Belle II Software  release-06-02-00
ARICHNtupleModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 // Own include
10 #include <arich/modules/arichNtuple/ARICHNtupleModule.h>
11 
12 // Hit classes
13 #include <framework/dataobjects/EventMetaData.h>
14 #include <arich/dataobjects/ARICHLikelihood.h>
15 #include <arich/dataobjects/ARICHAeroHit.h>
16 #include <arich/dataobjects/ARICHTrack.h>
17 #include <arich/dataobjects/ARICHPhoton.h>
18 #include <arich/dataobjects/ARICHHit.h>
19 #include <arich/dataobjects/ARICHInfo.h>
20 
21 // framework - DataStore
22 #include <framework/datastore/StoreObjPtr.h>
23 
24 #include <mdst/dataobjects/HitPatternCDC.h>
25 
26 #ifdef ALIGNMENT_USING_BHABHA
27 #include <mdst/dataobjects/ECLCluster.h>
28 #endif
29 
30 // framework aux
31 #include <framework/gearbox/Const.h>
32 #include <framework/logging/Logger.h>
33 
34 // ROOT
35 #include <TVector3.h>
36 #include <vector>
37 
38 using namespace std;
39 
40 namespace Belle2 {
46  //-----------------------------------------------------------------
47  // Register module
48  //-----------------------------------------------------------------
49 
50  REG_MODULE(ARICHNtuple)
51 
52  //-----------------------------------------------------------------
53  // Implementation
54  //-----------------------------------------------------------------
55 
57  m_file(0),
58  m_tree(0)
59  {
60  // set module description (e.g. insert text)
61  setDescription("The module saves variables needed for performance analysis, such as postion and momentum of the hit, likelihoods for hypotheses and number of photons.");
62 
63  // Add parameters
64  addParam("outputFile", m_outputFile, "ROOT output file name", string("ARICHNtuple.root"));
65  }
66 
67  ARICHNtupleModule::~ARICHNtupleModule()
68  {
69  }
70 
71  void ARICHNtupleModule::initialize()
72  {
73 
74  if (m_file) delete m_file;
75  m_file = new TFile(m_outputFile.c_str(), "RECREATE");
76  if (m_file->IsZombie()) {
77  B2FATAL("Couldn't open file '" << m_outputFile << "' for writing!");
78  return;
79  }
80 
81  m_tree = new TTree("arich", "ARICH validation ntuple");
82 
83  m_tree->Branch("evt", &m_arich.evt, "evt/I");
84  m_tree->Branch("run", &m_arich.run, "run/I");
85  m_tree->Branch("exp", &m_arich.exp, "exp/I");
86 
87  m_tree->Branch("charge", &m_arich.charge, "charge/S");
88  m_tree->Branch("pValue", &m_arich.pValue, "pValue/F");
89  m_tree->Branch("d0", &m_arich.z0, "pValue/F");
90  m_tree->Branch("z0", &m_arich.d0, "pValue/F");
91 
92 #ifdef ALIGNMENT_USING_BHABHA
93  m_tree->Branch("eop", &m_arich.eop, "eop/F");
94  m_tree->Branch("e9e21", &m_arich.e9e21, "e9e21/F");
95  m_tree->Branch("etot", &m_arich.etot, "etot/F");
96 #endif
97 
98  m_tree->Branch("PDG", &m_arich.PDG, "PDG/I");
99  m_tree->Branch("motherPDG", &m_arich.motherPDG, "motherPDG/I");
100  m_tree->Branch("primary", &m_arich.primary, "primary/S");
101  m_tree->Branch("seen", &m_arich.seen, "seen/S");
102  m_tree->Branch("scatter", &m_arich.scatter, "scatter/I");
103 
104  m_tree->Branch("rhoProd", &m_arich.rhoProd, "rhoProd/F");
105  m_tree->Branch("zProd", &m_arich.zProd, "zProd/F");
106  m_tree->Branch("phiProd", &m_arich.phiProd, "phiProd/F");
107  m_tree->Branch("rhoDec", &m_arich.rhoDec, "rhoDec/F");
108  m_tree->Branch("zDec", &m_arich.zDec, "zDec/F");
109  m_tree->Branch("phiDec", &m_arich.phiDec, "phiDec/F");
110  m_tree->Branch("status", &m_arich.status, "status/I");
111 
112  m_tree->Branch("detPhot", &m_arich.detPhot, "detPhot/I");
113  m_tree->Branch("numBkg", &m_arich.numBkg, "e/F:mu:pi:K:p:d");
114  m_tree->Branch("expPhot", &m_arich.expPhot, "e/F:mu:pi:K:p:d");
115  m_tree->Branch("logL", &m_arich.logL, "e/F:mu:pi:K:p:d");
116 
117  m_tree->Branch("recHit", &m_arich.recHit, "PDG/I:x/F:y:z:p:theta:phi");
118  m_tree->Branch("mcHit", &m_arich.mcHit, "PDG/I:x/F:y:z:p:theta:phi");
119  m_tree->Branch("winHit", &m_arich.winHit, "x/F:y");
120  m_tree->Branch("nrec", &m_arich.nRec, "nRec/I");
121  m_tree->Branch("nCDC", &m_arich.nCDC, "nCDC/I");
122  m_tree->Branch("inAcceptance", &m_arich.inAcc, "inAcc/O");
123  m_tree->Branch("photons", "std::vector<Belle2::ARICHPhoton>", &m_arich.photons);
124 
125  // required input
126  m_arichTracks.isRequired();
127  m_arichLikelihoods.isRequired();
128 
129  // optional input
130  m_tracks.isOptional();
131  m_arichMCPs.isOptional();
132  m_arichAeroHits.isOptional();
133  m_arichInfo.isOptional();
134 
135  StoreArray<ARICHHit> arichHits;
136  arichHits.isRequired();
137 
138  }
139 
140  void ARICHNtupleModule::beginRun()
141  {
142  }
143 
144  void ARICHNtupleModule::event()
145  {
146 
147  StoreObjPtr<EventMetaData> evtMetaData;
148  StoreArray<ARICHHit> arichHits;
149 
150  if (!m_arichTracks.isValid()) return;
151 
152 #ifdef ALIGNMENT_USING_BHABHA
153  double E_tot = 0.;
154  StoreArray<ECLCluster> eclClusters;
155  for (const auto& trk : m_tracks) {
156  const ECLCluster* eclTrack = trk.getRelated<ECLCluster>();
157  if (eclTrack) {
158  if (eclTrack->getEnergy() > 0.1) E_tot += eclTrack->getEnergy();
159  }
160  }
161  for (const auto& cluster : eclClusters) {
162  if (!cluster.isNeutral()) continue;
163  if (cluster.getHypothesisId() != 5 && cluster.getHypothesisId() != 1)
164  continue;
165  if (cluster.getEnergy() > 0.1) E_tot += cluster.getEnergy();
166  }
167 #endif
168 
169  for (const auto& arichTrack : m_arichTracks) {
170 
171  const ExtHit* extHit = arichTrack.getRelated<ExtHit>();
172 
173  const Track* mdstTrack = NULL;
174  if (extHit) mdstTrack = extHit->getRelated<Track>();
175 
176  const ARICHAeroHit* aeroHit = arichTrack.getRelated<ARICHAeroHit>();
177 
178  const ARICHLikelihood* lkh = NULL;
179  if (mdstTrack) lkh = mdstTrack->getRelated<ARICHLikelihood>();
180  else lkh = arichTrack.getRelated<ARICHLikelihood>();
181 
182  if (!lkh) continue;
183 
184  m_arich.clear();
185 
186  m_arich.evt = evtMetaData->getEvent();
187  m_arich.run = evtMetaData->getRun();
188  m_arich.exp = evtMetaData->getExperiment();
189 
190  // set hapd window hit if available
191  if (arichTrack.hitsWindow()) {
192  TVector2 winHit = arichTrack.windowHitPosition();
193  m_arich.winHit[0] = winHit.X();
194  m_arich.winHit[1] = winHit.Y();
195  }
196 
197  if (lkh->getFlag() == 1) m_arich.inAcc = 1;
198  m_arich.logL.e = lkh->getLogL(Const::electron);
199  m_arich.logL.mu = lkh->getLogL(Const::muon);
200  m_arich.logL.pi = lkh->getLogL(Const::pion);
201  m_arich.logL.K = lkh->getLogL(Const::kaon);
202  m_arich.logL.p = lkh->getLogL(Const::proton);
203  m_arich.logL.d = lkh->getLogL(Const::deuteron);
204 
205  m_arich.expPhot.e = lkh->getExpPhot(Const::electron);
206  m_arich.expPhot.mu = lkh->getExpPhot(Const::muon);
207  m_arich.expPhot.pi = lkh->getExpPhot(Const::pion);
208  m_arich.expPhot.K = lkh->getExpPhot(Const::kaon);
209  m_arich.expPhot.p = lkh->getExpPhot(Const::proton);
210  m_arich.expPhot.d = lkh->getExpPhot(Const::deuteron);
211 
212  m_arich.detPhot = lkh->getDetPhot();
213 
214 #ifdef ALIGNMENT_USING_BHABHA
215  m_arich.etot = E_tot;
216 #endif
217  m_arich.status = 1;
218 
219  m_arich.trgtype = 0;
220  if (m_arichInfo.getEntries() > 0) {
221  auto arichInfo = *m_arichInfo[0];
222  m_arich.trgtype = arichInfo.gettrgtype();
223  }
224 
225  const MCParticle* particle = 0;
226 
227  if (mdstTrack) {
228  const TrackFitResult* fitResult = mdstTrack->getTrackFitResultWithClosestMass(Const::pion);
229  if (fitResult) {
230  m_arich.pValue = fitResult->getPValue();
231  TVector3 trkPos = fitResult->getPosition();
232  m_arich.charge = fitResult->getChargeSign();
233  m_arich.z0 = trkPos.Z();
234  m_arich.d0 = (trkPos.XYvector()).Mod();
235  m_arich.nCDC = fitResult->getHitPatternCDC().getNHits();
236 #ifdef ALIGNMENT_USING_BHABHA
237  TVector3 trkMom = fitResult->getMomentum();
238  const ECLCluster* eclTrack = mdstTrack->getRelated<ECLCluster>();
239  if (eclTrack) {
240  m_arich.eop = eclTrack->getEnergy() / trkMom.Mag();
241  m_arich.e9e21 = eclTrack->getE9oE21();
242  }
243 #endif
244  }
245  m_arich.status += 10;
246  int fromARICHMCP = 0;
247  particle = mdstTrack->getRelated<MCParticle>();
248  if (!particle) { particle = mdstTrack->getRelated<MCParticle>("arichMCParticles"); fromARICHMCP = 1;}
249  if (particle) {
250  m_arich.PDG = particle->getPDG();
251  m_arich.primary = particle->getStatus(MCParticle::c_PrimaryParticle);
252  m_arich.seen = particle->hasSeenInDetector(Const::ARICH);
253  TVector3 prodVertex = particle->getProductionVertex();
254  m_arich.rhoProd = prodVertex.Perp();
255  m_arich.zProd = prodVertex.Z();
256  m_arich.phiProd = prodVertex.Phi();
257  TVector3 decVertex = particle->getDecayVertex();
258  m_arich.rhoDec = decVertex.Perp();
259  m_arich.zDec = decVertex.Z();
260  m_arich.phiDec = decVertex.Phi();
261 
262  if (!fromARICHMCP) {
263  MCParticle* mother = particle->getMother();
264  if (mother) m_arich.motherPDG = mother->getPDG();
265  std::vector<Belle2::MCParticle*> daughs = particle->getDaughters();
266  for (const auto daugh : daughs) {
267  if (daugh == NULL) continue;
268  if (daugh->getPDG() == particle->getPDG()) m_arich.scatter = 1;
269  }
270  }
271  }
272  }
273 
274 
275  // get reconstructed photons associated with track
276  const std::vector<ARICHPhoton>& photons = arichTrack.getPhotons();
277  m_arich.nRec = photons.size();
278  int nphot = 0;
279  for (auto it = photons.begin(); it != photons.end(); ++it) {
280  ARICHPhoton iph = *it;
281  if (iph.getHitID() < arichHits.getEntries()) {
282  auto ihit = *arichHits[iph.getHitID()];
283  int module = ihit.getModule();
284  int channel = ihit.getChannel();
285  int chid = 1 << 16 | module << 8 | channel;
286  iph.setHitID(chid);
287  }
288  m_arich.photons.push_back(iph);
289  nphot++;
290  }
291 
292  TVector3 recPos = arichTrack.getPosition();
293  m_arich.recHit.x = recPos.X();
294  m_arich.recHit.y = recPos.Y();
295  m_arich.recHit.z = recPos.Z();
296 
297  TVector3 recMom = arichTrack.getDirection() * arichTrack.getMomentum();
298  m_arich.recHit.p = recMom.Mag();
299  m_arich.recHit.theta = recMom.Theta();
300  m_arich.recHit.phi = recMom.Phi();
301 
302  if (aeroHit) {
303  TVector3 truePos = aeroHit->getPosition();
304  m_arich.mcHit.x = truePos.X();
305  m_arich.mcHit.y = truePos.Y();
306  m_arich.mcHit.z = truePos.Z();
307 
308  TVector3 trueMom = aeroHit->getMomentum();
309  m_arich.mcHit.p = trueMom.Mag();
310  m_arich.mcHit.theta = trueMom.Theta();
311  m_arich.mcHit.phi = trueMom.Phi();
312 
313  m_arich.mcHit.PDG = aeroHit->getPDG();
314  m_arich.status += 1000;
315 
316  if (!particle) {
317  particle = aeroHit->getRelated<MCParticle>();
318  if (particle) {
319  m_arich.PDG = particle->getPDG();
320  MCParticle* mother = particle->getMother();
321  if (mother) m_arich.motherPDG = mother->getPDG();
322  m_arich.primary = particle->getStatus(MCParticle::c_PrimaryParticle);
323  m_arich.seen = particle->hasSeenInDetector(Const::ARICH);
324  TVector3 prodVertex = particle->getProductionVertex();
325  m_arich.rhoProd = prodVertex.Perp();
326  m_arich.zProd = prodVertex.Z();
327  m_arich.phiProd = prodVertex.Phi();
328  TVector3 decVertex = particle->getDecayVertex();
329  m_arich.rhoDec = decVertex.Perp();
330  m_arich.zDec = decVertex.Z();
331  m_arich.phiDec = decVertex.Phi();
332 
333  std::vector<Belle2::MCParticle*> daughs = particle->getDaughters();
334  for (const auto daugh : daughs) {
335  if (daugh->getPDG() == particle->getPDG()) m_arich.scatter = 1;
336  }
337  }
338  }
339  }
340  m_tree->Fill();
341  }
342  }
343 
344 
345 
346  void ARICHNtupleModule::endRun()
347  {
348  }
349 
350  void ARICHNtupleModule::terminate()
351  {
352  m_file->cd();
353  m_tree->Write();
354  m_file->Close();
355  }
356 
358 } // end Belle2 namespace
359 
Datastore class that holds information on track parameters at the entrance in aerogel.
Definition: ARICHAeroHit.h:27
This is a class to store ARICH likelihoods in the datastore.
ARICH reconstruction efficiency test module.
Struct for ARICH reconstructed photon (hit related to track) information.
Definition: ARICHPhoton.h:25
void setHitID(int id)
Set ID of corresponding ARICHHit.
Definition: ARICHPhoton.h:48
int getHitID()
Get ID of corresponding ARICHHit.
Definition: ARICHPhoton.h:230
ECL cluster data.
Definition: ECLCluster.h:27
double getE9oE21() const
Return E9/E21 (shower shape variable).
Definition: ECLCluster.h:277
double getEnergy(EHypothesisBit hypothesis) const
Return Energy (GeV).
Definition: ECLCluster.cc:19
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:30
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
Base class for Modules.
Definition: Module.h:72
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
double getPValue() const
Getter for Chi2 Probability of the track fit.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
TVector3 getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
Class that bundles various TrackFitResults.
Definition: Track.h:25
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.