Belle II Software  release-06-01-15
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 /* Belle 2 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 std;
23 using namespace Belle2;
24 using namespace CLHEP;
25 
26 //-----------------------------------------------------------------
27 // Register the Module
28 //-----------------------------------------------------------------
29 REG_MODULE(BKLMAna)
30 
31 //-----------------------------------------------------------------
32 // Implementation
33 //-----------------------------------------------------------------
34 
36  m_file(nullptr),
37  m_extTree(nullptr),
38  m_run(0),
39  m_nExtHit(0),
40  m_hdistance(nullptr),
41  m_totalMom(nullptr),
42  m_passMom(nullptr),
43  m_effiMom(nullptr),
44  m_totalYX(nullptr),
45  m_passYX(nullptr),
46  m_effiYX(nullptr),
47  m_totalYZ(nullptr),
48  m_passYZ(nullptr),
49  m_effiYZ(nullptr)
50 {
51  for (int i = 0; i < 200; ++i) {
52  m_extx[i] = 0.0;
53  m_exty[i] = 0.0;
54  m_extz[i] = 0.0;
55  }
56  for (int i = 0; i < 15; ++i) {
57  m_totalThephi[i] = nullptr;
58  m_passThephi[i] = nullptr;
59  m_effiThephi[i] = nullptr;
60  m_totalTrkThephi[i] = nullptr;
61  m_passTrkThephi[i] = nullptr;
62  m_effiTrkThephi[i] = nullptr;
63  }
64  setDescription("analyze bklm efficiency associated to CDC, check performance of bklm et al.");
65  addParam("filename", m_filename, "Output root filename", string("bklmana.root"));
66 }
67 
68 BKLMAnaModule::~BKLMAnaModule()
69 {
70 }
71 
72 void BKLMAnaModule::initialize()
73 {
74  hits2D.isRequired();
75  extHits.isRequired();
76  tracks.isRequired();
77 
78 
79  m_file = new TFile(m_filename.c_str(), "RECREATE");
80 
81  m_extTree = new TTree("exthit", "ext hit");
82  m_extTree->Branch("run", &m_run, "m_run/I");
83  m_extTree->Branch("nExtHit", &m_nExtHit, "m_nExtHit/I");
84  m_extTree->Branch("x", &m_extx, "m_extx[m_nExtHit]/F");
85  m_extTree->Branch("y", &m_exty, "m_exty[m_nExtHit]/F");
86  m_extTree->Branch("z", &m_extz, "m_extz[m_nExtHit]/F");
87 
88  float phiBins = 36;
89  float phiMin = 0.0;
90  float phiMax = 360.0;
91  float thetaBins = 19;
92  float thetaMin = 35.0;
93  float thetaMax = 130.0;
94  char hname[50];
95  for (int iL = 0; iL < 15 ; iL ++) {
96  //based on theta phi of ExtHit position
97  sprintf(hname, "denominator_Layer%i", iL + 1);
98  m_totalThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
99  m_totalThephi[iL]->GetXaxis()->SetTitle("phi");
100  m_totalThephi[iL]->GetYaxis()->SetTitle("theta");
101  sprintf(hname, "numerator_Layer%i", iL + 1);
102  m_passThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
103  sprintf(hname, "effi_Layer%i", iL + 1);
104  m_passThephi[iL]->GetXaxis()->SetTitle("phi");
105  m_passThephi[iL]->GetYaxis()->SetTitle("theta");
106  m_effiThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
107  m_effiThephi[iL]->GetXaxis()->SetTitle("phi");
108  m_effiThephi[iL]->GetYaxis()->SetTitle("theta");
109  //based on theta phi of trk
110  sprintf(hname, "Denominator_Layer%i", iL + 1);
111  m_totalTrkThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
112  m_totalTrkThephi[iL]->GetXaxis()->SetTitle("trk.phi");
113  m_totalTrkThephi[iL]->GetYaxis()->SetTitle("trk.theta");
114  sprintf(hname, "Numerator_Layer%i", iL + 1);
115  m_passTrkThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
116  sprintf(hname, "Effi_Layer%i", iL + 1);
117  m_passTrkThephi[iL]->GetXaxis()->SetTitle("trk.phi");
118  m_passTrkThephi[iL]->GetYaxis()->SetTitle("trk.theta");
119  m_effiTrkThephi[iL] = new TH2F(hname, hname, phiBins, phiMin, phiMax, thetaBins, thetaMin, thetaMax);
120  m_effiTrkThephi[iL]->GetXaxis()->SetTitle("trk.phi");
121  m_effiTrkThephi[iL]->GetYaxis()->SetTitle("trk.theta");
122  }
123  float gmin = -350;
124  float gmax = 350;
125  int gNbin = 150;
126  float mommin = 0.;
127  float mommax = 15;
128  int mNbin = 30;
129 
130  m_hdistance = new TH1F("m_hdistance", " distance between mathced extHit and bklmHit2d ", 100, 0, 30);
131 
132  m_totalMom = new TH1F("m_totalMom", " denominator vs. p", mNbin, mommin, mommax);
133  m_passMom = new TH1F("m_passMom", " numerator vs. p", mNbin, mommin, mommax);
134  m_effiMom = new TH1F("m_effiMom", " effi. vs. p", mNbin, mommin, mommax);
135  m_effiMom->GetXaxis()->SetTitle("p (GeV)");
136  m_effiMom->GetYaxis()->SetTitle("Efficiency (GeV)");
137 
138  m_totalYX = new TH2F("m_totalYX", " denominator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
139  m_passYX = new TH2F("m_passYX", " numerator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
140  m_totalYZ = new TH2F("m_totalYZ", " denominator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
141  m_passYZ = new TH2F("m_passYZ", " numerator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
142  m_effiYX = new TH2F("m_effiYX", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
143  m_effiYZ = new TH2F("m_effiYZ", " effi. Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
144  m_effiYX->GetXaxis()->SetTitle("x (cm)");
145  m_effiYX->GetYaxis()->SetTitle("y (cm)");
146  m_effiYZ->GetXaxis()->SetTitle("z (cm)");
147  m_effiYZ->GetYaxis()->SetTitle("y (cm)");
148 
149 }
150 
151 void BKLMAnaModule::beginRun()
152 {
153 }
154 
155 void BKLMAnaModule::event()
156 {
157  StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
158  //unsigned long eventNumber = eventMetaData->getEvent();
159  unsigned long runNumber = eventMetaData->getRun();
160  //unsigned long expNumber = eventMetaData->getExperiment();
161  //StoreArray<RecoTrack> recoTracks;
162  //numRecoTrk = recoTracks.getEntries();
163  //StoreArray<BKLMTrack> bklmtracks;
164  //StoreArray<TrackFitResult> trackFitResults;
165  //StoreArray<MuidHit> muidHits;
166  //StoreArray<Muid> muids;
167  //StoreArray<MCParticle> mcParticles;
168 
169  //set<int> m_pointUsed;
170  //m_pointUsed.clear();
171  //check extHits
172  //all ExtHit in bklm scope in each event, should not be many
173  int nExtHit = 0;
174  for (int t = 0; t < extHits.getEntries(); t++) {
175  ExtHit* exthit = extHits[t];
176  if (exthit->getDetectorID() != Const::EDetector::BKLM)
177  continue;
178  m_extx[nExtHit] = exthit->getPosition()[0];
179  m_exty[nExtHit] = exthit->getPosition()[1];
180  m_extz[nExtHit] = exthit->getPosition()[2];
181  nExtHit++;
182  if (nExtHit > 199)
183  break;
184  }
185  m_run = runNumber;
186  m_nExtHit = nExtHit;
187 
188 //the second way, require muid
189  for (int k = 0; k < tracks.getEntries(); k++) {
190  Track* track = tracks[k];
191  // load the muon fit hypothesis or the hypothesis which is the clostes in mass to a muon
192  // the tracking will not always fit a muon hypothesis
193  const TrackFitResult* fitres = track->getTrackFitResultWithClosestMass(Belle2::Const::muon);
194  double mom = fitres->getMomentum().Mag();
195  //double pt = fitres->getTransverseMomentum();
196  TLorentzVector p4 = fitres->get4Momentum();
197  double trkphi = p4.Vect().Phi() * 180.0 / CLHEP::pi;
198  double trktheta = p4.Vect().Theta() * 180.0 / CLHEP::pi;
199  if (trkphi < 0)
200  trkphi = trkphi + 360.0;
201  RelationVector<BKLMHit2d> relatedHit2D = track->getRelationsTo<BKLMHit2d>();
202  RelationVector<ExtHit> relatedExtHit = track->getRelationsTo<ExtHit>();
203  for (unsigned int t = 0; t < relatedExtHit.size(); t++) {
204  ExtHit* exthit = relatedExtHit[t];
205  if (exthit->getDetectorID() != Const::EDetector::BKLM)
206  continue;
207  int module = exthit->getCopyID();
208  int section = BKLMElementNumbers::getSectionByModule(module);
209  int sector = BKLMElementNumbers::getSectorByModule(module);
210  int layer = BKLMElementNumbers::getLayerByModule(module);
211  bool crossed = false; // should be only once ?
212  KLMMuidLikelihood* muid = track->getRelatedTo<KLMMuidLikelihood>();
213  if (!muid)
214  continue;
215  int extPattern = muid->getExtLayerPattern();
216  if ((extPattern & (1 << (layer - 1))) != 0)
217  crossed = true;
218  if (!crossed)
219  continue;
220 
221  TVector3 extMom = exthit->getMomentum();
222  TVector3 extVec = exthit->getPosition();
223  bool matched = false;
224  m_totalYX->Fill(extVec[0], extVec[1]);
225  m_totalYZ->Fill(extVec[2], extVec[1]);
226  float phi = extVec.Phi() * 180.0 / CLHEP::pi;
227  float theta = extVec.Theta() * 180.0 / CLHEP::pi;
228  if (phi < 0)
229  phi = phi + 360.0;
230  m_totalThephi[layer - 1]->Fill(phi, theta);
231  m_totalTrkThephi[layer - 1]->Fill(trkphi, trktheta);
232  m_totalMom->Fill(mom);
233  //look for mateched BKLM2dHit
234  //for (unsigned int mHit = 0; mHit < relatedHit2D.size(); mHit++) {
235  // BKLMHit2d* hit = relatedHit2D[mHit];
236  for (int mHit = 0; mHit < hits2D.getEntries(); mHit++) {
237  BKLMHit2d* hit = hits2D[mHit];
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  TVector3 position = hit->getGlobalPosition();
246  TVector3 distance = extVec - position;
247  //on same track, same sector, same layer, we should believe extHit and BKLMHit2d are matched.
248  //let's record the distance to check, should be small
249  m_hdistance->Fill(distance.Mag());
250  if (distance.Mag() < 20)
251  matched = true;
252  //m_pointUsed.insert(m);
253  if (matched)
254  break;
255  }
256  if (matched) {
257  m_passYX->Fill(extVec[0], extVec[1]);
258  m_passYZ->Fill(extVec[2], extVec[1]);
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 
268 void BKLMAnaModule::endRun()
269 {
270 }
271 
272 void BKLMAnaModule::terminate()
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 
A class defining a module that perform efficiencies studies on bklm.
Definition: BKLMAnaModule.h:34
Store one BKLM strip hit as a ROOT object.
Definition: BKLMHit2d.h:32
static const ChargedStable muon
muon particle
Definition: Const.h:541
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:30
TVector3 getMomentum() const
Get momentum at this extrapolation hit.
Definition: ExtHit.h:147
int getCopyID() const
Get detector-element ID of sensitive element within detector.
Definition: ExtHit.h:124
TVector3 getPosition() const
Get position of this extrapolation hit.
Definition: ExtHit.h:143
Const::EDetector getDetectorID() const
Get detector ID of this extrapolation hit.
Definition: ExtHit.h:120
Class to store the likelihoods from KLM with additional informations related to the extrapolation.
Base class for Modules.
Definition: Module.h:72
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:95
Values of the result of a track fit with a given particle hypothesis.
TLorentzVector get4Momentum() const
Getter for the 4Momentum at the closest approach of the track in the r/phi projection.
TVector3 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
#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.