Belle II Software  release-08-01-10
BKLMAnaModule.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 header. */
10 #include <klm/bklm/modules/bklmAna/BKLMAnaModule.h>
11 
12 /* Basf2 headers. */
13 #include <framework/dataobjects/EventMetaData.h>
14 #include <framework/datastore/StoreObjPtr.h>
15 #include <mdst/dataobjects/Track.h>
16 #include <mdst/dataobjects/TrackFitResult.h>
17 #include <klm/dataobjects/KLMMuidLikelihood.h>
18 
19 /* CLHEP headers. */
20 #include <CLHEP/Units/SystemOfUnits.h>
21 
22 using namespace Belle2;
23 using namespace CLHEP;
24 
25 //-----------------------------------------------------------------
26 // Register the Module
27 //-----------------------------------------------------------------
28 REG_MODULE(BKLMAna);
29 
30 //-----------------------------------------------------------------
31 // Implementation
32 //-----------------------------------------------------------------
33 
35  m_file(nullptr),
36  m_extTree(nullptr),
37  m_run(0),
38  m_nExtHit(0),
39  m_hdistance(nullptr),
40  m_totalMom(nullptr),
41  m_passMom(nullptr),
42  m_effiMom(nullptr),
43  m_totalYX(nullptr),
44  m_passYX(nullptr),
45  m_effiYX(nullptr),
46  m_totalYZ(nullptr),
47  m_passYZ(nullptr),
48  m_effiYZ(nullptr)
49 {
50  for (int i = 0; i < 200; ++i) {
51  m_extx[i] = 0.0;
52  m_exty[i] = 0.0;
53  m_extz[i] = 0.0;
54  }
55  for (int i = 0; i < 15; ++i) {
56  m_totalThephi[i] = nullptr;
57  m_passThephi[i] = nullptr;
58  m_effiThephi[i] = nullptr;
59  m_totalTrkThephi[i] = nullptr;
60  m_passTrkThephi[i] = nullptr;
61  m_effiTrkThephi[i] = nullptr;
62  }
63  setDescription("analyze bklm efficiency associated to CDC, check performance of bklm et al.");
64  addParam("filename", m_filename, "Output root filename", std::string("bklmana.root"));
65 }
66 
68 {
69 }
70 
72 {
73  hits2D.isRequired();
74  extHits.isRequired();
75  tracks.isRequired();
76 
77 
78  m_file = new TFile(m_filename.c_str(), "RECREATE");
79 
80  m_extTree = new TTree("exthit", "ext hit");
81  m_extTree->Branch("run", &m_run, "m_run/I");
82  m_extTree->Branch("nExtHit", &m_nExtHit, "m_nExtHit/I");
83  m_extTree->Branch("x", &m_extx, "m_extx[m_nExtHit]/F");
84  m_extTree->Branch("y", &m_exty, "m_exty[m_nExtHit]/F");
85  m_extTree->Branch("z", &m_extz, "m_extz[m_nExtHit]/F");
86 
87  float phiBins = 36;
88  float phiMin = 0.0;
89  float phiMax = 360.0;
90  float thetaBins = 19;
91  float thetaMin = 35.0;
92  float thetaMax = 130.0;
93  char hname[50];
94  for (int iL = 0; iL < 15 ; iL ++) {
95  //based on theta phi of ExtHit position
96  sprintf(hname, "denominator_Layer%i", iL + 1);
97  m_totalThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
98  m_totalThephi[iL]->GetXaxis()->SetTitle("phi");
99  m_totalThephi[iL]->GetYaxis()->SetTitle("theta");
100  sprintf(hname, "numerator_Layer%i", iL + 1);
101  m_passThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
102  sprintf(hname, "effi_Layer%i", iL + 1);
103  m_passThephi[iL]->GetXaxis()->SetTitle("phi");
104  m_passThephi[iL]->GetYaxis()->SetTitle("theta");
105  m_effiThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
106  m_effiThephi[iL]->GetXaxis()->SetTitle("phi");
107  m_effiThephi[iL]->GetYaxis()->SetTitle("theta");
108  //based on theta phi of trk
109  sprintf(hname, "Denominator_Layer%i", iL + 1);
110  m_totalTrkThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
111  m_totalTrkThephi[iL]->GetXaxis()->SetTitle("trk.phi");
112  m_totalTrkThephi[iL]->GetYaxis()->SetTitle("trk.theta");
113  sprintf(hname, "Numerator_Layer%i", iL + 1);
114  m_passTrkThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
115  sprintf(hname, "Effi_Layer%i", iL + 1);
116  m_passTrkThephi[iL]->GetXaxis()->SetTitle("trk.phi");
117  m_passTrkThephi[iL]->GetYaxis()->SetTitle("trk.theta");
118  m_effiTrkThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
119  m_effiTrkThephi[iL]->GetXaxis()->SetTitle("trk.phi");
120  m_effiTrkThephi[iL]->GetYaxis()->SetTitle("trk.theta");
121  }
122  float gmin = -350;
123  float gmax = 350;
124  int gNbin = 150;
125  float mommin = 0.;
126  float mommax = 15;
127  int mNbin = 30;
128 
129  m_hdistance = new TH1F("m_hdistance", " distance between mathced extHit and bklmHit2d ", 100, 0, 30);
130 
131  m_totalMom = new TH1F("m_totalMom", " denominator vs. p", mNbin, mommin, mommax);
132  m_passMom = new TH1F("m_passMom", " numerator vs. p", mNbin, mommin, mommax);
133  m_effiMom = new TH1F("m_effiMom", " effi. vs. p", mNbin, mommin, mommax);
134  m_effiMom->GetXaxis()->SetTitle("p (GeV)");
135  m_effiMom->GetYaxis()->SetTitle("Efficiency (GeV)");
136 
137  m_totalYX = new TH2F("m_totalYX", " denominator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
138  m_passYX = new TH2F("m_passYX", " numerator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
139  m_totalYZ = new TH2F("m_totalYZ", " denominator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
140  m_passYZ = new TH2F("m_passYZ", " numerator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
141  m_effiYX = new TH2F("m_effiYX", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
142  m_effiYZ = new TH2F("m_effiYZ", " effi. Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
143  m_effiYX->GetXaxis()->SetTitle("x (cm)");
144  m_effiYX->GetYaxis()->SetTitle("y (cm)");
145  m_effiYZ->GetXaxis()->SetTitle("z (cm)");
146  m_effiYZ->GetYaxis()->SetTitle("y (cm)");
147 
148 }
149 
151 {
152 }
153 
155 {
156  StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
157  //unsigned long eventNumber = eventMetaData->getEvent();
158  unsigned long runNumber = eventMetaData->getRun();
159  //unsigned long expNumber = eventMetaData->getExperiment();
160  //StoreArray<RecoTrack> recoTracks;
161  //numRecoTrk = recoTracks.getEntries();
162  //StoreArray<BKLMTrack> bklmtracks;
163  //StoreArray<TrackFitResult> trackFitResults;
164  //StoreArray<MuidHit> muidHits;
165  //StoreArray<Muid> muids;
166  //StoreArray<MCParticle> mcParticles;
167 
168  //set<int> m_pointUsed;
169  //m_pointUsed.clear();
170  //check extHits
171  //all ExtHit in bklm scope in each event, should not be many
172  int nExtHit = 0;
173  for (int t = 0; t < extHits.getEntries(); t++) {
174  ExtHit* exthit = extHits[t];
175  if (exthit->getDetectorID() != Const::EDetector::BKLM)
176  continue;
177  m_extx[nExtHit] = exthit->getPosition().X();
178  m_exty[nExtHit] = exthit->getPosition().Y();
179  m_extz[nExtHit] = exthit->getPosition().Z();
180  nExtHit++;
181  if (nExtHit > 199)
182  break;
183  }
184  m_run = runNumber;
185  m_nExtHit = nExtHit;
186 
187 //the second way, require muid
188  for (int k = 0; k < tracks.getEntries(); k++) {
189  Track* track = tracks[k];
190  // load the muon fit hypothesis or the hypothesis which is the clostes in mass to a muon
191  // the tracking will not always fit a muon hypothesis
192  const TrackFitResult* fitres = track->getTrackFitResultWithClosestMass(Belle2::Const::muon);
193  double mom = fitres->getMomentum().R();
194  //double pt = fitres->getTransverseMomentum();
195  ROOT::Math::XYZVector p3 = fitres->getMomentum();
196  double trkphi = p3.Phi() * TMath::RadToDeg();
197  double trktheta = p3.Theta() * TMath::RadToDeg();
198  if (trkphi < 0)
199  trkphi = trkphi + 360.0;
200  RelationVector<KLMHit2d> relatedHit2D = track->getRelationsTo<KLMHit2d>();
201  RelationVector<ExtHit> relatedExtHit = track->getRelationsTo<ExtHit>();
202  for (unsigned int t = 0; t < relatedExtHit.size(); t++) {
203  ExtHit* exthit = relatedExtHit[t];
204  if (exthit->getDetectorID() != Const::EDetector::BKLM)
205  continue;
206  int module = exthit->getCopyID();
207  int section = BKLMElementNumbers::getSectionByModule(module);
208  int sector = BKLMElementNumbers::getSectorByModule(module);
209  int layer = BKLMElementNumbers::getLayerByModule(module);
210  bool crossed = false; // should be only once ?
211  KLMMuidLikelihood* muid = track->getRelatedTo<KLMMuidLikelihood>();
212  if (!muid)
213  continue;
214  int extPattern = muid->getExtLayerPattern();
215  if ((extPattern & (1 << (layer - 1))) != 0)
216  crossed = true;
217  if (!crossed)
218  continue;
219 
220  ROOT::Math::XYZVector extVec = exthit->getPosition();
221  bool matched = false;
222  m_totalYX->Fill(extVec.X(), extVec.Y());
223  m_totalYZ->Fill(extVec.Z(), extVec.Y());
224  float phi = extVec.Phi() * TMath::RadToDeg();
225  float theta = extVec.Theta() * TMath::RadToDeg();
226  if (phi < 0)
227  phi = phi + 360.0;
228  m_totalThephi[layer - 1]->Fill(phi, theta);
229  m_totalTrkThephi[layer - 1]->Fill(trkphi, trktheta);
230  m_totalMom->Fill(mom);
231  //look for mateched BKLM2dHit
232  //for (unsigned int mHit = 0; mHit < relatedHit2D.size(); mHit++) {
233  // KLMHit2d* hit = relatedHit2D[mHit];
234  for (int mHit = 0; mHit < hits2D.getEntries(); mHit++) {
235  KLMHit2d* hit = hits2D[mHit];
236  if (hit->getSubdetector() != KLMElementNumbers::c_BKLM)
237  continue;
238  //if(!hit->inRPC()) continue;
239  if (hit->getSection() != section)
240  continue;
241  if (hit->getSector() != sector)
242  continue;
243  if (hit->getLayer() != layer)
244  continue;
245  ROOT::Math::XYZVector position = hit->getPosition();
246  ROOT::Math::XYZVector distance = extVec - position;
247  //on same track, same sector, same layer, we should believe extHit and KLMHit2d are matched.
248  //let's record the distance to check, should be small
249  m_hdistance->Fill(distance.R());
250  if (distance.R() < 20)
251  matched = true;
252  //m_pointUsed.insert(m);
253  if (matched)
254  break;
255  }
256  if (matched) {
257  m_passYX->Fill(extVec.X(), extVec.Y());
258  m_passYZ->Fill(extVec.Z(), extVec.Y());
259  m_passTrkThephi[layer - 1]->Fill(trkphi, trktheta);
260  m_passThephi[layer - 1]->Fill(phi, theta);
261  m_passMom->Fill(mom);
262  }
263  }//end of loop ext hit
264  }//end of loop tracks
265  m_extTree->Fill();
266 }
267 
269 {
270 }
271 
273 {
274  for (int iL = 0; iL < 15; iL ++) {
275  for (int i = 0; i < m_totalThephi[iL]->GetNbinsX(); i++) {
276  for (int j = 0; j < m_totalThephi[iL]->GetNbinsY(); j++) {
277  float num = m_passThephi[iL]->GetBinContent(i + 1, j + 1);
278  float denom = m_totalThephi[iL]->GetBinContent(i + 1, j + 1);
279  if (num > 0) {
280  m_effiThephi[iL]->SetBinContent(i + 1, j + 1, num / denom);
281  m_effiThephi[iL]->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
282  } else {
283  m_effiThephi[iL]->SetBinContent(i + 1, j + 1, 0);
284  m_effiThephi[iL]->SetBinError(i + 1, j + 1, 0);
285  }
286  }
287  }
288  }
289 
290  for (int iL = 0; iL < 15; iL ++) {
291  for (int i = 0; i < m_totalTrkThephi[iL]->GetNbinsX(); i++) {
292  for (int j = 0; j < m_totalTrkThephi[iL]->GetNbinsY(); j++) {
293  float num = m_passTrkThephi[iL]->GetBinContent(i + 1, j + 1);
294  float denom = m_totalTrkThephi[iL]->GetBinContent(i + 1, j + 1);
295  if (num > 0) {
296  m_effiTrkThephi[iL]->SetBinContent(i + 1, j + 1, num / denom);
297  m_effiTrkThephi[iL]->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
298  } else {
299  m_effiTrkThephi[iL]->SetBinContent(i + 1, j + 1, 0);
300  m_effiTrkThephi[iL]->SetBinError(i + 1, j + 1, 0);
301  }
302  }
303  }
304  }
305 
306 
307  for (int i = 0; i < m_totalYX->GetNbinsX(); i++) {
308  for (int j = 0; j < m_totalYX->GetNbinsY(); j++) {
309  float num = m_passYX->GetBinContent(i + 1, j + 1);
310  float denom = m_totalYX->GetBinContent(i + 1, j + 1);
311  if (num > 0) {
312  m_effiYX->SetBinContent(i + 1, j + 1, num / denom);
313  m_effiYX->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
314  } else {
315  m_effiYX->SetBinContent(i + 1, j + 1, 0);
316  m_effiYX->SetBinError(i + 1, j + 1, 0);
317  }
318 
319  num = m_passYZ->GetBinContent(i + 1, j + 1);
320  denom = m_totalYZ->GetBinContent(i + 1, j + 1);
321  if (num > 0) {
322  m_effiYZ->SetBinContent(i + 1, j + 1, num / denom);
323  m_effiYZ->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
324  } else {
325  m_effiYZ->SetBinContent(i + 1, j + 1, 0);
326  m_effiYZ->SetBinError(i + 1, j + 1, 0);
327  }
328  }
329  }
330 
331  for (int i = 0; i < m_totalMom->GetNbinsX(); i++) {
332  float num = m_passMom->GetBinContent(i + 1);
333  float denom = m_totalMom->GetBinContent(i + 1);
334  if (num > 0) {
335  m_effiMom->SetBinContent(i + 1, num / denom);
336  m_effiMom->SetBinError(i + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
337  } else {
338  m_effiMom->SetBinContent(i + 1, 0);
339  m_effiMom->SetBinError(i + 1, 0);
340  }
341  }
342 
343  m_file->cd();
344  m_hdistance->Write();
345  m_totalYX->Write();
346  m_passYX->Write();
347  m_totalYZ->Write();
348  m_passYZ->Write();
349  m_effiYX->Write();
350  m_effiYZ->Write();
351  m_totalMom->Write();
352  m_passMom->Write();
353  m_effiMom->Write();
354  for (int i = 0; i < 15; i++) {
355  m_totalThephi[i]->Write();
356  m_passThephi[i]->Write();
357  m_effiThephi[i]->Write();
358  m_totalTrkThephi[i]->Write();
359  m_passTrkThephi[i]->Write();
360  m_effiTrkThephi[i]->Write();
361  }
362  m_extTree->Write();
363  m_file->Close();
364 }
365 
TH2F * m_passYZ
Y Z coordinate of extHit with matched bklmHit2d.
TH2F * m_effiYX
associated efficiencies vs glaoble y, x coordinate
float m_exty[200]
keep the global position of validated ExtHit in BKLM
Definition: BKLMAnaModule.h:83
void initialize() override
Initialize at start of job.
float m_extx[200]
keep the global position of validated ExtHit in BKLM
Definition: BKLMAnaModule.h:80
void event() override
This method is called for each event.
void endRun() override
Do any needed actions at the end of a simulation run.
void terminate() override
Terminate at the end of job.
StoreArray< ExtHit > extHits
extHits StoreArray
BKLMAnaModule()
Constructor.
TH2F * m_passThephi[15]
histogram of entries with matched bklmHit2d of extHit vs.theta phi angle of ExtHit
~BKLMAnaModule()
Destructor.
TH2F * m_totalTrkThephi[15]
histogram of total entries of extHit associating to ExtHit of each layer vs. theta phi angle of the t...
TH2F * m_totalYX
Y X coordinate of total entries of extHit that in klm region.
void beginRun() override
Do any needed actions at the start of a simulation run.
TH1F * m_hdistance
bklm GeometryPar
Definition: BKLMAnaModule.h:93
TH1F * m_totalMom
histogram of total entries of extHit vs. track momentum
Definition: BKLMAnaModule.h:96
StoreArray< KLMHit2d > hits2D
hits2D StoreArray
TFile * m_file
output TFile
Definition: BKLMAnaModule.h:62
int m_nExtHit
number of validated ExtHit in BKLM in one event
Definition: BKLMAnaModule.h:77
float m_extz[200]
keep the global position of validated ExtHit in BKLM
Definition: BKLMAnaModule.h:86
TH2F * m_totalThephi[15]
histogram of total entries of extHit associating to ExtHit of each layer vs.theta phi angle of ExtHit
TH2F * m_effiThephi[15]
efficiency associating to ExtHit of each layer vs.theta phi angle of ExtHit
std::string m_filename
name of the output TFile
Definition: BKLMAnaModule.h:65
TH2F * m_effiYZ
associated efficiencies vs glaoble y, z coordinate
TH1F * m_effiMom
efficiency associating to ExtHit vs. track momentum
StoreArray< Track > tracks
tracks StoreArray
TH2F * m_totalYZ
Y Z coordinate of total entries of extHit that in klm region.
TH2F * m_passTrkThephi[15]
histogram of extHit with matched bklmHit2d of each layer vs. theta phi angle of the track
TH2F * m_effiTrkThephi[15]
efficiency associating to ExtHit of each layer vs. theta phi angle of the track
TH2F * m_passYX
Y X coordinate of extHit with matched bklmHit2d.
int m_run
run number
Definition: BKLMAnaModule.h:74
TTree * m_extTree
TTree that holds several variable of interest.
Definition: BKLMAnaModule.h:68
TH1F * m_passMom
histogram of entries with matched bklmHit2d of extHit vs. track momentum
Definition: BKLMAnaModule.h:99
static int getSectorByModule(int module)
Get sector number by module identifier.
static int getLayerByModule(int module)
Get layer number by module identifier.
static int getSectionByModule(int module)
Get section number by module identifier.
static const ChargedStable muon
muon particle
Definition: Const.h:651
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:32
int getCopyID() const
Get detector-element ID of sensitive element within detector.
Definition: ExtHit.h:126
Const::EDetector getDetectorID() const
Get detector ID of this extrapolation hit.
Definition: ExtHit.h:122
ROOT::Math::XYZVector getPosition() const
Get position of this extrapolation hit.
Definition: ExtHit.h:145
KLM 2d hit.
Definition: KLMHit2d.h:33
Class to store the likelihoods from KLM with additional informations related to the extrapolation.
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Class for type safe access to objects that are referred to in relations.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Values of the result of a track fit with a given particle hypothesis.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Class that bundles various TrackFitResults.
Definition: Track.h:25
REG_MODULE(arichBtest)
Register the Module.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.