1 #if !defined(__CINT__) || defined(__MAKECINT__)
9 #include <TGraphErrors.h>
12 #include <TTreeFormula.h>
13 #include <TEventList.h>
14 #include <TObjString.h>
20 #include <TInterpreter.h>
37 TH1F*
ZoomFromTree(TH1F* h, TTree* atree, Int_t n,
const char* aVar, UShort_t aScaleFactor=1);
38 Double_t
GetTreeMinimum(TTree* aTree, Int_t n,
const char* columname);
39 Double_t
GetTreeMaximum(TTree* aTree, Int_t n,
const char* columname);
51 if(newMax) {nMax = newMax; N=0;
return 0;}
59 int PlotEMCALQATrendingTree(
const char* filename=
"trending.root",Bool_t SavePlots=0, TString expr =
"", TString fTrigger=
"")
62 QAPATH = TString(gSystem->Getenv(
"QAPATH"));
66 TFile* f = TFile::Open(filename);
68 TTree* tree = (TTree*)f->Get(
"trending");
69 if (! tree) {Error(
"PlotEMCALQATrendingTree",
"No Tree found!");
return -1;}
70 TFile* fout =
new TFile(Form(
"%s/trendingPlots.root",
QAPATH.Data()),
"RECREATE");
72 TList* TriggersList =
new TList();
75 tree->Draw(
"fTrigger",
"",
"goff");
77 for(Int_t i = 0 ; i < tree->GetSelectedRows() ; i++){
79 obj = tree->GetVar1()->PrintValue(0);
80 if(! TriggersList->FindObject(obj)) {TriggersList->Add(
new TObjString(obj));}
85 if(!fTrigger.Contains(
"QA")) {fTrigger =
"CaloQA_" + fTrigger;}
86 TriggersList->Add(
new TObjString(fTrigger.Data()));
88 TIter next(TriggersList);
90 while ((obj1 = next()))
104 TCut trig = Form(
"fTrigger==\"%s\"",Trig);
105 TCut NotZero = TCut(
"Nevent>0.");
108 if (Expr.Contains(
".C"))
110 Info(
"PlotEMCALQATrendingTree",Form(
"Additional selections from %s: ", Expr.Data()));
111 gInterpreter->ExecuteMacro(Expr.Data());
112 select = trig + expr;
115 if (! tree) {Error(
"PlotEMCALQATrendingTree",
"No Tree found!");
return -1;}
118 TString* fCalorimeter;
125 tree->SetBranchAddress(
"fDate",&dtime);
126 tree->SetBranchAddress(
"nSM",&CurN);
127 tree->SetBranchAddress(
"fCalorimeter",&fCalorimeter);
128 tree->SetBranchAddress(
"system",&system);
129 tree->SetBranchAddress(
"period",&period);
130 tree->SetBranchAddress(
"pass",&pass);
131 tree->SetBranchAddress(
"fTrigger",&fTrigger);
134 tree->SetEventList(0);
135 tree->Draw(
">>elist",select);
136 tree->Draw(
">>listNotZero",select+NotZero);
137 TEventList* listNotZero = (TEventList*)gDirectory->Get(
"listNotZero");
138 TEventList* elist = (TEventList*)gDirectory->Get(
"elist");
139 tree->SetEventList(elist);
140 if(! elist->GetN()) { Error(
"PlotEMCALQATrendingTree",
"The current selection doess not match any entry!");
return -2; }
141 CurN = tree->GetMinimum(
"nSM");
142 const Int_t n = CurN;
143 tree->GetEntry(elist->GetEntry(0));
145 TGraphErrors* AverNclustersSM[n];
146 TGraphErrors* AverNcellsPerClusterSM[n];
147 TGraphErrors* AverESM[n];
148 TGraphErrors* AverMeanSM[n];
149 TGraphErrors* AverWidthSM[n];
150 TGraphErrors* AverNpi0SM[n];
154 TString base =
QAPATH + period->Data() +
"_" + pass->Data() +
"_";
157 TString ClusterAverages ; ClusterAverages = base +
"ClAvNew" + (*fTrigger)(r) +
".png";
158 TString Entries; Entries = base +
"Nentries" + (*fTrigger)(r) +
".png";
159 TString ClusterAveragesEnergy; ClusterAveragesEnergy = base +
"ClAvEne" + (*fTrigger)(r) +
".png";
160 TString ClusterAveragesEnergy2; ClusterAveragesEnergy2 = base +
"ClAvEne" + (*fTrigger)(r) +
".pdf";
161 TString ClusterAveragesEntries; ClusterAveragesEntries = base +
"ClAvEnt" + (*fTrigger)(r) +
".png";
162 TString ClusterAveragesEntries2; ClusterAveragesEntries2 = base +
"ClAvEnt" + (*fTrigger)(r) +
".pdf";
163 TString ClusterAveragesCells; ClusterAveragesCells = base +
"ClAvCells" + (*fTrigger)(r) +
".png";
164 TString ClusterAveragesCells2; ClusterAveragesCells2 = base +
"ClAvCells" + (*fTrigger)(r) +
".pdf";
165 TString Pi0Entries; Pi0Entries = base +
"Pi0Entries" + (*fTrigger)(r) +
".png";
166 TString Pi0Entries2; Pi0Entries2 = base +
"Pi0Entries" + (*fTrigger)(r) +
".pdf";
167 TString Pi0Mass; Pi0Mass = base +
"Pi0Mass" + (*fTrigger)(r) +
".png";
168 TString Pi0Mass2; Pi0Mass2 = base +
"Pi0Mass" + (*fTrigger)(r) +
".pdf";
169 TString Pi0Width; Pi0Width = base +
"Pi0Width" + (*fTrigger)(r) +
".png";
170 TString Pi0Width2; Pi0Width2 = base +
"Pi0Width" + (*fTrigger)(r) +
".pdf";
173 int nEmptyRuns = tree->Draw(
"run",
"Nevent==0",
"goff");
174 if (nEmptyRuns && (nEmptyRuns != -1)) {
175 Info(
"PlotEMCALQATrendingTree",Form(
"The following %i runs are empty for trigger %s:",nEmptyRuns,Trig));
176 for(Int_t i = 0 ; i < nEmptyRuns ; i++){
177 cout<<tree->GetV1()[i]<<endl;
181 int nRun = tree->Draw(
"run",
"",
"goff");
183 TH1F* h1 =
new TH1F(
"h1",
"dummy", nRun, 0., nRun+0.5);
184 TGaxis::SetMaxDigits(3);
186 h1->SetStats(kFALSE) ;
187 h1->SetAxisRange(0, nRun,
"X") ;
188 h1->GetXaxis()->SetTitle(
"RUN Index");
189 h1->GetXaxis()->SetTitleOffset(1.86);
190 h1->GetXaxis()->SetTitleSize(0.03);
192 for(Int_t i = 0 ; i < nRun ; i++){
194 label+=tree->GetV1()[i];
195 h1->GetXaxis()->SetBinLabel(i+1,label.Data());
196 h1->GetXaxis()->LabelsOption(
"v");
200 TCanvas* c1 =
new TCanvas(
"Nevents",
"Nb of events", 1000, 500);
202 c1->SetBorderSize(0);
203 c1->SetFrameBorderMode(0);
204 gStyle->SetOptStat(0);
205 gPad->SetLeftMargin(0.08);
206 gPad->SetRightMargin(0.02);
208 tree->Draw(
"NextInt():Nevent",
"",
"goff");
209 h1->GetYaxis()->SetTitle(
"N_{events}");
211 if (h1->GetMinimum() > 0.) {c1->SetLogy();}
214 TGraph* Nevents =
new TGraph(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2());
215 Nevents->SetMarkerStyle(20);
216 Nevents->SetMarkerColor(1);
217 Nevents->SetLineColor(2);
218 Nevents->Draw(
"same lp") ;
221 if(SavePlots) c1->SaveAs(Entries);
223 TCanvas* c2 =
new TCanvas(
"ClusterAveragesEvents",
"Mean Nb of Cluster per Event", 1000, 500);
225 c2->SetBorderSize(0);
226 c2->SetFrameBorderMode(0);
229 gPad->SetLeftMargin(0.08);
230 gPad->SetRightMargin(0.02);
233 TH1F* h2 = (TH1F*)h1->Clone(
"");
234 h2->GetYaxis()->SetTitle(
"<N_{clusters}>/event");
236 h2->GetXaxis()->SetTitle(
"RUN Index");
237 h2->GetXaxis()->SetTitleOffset(1.86);
238 h2->GetXaxis()->SetTitleSize(0.03);
241 tree->Draw(
"NextInt():ClusterMean:xe:ClusterRMS",
"",
"goff");
242 TGraphErrors * AverNclusters =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(), tree->GetV4());
243 AverNclusters->SetMarkerStyle(20);
244 AverNclusters->SetMarkerColor(1);
245 AverNclusters->Draw(
"same P") ;
247 for(Int_t ism = 0 ; ism < n ; ism++){
248 tree->Draw(Form(
"NextInt():ClusterMeanSM[%i]:xe:ClusterRMSSM[%i]",ism,ism),
"",
"goff");
249 AverNclustersSM[ism] =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
250 if (ism !=8)AverNclustersSM[ism]->SetMarkerColor(ism<10?ism+2:ism+1);
else AverNclustersSM[ism]->SetMarkerColor(7);
251 AverNclustersSM[ism]->SetMarkerStyle(21+(ism<10 ? ism: ism-10));
253 AverNclustersSM[ism]->Draw(
"same P");
256 TLegend* l2 =
new TLegend(0.123, 0.744, 0.933, 0.894);
257 l2->SetNColumns((n+1)/2.);
259 l2->SetBorderSize(0);
260 l2->SetTextSize(0.04);
261 l2->SetHeader(Form(
"<# of clusters> in %s (period %s trigger %s)",fCalorimeter->Data(),period->Data(),((*fTrigger)(r)).Data()));
262 l2->AddEntry(AverNclusters,
"average",
"p");
263 for(Int_t ism = 0 ; ism < n ; ism++){
264 TString projname = Form(
"SM %d",ism);
265 l2->AddEntry(AverNclustersSM[ism],projname.Data(),
"p");
269 if(SavePlots) c2->SaveAs(ClusterAveragesEntries);
270 if(SavePlots==2) c2->SaveAs(ClusterAveragesEntries2);
272 TCanvas* c3 =
new TCanvas(
"ClusterAveragesEnergy",
"Mean Cluster Energy", 1000, 500);
274 c3->SetBorderSize(0);
275 c3->SetFrameBorderMode(0);
278 gPad->SetLeftMargin(0.08);
279 gPad->SetRightMargin(0.02);
282 TH1F* h3 = (TH1F*)h1->Clone(
"");
283 h3->GetYaxis()->SetTitle(
"<E> (GeV)");
285 h3->GetXaxis()->SetTitle(
"RUN Index");
286 h3->GetXaxis()->SetTitleOffset(1.86);
287 h3->GetXaxis()->SetTitleSize(0.03);
290 tree->Draw(
"NextInt():EtotalMean:xe:EtotalRMS",
"",
"goff");
291 TGraphErrors * AverE =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
292 AverE->SetMarkerStyle(20);
293 AverE->SetMarkerColor(1);
294 AverE->Draw(
"same P");
296 for(Int_t ism = 0 ; ism < n ; ism++){
298 tree->Draw(Form(
"NextInt():EtotalMeanSM[%i]:xe:EtotalRMSSM[%i]",ism,ism),
"",
"goff");
299 AverESM[ism] =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
300 if (ism !=8)AverESM[ism]->SetMarkerColor(ism<10?ism+2:ism+1);
else AverESM[ism]->SetMarkerColor(7);
301 AverESM[ism]->SetMarkerStyle(21+(ism<10 ? ism: ism-10));
302 AverESM[ism]->Draw(
"same P");
306 TLegend* l3 =
new TLegend(0.123, 0.744, 0.933, 0.894);
307 l3->SetNColumns((n+1)/2.);
309 l3->SetBorderSize(0);
310 l3->SetTextSize(0.04);
311 l3->SetHeader(Form(
"<E> in %s (period %s trigger %s)",fCalorimeter->Data(),period->Data(),((*fTrigger)(r)).Data()));
312 l3->AddEntry(AverE,
"average",
"p");
313 for(Int_t ism = 0 ; ism < n ; ism++){
314 TString projname = Form(
"SM %d",ism);
315 l3->AddEntry(AverESM[ism],projname.Data(),
"p");
319 if(SavePlots) c3->SaveAs(ClusterAveragesEnergy);
320 if(SavePlots==2) c3->SaveAs(ClusterAveragesEnergy2);
322 TCanvas* c4 =
new TCanvas(
"ClusterAveragesCells",
"Mean Nb of Cells per Cluster", 1000, 500);
324 c4->SetBorderSize(0);
325 c4->SetFrameBorderMode(0);
328 gPad->SetLeftMargin(0.08);
329 gPad->SetRightMargin(0.02);
332 TH1F* h4 = (TH1F*)h1->Clone(
"");
333 h4->GetYaxis()->SetTitle(
"<N_{CellsPerCluster}>");
335 h4->GetXaxis()->SetTitle(
"RUN Index");
336 h4->GetXaxis()->SetTitleOffset(1.86);
337 h4->GetXaxis()->SetTitleSize(0.03);
341 tree->Draw(
"NextInt():CellPerClusterMean:xe:CellPerClusterRMS",
"",
"goff");
342 TGraphErrors * AverCellPerCluster =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
343 AverCellPerCluster->SetMarkerStyle(20);
344 AverCellPerCluster->SetMarkerColor(1);
346 for(Int_t ism = 0 ; ism < n ; ism++){
347 tree->Draw(Form(
"NextInt():CellPerClusterMeanSM[%i]:xe:CellPerClusterRMSSM[%i]",ism,ism),
"",
"goff");
348 AverNcellsPerClusterSM[ism] =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
349 if (ism !=8)AverNcellsPerClusterSM[ism]->SetMarkerColor(ism<10?ism+2:ism+1);
else AverNcellsPerClusterSM[ism]->SetMarkerColor(7);
350 AverNcellsPerClusterSM[ism]->SetMarkerStyle(21+(ism<10 ? ism: ism-10));
351 AverNcellsPerClusterSM[ism]->Draw(
"same P");
355 TLegend* l4 =
new TLegend(0.123, 0.744, 0.933, 0.894);
356 l4->SetNColumns((n+1)/2.);
358 l4->SetBorderSize(0);
359 l4->SetTextSize(0.04);
360 l4->SetHeader(Form(
"<# of cells per cluster> in %s (period %s trigger %s)",fCalorimeter->Data(),period->Data(),((*fTrigger)(r)).Data()));
361 l4->AddEntry(AverCellPerCluster,
"average",
"p");
362 for(Int_t ism = 0 ; ism < n ; ism++){
363 TString projname = Form(
"SM %d",ism);
364 l4->AddEntry(AverNcellsPerClusterSM[ism],projname.Data(),
"p");
368 if(SavePlots) c4->SaveAs(ClusterAveragesCells);
370 TCanvas* c5 =
new TCanvas(
"Pi0Position",
"Mean Pi0 Mass", 1000, 500);
372 c5->SetBorderSize(0);
373 c5->SetFrameBorderMode(0);
376 gStyle->SetOptStat(0);
378 gPad->SetLeftMargin(0.08);
379 gPad->SetRightMargin(0.02);
382 TH1F * h5 = (TH1F*)h1->Clone(
"");
384 h5->GetXaxis()->SetTitle(
"RUN Index");
385 h5->GetXaxis()->SetTitleOffset(1.86);
386 h5->GetXaxis()->SetTitleSize(0.03);
387 h5->GetYaxis()->SetTitle(
"Mean_{#pi^{0}}");
391 tree->Draw(
"NextInt():MeanPos:xe:MeanPosErr",
"",
"goff");
392 TGraphErrors * AverMean =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
393 AverMean->SetMarkerStyle(20);
394 AverMean->SetMarkerColor(1);
395 AverMean->Draw(
"same P");
397 for(Int_t ism = 0 ; ism < n ; ism++){
399 tree->Draw(Form(
"NextInt():MeanPosSM[%i]:xe:MeanPosErrSM[%i]",ism,ism),
"",
"goff");
400 AverMeanSM[ism] =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
401 if (ism !=8)AverMeanSM[ism]->SetMarkerColor(ism<10?ism+2:ism+1);
else AverMeanSM[ism]->SetMarkerColor(7);
402 AverMeanSM[ism]->SetMarkerStyle(21+(ism<10 ? ism: ism-10));
403 AverMeanSM[ism]->Draw(
"same P");
407 TLegend* l5 =
new TLegend(0.123, 0.744, 0.933, 0.894);
408 l5->SetNColumns((n+1)/2.);
410 l5->SetBorderSize(0);
411 l5->SetTextSize(0.04);
412 l5->SetHeader(Form(
"<M_{#pi^{0}}> (MeV) in %s (period %s trigger %s)",fCalorimeter->Data(),period->Data(),((*fTrigger)(r)).Data()));
413 l5->AddEntry(AverMean,
"average",
"p");
414 for(Int_t ism = 0 ; ism < n ; ism++){
415 TString projname = Form(
"SM %d",ism);
416 l5->AddEntry(AverMeanSM[ism],projname.Data(),
"p");
421 if(SavePlots) c5->SaveAs(Pi0Mass);
424 TCanvas* c6 =
new TCanvas(
"Pi0Width",
"Mean Pi0 Width", 1000, 500);
426 c6->SetBorderSize(0);
427 c6->SetFrameBorderMode(0);
430 gPad->SetLeftMargin(0.08);
431 gPad->SetRightMargin(0.02);
434 TH1F* h6 = (TH1F*)h1->Clone(
"");
436 h6->GetXaxis()->SetTitle(
"RUN Index");
437 h6->GetXaxis()->SetTitleOffset(1.86);
438 h6->GetXaxis()->SetTitleSize(0.03);
439 h6->GetYaxis()->SetTitle(
"#sigma_{#pi^{0}}");
442 tree->Draw(
"NextInt():Width:xe:WidthErr",
"",
"goff");
443 TGraphErrors * AverWidth =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
444 AverWidth->SetMarkerStyle(20);
445 AverWidth->SetMarkerColor(1);
446 AverWidth->Draw(
"same P");
448 for(Int_t ism = 0 ; ism < n ; ism++){
449 tree->Draw(Form(
"NextInt():WidthSM[%i]:xe:WidthErrSM[%i]",ism,ism),
"",
"goff");
450 AverWidthSM[ism] =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
451 if (ism !=8)AverWidthSM[ism]->SetMarkerColor(ism<10?ism+2:ism+1);
else AverWidthSM[ism]->SetMarkerColor(7);
452 AverWidthSM[ism]->SetMarkerStyle(21+(ism<10 ? ism: ism-10));
453 AverWidthSM[ism]->Draw(
"same P");
457 TLegend* l6 =
new TLegend(0.123, 0.744, 0.933, 0.894);
458 l6->SetNColumns((n+1)/2.);
460 l6->SetBorderSize(0);
461 l6->SetTextSize(0.04);
462 l6->SetHeader(Form(
"#sigma_{#pi^{0}} in %s (period %s trigger %s)",fCalorimeter->Data(),period->Data(),((*fTrigger)(r)).Data()));
463 l6->AddEntry(AverWidth,
"total",
"p");
464 for(Int_t ism = 0 ; ism < n ; ism++){
465 TString projname = Form(
"SM %d",ism);
466 l6->AddEntry(AverWidthSM[ism],projname.Data(),
"p");
470 if(SavePlots) c6->SaveAs(Pi0Width);
472 TCanvas* c7 =
new TCanvas(
"Npi0",
"Mean Nb of Pi0", 1000, 500);
474 c7->SetBorderSize(0);
475 c7->SetFrameBorderMode(0);
478 gPad->SetLeftMargin(0.08);
479 gPad->SetRightMargin(0.02);
482 TH1F* h7 = (TH1F*)h1->Clone(
"");
484 if (h7->GetMinimum() > 0.) {c7->SetLogy();}
485 h7->GetXaxis()->SetTitle(
"RUN Index");
486 h7->GetXaxis()->SetTitleOffset(1.86);
487 h7->GetXaxis()->SetTitleSize(0.03);
488 h7->GetYaxis()->SetTitle(
"<N_{#pi^{0}}>/event");
491 tree->Draw(
"NextInt():Npi0:xe:Npi0Err",
"",
"goff");
492 if (tree->GetMinimum(
"Npi0") > 1) c4->SetLogy();
493 TGraphErrors * AverNpi0 =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
494 AverNpi0->SetMarkerStyle(20);
495 AverNpi0->SetMarkerColor(1);
496 AverNpi0->Draw(
"same P");
498 for(Int_t ism = 0 ; ism < n ; ism++){
499 tree->Draw(Form(
"NextInt():Npi0SM[%i]:xe:Npi0ErrSM[%i]",ism,ism),
"",
"goff");
500 AverNpi0SM[ism] =
new TGraphErrors(tree->GetSelectedRows(), tree->GetV1(), tree->GetV2(),tree->GetV3(),tree->GetV4());
501 if (ism !=8)AverNpi0SM[ism]->SetMarkerColor(ism<10?ism+2:ism+1);
else AverNpi0SM[ism]->SetMarkerColor(7);
502 AverNpi0SM[ism]->SetMarkerStyle(21+(ism<10 ? ism: ism-10));
503 AverNpi0SM[ism]->Draw(
"same P");
506 TLegend* l7 =
new TLegend(0.123, 0.744, 0.933, 0.894);
507 l7->SetNColumns((n+1)/2.);
509 l7->SetBorderSize(0);
510 l7->SetTextSize(0.04);
511 l7->SetHeader(Form(
"<N_{#pi^{0}}>/event in %s (period %s trigger %s)",fCalorimeter->Data(),period->Data(),((*fTrigger)(r)).Data()));
512 l7->AddEntry(AverNpi0,
"total",
"p");
513 for(Int_t ism = 0 ; ism < n ; ism++){
514 TString projname = Form(
"SM %d",ism);
515 l7->AddEntry(AverNpi0SM[ism],projname.Data(),
"p");
519 if(SavePlots) c7->SaveAs(Pi0Entries);
520 if(SavePlots==2) c7->SaveAs(Pi0Entries2);
522 fout->mkdir(Form(
"%s/%s/%s/%s",period->Data(),pass->Data(),
"TrendingQA",fTrigger->Data()));
524 fout->Cd(Form(
"%s/%s/%s/%s",period->Data(),pass->Data(),
"TrendingQA",fTrigger->Data()));
526 gROOT->GetListOfCanvases()->Write();
527 gROOT->GetListOfCanvases()->Delete();
529 if((!Expr.IsNull()) && (!Expr.EndsWith(
".root"))) elist->Write();
530 if(listNotZero) {listNotZero->Reset();}
531 if(elist) {elist->Reset();}
539 TH1F*
ZoomFromTree(TH1F* h, TTree* atree, Int_t n,
const char* aVar, UShort_t aScaleFactor)
542 atree->SetEventList(0);
543 TEventList *listNotZero = (TEventList*)gDirectory->Get(
"listNotZero");
544 atree->SetEventList(listNotZero);
548 double offset = 30*((treeMax - treeMin)/(100.*aScaleFactor));
550 if(treeMin != -treeMax){
551 h->SetMinimum(TMath::Max(0.,treeMin-offset));
552 h->SetMaximum(treeMax+2*offset);
555 atree->SetEventList(0);
556 TEventList *elist = (TEventList*)gDirectory->Get(
"elist");
557 atree->SetEventList(elist);
567 TLeaf* leaf = aTree->GetLeaf(columname);
571 TBranch* branch = leaf->GetBranch();
572 Double_t cmin = 3.40282e+38;
573 for (Long64_t i = 0; i < aTree->GetEntries(); ++i) {
574 Long64_t entryNumber = aTree->GetEntryNumber(i);
575 if (entryNumber < 0)
break;
576 branch->GetEntry(entryNumber);
577 for (Int_t j = 0;j < TMath::Min(leaf->GetLen(),n); ++j) {
578 Double_t val = leaf->GetValue(j);
593 TLeaf* leaf = aTree->GetLeaf(columname);
597 TBranch* branch = leaf->GetBranch();
598 Double_t cmax = - 3.40282e+38;
599 for (Long64_t i = 0; i < aTree->GetEntries(); ++i) {
600 Long64_t entryNumber = aTree->GetEntryNumber(i);
601 if (entryNumber < 0)
break;
602 branch->GetEntry(entryNumber);
603 for (Int_t j = 0; j < TMath::Min(leaf->GetLen(),n); ++j) {
604 Double_t val = leaf->GetValue(j);
Double_t GetTreeMinimum(TTree *aTree, Int_t n, const char *columname)
Double_t GetTreeMaximum(TTree *aTree, Int_t n, const char *columname)
int NextInt(int newMax=0)
TH1F * ZoomFromTree(TH1F *h, TTree *atree, Int_t n, const char *aVar, UShort_t aScaleFactor=1)
int PlotEMCALQATrendingTree(TTree *tree, const char *trig, TFile *fout, Bool_t SavePlots, TString expr)