Belle II Software  release-05-01-25
HitLevelInfoWriter.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2012 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Jake Bennett, Jitendra Kumar
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <reconstruction/modules/HitLevelInfoWriter/HitLevelInfoWriter.h>
12 #include <reconstruction/dataobjects/DedxConstants.h>
13 #include <mdst/dataobjects/PIDLikelihood.h>
14 
15 using namespace Belle2;
16 using namespace Dedx;
17 
18 REG_MODULE(HitLevelInfoWriter)
19 
21 {
22 
23  setDescription("Extract dE/dx information for calibration development.");
24 
25  addParam("outputBaseName", m_strOutputBaseName, "Suffix for output file name", std::string("HLInfo.root"));
26  addParam("particleLists", m_strParticleList, "Vector of ParticleLists to save", std::vector<std::string>());
27  addParam("enableHitLevel", enableHitLevel, "True or False for Hit level variables", false);
28  addParam("enableExtraVar", enableExtraVar, "True or False for extra track/hit level variables", false);
29  addParam("nodeadwire", nodeadwire, "True or False for deadwire hit variables", true);
30 
31 }
32 
34 
36 {
37 
38  B2INFO("Creating a ROOT file for the hit level information...");
39 
40  // required inputs
41  m_dedxTracks.isRequired();
42  m_tracks.isRequired();
43  m_trackFitResults.isRequired();
44  m_eclClusters.isOptional();
45  m_klmClusters.isOptional();
46 
47  // build a map to relate input strings to the right particle type
48  std::map<std::string, std::string> pdgMap = {{"pi+", "211"}, {"K+", "321"}, {"mu+", "13"}, {"e+", "11"}, {"p+", "2212"}, {"deuteron", "1000010020"}};
49 
50  // if no particle lists are given, write out all tracks
51  if (m_strParticleList.size() == 0) bookOutput(m_strOutputBaseName);
52 
53  // create a new output file for each particle list specified
54  for (unsigned int i = 0; i < m_strParticleList.size(); i++) {
55  // strip the name of the particle lists to make this work
56  std::string pdg = pdgMap[m_strParticleList[i].substr(0, m_strParticleList[i].find(":"))];
57  std::string filename = std::string(m_strOutputBaseName + "_PID" + pdg + ".root");
58  bookOutput(filename);
59  }
60 }
61 
63 {
64 
65  StoreArray<CDCDedxTrack> dedxTracks;
66  int nParticleList = m_strParticleList.size();
67 
68  // **************************************************
69  // LOOP OVER dE/dx measurements for all tracks if
70  // no particle list is specified.
71  // **************************************************
72  if (nParticleList == 0) {
73  for (int idedx = 0; idedx < m_dedxTracks.getEntries(); idedx++) {
74 
75  CDCDedxTrack* dedxTrack = m_dedxTracks[idedx];
76  if (!dedxTrack) {
77  B2WARNING("No dedx related track...");
78  continue;
79  }
80 
81  const Track* track = dedxTrack->getRelatedFrom<Track>();
82  if (!track) {
83  B2WARNING("No related track...");
84  continue;
85  }
86 
87  const TrackFitResult* fitResult = track->getTrackFitResultWithClosestMass(Const::pion);
88  if (!fitResult) {
89  B2WARNING("No related fit for this track...");
90  continue;
91  }
92 
93  if (dedxTrack->size() == 0 || dedxTrack->size() > 200) continue;
94  if (dedxTrack->getCosTheta() < -1.0 || dedxTrack->getCosTheta() > 1.0) continue;
95 
96  // fill the event meta data
97  StoreObjPtr<EventMetaData> evtMetaData;
98  m_expID = evtMetaData->getExperiment();
99  m_runID = evtMetaData->getRun();
100  m_eventID = evtMetaData->getEvent();
101 
102  // fill the E/P
103  const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
104  if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
105  m_eop = (eclCluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons)) / (fitResult->getMomentum().Mag());
107  m_e1_9 = eclCluster->getE1oE9();
108  m_e9_21 = eclCluster->getE9oE21();
109  // fill the muon depth
110  const KLMCluster* klmCluster = eclCluster->getRelated<KLMCluster>();
111  if (klmCluster) m_klmLayers = klmCluster->getLayers();
112  }
113 
114  // fill the TTree with the Track information
115  fillTrack(fitResult);
116 
117  // fill the TTree with the CDCDedxTrack information
118  fillDedx(dedxTrack);
119 
120  // fill the TTree
121  m_tree[0]->Fill();
122  }
123  }
124 
125  // **************************************************
126  //
127  // LOOP OVER particles in the given particle lists
128  //
129  // **************************************************
130 
131  for (int iList = 0; iList < nParticleList; iList++) {
132  // make sure the list exists and is not empty
133  StoreObjPtr<ParticleList> particlelist(m_strParticleList[iList]);
134  if (!particlelist or particlelist->getListSize(true) == 0) {
135  //B2WARNING("ParticleList " << m_strParticleList[iList] << " not found or empty, skipping");
136  continue;
137  }
138 
139  // loop over the particles in the list and follow the links to the
140  // dE/dx information (Particle -> PIDLikelihood -> Track -> CDCDedxTrack)
141  for (unsigned int iPart = 0; iPart < particlelist->getListSize(true); iPart++) {
142  Particle* part = particlelist->getParticle(iPart, true);
143  if (!part) {
144  B2WARNING("No particles...");
145  continue;
146  }
147  PIDLikelihood* pid = part->getRelatedTo<PIDLikelihood>();
148  if (!pid) {
149  B2WARNING("No related PID likelihood...");
150  continue;
151  }
152  Track* track = pid->getRelatedFrom<Track>();
153  if (!track) {
154  B2WARNING("No related track...");
155  continue;
156  }
157  CDCDedxTrack* dedxTrack = track->getRelatedTo<CDCDedxTrack>();
158  if (!dedxTrack) {
159  B2WARNING("No related CDCDedxTrack...");
160  continue;
161  }
162  std::string ptype = m_strParticleList[iList].substr(0, m_strParticleList[iList].find(":"));
163  const TrackFitResult* fitResult = track->getTrackFitResult(Const::pion);
164  if (ptype != "pi+") {
165  if (ptype == "K+") fitResult = track->getTrackFitResultWithClosestMass(Const::kaon);
166  else if (ptype == "p+") fitResult = track->getTrackFitResultWithClosestMass(Const::proton);
167  else if (ptype == "deuteron") fitResult = track->getTrackFitResultWithClosestMass(Const::deuteron);
168  else if (ptype == "mu+") fitResult = track->getTrackFitResultWithClosestMass(Const::muon);
169  else if (ptype == "e+") fitResult = track->getTrackFitResultWithClosestMass(Const::electron);
170  }
171  if (!fitResult) {
172  B2WARNING("No related fit for this track...");
173  continue;
174  }
175 
176  if (dedxTrack->size() == 0 || dedxTrack->size() > 200) continue;
177  //if out CDC (dont add as we dont correct via correctionmodules)
178  if (dedxTrack->getCosTheta() < TMath::Cos(150.0 * TMath::DegToRad()))continue; //-0.866
179  if (dedxTrack->getCosTheta() > TMath::Cos(17.0 * TMath::DegToRad())) continue; //0.95
180 
181  // fill the event meta data
182  StoreObjPtr<EventMetaData> evtMetaData;
183  m_expID = evtMetaData->getExperiment();
184  m_runID = evtMetaData->getRun();
185  m_eventID = evtMetaData->getEvent();
186 
187  // fill the E/P
188  const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
189  if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
190  m_eop = (eclCluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons)) / (fitResult->getMomentum().Mag());
191 
192  // fill the muon depth
193  const KLMCluster* klmCluster = eclCluster->getRelated<KLMCluster>();
194  if (klmCluster) m_klmLayers = klmCluster->getLayers();
195  }
196 
197  // fill the TTree with the Track information
198  fillTrack(fitResult);
199 
200  // fill the TTree with the CDCDedxTrack information
201  fillDedx(dedxTrack);
202 
203  // fill the TTree
204  m_tree[iList]->Fill();
205  }
206  }
207 }
208 
210 {
211 
212  for (unsigned int i = 0; i < m_file.size(); i++) {
213  B2INFO("Done writing out the hit level information...\t" << m_tree[i]->GetEntries() << " tracks");
214 
215  // write the ttree to a root file
216  m_file[i]->cd();
217  m_tree[i]->Write();
218  m_file[i]->Close();
219  }
220 }
221 
222 void
224 {
225  TVector3 trackMom = fitResult->getMomentum();
226  m_p = trackMom.Mag();
227  m_pt = trackMom.Pt();
228  m_phi = trackMom.Phi();
229 
230  m_theta = trackMom.Theta() * 180. / TMath::Pi(); //in degree
231  if (m_theta > 17. && m_theta < 150.)m_inCDC = 1;
232  else m_inCDC = 0;
233 
234  if (fitResult->getChargeSign() < 0) {
235  m_p *= -1;
236  m_pt *= -1;
237  }
238 
239  TVector3 trackPos = fitResult->getPosition();
240  m_vx0 = trackPos.x();
241  m_vy0 = trackPos.y();
242  m_vz0 = trackPos.z();
243 
244  m_d0 = fitResult->getD0();
245  m_z0 = fitResult->getZ0();
246  m_chi2 = fitResult->getPValue();
247  m_tanlambda = fitResult->getTanLambda();
248  m_phi0 = fitResult->getPhi0();
249  m_nCDChits = fitResult->getHitPatternCDC().getNHits();
250 
251  static DBObjPtr<BeamSpot> beamSpotDB;
252  const auto& frame = ReferenceFrame::GetCurrent();
253  UncertainHelix helix = fitResult->getUncertainHelix();
254  helix.passiveMoveBy(beamSpotDB->getIPPosition());
255  m_dr = frame.getVertex(helix.getPerigee()).Perp();
256  m_dphi = frame.getVertex(helix.getPerigee()).Phi();
257  m_dz = frame.getVertex(helix.getPerigee()).Z();
258 }
259 
260 void
262 {
263  // clear the containers first
264  clearEntries();
265 
266  m_trackID = dedxTrack->trackID();
267  m_length = dedxTrack->getLength();
268  m_charge = dedxTrack->getCharge();
269  m_cosTheta = dedxTrack->getCosTheta();
270  m_PDG = dedxTrack->getPDG();
271 
272  m_pCDC = dedxTrack->getMomentum();
273  if (m_charge < 0) {
274  m_pCDC *= -1;
275  }
276 
277  h_nhits = dedxTrack->size();
278  l_nhits = dedxTrack->getNLayerHits();
279  l_nhitsused = dedxTrack->getNLayerHitsUsed();
280 
281  m_mean = dedxTrack->getDedxMean();
282  m_trunc = dedxTrack->getDedx();
283  m_truncNoSat = dedxTrack->getDedxNoSat();
284  m_error = dedxTrack->getDedxError();
285 
286  // Get the calibration constants
287  m_scale = m_DBScaleFactor->getScaleFactor();
288  m_runGain = m_DBRunGain->getRunGain();
289  m_cosCor = m_DBCosineCor->getMean(m_cosTheta);
290 
291  if (m_cosTheta <= -0.850 || m_cosTheta >= 0.950) {
292  m_cosEdgeCor = m_DBCosEdgeCor->getMean(m_cosTheta);
293  } else {
294  m_cosEdgeCor = 1.0;
295  }
296 
297 
298  m_chie = dedxTrack->getChi(0);
299  m_chimu = dedxTrack->getChi(1);
300  m_chipi = dedxTrack->getChi(2);
301  m_chik = dedxTrack->getChi(3);
302  m_chip = dedxTrack->getChi(4);
303  m_chid = dedxTrack->getChi(5);
304 
305  m_pmeane = dedxTrack->getPmean(0);
306  m_pmeanmu = dedxTrack->getPmean(1);
307  m_pmeanpi = dedxTrack->getPmean(2);
308  m_pmeank = dedxTrack->getPmean(3);
309  m_pmeanp = dedxTrack->getPmean(4);
310  m_pmeand = dedxTrack->getPmean(5);
311 
312  m_prese = dedxTrack->getPreso(0);
313  m_presmu = dedxTrack->getPreso(1);
314  m_prespi = dedxTrack->getPreso(2);
315  m_presk = dedxTrack->getPreso(3);
316  m_presp = dedxTrack->getPreso(4);
317  m_presd = dedxTrack->getPreso(5);
318 
319  // Get the vector of dE/dx values for all layers
320  double lout = 0, lin = 0, increment = 0;
321  int lastlayer = 0;
322  for (int il = 0; il < l_nhits; ++il) {
323  l_nhitscombined[il] = dedxTrack->getNHitsCombined(il);
324  l_wirelongesthit[il] = dedxTrack->getWireLongestHit(il);
325  l_layer[il] = dedxTrack->getLayer(il);
326  l_path[il] = dedxTrack->getLayerPath(il);
327  l_dedx[il] = dedxTrack->getLayerDedx(il);
328 
329  if (l_layer[il] > lastlayer) lout++;
330  else if (l_layer[il] < lastlayer) lin++;
331  else continue;
332 
333  lastlayer = l_layer[il];
334  increment++;
335  }
336  m_ioasym = (lout - lin) / increment;
337 
338  // Get the vector of dE/dx values for all hits
339  if (enableHitLevel) {
340  for (int ihit = 0; ihit < h_nhits; ++ihit) {
341 
342  if (nodeadwire && m_DBWireGains->getWireGain(h_wire[ihit]) == 0)continue;
343  h_lwire[ihit] = dedxTrack->getWireInLayer(ihit);
344  h_wire[ihit] = dedxTrack->getWire(ihit);
345  h_layer[ihit] = dedxTrack->getHitLayer(ihit);
346  h_path[ihit] = dedxTrack->getPath(ihit);
347  h_dedx[ihit] = dedxTrack->getDedx(ihit);
348  h_adcraw[ihit] = dedxTrack->getADCBaseCount(ihit);
349  h_adccorr[ihit] = dedxTrack->getADCCount(ihit);
350  h_doca[ihit] = dedxTrack->getDoca(ihit);
351  h_ndoca[ihit] = h_doca[ihit] / dedxTrack->getCellHalfWidth(ihit);
352  h_ndocaRS[ihit] = dedxTrack->getDocaRS(ihit) / dedxTrack->getCellHalfWidth(ihit);
353  h_enta[ihit] = dedxTrack->getEnta(ihit);
354  h_entaRS[ihit] = dedxTrack->getEntaRS(ihit);
355  h_driftT[ihit] = dedxTrack->getDriftT(ihit);
356  h_driftD[ihit] = dedxTrack->getDriftD(ihit);
357 
358  // Get extra variable from tracking
359  if (enableExtraVar) {
360  h_WeightPionHypo[ihit] = dedxTrack->getWeightPionHypo(ihit);
361  h_WeightKaonHypo[ihit] = dedxTrack->getWeightKaonHypo(ihit);
362  h_WeightProtonHypo[ihit] = dedxTrack->getWeightProtonHypo(ihit);
363  h_foundByTrackFinder[ihit] = dedxTrack->getFoundByTrackFinder(ihit);
364  }
365  // Get the calibration constants
366  h_facnladc[ihit] = dedxTrack->getNonLADCCorrection(ihit);
367  h_wireGain[ihit] = m_DBWireGains->getWireGain(h_wire[ihit]);
368  h_twodCor[ihit] = m_DB2DCell->getMean(h_layer[ihit], h_ndocaRS[ihit], h_entaRS[ihit]);
369  h_onedCor[ihit] = m_DB1DCell->getMean(h_layer[ihit], h_entaRS[ihit]);
370  }
371  }
372 }
373 
374 void
376 {
377  for (int il = 0; il < 200; ++il) {
378  l_nhitscombined[il] = 0;
379  l_wirelongesthit[il] = 0;
380  l_layer[il] = 0;
381  l_path[il] = 0;
382  l_dedx[il] = 0;
383  }
384 
385  if (enableHitLevel) {
386  for (int ihit = 0; ihit < 200; ++ihit) {
387  h_lwire[ihit] = 0;
388  h_wire[ihit] = 0;
389  h_layer[ihit] = 0;
390  h_path[ihit] = 0;
391  h_dedx[ihit] = 0;
392  h_adcraw[ihit] = 0;
393  h_adccorr[ihit] = 0;
394  h_doca[ihit] = 0;
395  h_ndoca[ihit] = 0;
396  h_ndocaRS[ihit] = 0;
397  h_enta[ihit] = 0;
398  h_entaRS[ihit] = 0;
399  h_driftT[ihit] = 0;
400  // h_driftD[ihit] = 0;
401  h_facnladc[ihit] = 0;
402  h_wireGain[ihit] = 0;
403  h_twodCor[ihit] = 0;
404  h_onedCor[ihit] = 0;
405  h_WeightPionHypo[ihit] = 0;
406  h_WeightKaonHypo[ihit] = 0;
407  h_WeightProtonHypo[ihit] = 0;
408  h_foundByTrackFinder[ihit] = 0;
409  }
410  }
411 }
412 
413 void
415 {
416  // register output root file
417  m_file.push_back(new TFile(filename.c_str(), "RECREATE"));
418  m_tree.push_back(new TTree("track", "dE/dx information"));
419 
420  int i = m_tree.size() - 1;
421  m_tree[i]->SetDirectory(0);
422 
423  // event level information (from emd)
424  m_tree[i]->Branch("exp", &m_expID, "exp/I");
425  m_tree[i]->Branch("run", &m_runID, "run/I");
426  m_tree[i]->Branch("event", &m_eventID, "event/I");
427 
428  // track level information (from tfr)
429  m_tree[i]->Branch("d0", &m_d0, "d0/D");
430  m_tree[i]->Branch("z0", &m_z0, "z0/D");
431  m_tree[i]->Branch("dz", &m_dz, "dz/D");
432  m_tree[i]->Branch("dr", &m_dr, "dr/D");
433  m_tree[i]->Branch("dphi", &m_dphi, "dphi/D");
434  m_tree[i]->Branch("vx0", &m_vx0, "vx0/D");
435  m_tree[i]->Branch("vy0", &m_vy0, "vy0/D");
436  m_tree[i]->Branch("vz0", &m_vz0, "vz0/D");
437  m_tree[i]->Branch("tanlambda", &m_tanlambda, "tanlambda/D");
438  m_tree[i]->Branch("phi0", &m_phi0, "phi0/D");
439  m_tree[i]->Branch("chi2", &m_chi2, "chi2/D");
440 
441  // track level information (from cdt)
442  m_tree[i]->Branch("nCDChits", &m_nCDChits, "nCDChits/D");
443  m_tree[i]->Branch("inCDC", &m_inCDC, "inCDC/I");
444  m_tree[i]->Branch("track", &m_trackID, "track/I");
445  m_tree[i]->Branch("length", &m_length, "length/D");
446  m_tree[i]->Branch("charge", &m_charge, "charge/I");
447  m_tree[i]->Branch("costh", &m_cosTheta, "costh/D");
448  m_tree[i]->Branch("pF", &m_pCDC, "pF/D");
449  m_tree[i]->Branch("p", &m_p, "p/D");
450  m_tree[i]->Branch("pt", &m_pt, "pt/D");
451  m_tree[i]->Branch("ioasym", &m_ioasym, "ioasym/D");
452  m_tree[i]->Branch("phi", &m_phi, "phi/D");
453  m_tree[i]->Branch("pdg", &m_PDG, "pdg/D");
454  // track level dE/dx measurements
455  m_tree[i]->Branch("mean", &m_mean, "mean/D");
456  m_tree[i]->Branch("dedx", &m_trunc, "dedx/D");
457  m_tree[i]->Branch("dedxnosat", &m_truncNoSat, "dedxnosat/D");
458  m_tree[i]->Branch("dedxerr", &m_error, "dedxerr/D");
459 
460 
461  // other track level information
462  m_tree[i]->Branch("eop", &m_eop, "eop/D");
463  m_tree[i]->Branch("e", &m_e, "e/D");
464  m_tree[i]->Branch("e1_9", &m_e1_9, "e1_9/D");
465  m_tree[i]->Branch("e9_21", &m_e9_21, "e9_21/D");
466  m_tree[i]->Branch("klmLayers", &m_klmLayers, "klmLayers/I");
467 
468  // calibration constants
469  m_tree[i]->Branch("scale", &m_scale, "scale/D");
470  m_tree[i]->Branch("coscor", &m_cosCor, "coscor/D");
471  m_tree[i]->Branch("cosedgecor", &m_cosEdgeCor, "cosedgecor/D");
472  m_tree[i]->Branch("rungain", &m_runGain, "rungain/D");
473 
474  // PID values
475  m_tree[i]->Branch("chiE", &m_chie, "chiE/D");
476  m_tree[i]->Branch("chiMu", &m_chimu, "chiMu/D");
477  m_tree[i]->Branch("chiPi", &m_chipi, "chiPi/D");
478  m_tree[i]->Branch("chiK", &m_chik, "chiK/D");
479  m_tree[i]->Branch("chiP", &m_chip, "chiP/D");
480  m_tree[i]->Branch("chiD", &m_chid, "chiD/D");
481 
482  m_tree[i]->Branch("pmeanE", &m_pmeane, "pmeanE/D");
483  m_tree[i]->Branch("pmeanMu", &m_pmeanmu, "pmeanMu/D");
484  m_tree[i]->Branch("pmeanPi", &m_pmeanpi, "pmeanPi/D");
485  m_tree[i]->Branch("pmeanK", &m_pmeank, "pmeanK/D");
486  m_tree[i]->Branch("pmeanP", &m_pmeanp, "pmeanP/D");
487  m_tree[i]->Branch("pmeanD", &m_pmeand, "pmeanD/D");
488 
489  m_tree[i]->Branch("presE", &m_prese, "presE/D");
490  m_tree[i]->Branch("presMu", &m_presmu, "presMu/D");
491  m_tree[i]->Branch("presPi", &m_prespi, "presPi/D");
492  m_tree[i]->Branch("presK", &m_presk, "presK/D");
493  m_tree[i]->Branch("presP", &m_presp, "presP/D");
494  m_tree[i]->Branch("presD", &m_presd, "presD/D");
495 
496  // layer level information
497  m_tree[i]->Branch("lNHits", &l_nhits, "lNHits/I");
498  m_tree[i]->Branch("lNHitsUsed", &l_nhitsused, "lNHitsUsed/I");
499  m_tree[i]->Branch("lNHitsCombined", l_nhitscombined, "lNHitsCombined[lNHits]/I");
500  m_tree[i]->Branch("lWireLongestHit", l_wirelongesthit, "lWireLongestHit[lNHits]/I");
501  m_tree[i]->Branch("lLayer", l_layer, "lLayer[lNHits]/I");
502  m_tree[i]->Branch("lPath", l_path, "lPath[lNHits]/D");
503  m_tree[i]->Branch("lDedx", l_dedx, "lDedx[lNHits]/D");
504 
505  // hit level information
506  if (enableHitLevel) {
507  m_tree[i]->Branch("hNHits", &h_nhits, "hNHits/I");
508  m_tree[i]->Branch("hLWire", h_lwire, "hLWire[hNHits]/I");
509  m_tree[i]->Branch("hWire", h_wire, "hWire[hNHits]/I");
510  m_tree[i]->Branch("hLayer", h_layer, "hLayer[hNHits]/I");
511  m_tree[i]->Branch("hPath", h_path, "hPath[hNHits]/D");
512  m_tree[i]->Branch("hDedx", h_dedx, "hDedx[hNHits]/D");
513  m_tree[i]->Branch("hADCRaw", h_adcraw, "hADCRaw[hNHits]/D");
514  m_tree[i]->Branch("hADCCorr", h_adccorr, "hADCCorr[hNHits]/D");
515  m_tree[i]->Branch("hDoca", h_doca, "hDoca[hNHits]/D");
516  m_tree[i]->Branch("hNDoca", h_ndoca, "hNDoca[hNHits]/D");
517  m_tree[i]->Branch("hNDocaRS", h_ndocaRS, "hNDocaRS[hNHits]/D");
518  m_tree[i]->Branch("hEnta", h_enta, "hEnta[hNHits]/D");
519  m_tree[i]->Branch("hEntaRS", h_entaRS, "hEntaRS[hNHits]/D");
520  m_tree[i]->Branch("hDriftT", h_driftT, "hDriftT[hNHits]/D");
521  m_tree[i]->Branch("hDriftD", h_driftD, "hDriftD[hNHits]/D");
522  m_tree[i]->Branch("hFacnlADC", h_facnladc, "hFacnlADC[hNHits]/D");
523  m_tree[i]->Branch("hWireGain", h_wireGain, "hWireGain[hNHits]/D");
524  m_tree[i]->Branch("hTwodcor", h_twodCor, "hTwodcor[hNHits]/D");
525  m_tree[i]->Branch("hOnedcor", h_onedCor, "hOnedcor[hNHits]/D");
526 
527  if (enableExtraVar) {
528  m_tree[i]->Branch("hWeightPionHypo", h_WeightPionHypo, "hWeightPionHypo[hNHits]/D");
529  m_tree[i]->Branch("hWeightKaonHypo", h_WeightKaonHypo, "hWeightKaonHypo[hNHits]/D");
530  m_tree[i]->Branch("hWeightProtonHypo", h_WeightProtonHypo, "hWeightProtonHypo[hNHits]/D");
531  m_tree[i]->Branch("hFoundByTrackFinder", h_foundByTrackFinder, "hFoundByTrackFinder[hNHits]/I");
532  }
533  }
534 }
Belle2::CDCDedxTrack::getPDG
double getPDG() const
Get the identity of the particle.
Definition: CDCDedxTrack.h:111
Belle2::CDCDedxTrack::getDedxNoSat
double getDedxNoSat() const
Get dE/dx truncated mean without the saturation correction for this track.
Definition: CDCDedxTrack.h:117
Belle2::CDCDedxTrack::getHitLayer
int getHitLayer(int i) const
Return the (global) layer number for a hit.
Definition: CDCDedxTrack.h:205
Belle2::UncertainHelix
This class represents an ideal helix in perigee parameterization including the covariance matrix of t...
Definition: UncertainHelix.h:40
Belle2::ECLCluster::getE9oE21
double getE9oE21() const
Return E9/E21 (shower shape variable).
Definition: ECLCluster.h:295
Belle2::HitLevelInfoWriterModule::event
virtual void event() override
This method is called for each event.
Definition: HitLevelInfoWriter.cc:62
Belle2::CDCDedxTrack::getDedx
double getDedx() const
Get dE/dx truncated mean for this track.
Definition: CDCDedxTrack.h:114
Belle2::CDCDedxTrack::getEntaRS
double getEntaRS(int i) const
Return rescaled enta value for cell height=width assumption.
Definition: CDCDedxTrack.h:232
Belle2::KLMCluster::getLayers
int getLayers() const
Get number of layers with hits.
Definition: KLMCluster.h:76
Belle2::CDCDedxTrack::getChi
double getChi(int i) const
Return the PID (chi) value.
Definition: CDCDedxTrack.h:275
Belle2::TrackFitResult::getMomentum
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Definition: TrackFitResult.h:116
Belle2::ECLCluster
ECL cluster data.
Definition: ECLCluster.h:39
Belle2::CDCDedxTrack::getWire
int getWire(int i) const
Return the sensor ID for this hit: wire number for CDC (0-14336)
Definition: CDCDedxTrack.h:202
Belle2::Const::electron
static const ChargedStable electron
electron particle
Definition: Const.h:533
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::HitLevelInfoWriterModule
Extracts dE/dx information for calibration testing.
Definition: HitLevelInfoWriter.h:61
Belle2::CDCDedxTrack::getLayer
int getLayer(int i) const
Return the (global) layer number for a layer hit.
Definition: CDCDedxTrack.h:182
Belle2::HitLevelInfoWriterModule::fillDedx
void fillDedx(CDCDedxTrack *dedxTrack)
Fill the TTree with the information from a CDCDedxTrack object.
Definition: HitLevelInfoWriter.cc:261
Belle2::TrackFitResult::getPValue
double getPValue() const
Getter for Chi2 Probability of the track fit.
Definition: TrackFitResult.h:163
Belle2::CDCDedxTrack::getMomentum
double getMomentum() const
Return the track momentum valid in the CDC.
Definition: CDCDedxTrack.h:129
Belle2::ECLCluster::EHypothesisBit::c_nPhotons
@ c_nPhotons
CR is split into n photons (N1)
Belle2::RelationsInterface::getRelated
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
Definition: RelationsObject.h:280
Belle2::HitLevelInfoWriterModule::fillTrack
void fillTrack(const TrackFitResult *fitResult)
Fill the TTree with the information from the track fit.
Definition: HitLevelInfoWriter.cc:223
Belle2::PIDLikelihood
Class to collect log likelihoods from TOP, ARICH, dEdx, ECL and KLM aimed for output to mdst includes...
Definition: PIDLikelihood.h:37
Belle2::TrackFitResult::getTanLambda
double getTanLambda() const
Getter for tanLambda.
Definition: TrackFitResult.h:206
Belle2::CDCDedxTrack::getADCCount
int getADCCount(int i) const
Return the adcCount for this hit.
Definition: CDCDedxTrack.h:214
Belle2::UncertainHelix::passiveMoveBy
double passiveMoveBy(const TVector3 &by)
Moves origin of the coordinate system (passive transformation) by the given vector.
Definition: UncertainHelix.h:121
Belle2::CDCDedxTrack::getLength
double getLength() const
Return the total path length for this track.
Definition: CDCDedxTrack.h:138
Belle2::ECLCluster::hasHypothesis
bool hasHypothesis(EHypothesisBit bitmask) const
Return if specific hypothesis bit is set.
Definition: ECLCluster.h:374
Belle2::CDCDedxTrack::getEnta
double getEnta(int i) const
Return the entrance angle in the CDC cell for this hit.
Definition: CDCDedxTrack.h:226
Belle2::CDCDedxTrack::size
int size() const
Return the number of hits for this track.
Definition: CDCDedxTrack.h:196
Belle2::TrackFitResult
Values of the result of a track fit with a given particle hypothesis.
Definition: TrackFitResult.h:59
Belle2::CDCDedxTrack::getLayerDedx
double getLayerDedx(int i) const
Return the total dE/dx for this layer.
Definition: CDCDedxTrack.h:188
Belle2::TrackFitResult::getPosition
TVector3 getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
Definition: TrackFitResult.h:109
Belle2::CDCDedxTrack::getCharge
int getCharge() const
Return the charge for this track.
Definition: CDCDedxTrack.h:135
Belle2::HitLevelInfoWriterModule::bookOutput
void bookOutput(std::string filename)
Create the output TFiles and TTrees.
Definition: HitLevelInfoWriter.cc:414
Belle2::TrackFitResult::getHitPatternCDC
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
Definition: TrackFitResult.cc:120
Belle2::TrackFitResult::getZ0
double getZ0() const
Getter for z0.
Definition: TrackFitResult.h:200
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::Const::kaon
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:536
Belle2::CDCDedxTrack::getDriftD
double getDriftD(int i) const
Return the drift distance for this hit.
Definition: CDCDedxTrack.h:247
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::CDCDedxTrack::getFoundByTrackFinder
int getFoundByTrackFinder(int i) const
Return the TrackFinder which added the given hit to track.
Definition: CDCDedxTrack.h:263
Belle2::ECLCluster::getEnergy
double getEnergy(const EHypothesisBit &hypothesis) const
Return Energy (GeV).
Definition: ECLCluster.cc:21
Belle2::Const::pion
static const ChargedStable pion
charged pion particle
Definition: Const.h:535
Belle2::CDCDedxTrack::getPath
double getPath(int i) const
Return the path length through the cell for this hit.
Definition: CDCDedxTrack.h:208
Belle2::CDCDedxTrack::getLayerPath
double getLayerPath(int i) const
Return the distance travelled in this layer.
Definition: CDCDedxTrack.h:185
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::CDCDedxTrack::getWireLongestHit
int getWireLongestHit(int i) const
Return the wire number of the longest hit per layer.
Definition: CDCDedxTrack.h:179
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::Const::deuteron
static const ChargedStable deuteron
deuteron particle
Definition: Const.h:538
Belle2::CDCDedxTrack::getDedxMean
double getDedxMean() const
Get the dE/dx mean for this track.
Definition: CDCDedxTrack.h:123
Belle2::CDCDedxTrack::getWireInLayer
int getWireInLayer(int i) const
Return the sensor ID for this hit: wire number in the layer.
Definition: CDCDedxTrack.h:199
Belle2::CDCDedxTrack::getDoca
double getDoca(int i) const
Return the distance of closest approach to the sense wire for this hit.
Definition: CDCDedxTrack.h:223
Belle2::CDCDedxTrack::getNHitsCombined
int getNHitsCombined(int i) const
Return the number of hits combined per layer.
Definition: CDCDedxTrack.h:176
Belle2::TrackFitResult::getPhi0
double getPhi0() const
Getter for phi0.
Definition: TrackFitResult.h:184
Belle2::CDCDedxTrack::getDriftT
int getDriftT(int i) const
Return the drift time for this hit.
Definition: CDCDedxTrack.h:235
Belle2::CDCDedxTrack::getNLayerHitsUsed
double getNLayerHitsUsed() const
Return the number of hits used to determine the truncated mean.
Definition: CDCDedxTrack.h:173
Belle2::CDCDedxTrack::getDedxError
double getDedxError() const
Get the error on the dE/dx truncated mean for this track.
Definition: CDCDedxTrack.h:120
Belle2::KLMCluster
KLM cluster data.
Definition: KLMCluster.h:38
Belle2::CDCDedxTrack::getCellHalfWidth
double getCellHalfWidth(int i) const
Return the half-width of the CDC cell.
Definition: CDCDedxTrack.h:244
Belle2::CDCDedxTrack::getADCBaseCount
int getADCBaseCount(int i) const
Return the base adcCount (no non-linearity) for this hit.
Definition: CDCDedxTrack.h:217
Belle2::HitLevelInfoWriterModule::~HitLevelInfoWriterModule
virtual ~HitLevelInfoWriterModule()
Destructor.
Definition: HitLevelInfoWriter.cc:33
Belle2::Particle
Class to store reconstructed particles.
Definition: Particle.h:77
Belle2::CDCDedxTrack::getPreso
double getPreso(int i) const
Return the PID (predicted reso) value.
Definition: CDCDedxTrack.h:281
Belle2::CDCDedxTrack::getWeightPionHypo
double getWeightPionHypo(int i) const
Return the max weights from KalmanFitterInfo using pion hypothesis.
Definition: CDCDedxTrack.h:266
Belle2::CDCDedxTrack::getWeightKaonHypo
double getWeightKaonHypo(int i) const
Return the max weights from KalmanFitterInfo using kaon hypothesis.
Definition: CDCDedxTrack.h:269
Belle2::TrackFitResult::getUncertainHelix
UncertainHelix getUncertainHelix() const
Conversion to framework Uncertain Helix (i.e., with covariance).
Definition: TrackFitResult.h:234
Belle2::ECLCluster::getE1oE9
double getE1oE9() const
Return E1/E9 (shower shape variable).
Definition: ECLCluster.h:292
Belle2::ReferenceFrame::GetCurrent
static const ReferenceFrame & GetCurrent()
Get current rest frame.
Definition: ReferenceFrame.cc:28
Belle2::Const::proton
static const ChargedStable proton
proton particle
Definition: Const.h:537
Belle2::CDCDedxTrack::getPmean
double getPmean(int i) const
Return the PID (predicted mean) value.
Definition: CDCDedxTrack.h:278
Belle2::CDCDedxTrack::trackID
double trackID() const
Return the track ID.
Definition: CDCDedxTrack.h:108
Belle2::CDCDedxTrack::getCosTheta
double getCosTheta() const
Return cos(theta) for this track.
Definition: CDCDedxTrack.h:132
Belle2::HitPatternCDC::getNHits
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
Definition: HitPatternCDC.cc:49
Belle2::Const::muon
static const ChargedStable muon
muon particle
Definition: Const.h:534
Belle2::Track
Class that bundles various TrackFitResults.
Definition: Track.h:35
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::CDCDedxTrack
Debug output for CDCDedxPID module.
Definition: CDCDedxTrack.h:36
Belle2::CDCDedxTrack::getNLayerHits
int getNLayerHits() const
Return the number of layer hits for this track.
Definition: CDCDedxTrack.h:170
Belle2::CDCDedxTrack::getDocaRS
double getDocaRS(int i) const
Return rescaled doca value for cell height=width assumption.
Definition: CDCDedxTrack.h:229
Belle2::HitLevelInfoWriterModule::initialize
virtual void initialize() override
Initialize the module.
Definition: HitLevelInfoWriter.cc:35
Belle2::CDCDedxTrack::getWeightProtonHypo
double getWeightProtonHypo(int i) const
Return the max weights from KalmanFitterInfo using proton hypothesis.
Definition: CDCDedxTrack.h:272
Belle2::CDCDedxTrack::getNonLADCCorrection
double getNonLADCCorrection(int i) const
Return the factor introduce for adcCount (non-linearity) correction.
Definition: CDCDedxTrack.h:220
Belle2::TrackFitResult::getChargeSign
short getChargeSign() const
Return track charge (1 or -1).
Definition: TrackFitResult.h:160
Belle2::HitLevelInfoWriterModule::clearEntries
void clearEntries()
Clear the arrays before filling an event.
Definition: HitLevelInfoWriter.cc:375
Belle2::RelationsInterface::getRelatedFrom
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
Definition: RelationsObject.h:265
Belle2::HitLevelInfoWriterModule::terminate
virtual void terminate() override
End of the event processing.
Definition: HitLevelInfoWriter.cc:209
Belle2::TrackFitResult::getD0
double getD0() const
Getter for d0.
Definition: TrackFitResult.h:178