BmnRoot
Loading...
Searching...
No Matches
BmnEfficiencyTools.cxx
Go to the documentation of this file.
2
4#include "BmnEfficiency.h"
5#include "TFile.h"
6
8 : fPeriod(7)
9{
10 fRunAna->Print();
11}
12
14 : fAna(nullptr)
15 , fInnTracker(nullptr)
16 , fPeriod(7)
17 , fRunTrigInfo(nullptr)
18 , fDstName(dst)
19 , fGeomFile("")
20 , fOutFile("")
21 , isMc(kTRUE)
22 , isL1(kFALSE)
23 , effGem(nullptr)
24 , effSilicon(nullptr)
25{
26 // Inner tracker config ...
27 fInnTracker = new BmnInnerTrackerGeometryDraw();
28
29 fAna = new FairRunAna();
30}
31
33 : fAna(nullptr)
34 , fInnTracker(nullptr)
35 , fPeriod(7)
36 , startRun(3589)
37 , finishRun(5186)
38 , fRunTrigInfo(nullptr)
39 , fDataPath("")
40 , fDstName("")
41 , fGeomFile("")
42 , fOutFile("")
43 , isMc(kFALSE)
44 , isL1(kFALSE)
45 , effGem(nullptr)
46 , effSilicon(nullptr)
47{
48 // Setting zero content vectors for beams, targets and triggers ...
49 // To be checked when processing a file list created manually ...
50 fBeams.resize(0);
51 fTargets.resize(0);
52 fTriggers.resize(0);
53
54 fManualList.resize(0);
55
56 // Initializing trigg. info ...
57 fRunTrigInfo = new BmnDataTriggerInfo();
58
59 // Inner tracker config ...
60 fInnTracker = new BmnInnerTrackerGeometryDraw();
61
62 fAna = new FairRunAna();
63}
64
66{
67
68 // Setting inn. tracker geometry ...
69 if (fGeomFile.IsNull())
70 Fatal("BmnEfficiencyTools::Process()", "No geometry file passed!!!");
71 else
72 TGeoManager::Import(fGeomFile.Data());
73
74 // Creating a list of dst files corresponding to the criteria established ...
75
76 vector<Int_t> selectedRuns;
77 vector<TString> selectedDst;
78
79 if (!isMc && !fManualList.size()) {
80 // Loop over all dst files recorded ...
81 for (Int_t iDst = startRun; iDst < finishRun; iDst++) {
82
83 // Getting curr. trigger ...
84 TString trigger = fRunTrigInfo->GetTrigger(iDst);
85 if (trigger.IsNull())
86 continue;
87
88 // Getting curr. target ...
89 UniRun* pCurrentRun = UniRun::GetRun(fPeriod, iDst);
90 if (!pCurrentRun)
91 continue;
92
93 TString targ = pCurrentRun->GetTargetParticle();
94 if (targ == "")
95 continue;
96
97 // Getting beam ...
98 TString beam = pCurrentRun->GetBeamParticle();
99
100 auto triggerSearch = find(fTriggers.begin(), fTriggers.end(), trigger);
101 auto targetSearch = find(fTargets.begin(), fTargets.end(), *targ);
102 auto beamSearch = find(fBeams.begin(), fBeams.end(), beam);
103
104 if (targetSearch != fTargets.end() && triggerSearch != fTriggers.end() && beamSearch != fBeams.end())
105 selectedRuns.push_back(iDst);
106 }
107
108 for (auto run : selectedRuns) {
109 TString dstFile =
110 fDataPath + "/"
111 + ((!isL1) ? TString::Format("bmndst_%d.root", run) : (TString::Format("dst-cbm2bmn-%d.root", run)));
112 selectedDst.push_back(dstFile);
113 }
114 } else if (!isMc && fManualList.size())
115 for (auto it : fManualList)
116 selectedDst.push_back(fDataPath + "/" + TString::Format("dst-cbm2bmn-%d.root", it));
117 else
118 selectedDst.push_back(fDstName);
119
120 const Int_t nStatsGem = fInnTracker->GetGemGeometry()->GetNStations();
121 const Int_t nStatsSil = fInnTracker->GetSiliconGeometry()->GetNStations();
122
123 // Preparing containers to store 2d-efficiency ...
124 effGem = new TClonesArray("EffStore2D");
125 for (Int_t iStat = 0; iStat < nStatsGem; iStat++)
126 new ((*effGem)[effGem->GetEntriesFast()]) EffStore2D("GEM", iStat, fYRangesGem);
127
128 effSilicon = new TClonesArray("EffStore2D");
129 for (Int_t iStat = 0; iStat < nStatsSil; iStat++)
130 new ((*effSilicon)[effSilicon->GetEntriesFast()]) EffStore2D("SILICON", iStat, fYRangesSilicon);
131
132 for (auto it = selectedDst.begin(); it != selectedDst.end(); it++) {
133 BmnEfficiency* eff = new BmnEfficiency(fAna, fInnTracker, *it);
134 eff->Efficiency(effGem, effSilicon, fYRangesGem, fYRangesSilicon);
135 delete eff;
136 }
137
138 TFile* fOut = nullptr;
139
140 if (fOutFile.IsNull())
141 fOut = new TFile("outFile.root", "recreate");
142 else
143 fOut = new TFile(Form("%s", fOutFile.Data()), "recreate");
144
145 // GEM ...
146 const Int_t nQp = 2;
147
148 for (Int_t iEff = 0; iEff < effGem->GetEntriesFast(); iEff++) {
149 EffStore2D* epsilon = (EffStore2D*)effGem->UncheckedAt(iEff);
150
151 Int_t stat = epsilon->Station();
152
153 for (Int_t iGemStat = 0; iGemStat < nStatsGem; iGemStat++) {
154 if (iGemStat != stat)
155 continue;
156
157 for (Int_t iQp = 0; iQp < nQp; iQp++) {
158 epsilon->efficiency(iQp)->Write();
159
160 for (size_t iRange = 0; iRange < fYRangesGem.find(iGemStat)->second.size(); iRange++)
161 epsilon->efficiency4range(iRange, iQp)->Write();
162 }
163 }
164 }
165
166 // SILICON ...
167 for (Int_t iEff = 0; iEff < effSilicon->GetEntriesFast(); iEff++) {
168 EffStore2D* epsilon = (EffStore2D*)effSilicon->UncheckedAt(iEff);
169
170 Int_t stat = epsilon->Station();
171
172 for (Int_t iSilStat = 0; iSilStat < nStatsSil; iSilStat++) {
173
174 if (iSilStat != stat)
175 continue;
176
177 for (Int_t iQp = 0; iQp < nQp; iQp++) {
178 epsilon->efficiency(iQp, "", "SILICON")->Write();
179 TH2F* hPassed = (TH2F*)epsilon->efficiency(iQp, "xP", "SILICON")->GetPassedHistogram();
180 if (!iSilStat) {
181 epsilon->efficiency(iQp, "xP", "SILICON")->Write();
182 hPassed->Write("#varepsilon (x, P)");
183 }
184
185 TString tmp = !iQp ? ">" : "<";
186 if (!iSilStat) {
187 TH1D* hMomPassed = hPassed->ProjectionY(Form("Momentum (passed), Q_{p} %s 0", tmp.Data()), 1,
188 hPassed->GetNbinsX());
189 DoNormalization(hMomPassed);
190 hMomPassed->Write();
191 }
192
193 for (size_t iRange = 0; iRange < fYRangesSilicon.find(iSilStat)->second.size(); iRange++)
194 epsilon->efficiency4range(iRange, iQp, "SILICON")->Write();
195 }
196 }
197 }
198
199 delete fOut;
200}
201
203{
204 if (h->GetEntries() == 0)
205 return;
206
207 // Collecting all bin contents ...
208 Double_t contentAll = 0.;
209 for (Int_t iBin = 1; iBin < h->GetNbinsX() + 1; iBin++)
210 contentAll += h->GetBinContent(iBin);
211 // Normalizing histo ...
212 for (Int_t iBin = 1; iBin < h->GetNbinsX() + 1; iBin++) {
213 h->SetBinContent(iBin, h->GetBinContent(iBin) / contentAll);
214 h->SetBinError(iBin, 0.);
215 }
216 h->GetYaxis()->SetRangeUser(0., 1.1 * h->GetMaximum());
217}
void Efficiency(TClonesArray *, TClonesArray *effSilicon, map< Int_t, vector< pair< Double_t, Double_t > > >, map< Int_t, vector< pair< Double_t, Double_t > > >)
BmnGemStripStationSet * GetGemGeometry()
BmnSiliconStationSet * GetSiliconGeometry()
TEfficiency * efficiency(Int_t qpBin=0, TString hType="", TString det="GEM")
TEfficiency * efficiency4range(Int_t range, Int_t qpBin=0, TString det="GEM")
Int_t Station()
static UniRun * GetRun(int period_number, int run_number)
get run from the database
Definition UniRun.cxx:177
TString GetBeamParticle()
get beam particle of the current run
Definition UniRun.h:113
TString GetTargetParticle()
get target particle of the current run
Definition UniRun.h:115
-clang-format