10 #include <klm/bklm/modules/bklmAna/BKLMAnaModule.h>
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>
20 #include <CLHEP/Units/SystemOfUnits.h>
24 using namespace CLHEP;
51 for (
int i = 0; i < 200; ++i) {
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;
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"));
68 BKLMAnaModule::~BKLMAnaModule()
72 void BKLMAnaModule::initialize()
79 m_file =
new TFile(m_filename.c_str(),
"RECREATE");
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");
92 float thetaMin = 35.0;
93 float thetaMax = 130.0;
95 for (
int iL = 0; iL < 15 ; iL ++) {
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");
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");
130 m_hdistance =
new TH1F(
"m_hdistance",
" distance between mathced extHit and bklmHit2d ", 100, 0, 30);
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)");
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)");
151 void BKLMAnaModule::beginRun()
155 void BKLMAnaModule::event()
159 unsigned long runNumber = eventMetaData->getRun();
174 for (
int t = 0; t < extHits.getEntries(); t++) {
175 ExtHit* exthit = extHits[t];
189 for (
int k = 0; k < tracks.getEntries(); k++) {
190 Track* track = tracks[k];
197 double trkphi = p4.Vect().Phi() * 180.0 / CLHEP::pi;
198 double trktheta = p4.Vect().Theta() * 180.0 / CLHEP::pi;
200 trkphi = trkphi + 360.0;
203 for (
unsigned int t = 0; t < relatedExtHit.size(); t++) {
204 ExtHit* exthit = relatedExtHit[t];
208 int section = BKLMElementNumbers::getSectionByModule(module);
209 int sector = BKLMElementNumbers::getSectorByModule(module);
210 int layer = BKLMElementNumbers::getLayerByModule(module);
211 bool crossed =
false;
215 int extPattern = muid->getExtLayerPattern();
216 if ((extPattern & (1 << (layer - 1))) != 0)
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;
230 m_totalThephi[layer - 1]->Fill(phi, theta);
231 m_totalTrkThephi[layer - 1]->Fill(trkphi, trktheta);
232 m_totalMom->Fill(mom);
236 for (
int mHit = 0; mHit < hits2D.getEntries(); mHit++) {
239 if (hit->getSection() != section)
241 if (hit->getSector() != sector)
243 if (hit->getLayer() != layer)
245 TVector3 position = hit->getGlobalPosition();
246 TVector3 distance = extVec - position;
249 m_hdistance->Fill(distance.Mag());
250 if (distance.Mag() < 20)
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);
268 void BKLMAnaModule::endRun()
272 void BKLMAnaModule::terminate()
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);
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)));
283 m_effiThephi[iL]->SetBinContent(i + 1, j + 1, 0);
284 m_effiThephi[iL]->SetBinError(i + 1, j + 1, 0);
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);
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)));
299 m_effiTrkThephi[iL]->SetBinContent(i + 1, j + 1, 0);
300 m_effiTrkThephi[iL]->SetBinError(i + 1, j + 1, 0);
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);
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)));
315 m_effiYX->SetBinContent(i + 1, j + 1, 0);
316 m_effiYX->SetBinError(i + 1, j + 1, 0);
319 num = m_passYZ->GetBinContent(i + 1, j + 1);
320 denom = m_totalYZ->GetBinContent(i + 1, j + 1);
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)));
325 m_effiYZ->SetBinContent(i + 1, j + 1, 0);
326 m_effiYZ->SetBinError(i + 1, j + 1, 0);
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);
335 m_effiMom->SetBinContent(i + 1, num / denom);
336 m_effiMom->SetBinError(i + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
338 m_effiMom->SetBinContent(i + 1, 0);
339 m_effiMom->SetBinError(i + 1, 0);
344 m_hdistance->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();
A class defining a module that perform efficiencies studies on bklm.
Store one BKLM strip hit as a ROOT object.
static const ChargedStable muon
muon particle
Store one Ext hit as a ROOT object.
TVector3 getMomentum() const
Get momentum at this extrapolation hit.
int getCopyID() const
Get detector-element ID of sensitive element within detector.
TVector3 getPosition() const
Get position of this extrapolation hit.
Const::EDetector getDetectorID() const
Get detector ID of this extrapolation hit.
Class to store the likelihoods from KLM with additional informations related to the extrapolation.
Class for type safe access to objects that are referred to in relations.
Type-safe access to single objects in the data store.
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.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.