ATLAS Offline Software
Loading...
Searching...
No Matches
MonitoringFile_MuonTrackPostProcess.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4
5
6
8//
9//Post processing algorithm for Muon Track Physics Monitoring
10//
12
14
15#include <iostream>
16#include <algorithm>
17#include <fstream>
18#include <cmath>
19#include <cstdlib>
20#include <sstream>
21#include <stdio.h>
22#include <map>
23#include <iomanip>
24#include <set>
25
26#include "TH1F.h"
27#include "TH2F.h"
28#include "TFile.h"
29#include "TClass.h"
30#include "TKey.h"
31#include "TMath.h"
32#include "TF1.h"
33#include "TTree.h"
34#include "TBranch.h"
35#include "TList.h"
36#include "TProfile.h"
37#include "TMinuit.h"
38
39const char* const SegStationName[17] = {
40 "BIS", "BIL", "BMS", "BML", "BOS", "BOL", "BEE",
41 "EIS", "EIL", "EMS", "EML", "EOS", "EOL", "EES",
42 "EEL", "CSS", "CSL"
43};//For filling in monitoring plots
44
45double breitgausfun(double* x, double* par) {
46 //Fit parameters:
47 //par[0]=Width (scale) Breit-Wigner
48 //par[1]=Most Probable (MP, location) Breit mean
49 //par[2]=Total area (integral -inf to inf, normalization constant)
50 //par[3]=Width (sigma) of convoluted Gaussian function
51
52 // Numeric constants
53 double invsq2pi = 0.3989422804014; // (2 pi)^(-1/2)
54 // Control constants
55 double np = 100.0; // number of convolution steps
56 double sc = 4; // convolution extends to +-sc Gaussian sigmas
57
58 // Variables
59 double xx;
60 double fland;
61 double sum = 0.0;
62 double xlow, xupp;
63 double step;
64 double i;
65
66 // Range of convolution integral
67 xlow = x[0] - sc * par[3];
68 xupp = x[0] + sc * par[3];
69 step = (xupp - xlow) / np;
70
71 // Convolution integral of Breit and Gaussian by sum
72 for (i = 1.0; i <= np / 2; i++) {
73 xx = xlow + (i - .5) * step;
74 fland = TMath::BreitWigner(xx, par[1], par[0]);
75 sum += fland * TMath::Gaus(x[0], xx, par[3]);
76
77 xx = xupp - (i - .5) * step;
78 fland = TMath::BreitWigner(xx, par[1], par[0]);
79 sum += fland * TMath::Gaus(x[0], xx, par[3]);
80 }
81 return(par[2] * step * sum * invsq2pi / par[3]);
82}
83
84void TwoDto2D_Eff(TH2* Numerator, TH2* Denominator, TH2* Efficiency, bool rebin2d = false) {
85 //the input histograms must have the same dimension!
86 if (Numerator == NULL || Denominator == NULL || Efficiency == NULL) {
87 return;
88 }
89 if (rebin2d) {
90 Numerator->Rebin2D();//here change the default binnning of eta-phi
91 Efficiency->Rebin2D();//here change the default binnning of eta-phi
92 }
93 //then check the dimensions
94 int n_xbins = Numerator->GetNbinsX();
95 if (Denominator->GetNbinsX() != n_xbins || Efficiency->GetNbinsX() != n_xbins) {
96 return;
97 }
98 int n_ybins = Numerator->GetNbinsY();
99 if (Denominator->GetNbinsY() != n_ybins || Efficiency->GetNbinsY() != n_ybins) {
100 return;
101 }
102
103 //after protection
104 for (int bin_itrX = 1; bin_itrX < Efficiency->GetNbinsX() + 1; bin_itrX++) {
105 for (int bin_itrY = 1; bin_itrY < Efficiency->GetNbinsY() + 1; bin_itrY++) {
106 if (Denominator->GetBinContent(bin_itrX, bin_itrY) == 0) continue;
107 Efficiency->SetBinContent(bin_itrX, bin_itrY,
108 Numerator->GetBinContent(bin_itrX,
109 bin_itrY) / Denominator->GetBinContent(bin_itrX, bin_itrY));
110 }
111 }
112 Efficiency->Write("", TObject::kOverwrite);
113 return;
114}
115
116void TwoDto1D_Mean(TH2* h_parent, TH1* h_child, int rebinning = 1) {
117 //the input histograms must have the same dimension!
118 if (h_parent == NULL || h_child == NULL) {
119 return;
120 }
121 if (h_parent->GetNbinsX() != h_child->GetNbinsX()) {
122 return;
123 }
124 //after protection
125 for (int bin_itrX = 1; bin_itrX < h_parent->GetNbinsX() + 1; bin_itrX++) {
126 double parent_event = 0;
127 double parent_sum = 0;
128 for (int bin_itrY = 1; bin_itrY < h_parent->GetNbinsY() + 1; bin_itrY++) {
129 parent_event += h_parent->GetBinContent(bin_itrX, bin_itrY);
130 parent_sum += h_parent->GetBinContent(bin_itrX, bin_itrY) * h_parent->GetYaxis()->GetBinCenter(bin_itrY);
131 }
132 if (parent_event == 0) {
133 continue;
134 }
135 h_child->SetBinContent(bin_itrX, parent_sum / parent_event);
136 }
137 std::string sHistTitle = h_child->GetTitle();
138 h_child->Rebin(rebinning);
139 h_child->SetTitle((sHistTitle + " per event").c_str());
140 h_child->Write("", TObject::kOverwrite);
141 return;
142}
143
144void TwoDto1D_Sum(TH2* h_parent, TH1* h_child, int rebinning = 2) {
145 //the input histograms must have the same dimension!
146 if (h_parent == NULL || h_child == NULL) {
147 return;
148 }
149 if (h_parent->GetNbinsX() != h_child->GetNbinsX()) {
150 return;
151 }
152 //after protection
153 for (int bin_itrX = 1; bin_itrX < h_parent->GetNbinsX() + 1; bin_itrX++) {
154 double parent_sum = 0;
155 for (int bin_itrY = 1; bin_itrY < h_parent->GetNbinsY() + 1; bin_itrY++) {
156 parent_sum += h_parent->GetBinContent(bin_itrX, bin_itrY) * h_parent->GetYaxis()->GetBinCenter(bin_itrY);
157 }
158 if (parent_sum == 0) {
159 continue;
160 }
161 h_child->SetBinContent(bin_itrX, parent_sum);
162 }
163 h_child->Rebin(rebinning);
164 h_child->Write("", TObject::kOverwrite);
165 return;
166}
167
168void SetMassInfo(int iBin, TH1* InputHist, TH1* OutMean, TH1* OutSigma, const std::string& recalg_path) {
169 if (InputHist == NULL || OutMean == NULL || OutSigma == NULL) {
170 return;
171 }
172 if (InputHist->GetMaximum() > 10.0) {
173 if (recalg_path == "Z") {//only for Z
174 TF1* fit1 = new TF1("fit1", breitgausfun, 76, 106, 4);
175 fit1->SetLineColor(kRed);
176 fit1->SetParameter(0, 3.0);//par[0]=Width (scale) Breit-Wigner
177 fit1->SetParameter(1, 91.2);//par[1]=Most Probable (MP, location) Breit mean
178 fit1->SetParameter(2, InputHist->GetEntries());//par[2]=Total area (integral -inf to inf, normalization constant)
179 fit1->SetParameter(3, 1.0);//par[3]=Width (sigma) of convoluted Gaussian function
180 InputHist->Fit("fit1", "q", "", 77, 105);
181 fit1->Draw();
182 OutMean->SetBinContent(iBin, fit1->GetParameter(1));
183 OutMean->SetBinError(iBin, fit1->GetParError(1));
184 OutSigma->SetBinContent(iBin, fit1->GetParameter(0));
185 OutSigma->SetBinError(iBin, fit1->GetParError(0));
186 }
187 if (recalg_path == "Jpsi") {
188 InputHist->Fit("gaus", "q", "", 2.95, 3.25);
189 TF1* fit1 = (TF1*) InputHist->GetFunction("gaus");
190 fit1->Draw();
191 OutMean->SetBinContent(iBin, fit1->GetParameter(1));
192 OutMean->SetBinError(iBin, fit1->GetParError(1));
193 OutSigma->SetBinContent(iBin, fit1->GetParameter(2));
194 OutSigma->SetBinError(iBin, fit1->GetParError(2));
195 }
196 } else {
197 OutMean->SetBinContent(iBin, InputHist->GetMean(1));
198 OutMean->SetBinError(iBin, InputHist->GetMeanError(1));
199 OutSigma->SetBinContent(iBin, InputHist->GetRMS(1));
200 OutSigma->SetBinError(iBin, InputHist->GetRMSError(1));
201 }
202 return;
203}
204
205namespace dqutils {
206 //main function
207 void MonitoringFile::MuonTrackPostProcess(const std::string& inFilename, bool isIncremental) {
208 if (isIncremental) {
209 return;
210 }
211 MonitoringFile::MuonTrack_Main(inFilename, "");
212 MonitoringFile::MuonTrack_Main(inFilename, "NoTrig/");
213 return;
214 }
215
216 //subfunctions
217 void MonitoringFile::MuonTrack_Main(const std::string& inFilename, const std::string& dirname_in) {
218 TString stmp = dirname_in;
219 stmp.ReplaceAll("/", "_");//name clean
220 std::string plotdirname = stmp.Data();//set the plottting dir name
221
222 std::string dirname = "MuonPhysics/" + dirname_in;//give it the full path
223 TFile* f = TFile::Open(inFilename.c_str(), "UPDATE");
224
225 if (f == 0) {
226 std::cerr << "MuonTrackMonitoring(): " << "Input file not opened \n";
227 return;
228 }
229 if (f->GetSize() < 1000.) {
230 std::cerr << "MuonTrackMonitoring(): " << "Input file empty \n";
231 return;
232 }
233 // get run directory name
234 //Seemingly unnecessary lines are necessary
235 TIter nextcd0(gDirectory->GetListOfKeys());
236 TKey* key0 = (TKey*) nextcd0();
237 if (key0 == 0) return;
238
239 TDirectory* dir0 = dynamic_cast<TDirectory*> (key0->ReadObj());
240 if (dir0 == 0) return;
241
242 dir0->cd();
243 std::string runNumber = dir0->GetName();
244 std::string motherDir = runNumber + "/" + dirname;
246
247 //Do the segment part
248 std::string mDir = motherDir + "Segments/";
249 if (!f->cd(mDir.c_str())) return;
250
251 TIter nextcd1(gDirectory->GetListOfKeys());
252 while (TKey* key1 = dynamic_cast<TKey*>(nextcd1())) {
253 //While in the segments
254 std::string recalg_path = key1->GetName();
255 std::string recalg_fullStr = mDir + key1->GetName();
256 TDirectory* dir1 = f->GetDirectory(recalg_fullStr.c_str());
257 if (!dir1) continue;
258 dir1->cd();
259
260 // Divide the efficiency histograms
261 TH2F* h_EffNumerator =
262 (TH2F*) dir1->Get(Form("%sSegments_%s_eff_chamberIndex_perSector_numerator", plotdirname.c_str(),
263 recalg_path.c_str()));
264 TH2F* h_EffDenominator =
265 (TH2F*) dir1->Get(Form("%sSegments_%s_eff_chamberIndex_perSector_denominator", plotdirname.c_str(),
266 recalg_path.c_str()));
267 TH2F* h_Efficiency =
268 (TH2F*) dir1->Get(Form("%sSegments_%s_eff_chamberIndex_perSector", plotdirname.c_str(), recalg_path.c_str()));
269
270 TwoDto2D_Eff(h_EffNumerator, h_EffDenominator, h_Efficiency);
271
272 //add the efficiency for precision
273 for (int i = 0; i < 17; i++) {
274 TH2F* seg_prec_EffNumerator =
275 (TH2F*) dir1->Get(Form("%sSegments_%s_%s_etastation_nPrechit", plotdirname.c_str(), recalg_path.c_str(),
276 SegStationName[i]));
277 TH2F* seg_prec_EffDenominator =
278 (TH2F*) dir1->Get(Form("%sSegments_%s_eff_%s_etastation_nPrechit", plotdirname.c_str(), recalg_path.c_str(),
279 SegStationName[i]));
280 TH2F* seg_prec_Efficiency =
281 (TH2F*) dir1->Get(Form("%sSegments_%s_eff_%s_etastation_nPrechit", plotdirname.c_str(), recalg_path.c_str(),
282 SegStationName[i]));
283
284 TwoDto2D_Eff(seg_prec_EffNumerator, seg_prec_EffDenominator, seg_prec_Efficiency);
285
286 TH2F* seg_trig_EffNumerator =
287 (TH2F*) dir1->Get(Form("%sSegments_%s_%s_etastation_nTrighit", plotdirname.c_str(), recalg_path.c_str(),
288 SegStationName[i]));
289 TH2F* seg_trig_EffDenominator =
290 (TH2F*) dir1->Get(Form("%sSegments_%s_eff_%s_etastation_nTrighit", plotdirname.c_str(), recalg_path.c_str(),
291 SegStationName[i]));
292 TH2F* seg_trig_Efficiency =
293 (TH2F*) dir1->Get(Form("%sSegments_%s_eff_%s_etastation_nTrighit", plotdirname.c_str(), recalg_path.c_str(),
294 SegStationName[i]));
295
296 TwoDto2D_Eff(seg_trig_EffNumerator, seg_trig_EffDenominator, seg_trig_Efficiency);
297 }
298 }//ends different subfolder for segment efficiency
299
300 //Do the muon part
301 std::string mDir_muons = motherDir + "Muons/";
302 if (!f->cd(mDir_muons.c_str())) return;
303
304 TIter nextcd_muons(gDirectory->GetListOfKeys());
305 while (TKey* key1 = dynamic_cast<TKey*>(nextcd_muons())) {
306 //While in the segments
307 std::string recalg_path = key1->GetName();
308 std::string recalg_fullStr = mDir_muons + key1->GetName();
309 TDirectory* dir1 = f->GetDirectory(recalg_fullStr.c_str());
310 if (!dir1) continue;
311 dir1->cd();
312
313 std::string muonqualstr[4] = {
314 "Tight", "Medium", "Loose", "Veryloose"
315 };
316 // Divide the efficiency histograms
317 TH2F* h_EffDenominator =
318 (TH2F*) dir1->Get(Form("%sMuons_%s_Origin_eta_phi", plotdirname.c_str(), recalg_path.c_str()));
319 //m_EffDenominator->Rebin2D();//here change the default binnning of eta-phi! disabled once we are in 64 bins
320 for (int i = 0; i < 4; i++) {
321 TH2F* h_EffNumerator =
322 (TH2F*) dir1->Get(Form("%sMuons_%s_%s_eta_phi", plotdirname.c_str(), recalg_path.c_str(),
323 muonqualstr[i].c_str()));
324 TH2F* h_Efficiency =
325 (TH2F*) dir1->Get(Form("%sMuons_%s_%s_eff", plotdirname.c_str(), recalg_path.c_str(), muonqualstr[i].c_str()));
326 TwoDto2D_Eff(h_EffNumerator, h_EffDenominator, h_Efficiency);//here change the default binnning of eta-phi
327 }
328
329 TH2F* eff_nPrec = (TH2F*) dir1->Get(Form("%sMuons_%s_eff_nPrec", plotdirname.c_str(), recalg_path.c_str()));
330 TH2F* eff_nPhi = (TH2F*) dir1->Get(Form("%sMuons_%s_eff_nPhi", plotdirname.c_str(), recalg_path.c_str()));
331 TH2F* eff_nTrigEta = (TH2F*) dir1->Get(Form("%sMuons_%s_eff_nTrigEta", plotdirname.c_str(), recalg_path.c_str()));
332 TH2F* eff_ndof = (TH2F*) dir1->Get(Form("%sMuons_%s_eff_ndof", plotdirname.c_str(), recalg_path.c_str()));
333 TH2F* eff_chi2 = (TH2F*) dir1->Get(Form("%sMuons_%s_eff_chi2", plotdirname.c_str(), recalg_path.c_str()));
334 TH2F* ID_eff_ndof = (TH2F*) dir1->Get(Form("%sMuons_%s_ID_eff_ndof", plotdirname.c_str(), recalg_path.c_str()));
335 TH2F* ID_eff_chi2 = (TH2F*) dir1->Get(Form("%sMuons_%s_ID_eff_chi2", plotdirname.c_str(), recalg_path.c_str()));
336 TH2F* MS_eff_ndof = (TH2F*) dir1->Get(Form("%sMuons_%s_MS_eff_ndof", plotdirname.c_str(), recalg_path.c_str()));
337 TH2F* MS_eff_chi2 = (TH2F*) dir1->Get(Form("%sMuons_%s_MS_eff_chi2", plotdirname.c_str(), recalg_path.c_str()));
338
339 TH2F* avg_hits_precision_inner =
340 (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_precision_inner", plotdirname.c_str(), recalg_path.c_str()));
341 TH2F* avg_hits_precision_middle =
342 (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_precision_middle", plotdirname.c_str(), recalg_path.c_str()));
343 TH2F* avg_hits_precision_outer =
344 (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_precision_outer", plotdirname.c_str(), recalg_path.c_str()));
345 TH2F* avg_hits_precision_extended =
346 (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_precision_extended", plotdirname.c_str(), recalg_path.c_str()));
347
348 TH2F* avg_hits_trigger_layer1 =
349 (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_trigger_layer1", plotdirname.c_str(), recalg_path.c_str()));
350 TH2F* avg_hits_trigger_layer2 =
351 (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_trigger_layer2", plotdirname.c_str(), recalg_path.c_str()));
352 TH2F* avg_hits_trigger_layer3 =
353 (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_trigger_layer3", plotdirname.c_str(), recalg_path.c_str()));
354 TH2F* avg_hits_trigger_layer4 =
355 (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_trigger_layer4", plotdirname.c_str(), recalg_path.c_str()));
356
357 TH2F* avg_hits_ibl = (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_ibl", plotdirname.c_str(), recalg_path.c_str()));
358 TH2F* avg_hits_pix = (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_pix", plotdirname.c_str(), recalg_path.c_str()));
359 TH2F* avg_hits_sct = (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_sct", plotdirname.c_str(), recalg_path.c_str()));
360 TH2F* avg_hits_trt = (TH2F*) dir1->Get(Form("%sMuons_%s_avg_hits_trt", plotdirname.c_str(), recalg_path.c_str()));
361
362 TH2F* avg_ddpt_idme = (TH2F*) dir1->Get(Form("%sMuons_%s_avg_ddpt_idme", plotdirname.c_str(), recalg_path.c_str()));
363 TH2F* avg_dptsignif = (TH2F*) dir1->Get(Form("%sMuons_%s_avg_dptsignif", plotdirname.c_str(), recalg_path.c_str()));
364
365 TwoDto2D_Eff(eff_nPrec, h_EffDenominator, eff_nPrec);
366 TwoDto2D_Eff(eff_nPhi, h_EffDenominator, eff_nPhi);
367 TwoDto2D_Eff(eff_nTrigEta, h_EffDenominator, eff_nTrigEta);
368 TwoDto2D_Eff(eff_ndof, h_EffDenominator, eff_ndof);
369 TwoDto2D_Eff(eff_chi2, h_EffDenominator, eff_chi2);
370 TwoDto2D_Eff(ID_eff_ndof, h_EffDenominator, ID_eff_ndof);
371 TwoDto2D_Eff(ID_eff_chi2, h_EffDenominator, ID_eff_chi2);
372 TwoDto2D_Eff(MS_eff_ndof, h_EffDenominator, MS_eff_ndof);
373 TwoDto2D_Eff(MS_eff_chi2, h_EffDenominator, MS_eff_chi2);
374
375 TwoDto2D_Eff(avg_hits_precision_inner, h_EffDenominator, avg_hits_precision_inner);
376 TwoDto2D_Eff(avg_hits_precision_middle, h_EffDenominator, avg_hits_precision_middle);
377 TwoDto2D_Eff(avg_hits_precision_outer, h_EffDenominator, avg_hits_precision_outer);
378 TwoDto2D_Eff(avg_hits_precision_extended, h_EffDenominator, avg_hits_precision_extended);
379
380 TwoDto2D_Eff(avg_hits_trigger_layer1, h_EffDenominator, avg_hits_trigger_layer1);
381 TwoDto2D_Eff(avg_hits_trigger_layer2, h_EffDenominator, avg_hits_trigger_layer2);
382 TwoDto2D_Eff(avg_hits_trigger_layer3, h_EffDenominator, avg_hits_trigger_layer3);
383 TwoDto2D_Eff(avg_hits_trigger_layer4, h_EffDenominator, avg_hits_trigger_layer4);
384
385 TwoDto2D_Eff(avg_hits_ibl, h_EffDenominator, avg_hits_ibl);
386 TwoDto2D_Eff(avg_hits_pix, h_EffDenominator, avg_hits_pix);
387 TwoDto2D_Eff(avg_hits_sct, h_EffDenominator, avg_hits_sct);
388 TwoDto2D_Eff(avg_hits_trt, h_EffDenominator, avg_hits_trt);
389
390 TwoDto2D_Eff(avg_ddpt_idme, h_EffDenominator, avg_ddpt_idme);
391 TwoDto2D_Eff(avg_dptsignif, h_EffDenominator, avg_dptsignif);
392 }//ends different subfolder for muon efficiency
393
394
395 //Do the luminoisty part
396 std::string mDir_lb = motherDir + "Overview/";
397 if (!f->cd(mDir_lb.c_str())) return;
398
399 TIter nextcd_lb(gDirectory->GetListOfKeys());
400 while (TKey* key1 = dynamic_cast<TKey*>(nextcd_lb())) {
401 //While in the segments
402 std::string recalg_path = key1->GetName();
403 std::string recalg_fullStr = mDir_lb + key1->GetName();
404 TDirectory* dir1 = f->GetDirectory(recalg_fullStr.c_str());
405 if (!dir1) continue;
406 dir1->cd();
407
408 std::string montype[3] = {
409 "Segment", "MuonTrack", "Muon"
410 };
411 // Divide the efficiency histograms
412 for (int i = 0; i < 3; i++) {
413 TH2F* h_parent_lb =
414 (TH2F*) dir1->Get(Form("%sOverview_%s_n%s_LB_2D", plotdirname.c_str(), recalg_path.c_str(), montype[i].c_str()));
415 TH1F* h_child_lb =
416 (TH1F*) dir1->Get(Form("%sOverview_%s_n%s_LB", plotdirname.c_str(), recalg_path.c_str(), montype[i].c_str()));
417 //TH2F* h_parent_inst = (TH2F*)dir1->Get(Form("%sOverview_%s_n%s_Inst_2D", plotdirname.c_str(),
418 // recalg_path.c_str(), montype[i].c_str()));
419 //TH1F* h_child_inst = (TH1F*)dir1->Get(Form("%sOverview_%s_n%s_Inst", plotdirname.c_str(), recalg_path.c_str(),
420 // montype[i].c_str()));
421 //TH2F* h_parent_intlumi = (TH2F*)dir1->Get(Form("%sOverview_%s_n%s_IntLumi_2D", plotdirname.c_str(),
422 // recalg_path.c_str(), montype[i].c_str()));
423 //TH1F* h_child_intlumi = (TH1F*)dir1->Get(Form("%sOverview_%s_n%s_IntLumi", plotdirname.c_str(),
424 // recalg_path.c_str(), montype[i].c_str()));
425 TwoDto1D_Mean(h_parent_lb, h_child_lb);
426 //TwoDto1D_Mean(h_parent_inst, h_child_inst);
427 //TwoDto1D_Mean(h_parent_intlumi, h_child_intlumi);
428 }
429 std::string resonance[2] = {
430 "Z", "Jpsi"
431 };
432 // Divide the efficiency histograms
433 for (int i = 0; i < 2; i++) {
434 TH2F* h_parent_lb =
435 (TH2F*) dir1->Get(Form("%sOverview_%s_n%s_LB_2D", plotdirname.c_str(), recalg_path.c_str(),
436 resonance[i].c_str()));
437 TH1F* h_child_lb =
438 (TH1F*) dir1->Get(Form("%sOverview_%s_n%s_LB", plotdirname.c_str(), recalg_path.c_str(), resonance[i].c_str()));
439 //TH2F* h_parent_inst = (TH2F*)dir1->Get(Form("%sOverview_%s_n%s_Inst_2D", plotdirname.c_str(),
440 // recalg_path.c_str(), resonance[i].c_str()));
441 //TH1F* h_child_inst = (TH1F*)dir1->Get(Form("%sOverview_%s_n%s_Inst", plotdirname.c_str(), recalg_path.c_str(),
442 // resonance[i].c_str()));
443 //TH2F* h_parent_intlumi = (TH2F*)dir1->Get(Form("%sOverview_%s_n%s_IntLumi_2D", plotdirname.c_str(),
444 // recalg_path.c_str(), resonance[i].c_str()));
445 //TH1F* h_child_intlumi = (TH1F*)dir1->Get(Form("%sOverview_%s_n%s_IntLumi", plotdirname.c_str(),
446 // recalg_path.c_str(), resonance[i].c_str()));
447 TwoDto1D_Sum(h_parent_lb, h_child_lb);
448 //TwoDto1D_Sum(h_parent_inst, h_child_inst);
449 //TwoDto1D_Sum(h_parent_intlumi, h_child_intlumi);
450 }
451 }//ends different subfolder for luminosity
452
453
454 //Do the muon part; only for the main directory!
455 if (dirname.find("NoTrig") == std::string::npos) {
456 //std::cout << "get to trackphys " << std::endl;
457 std::string mDir_phys = motherDir + "MuonTrkPhys/";
458 if (!f->cd(mDir_phys.c_str())) return;
459
460 TIter nextcd_phys(gDirectory->GetListOfKeys());
461 while (TKey* key1 = dynamic_cast<TKey*>(nextcd_phys())) {
462 //While in the segments
463 std::string recalg_path = key1->GetName();
464 std::string recalg_fullStr = mDir_phys + key1->GetName();
465 TDirectory* dir1 = f->GetDirectory(recalg_fullStr.c_str());
466 if (!dir1) continue;
467 dir1->cd();
468
469 TH1* h_Mass_Mean = (TH1F*) dir1->Get(Form("m_%s_M_Mean", recalg_path.c_str()));
470 TH1* h_Mass_Sigma = (TH1F*) dir1->Get(Form("m_%s_M_Sigma", recalg_path.c_str()));
471 // Get each of the mass histograms
472 std::string det_region[4] = {
473 "EC", "BC", "BA", "EA"
474 };
475 for (int i = 0; i < 4; i++) {
476 for (int j = 0; j < 4; j++) {
477 TH1* h_Mass_region = (TH1F*) dir1->Get(Form("m_%s_M_%s_%s", recalg_path.c_str(),
478 det_region[i].c_str(), det_region[j].c_str()));
479 //std::cout << " bin " << i * 4 + (j + 1) << " content " << det_region[i] << " " << det_region[j] <<
480 // std::endl;
481 SetMassInfo(i * 4 + (j + 1), h_Mass_region, h_Mass_Mean, h_Mass_Sigma, recalg_path);
482 if (h_Mass_region != NULL) h_Mass_region->Write("", TObject::kOverwrite);
483 }
484 }
485 if (h_Mass_Mean != NULL) h_Mass_Mean->Write("", TObject::kOverwrite);
486 if (h_Mass_Sigma != NULL) h_Mass_Sigma->Write("", TObject::kOverwrite);
487 }
488 }//ends different subfolder for muon efficiency
489
490 f->Close();
491 delete f;
492 return;
493 } //ends function
494} //namespace
static Double_t sc
void SetMassInfo(int iBin, TH1 *InputHist, TH1 *OutMean, TH1 *OutSigma, const std::string &recalg_path)
void TwoDto1D_Mean(TH2 *h_parent, TH1 *h_child, int rebinning=1)
void TwoDto2D_Eff(TH2 *Numerator, TH2 *Denominator, TH2 *Efficiency, bool rebin2d=false)
double breitgausfun(double *x, double *par)
void TwoDto1D_Sum(TH2 *h_parent, TH1 *h_child, int rebinning=2)
const char *const SegStationName[17]
#define x
static void MuonTrackPostProcess(const std::string &inFileName, bool isIncremental=false)
static void MuonTrack_Main(const std::string &inFileName, const std::string &dirname_in)
a structure to hold an efficiency together with a variable number of uncertainties
std::string dirname(std::string name)
Definition utils.cxx:200