|
| 1 | +// Copyright 2019-2020 CERN and copyright holders of ALICE O2. |
| 2 | +// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders. |
| 3 | +// All rights not expressly granted are reserved. |
| 4 | +// |
| 5 | +// This software is distributed under the terms of the GNU General Public |
| 6 | +// License v3 (GPL Version 3), copied verbatim in the file "COPYING". |
| 7 | +// |
| 8 | +// In applying this license CERN does not waive the privileges and immunities |
| 9 | +// granted to it by virtue of its status as an Intergovernmental Organization |
| 10 | +// or submit itself to any jurisdiction. |
| 11 | + |
| 12 | +/// \file upcPolarisationJpsiIncoh.cxx |
| 13 | +/// \brief Workflow to analyse UPC forward events and perform J/psi polarization selections |
| 14 | +/// \author Niveditha Ram, IP2I <niv.ram@cern.ch> |
| 15 | +/// \ingroup PWGUD |
| 16 | +/// executable name: o2-analysis-ud-upc-polarisation-jpsiincoh |
| 17 | + |
| 18 | +#include "PWGUD/DataModel/UDTables.h" |
| 19 | + |
| 20 | +#include "Common/Core/RecoDecay.h" |
| 21 | + |
| 22 | +#include "CCDB/BasicCCDBManager.h" |
| 23 | +#include "CommonConstants/PhysicsConstants.h" |
| 24 | +#include "DataFormatsParameters/GRPECSObject.h" |
| 25 | +#include "DataFormatsParameters/GRPLHCIFData.h" |
| 26 | +#include "Framework/AnalysisDataModel.h" |
| 27 | +#include "Framework/AnalysisTask.h" |
| 28 | +#include "Framework/runDataProcessing.h" |
| 29 | + |
| 30 | +#include "Math/Vector4D.h" |
| 31 | +#include "TMath.h" |
| 32 | +#include "TRandom3.h" |
| 33 | +#include "TSystem.h" |
| 34 | + |
| 35 | +#include <unordered_map> |
| 36 | +#include <vector> |
| 37 | + |
| 38 | +using namespace ROOT::Math; |
| 39 | + |
| 40 | +// table for saving tree with info on data |
| 41 | +namespace dimu |
| 42 | +{ |
| 43 | +// dimuon |
| 44 | +DECLARE_SOA_COLUMN(RunNumber, runNumber, int); |
| 45 | +DECLARE_SOA_COLUMN(M, m, float); |
| 46 | +DECLARE_SOA_COLUMN(Pt, pt, float); |
| 47 | +DECLARE_SOA_COLUMN(Rap, rap, float); |
| 48 | +DECLARE_SOA_COLUMN(Phi, phi, float); |
| 49 | +} // namespace dimu |
| 50 | + |
| 51 | +namespace o2::aod |
| 52 | +{ |
| 53 | +DECLARE_SOA_TABLE(DiMu, "AOD", "DIMU", |
| 54 | + dimu::RunNumber, |
| 55 | + dimu::M, dimu::Pt, dimu::Rap, dimu::Phi); |
| 56 | +} // namespace o2::aod |
| 57 | +using namespace o2; |
| 58 | +using namespace o2::framework; |
| 59 | +using namespace o2::framework::expressions; |
| 60 | + |
| 61 | +// constants used in the track selection |
| 62 | +const float kRAbsMin = 17.6; |
| 63 | +const float kRAbsMax = 89.5; |
| 64 | +const float kPDca = 200.; |
| 65 | +float kEtaMin = -4.0; |
| 66 | +float kEtaMax = -2.5; |
| 67 | +const float kPtMin = 0.; |
| 68 | +const float kMaxAmpV0A = 100.; |
| 69 | +const int kReqMatchMIDTracks = 2; |
| 70 | +const int kReqMatchMFTTracks = 2; |
| 71 | +const int kMaxChi2MFTMatch = 30; |
| 72 | +struct UpcPolarisationJpsiIncoh { |
| 73 | + |
| 74 | + using CandidatesFwd = soa::Join<o2::aod::UDCollisions, o2::aod::UDCollisionsSelsFwd>; |
| 75 | + using ForwardTracks = soa::Join<o2::aod::UDFwdTracks, o2::aod::UDFwdTracksExtra>; |
| 76 | + using CompleteFwdTracks = soa::Join<ForwardTracks, o2::aod::UDMcFwdTrackLabels>; |
| 77 | + |
| 78 | + Produces<o2::aod::DiMu> dimuSel; |
| 79 | + // defining histograms using histogram registry: different histos for the different process functions |
| 80 | + HistogramRegistry registry{"registry", {}, OutputObjHandlingPolicy::AnalysisObject, true, true}; |
| 81 | + // CONFIGURABLES |
| 82 | + static constexpr double Pi = o2::constants::math::PI; |
| 83 | + // pT of muon pairs |
| 84 | + Configurable<int> nBinsPt{"nBinsPt", 250, "N bins in pT histo"}; |
| 85 | + Configurable<float> lowPt{"lowPt", 0., "lower limit in pT histo"}; |
| 86 | + Configurable<float> highPt{"highPt", 2, "upper limit in pT histo"}; |
| 87 | + // mass of muon pairs |
| 88 | + Configurable<int> nBinsMass{"nBinsMass", 500, "N bins in mass histo"}; |
| 89 | + Configurable<float> lowMass{"lowMass", 0., "lower limit in mass histo"}; |
| 90 | + Configurable<float> highMass{"highMass", 10., "upper limit in mass histo"}; |
| 91 | + // eta of muon pairs |
| 92 | + Configurable<int> nBinsEta{"nBinsEta", 600, "N bins in eta histo"}; |
| 93 | + Configurable<float> lowEta{"lowEta", -10., "lower limit in eta histo"}; |
| 94 | + Configurable<float> highEta{"highEta", -2., "upper limit in eta histo"}; |
| 95 | + // rapidity of muon pairs |
| 96 | + Configurable<int> nBinsRapidity{"nBinsRapidity", 250, "N bins in rapidity histo"}; |
| 97 | + Configurable<float> lowRapidity{"lowRapidity", -4.5, "lower limit in rapidity histo"}; |
| 98 | + Configurable<float> highRapidity{"highRapidity", -2., "upper limit in rapidity histo"}; |
| 99 | + // phi of muon pairs |
| 100 | + Configurable<int> nBinsPhi{"nBinsPhi", 600, "N bins in phi histo"}; |
| 101 | + Configurable<float> lowPhi{"lowPhi", -Pi, "lower limit in phi histo"}; |
| 102 | + Configurable<float> highPhi{"highPhi", Pi, "upper limit in phi histo"}; |
| 103 | + // Analysis cuts |
| 104 | + Configurable<float> maxJpsiMass{"maxJpsiMass", 3.18, "Maximum of the jpsi peak for peak cut"}; |
| 105 | + Configurable<float> minJpsiMass{"minJpsiMass", 3.0, "Minimum of the jpsi peak for peak cut"}; |
| 106 | + // my track type |
| 107 | + // 0 = MCH-MID-MFT |
| 108 | + // 1 = MCH-MID |
| 109 | + Configurable<int> myTrackType{"myTrackType", 1, "My track type"}; |
| 110 | + |
| 111 | + void init(InitContext&) |
| 112 | + { |
| 113 | + // axis |
| 114 | + const AxisSpec axisPt{nBinsPt, lowPt, highPt, "#it{p}_{T} GeV/#it{c}"}; |
| 115 | + const AxisSpec axisMass{nBinsMass, lowMass, highMass, "m_{#mu#mu} GeV/#it{c}^{2}"}; |
| 116 | + const AxisSpec axisEta{nBinsEta, lowEta, highEta, "#eta"}; |
| 117 | + const AxisSpec axisRapidity{nBinsRapidity, lowRapidity, highRapidity, "Rapidity"}; |
| 118 | + const AxisSpec axisPhi{nBinsPhi, lowPhi, highPhi, "#varphi"}; |
| 119 | + // histos |
| 120 | + // data and reco MC |
| 121 | + registry.add("hMass", "Invariant mass of muon pairs;;#counts", kTH1D, {axisMass}); |
| 122 | + registry.add("hPt", "Transverse momentum of muon pairs;;#counts", kTH1D, {axisPt}); |
| 123 | + registry.add("hEta", "Pseudorapidty of muon pairs;;#counts", kTH1D, {axisEta}); |
| 124 | + registry.add("hRapidity", "Rapidty of muon pairs;;#counts", kTH1D, {axisRapidity}); |
| 125 | + registry.add("hPhi", "#varphi of muon pairs;;#counts", kTH1D, {axisPhi}); |
| 126 | + } |
| 127 | + |
| 128 | + // template function that fills a map with the collision id of each udcollision as key |
| 129 | + // and a vector with the tracks |
| 130 | + // map == (key, element) == (udCollisionId, vector of trks) |
| 131 | + template <typename TTracks> |
| 132 | + void collectCandIDs(std::unordered_map<int32_t, std::vector<int32_t>>& tracksPerCand, TTracks& tracks) |
| 133 | + { |
| 134 | + for (const auto& tr : tracks) { |
| 135 | + int32_t candId = tr.udCollisionId(); |
| 136 | + if (candId < 0) { |
| 137 | + continue; |
| 138 | + } |
| 139 | + tracksPerCand[candId].push_back(tr.globalIndex()); |
| 140 | + } |
| 141 | + } |
| 142 | + |
| 143 | + // template function that fills a map with the collision id of each udmccollision as key |
| 144 | + // and a vector with the tracks |
| 145 | + // map == (key, element) == (udMcCollisionId, vector of mc particles) |
| 146 | + template <typename TTracks> |
| 147 | + void collectMcCandIDs(std::unordered_map<int32_t, std::vector<int32_t>>& tracksPerCand, TTracks& tracks) |
| 148 | + { |
| 149 | + for (const auto& tr : tracks) { |
| 150 | + int32_t candId = tr.udMcCollisionId(); |
| 151 | + if (candId < 0) { |
| 152 | + continue; |
| 153 | + } |
| 154 | + tracksPerCand[candId].push_back(tr.globalIndex()); |
| 155 | + } |
| 156 | + } |
| 157 | + |
| 158 | + // struct used to store the ZDC info in a map |
| 159 | + struct ZDCinfo { |
| 160 | + float timeA; |
| 161 | + float timeC; |
| 162 | + float enA; |
| 163 | + float enC; |
| 164 | + int32_t id; |
| 165 | + }; |
| 166 | + |
| 167 | + // function that fills a map with the collision id of each udcollision as key |
| 168 | + // and a ZDCinfo struct with the ZDC information |
| 169 | + void collectCandZDCInfo(std::unordered_map<int32_t, ZDCinfo>& zdcPerCand, o2::aod::UDZdcsReduced const& ZDCs) |
| 170 | + { |
| 171 | + |
| 172 | + for (const auto& zdc : ZDCs) { |
| 173 | + int32_t candId = zdc.udCollisionId(); |
| 174 | + if (candId < 0) { |
| 175 | + continue; |
| 176 | + } |
| 177 | + |
| 178 | + zdcPerCand[candId].timeA = zdc.timeZNA(); |
| 179 | + zdcPerCand[candId].timeC = zdc.timeZNC(); |
| 180 | + zdcPerCand[candId].enA = zdc.energyCommonZNA(); |
| 181 | + zdcPerCand[candId].enC = zdc.energyCommonZNC(); |
| 182 | + |
| 183 | + // take care of the infinity |
| 184 | + if (std::isinf(zdcPerCand[candId].timeA)) |
| 185 | + zdcPerCand[candId].timeA = -999; |
| 186 | + if (std::isinf(zdcPerCand[candId].timeC)) |
| 187 | + zdcPerCand[candId].timeC = -999; |
| 188 | + if (std::isinf(zdcPerCand[candId].enA)) |
| 189 | + zdcPerCand[candId].enA = -999; |
| 190 | + if (std::isinf(zdcPerCand[candId].enC)) |
| 191 | + zdcPerCand[candId].enC = -999; |
| 192 | + } |
| 193 | + } |
| 194 | + |
| 195 | + // function to select muon tracks |
| 196 | + template <typename TTracks> |
| 197 | + bool isMuonSelected(const TTracks& fwdTrack) |
| 198 | + { |
| 199 | + float rAbs = fwdTrack.rAtAbsorberEnd(); |
| 200 | + float pDca = fwdTrack.pDca(); |
| 201 | + float pt = RecoDecay::pt(fwdTrack.px(), fwdTrack.py()); |
| 202 | + float eta = RecoDecay::eta(std::array{fwdTrack.px(), fwdTrack.py(), fwdTrack.pz()}); |
| 203 | + if (eta < kEtaMin || eta > kEtaMax) |
| 204 | + return false; |
| 205 | + if (pt < kPtMin) |
| 206 | + return false; |
| 207 | + if (rAbs < kRAbsMin || rAbs > kRAbsMax) |
| 208 | + return false; |
| 209 | + if (pDca > kPDca) |
| 210 | + return false; |
| 211 | + return true; |
| 212 | + } |
| 213 | + |
| 214 | + // function that processes the candidates: |
| 215 | + // it applies V0 selection, trk selection, kine selection, and fills the histograms |
| 216 | + // it also divides the data in neutron classes |
| 217 | + // used for real data |
| 218 | + void processCand(CandidatesFwd::iterator const& cand, |
| 219 | + ForwardTracks::iterator const& tr1, ForwardTracks::iterator const& tr2) |
| 220 | + { |
| 221 | + // V0 selection |
| 222 | + const auto& ampsV0A = cand.amplitudesV0A(); |
| 223 | + const auto& ampsRelBCsV0A = cand.ampRelBCsV0A(); |
| 224 | + for (unsigned int i = 0; i < ampsV0A.size(); ++i) { |
| 225 | + if (std::abs(ampsRelBCsV0A[i]) <= 1) { |
| 226 | + if (ampsV0A[i] > kMaxAmpV0A) |
| 227 | + return; |
| 228 | + } |
| 229 | + } |
| 230 | + |
| 231 | + // MCH-MID match selection |
| 232 | + int nMIDs = 0; |
| 233 | + if (tr1.chi2MatchMCHMID() > 0) |
| 234 | + nMIDs++; |
| 235 | + if (tr2.chi2MatchMCHMID() > 0) |
| 236 | + nMIDs++; |
| 237 | + if (nMIDs != kReqMatchMIDTracks) |
| 238 | + return; |
| 239 | + // MFT-MID match selection (if MFT is requested by the trackType) |
| 240 | + if (myTrackType == 0) { |
| 241 | + // if MFT is requested check that the tracks is inside the MFT acceptance |
| 242 | + kEtaMin = -3.6; |
| 243 | + kEtaMax = -2.5; |
| 244 | + |
| 245 | + int nMFT = 0; |
| 246 | + if (tr1.chi2MatchMCHMFT() > 0 && tr1.chi2MatchMCHMFT() < kMaxChi2MFTMatch) |
| 247 | + nMFT++; |
| 248 | + if (tr2.chi2MatchMCHMFT() > 0 && tr2.chi2MatchMCHMFT() < kMaxChi2MFTMatch) |
| 249 | + nMFT++; |
| 250 | + if (nMFT != kReqMatchMFTTracks) |
| 251 | + return; |
| 252 | + } |
| 253 | + // track selection |
| 254 | + if (!isMuonSelected(*tr1)) |
| 255 | + return; |
| 256 | + if (!isMuonSelected(*tr2)) |
| 257 | + return; |
| 258 | + |
| 259 | + // form Lorentz vectors |
| 260 | + auto mMu = o2::constants::physics::MassMuonMinus; |
| 261 | + LorentzVector<PxPyPzM4D<float>> p1(tr1.px(), tr1.py(), tr1.pz(), mMu); |
| 262 | + LorentzVector<PxPyPzM4D<float>> p2(tr2.px(), tr2.py(), tr2.pz(), mMu); |
| 263 | + LorentzVector p = p1 + p2; |
| 264 | + |
| 265 | + // cut on pair kinematics |
| 266 | + // select mass |
| 267 | + if (p.M() < lowMass) |
| 268 | + return; |
| 269 | + if (p.M() > highMass) |
| 270 | + return; |
| 271 | + // select pt |
| 272 | + if (p.Pt() < lowPt) |
| 273 | + return; |
| 274 | + if (p.Pt() > highPt) |
| 275 | + return; |
| 276 | + // select rapidity |
| 277 | + if (p.Rapidity() < lowRapidity) |
| 278 | + return; |
| 279 | + if (p.Rapidity() > highRapidity) |
| 280 | + return; |
| 281 | + // fill the histos without looking at neutron emission |
| 282 | + registry.fill(HIST("hMass"), p.M()); |
| 283 | + registry.fill(HIST("hPt"), p.Pt()); |
| 284 | + registry.fill(HIST("hEta"), p.Eta()); |
| 285 | + registry.fill(HIST("hRapidity"), p.Rapidity()); |
| 286 | + registry.fill(HIST("hPhi"), p.Phi()); |
| 287 | + |
| 288 | + dimuSel(cand.runNumber(), p.M(), p.Pt(), p.Rapidity(), p.Phi()); |
| 289 | + } |
| 290 | + // PROCESS FUNCTION |
| 291 | + void processData(CandidatesFwd const& eventCandidates, |
| 292 | + o2::aod::UDZdcsReduced const& ZDCs, |
| 293 | + ForwardTracks const& fwdTracks) |
| 294 | + { |
| 295 | + |
| 296 | + // map with the tracks |
| 297 | + std::unordered_map<int32_t, std::vector<int32_t>> tracksPerCand; |
| 298 | + // takes a tracks table with a coloumn of collision ID and makes it into a map of collision ID to each track. |
| 299 | + collectCandIDs(tracksPerCand, fwdTracks); |
| 300 | + |
| 301 | + // map with the ZDC info |
| 302 | + std::unordered_map<int32_t, ZDCinfo> zdcPerCand; |
| 303 | + collectCandZDCInfo(zdcPerCand, ZDCs); |
| 304 | + |
| 305 | + // loop over the candidates |
| 306 | + for (const auto& item : tracksPerCand) { |
| 307 | + int32_t trId1 = item.second[0]; |
| 308 | + int32_t trId2 = item.second[1]; |
| 309 | + int32_t candID = item.first; |
| 310 | + auto cand = eventCandidates.iteratorAt(candID); |
| 311 | + auto tr1 = fwdTracks.iteratorAt(trId1); |
| 312 | + auto tr2 = fwdTracks.iteratorAt(trId2); |
| 313 | + processCand(cand, tr1, tr2); |
| 314 | + } |
| 315 | + } |
| 316 | + |
| 317 | + PROCESS_SWITCH(UpcPolarisationJpsiIncoh, processData, "", true); |
| 318 | +}; |
| 319 | + |
| 320 | +WorkflowSpec defineDataProcessing(ConfigContext const& cfgc) |
| 321 | +{ |
| 322 | + return WorkflowSpec{ |
| 323 | + adaptAnalysisTask<UpcPolarisationJpsiIncoh>(cfgc), |
| 324 | + }; |
| 325 | +} |
0 commit comments