Belle II Software  release-06-02-00
CDCCosmicAnalysisModule.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 #include "cdc/modules/cdcCosmicAnalysis/CDCCosmicAnalysisModule.h"
10 #include <framework/geometry/BFieldManager.h>
11 #include <framework/gearbox/Const.h>
12 #include <framework/datastore/RelationArray.h>
13 
14 #include <Math/ProbFuncMathCore.h>
15 #include "algorithm"
16 
17 using namespace std;
18 using namespace Belle2;
19 using namespace CDC;
20 
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(CDCCosmicAnalysis)
26 
27 
28 //-----------------------------------------------------------------
29 // Implementation
30 //-----------------------------------------------------------------
31 
33 {
34  setDescription("Module for harvesting parameters of the two half (up/down) of a cosmic track for performance study");
35  setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
36  addParam("RecoTracksColName", m_recoTrackArrayName, "Name of collectrion hold RecoTracks", std::string(""));
37  addParam("Output", m_outputFileName, "output file name", string("twotracks.root"));
38  addParam("EventT0Extraction", m_eventT0Extraction, "use event t0 extract t0 or not", true);
39  addParam("treeName", m_treeName, "Output tree name", string("tree"));
40  addParam("phi0InRad", m_phi0InRad, "Phi0 in unit of radian, true: rad, false: deg", true);
41  addParam("StoreTrackParErrors", m_storeTrackParErrors,
42  "Store Track Parameter errors (true) or not (false)", false);
43 }
44 
45 CDCCosmicAnalysisModule::~CDCCosmicAnalysisModule()
46 {
47 }
48 void CDCCosmicAnalysisModule::initialize()
49 {
50  // Register histograms (calls back defineHisto)
51  m_Tracks.isRequired(m_trackArrayName);
52  m_RecoTracks.isRequired(m_recoTrackArrayName);
53  m_TrackFitResults.isRequired(m_trackFitResultArrayName);
54  RelationArray relRecoTrackTrack(m_RecoTracks, m_Tracks, m_relRecoTrackTrackName);
55 
56  m_relRecoTrackTrackName = relRecoTrackTrack.getName();
57 
58  tfile = new TFile(m_outputFileName.c_str(), "RECREATE");
59  // tree = new TTree("treeTrk", "treeTrk");
60  tree = new TTree(m_treeName.c_str(), m_treeName.c_str());
61  tree->Branch("run", &run, "run/I");
62  tree->Branch("evtT0", &evtT0, "evtT0/D");
63  tree->Branch("charge", &charge, "charge/S");
64 
65  tree->Branch("Pval1", &Pval1, "Pval1/D");
66  tree->Branch("ndf1", &ndf1, "ndf1/D");
67  tree->Branch("Pt1", &Pt1, "Pt1/D");
68  tree->Branch("D01", &D01, "D01/D");
69  tree->Branch("Phi01", &Phi01, "Phi01/D");
70  // tree->Branch("Om1", &Om1, "Om1/D");
71  tree->Branch("Z01", &Z01, "Z01/D");
72  tree->Branch("tanLambda1", &tanLambda1, "tanLambda1/D");
73  tree->Branch("posSeed1", "TVector3", &posSeed1);
74  tree->Branch("Omega1", &Omega1, "Omega1/D");
75  tree->Branch("Mom1", "TVector3", &Mom1);
76 
77  tree->Branch("Pval2", &Pval2, "Pval2/D");
78  tree->Branch("ndf2", &ndf2, "ndf2/D");
79  tree->Branch("Pt2", &Pt2, "Pt2/D");
80  tree->Branch("D02", &D02, "D02/D");
81  tree->Branch("Phi02", &Phi02, "Phi02/D");
82  // tree->Branch("Om2", &Om2, "Om2/D");
83  tree->Branch("Z02", &Z02, "Z02/D");
84  tree->Branch("tanLambda2", &tanLambda2, "tanLambda2/D");
85  tree->Branch("posSeed2", "TVector3", &posSeed2);
86  tree->Branch("Omega2", &Omega2, "Omega2/D");
87  tree->Branch("Mom2", "TVector3", &Mom2);
88 
89  if (m_storeTrackParErrors) {
90  tree->Branch("eD01", &eD01, "eD01/D");
91  tree->Branch("ePhi01", &ePhi01, "ePhi01/D");
92  tree->Branch("eOm1", &eOm1, "eOm1/D");
93  tree->Branch("eZ01", &eZ01, "eZ01/D");
94  tree->Branch("etanL1", &etanL1, "etanL1/D");
95 
96  tree->Branch("eD02", &eD02, "eD02/D");
97  tree->Branch("ePhi02", &ePhi02, "ePhi02/D");
98  tree->Branch("eOm2", &eOm2, "eOm2/D");
99  tree->Branch("eZ02", &eZ02, "eZ02/D");
100  tree->Branch("etanL2", &etanL2, "etanL2/D");
101  }
102 
103 }
104 
105 void CDCCosmicAnalysisModule::beginRun()
106 {
107  B2Vector3D pos(0, 0, 0);
108  B2Vector3D bfield = BFieldManager::getFieldInTesla(pos);
109  if (bfield.Z() > 0.5) {
110  m_bField = true;
111  B2INFO("CDCCosmicAnalysis: Magnetic field is ON");
112  } else {
113  m_bField = false;
114  B2INFO("CDCCosmicAnalysis: Magnetic field is OFF");
115  }
116  B2INFO("CDCCosmicAnalysis: BField at (0,0,0) = " << bfield.Mag());
117 }
118 
119 void CDCCosmicAnalysisModule::event()
120 {
121  if (m_EventMetaData.isValid())
122  run = m_EventMetaData->getRun();
123 
124  evtT0 = 0;
125  if (m_eventT0Extraction) {
126  // event with is fail to extract event-t0 will be excluded
127  if (m_eventTimeStoreObject.isValid() && m_eventTimeStoreObject->hasEventT0()) {
128  evtT0 = m_eventTimeStoreObject->getEventT0();
129  } else {
130  return;
131  }
132  }
133 
134  // Loop over Tracks
135  int nTr = m_Tracks.getEntries();
136 
137  short charge2 = 0;
138  short charge1 = 0;
139  bool up(false), down(false);
140  for (int i = 0; i < nTr; ++i) {
141  const Belle2::Track* b2track = m_Tracks[i];
142  const Belle2::TrackFitResult* fitresult;
143  fitresult = b2track->getTrackFitResult(Const::muon);
144  if (!fitresult) {
145  B2WARNING("There is no track fit result for muon hypothesis; try with the closest mass hypothesis...");
146  fitresult = b2track->getTrackFitResultWithClosestMass(Const::muon);
147  if (!fitresult) {
148  B2WARNING("There is also no track fit reslt for the other mass hypothesis");
149  continue;
150  }
151  }
152  const Belle2::RecoTrack* recoTrack = b2track->getRelatedTo<Belle2::RecoTrack>(m_recoTrackArrayName);
153  if (!recoTrack) {
154  B2WARNING("Can not access RecoTrack of this Belle2::Track");
155  continue;
156  }
157 
158  TVector3 posSeed = recoTrack->getPositionSeed();
159  const genfit::FitStatus* fs = recoTrack->getTrackFitStatus();
160 
161  double ndf = fs->getNdf();
162  if (!m_bField) // in case of no magnetic field, #track par=4 instead of 5.
163  ndf += 1;
164 
165  double Chi2 = fs->getChi2();
166  double TrPval = std::max(0., ROOT::Math::chisquared_cdf_c(Chi2, ndf));
167  double Phi0 = fitresult->getPhi0();
168  if (m_phi0InRad == false) { // unit of degrees.
169  Phi0 *= 180 / M_PI;
170  }
171 
172  /*** Two track case.***/
173  if ((posSeed.Y() > 0 && !up) ||
174  (up && posSeed.Y()*posSeed1.Y() > 0 && ndf > ndf1)) {
175  charge1 = fitresult->getChargeSign();
176  posSeed1 = posSeed;
177  D01 = fitresult->getD0();
178  eD01 = sqrt(fitresult->getCovariance5()[0][0]);
179  Phi01 = Phi0;
180  ePhi01 = sqrt(fitresult->getCovariance5()[1][1]);
181  if (m_phi0InRad == false) { // unit of degrees.
182  ePhi01 *= 180 / M_PI;
183  }
184 
185  Omega1 = fitresult->getOmega();
186  Mom1 = fitresult->getMomentum();
187  eOm1 = sqrt(fitresult->getCovariance5()[2][2]);
188  Z01 = fitresult->getZ0();
189  eZ01 = sqrt(fitresult->getCovariance5()[3][3]);
190  tanLambda1 = fitresult->getTanLambda();
191  etanL1 = sqrt(fitresult->getCovariance5()[4][4]);
192  Pt1 = fitresult->getTransverseMomentum();
193  ndf1 = ndf;
194  Pval1 = TrPval;
195  up = true;
196  }
197 
198  if ((posSeed.Y() < 0 && !down) ||
199  (down && posSeed.Y()*posSeed2.Y() > 0 && ndf > ndf2)) {
200  charge2 = fitresult->getChargeSign();
201  posSeed2 = posSeed;
202  D02 = fitresult->getD0();
203  eD02 = sqrt(fitresult->getCovariance5()[0][0]);
204  Phi02 = Phi0;
205  ePhi02 = sqrt(fitresult->getCovariance5()[1][1]);
206  if (m_phi0InRad == false) { // unit of degrees.
207  ePhi02 *= 180 / M_PI;
208  }
209  Omega2 = fitresult->getOmega();
210  Mom2 = fitresult->getMomentum();
211  eOm2 = sqrt(fitresult->getCovariance5()[2][2]);
212  Z02 = fitresult->getZ0();
213  eZ02 = sqrt(fitresult->getCovariance5()[3][3]);
214  tanLambda2 = fitresult->getTanLambda();
215  etanL2 = sqrt(fitresult->getCovariance5()[4][4]);
216  Pt2 = fitresult->getTransverseMomentum();
217  ndf2 = ndf;
218  Pval2 = TrPval;
219  down = true;
220  }
221  }
222  if (m_bField && charge1 * charge2 == 0) return;
223  if (charge1 * charge2 >= 0 && up && down) {
224  charge = charge1;
225  tree->Fill();
226  }
227 
228 }
229 
230 void CDCCosmicAnalysisModule::endRun()
231 {
232 }
233 
234 void CDCCosmicAnalysisModule::terminate()
235 {
236  tfile->cd();
237  tree->Write();
238  tfile->Close();
239 }
240 
241 
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:420
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:144
Analysis module for CDC CR data.
Base class for Modules.
Definition: Module.h:72
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:76
const genfit::FitStatus * getTrackFitStatus(const genfit::AbsTrackRep *representation=nullptr) const
Return the track fit status for the given representation or for the cardinal one. You are not allowed...
Definition: RecoTrack.h:536
TVector3 getPositionSeed() const
Return the position seed stored in the reco track. ATTENTION: This is not the fitted position.
Definition: RecoTrack.h:476
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
const std::string & getName() const
Return name under which the object is saved in the DataStore.
Values of the result of a track fit with a given particle hypothesis.
TMatrixDSym getCovariance5() const
Getter for covariance matrix of perigee parameters in matrix form.
short getChargeSign() const
Return track charge (1 or -1).
double getOmega() const
Getter for omega.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
double getD0() const
Getter for d0.
double getTransverseMomentum() const
Getter for the absolute value of the transverse momentum at the perigee.
double getTanLambda() const
Getter for tanLambda.
double getZ0() const
Getter for z0.
double getPhi0() const
Getter for phi0.
Class that bundles various TrackFitResults.
Definition: Track.h:25
const TrackFitResult * getTrackFitResult(const Const::ChargedStable &chargedStable) const
Access to TrackFitResults.
Definition: Track.cc:17
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for a fit hypothesis with the closest mass.
Definition: Track.cc:68
Class where important numbers and properties of a fit can be stored.
Definition: FitStatus.h:80
double getChi2() const
Get chi^2 of the fit.
Definition: FitStatus.h:120
double getNdf() const
Get the degrees of freedom of the fit.
Definition: FitStatus.h:122
#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.