10 # include <TFitResult.h>
13 # include <THashList.h>
14 # include <TParameter.h>
36 const Color_t
cc[] = { kMagenta+2,
62 TObject*
GetO(TSeqCollection* l,
const char* name, TClass* cls)
65 Warning(
"GetO",
"No list passed");
68 TObject* o = l->FindObject(name);
70 Warning(
"GetO",
"No object named %s found in list %s",
76 if (!o->IsA()->InheritsFrom(cls)) {
77 Warning(
"GetO",
"Object %s read from %s is not a %s but a %s",
78 name, l->GetName(), cls->GetName(), o->ClassName());
84 TObject*
CopyO(TSeqCollection* l,
const char* name,
const char* newName,
89 TObject* copy = orig->Clone(newName ? newName : name);
103 return static_cast<TH1*
>(
GetO(l,name,TH1::Class()));
116 return static_cast<TH2*
>(
GetO(l,name,TH2::Class()));
128 TH1*
CopyH1(TSeqCollection* l,
const char* name,
const char* newName=0)
130 TH1* copy =
static_cast<TH1*
>(
CopyO(l,name,newName,TH1::Class()));
131 if (copy) copy->SetDirectory(0);
144 TH2*
CopyH2(TSeqCollection* l,
const char* name,
const char* newName=0)
146 TH2* copy =
static_cast<TH2*
>(
CopyO(l,name,newName,TH2::Class()));
147 if (copy) copy->SetDirectory(0);
170 if (a1->GetNbins() != a2->GetNbins()) {
171 ::Warning(
"CheckAxisNBins",
"Incompatible number %s bins: %d vs %d",
172 which, a1->GetNbins(), a2->GetNbins());
191 if (!TMath::AreEqualRel(a1->GetXmin(), a2->GetXmin(),1.E-12) ||
192 !TMath::AreEqualRel(a1->GetXmax(), a2->GetXmax(),1.E-12)) {
193 Warning(
"CheckAxisLimits",
194 "Limits of %s axis incompatible [%f,%f] vs [%f,%f]", which,
195 a1->GetXmin(), a1->GetXmax(), a2->GetXmin(), a2->GetXmax());
214 const TArrayD * h1Array = a1->GetXbins();
215 const TArrayD * h2Array = a2->GetXbins();
216 Int_t fN = h1Array->fN;
217 if ( fN == 0 )
return true;
218 if (h2Array->fN != fN) {
220 Warning(
"CheckAxisBins",
"Not equal number of %s bin limits: %d vs %d",
221 which, fN, h2Array->fN);
225 for (
int i = 0; i < fN; ++i) {
226 if (!TMath::AreEqualRel(h1Array->GetAt(i),h2Array->GetAt(i),1E-10)) {
227 Warning(
"CheckAxisBins",
228 "%s limit # %3d incompatible: %f vs %f",
229 which, i, h1Array->GetAt(i),h2Array->GetAt(i));
253 THashList *l1 = (
const_cast<TAxis*
>(a1))->GetLabels();
254 THashList *l2 = (
const_cast<TAxis*
>(a2))->GetLabels();
256 if (!l1 && !l2)
return true;
258 Warning(
"CheckAxisLabels",
"Difference in %s labels: %p vs %p",
263 if (l1->GetSize() != l2->GetSize()) {
264 Warning(
"CheckAxisLabels",
"Different number of %s labels: %d vs %d",
265 which, l1->GetSize(), l2->GetSize());
268 for (
int i = 1; i <= a1->GetNbins(); ++i) {
269 TString label1 = a1->GetBinLabel(i);
270 TString label2 = a2->GetBinLabel(i);
271 if (label1 != label2) {
272 Warning(
"CheckAxisLabels",
"%s label # %d not the same: '%s' vs '%s'",
273 which, i, label1.Data(), label2.Data());
315 if (h1 == h2)
return true;
317 if (h1->GetDimension() != h2->GetDimension()) {
318 Warning(
"CheckConsistency",
319 "%s and %s have different dimensions %d vs %d",
320 h1->GetName(), h2->GetName(),
321 h1->GetDimension(), h2->GetDimension());
324 Int_t dim = h1->GetDimension();
327 Bool_t alsoLbls = (h1->GetEntries() != 0 && h2->GetEntries() != 0);
328 if (!
CheckAxis(
"X", h1->GetXaxis(), h2->GetXaxis(), alsoLbls)) ret =
false;
330 !
CheckAxis(
"Y", h1->GetYaxis(), h2->GetYaxis(), alsoLbls)) ret =
false;
332 !
CheckAxis(
"Z", h1->GetZaxis(), h2->GetZaxis(), alsoLbls)) ret =
false;
349 TH1* ret =
static_cast<TH1*
>(templ->Clone(src->GetName()));
350 ret->SetDirectory(0);
352 for (
Int_t i = 1; i <= ret->GetNbinsX(); i++) {
353 Double_t x = ret->GetXaxis()->GetBinCenter(i);
354 Int_t b = src->GetXaxis()->FindBin(x);
355 if (b <= 0 || b >= src->GetNbinsX())
continue;
358 ret->SetBinContent(i, c);
359 ret->SetBinError (i, e);
376 TH2* ret =
static_cast<TH2*
>(templ->Clone(src->GetName()));
377 ret->SetDirectory(0);
379 for (
Int_t i = 1; i <= ret->GetNbinsX(); i++) {
380 Double_t x = ret->GetXaxis()->GetBinCenter(i);
381 Int_t bx = src->GetXaxis()->FindBin(x);
382 if (bx <= 0 || bx >= src->GetNbinsX())
continue;
383 for (
Int_t j = 1; j <= ret->GetNbinsY(); j++) {
384 Double_t y = ret->GetYaxis()->GetBinCenter(j);
385 Int_t by = src->GetYaxis()->FindBin(y);
387 Double_t e = src->GetBinError(bx,by);
388 ret->SetBinContent(i, j, c);
389 ret->SetBinError (i, j, e);
423 h->SetMarkerColor(color);
424 h->SetMarkerStyle(marker);
425 h->SetMarkerSize(size);
426 h->SetFillColor(color);
427 h->SetFillStyle(fill);
428 h->SetLineColor(color);
429 h->SetLineStyle(line);
430 h->SetLineWidth(width);
431 h->GetXaxis()->SetNdivisions(210);
432 h->GetYaxis()->SetNdivisions(210);
444 Printf(
"Content of %s - %s", h->GetName(), h->GetTitle());
445 for (
Int_t i = 1; i <= h->GetNbinsX(); i++) {
447 for (
Int_t j = 1; j <= h->GetNbinsY(); j++) {
450 if (TMath::IsNaN(c) || TMath::IsNaN(e))
451 printf(
"*** NAN ***");
453 printf(
"%.*f+/-%.*f ", prec, c, prec, e);
467 Printf(
"Content of %s - %s", h->GetName(), h->GetTitle());
468 for (
Int_t i = 1; i <= h->GetNbinsX(); i++) {
469 if (h->GetBinContent(i) <= 1e-6)
continue;
470 printf(
"%3d (%+6.3f): %.*f+/-%.*f\n", i,
471 h->GetXaxis()->GetBinCenter(i),
472 prec, h->GetBinContent(i),
473 prec, h->GetBinError(i));
502 Int_t bMin = h->FindBin(min+epsilon);
503 Int_t bMax = h->FindBin(max-epsilon);
504 if (bMin < 1) bMin = 0;
505 Double_t val = h->IntegralAndError(bMin, bMax, err);
507 if (TMath::Abs(h->GetXaxis()->GetXmin()+h->GetXaxis()->GetXmax())>=epsilon)
512 bMin = h->FindBin(-min+epsilon);
513 bMax = h->FindBin(-max-epsilon);
514 val += h->IntegralAndError(bMin, bMax, err2);
515 err = TMath::Sqrt(err*err+err2*err2);
540 if (TMath::Abs(n) < 1e-16 || TMath::Abs(d) < 1e-16)
return 0;
542 er = TMath::Sqrt(en*en/n/n + ed*ed/d/d);
570 for (
Int_t i = 1; i <= h->GetNbinsX(); i++) {
573 Double_t dr = (dc > 1e-6 ? 1/dc : 0);
575 for (
Int_t j = 1; j <= h->GetNbinsY(); j++) {
576 Double_t nc = h->GetBinContent(i,j);
580 Double_t se = sc*TMath::Sqrt(ns*ns+dq*dq);
582 h->SetBinContent(i,j,sc);
583 h->SetBinError (i,j,se);
608 for (
Int_t i = 1; i <= h->GetNbinsX(); i++) {
613 Double_t se = sc*TMath::Sqrt(es*es+rr*rr);
615 h->SetBinContent(i,sc);
616 h->SetBinError (i,se);
641 for (
Int_t i = 1; i <= h->GetNbinsX(); i++) {
642 for (
Int_t j = 1; j <= h->GetNbinsY(); j++) {
647 Double_t se = sc*TMath::Sqrt(es*es+rr*rr);
649 h->SetBinContent(i,j,sc);
650 h->SetBinError (i,j,se);
668 for (
Int_t iy = 1; iy <= in->GetNbinsY(); iy++) {
669 Double_t z = in->GetYaxis()->GetBinCenter(iy);
670 Int_t bin = ipz->GetXaxis()->FindBin(z);
671 Double_t nEv = ipz->GetBinContent(bin);
672 Double_t eEv = ipz->GetBinError (bin);
673 Double_t esc = (nEv > 0 ? 1./nEv : 0);
676 for (
Int_t ix = 1; ix <= in->GetNbinsX(); ix++) {
678 Double_t e = in->GetBinError (ix,iy);
683 if (full)se = sc * TMath::Sqrt(r*r+rE2);
685 Double_t scl = 1 / in->GetXaxis()->GetBinWidth(ix);
688 in->SetBinContent(ix, iy, scl*sc);
689 in->SetBinError (ix, iy, scl*se);
716 TH2* mask = other ? other : h;
719 Int_t nIPz = h->GetNbinsY();
720 Int_t nEta = h->GetNbinsX();
721 TH1* p = h->ProjectionX(name,1,nIPz,
"e");
725 p->SetYTitle(Form(
"#LT(%s)#GT", h->GetYaxis()->GetTitle()));
727 for (
Int_t etaBin = 1; etaBin <= nEta; etaBin++) {
733 for (
Int_t ipzBin = 1; ipzBin <= nIPz; ipzBin++) {
734 Double_t bc = mask->GetBinContent(etaBin, ipzBin);
735 if (bc < 1e-9)
continue;
736 Double_t be = mask->GetBinError (etaBin, ipzBin);
737 if (TMath::IsNaN(be))
continue;
747 TMath::Sort(j, hr.fArray, idx.fArray,
false);
756 for (k = 0; k < j; k++) {
758 Int_t ipzBin = hb[l];
762 Double_t x = TMath::Sqrt(nsume+hee*hee)/(nsum+hvv);
766 Double_t by = h->GetYaxis()->GetBinCenter(ipzBin);
767 Int_t ib = ipz ? ipz->FindBin(by) : 0;
769 nsum += h->GetBinContent(etaBin, ipzBin);
770 nsume += TMath::Power(h->GetBinError(etaBin, ipzBin), 2);
773 dsum += !ipz ? 1 : ipz->GetBinContent(ib);
774 dsume += !ipz ? 0 : TMath::Power(ipz->GetBinError(ib),2);
777 if (k == 0 || n == 0) {
782 Double_t norm = (mode > 0 ? n : dsum);
787 if (mode==2) ave = TMath::Sqrt(nsume)/n;
788 else ave = av*TMath::Sqrt(rne+rde);
790 Printf(
"%10s - bin %3d (%+6.3f) count=%2d n=%9.3f+/-%9.3f d=%9.3f+/-%9.3f "
791 "norm=%9.3f -> %9.3f +/- %9.3f",
792 h->GetName(), etaBin, h->GetXaxis()->GetBinCenter(etaBin), n,
793 nsum, TMath::Sqrt(nsume)/(mode == 2 ? n : nsum),
794 dsum, TMath::Sqrt(dsume)/(mode == 2 ? n : dsum),
798 p->SetBinContent(etaBin, av);
799 p->SetBinError (etaBin, ave);
801 if (mode == 0) p->Scale(1,
"width");
821 TH2* tmp =
static_cast<TH2*
>(h->Clone(
"tmp"));
822 tmp->SetDirectory(0);
825 TH1* ret =
Avg(tmp, 0, 2, name);
866 if (!deltaRec || !deltaBg) {
867 Error(
"GetDeltaScale",
"Missing input histogram(s): %p %p",
871 deltaRec->Scale(1./nEvRec);
872 deltaBg ->Scale(1./nEvBg);
873 Double_t top = TMath::Min(lim,deltaRec->GetXaxis()->GetXmax());
881 Printf(
"%20s: data=%9g+/-%9g inj=%9g+/-%9g -> %9g+/-%9g",
882 list->GetName(), iRec, eRec, iBg, eBg, r, eR);
906 TH1* deltaRec =
CopyH1(list,Form(
"b%d_TrData_WDist", b));
907 TH1* deltaBg =
CopyH1(list,Form(
"b%d_TrInj_WDist", b));
930 TSeqCollection* simList,
937 TH1* deltaReal =
CopyH1(realList,Form(
"b%d_TrData_WDist", b));
938 TH1* deltaSim =
CopyH1(simList, Form(
"b%d_TrData_WDist", b));
965 Double_t top = deltaReal->GetYaxis()->GetXmax();
966 TH1* ret = deltaReal->ProjectionX(
"scalar");
968 ret->SetTitle(
"Background scale (#it{k})");
969 ret->SetYTitle(
"#delta scale");
970 ret->SetXTitle(
"#it{#eta}");
971 ret->SetDirectory(0);
972 for (
Int_t i = 1; i <= ret->GetNbinsX(); i++) {
973 TH1* tmpReal = deltaReal->ProjectionY(
"tmpReal",i,i);
974 TH1* tmpSim = deltaSim ->ProjectionY(
"tmpSim", i,i);
976 ret->SetBinContent(i, r);
977 ret->SetBinError (i, eR);
1000 TH2* deltaRec =
GetH2(list,Form(
"b%d_TrData_WDvsEta", b));
1001 TH2* deltaBg =
GetH2(list,Form(
"b%d_TrInj_WDvsEta", b));
1021 TSeqCollection* simList,
1027 TH2* deltaReal =
CopyH2(realList,Form(
"b%d_TrData_WDvsEta", b));
1028 TH2* deltaSim =
CopyH2(simList, Form(
"b%d_TrData_WDvsEta", b));
1029 return GetDeltaDist(deltaReal, deltaSim, nEvReal, nEvSim, lim);
1054 TSeqCollection* list,
1060 TH2* in =
GetH2(list,Form(
"b%d_Tr%s_ZvEtaCutT",b,comb ?
"Comb" :
"Inj"));
1061 TH2* ret =
static_cast<TH2*
>(in->Clone(
"bg"));
1062 in ->SetDirectory(0);
1063 ret->SetDirectory(0);
1065 if (comb)
return ret;
1072 return Scale(ret, r, eR);
1110 TH2* beta =
static_cast<TH2*
>(bg->Clone(
"beta"));
1111 beta->Divide(simMeas);
1112 beta->SetDirectory(0);
1113 TH2* one =
static_cast<TH2*
>(bg->Clone(
"one"));
1115 one->SetDirectory(0);
1116 for (
Int_t i = 1; i <= one->GetNbinsX(); i++)
1117 for (
Int_t j = 1; j <= one->GetNbinsY(); j++)
1118 one->SetBinContent(i,j,1);
1120 meas->Multiply(one);
1134 TH1* h =
CopyH1(l, Form(
"b%d_Tr%s_WDist", b, which));
1137 h->SetYTitle(
"tracklets/event");
1138 h->SetXTitle(
"#Delta");
1140 Bool_t isData = (scale < 1e-6);
1141 Color_t
color = (isData ? kGreen+2 : kBlue+2);
1143 if (w.EqualTo(
"Data"))
SetAttr(h, color, isData ? 21 : 25);
1144 if (w.EqualTo(
"Inj"))
SetAttr(h, color, isData ? 23 : 32);
1145 if (w.EqualTo(
"Comb"))
SetAttr(h, color, 26);
1146 if (w.EqualTo(
"Prim"))
SetAttr(h, kMagenta+2, 28);
1147 if (w.EqualTo(
"Sec"))
SetAttr(h, kMagenta+2, 27);
1149 if (w.EqualTo(
"Data")) h->SetTitle(
"Measured");
1150 if (w.EqualTo(
"Inj")) h->SetTitle(
"Injection");
1151 if (w.EqualTo(
"Comb")) h->SetTitle(
"Combinatorial");
1152 if (w.EqualTo(
"Prim")) h->SetTitle(
"Primaries");
1153 if (w.EqualTo(
"Sec")) h->SetTitle(
"Secondaries");
1155 if (w.EqualTo(
"Inj")) {
1158 h->SetTitle(Form(
"%5.3f#pm%5.3f #times %s",
1159 r, eR, h->GetTitle()));
1163 Scale(h, scale, scaleE);
1164 h->SetTitle(Form(
"%5.3f#pm%5.3f #times %s",
1165 scale, scaleE, h->GetTitle()));
1257 TSeqCollection* realList,
1258 TSeqCollection* simList,
1267 const Double_t k = (realList == simList ? 1 : 1.3);
1268 realList->SetName(
"realList");
1269 simList ->SetName(
"simList");
1270 TH2* realMeas =
GetH2(realList, Form(
"b%d_TrData_ZvEtaCutT", b));
1271 TH2* realIPzC =
GetH2(realList,
"zv");
1272 TH1* realIPz = realIPzC->ProjectionX(
"realIPz", b+1,b+1,
"e");
1274 GetH2(simList, Form(
"b%d_TrData_ZvEtaCutT", b)));
1276 TH1* simIPz = simIPzC->ProjectionX(
"simIPz", b+1,b+1,
"e");
1278 GetH2(simList, Form(
"b%d_zvEtaPrimMC", b)));
1280 TH1* trueIPz = simIPzC->ProjectionX(
"trueIPz", b+1,b+1,
"e");
1281 TH2* simBg =
Coerce(realMeas,
GetBg(simComb, simList, b, 1, lim));
1283 (realComb ?simList : realList),
1291 if (scaleMode == 0) {
1294 Scale(realBg, scale, 0);
1296 else if (scaleMode == 1) {
1301 realIPz->GetEntries(), lim);
1303 simIPz->GetEntries(), lim);
1304 scale =
RatioE(realR, realRe, simR, simRe, scaleE);
1307 realIPz->GetEntries(),
1308 simIPz->GetEntries(),
1311 Scale(realBg, scale, scaleE);
1313 else if (scaleMode == 2) {
1318 tmpReal->Divide(tmpSim);
1319 tmpReal->SetName(
"k");
1320 tmpReal->SetTitle(
"Background scale");
1321 tmpReal->SetMarkerColor(tmpReal->GetMarkerColor());
1322 tmpReal->SetMarkerStyle(tmpReal->GetMarkerStyle()+4);
1323 tmpReal->SetLineColor(tmpReal->GetLineColor());
1325 Scale(realBg, tmpReal);
1328 realIPz->GetEntries(),
1329 simIPz ->GetEntries(), lim);
1330 TF1* kfit =
new TF1(
"kfit",
"pol0");
1331 histK->Fit(kfit,
"Q0+");
1332 scale = kfit->GetParameter(0);
1333 scaleE = kfit->GetParError(0);
1334 Scale(realBg, histK);
1343 Warning(
"CorrectOneBin",
"Real data (%s) and simulated data (%s) done "
1344 "with different binning", realList->GetName(), simList->GetName());
1347 realMeas->SetName(
"realMeas");
1348 simMeas ->SetName(
"simMeas");
1349 realBg ->SetName(
"realBg");
1350 simBg ->SetName(
"simBg");
1351 trueGen ->SetName(
"trueGen");
1352 realMeas->SetTitle(
"Measured (real)");
1353 realBg ->SetTitle(
"Background (real)");
1354 simMeas ->SetTitle(
"Measured (simulated)");
1355 simBg ->SetTitle(
"Background (simulated)");
1356 trueGen ->SetTitle(
"Generated");
1357 realMeas->SetDirectory(0);
1358 simMeas ->SetDirectory(0);
1359 trueGen ->SetDirectory(0);
1360 realIPz ->SetDirectory(0);
1361 simIPz ->SetDirectory(0);
1362 trueIPz ->SetDirectory(0);
1366 Double_t ipzMin = realMeas->GetYaxis()->GetXmin();
1367 Double_t ipzMax = realMeas->GetYaxis()->GetXmax();
1371 realIPz->SetMarkerColor(kRed+2); realIPz->SetMarkerStyle(20);
1372 simIPz ->SetMarkerColor(kBlue+2); simIPz ->SetMarkerStyle(21);
1373 trueIPz->SetMarkerColor(kBlack); trueIPz->SetMarkerStyle(24);
1376 for (
Int_t i = 1; i <= realIPz->GetNbinsX(); i++) {
1377 if (simIPz->GetBinContent(i) < 1e-6) {
1378 realIPz->SetBinContent(i,0);
1379 realIPz->SetBinError (i,0);
1395 realMeas->Scale(1/realNev);
1396 realBg ->Scale(1/(realComb ? simNev : realNev));
1397 realIPz ->Scale(1/realNev);
1398 simMeas ->Scale(1/simNev);
1399 simBg ->Scale(1/simNev);
1400 simIPz ->Scale(1/simNev);
1401 trueGen ->Scale(1/trueNev);
1402 trueIPz ->Scale(1/trueNev);
1406 TH2* realSig =
static_cast<TH2*
>(realMeas->Clone(
"realSig"));
1407 realSig->SetTitle(
"Signal (real)");
1408 realSig->SetDirectory(0);
1412 TH2* simSig =
static_cast<TH2*
>(simMeas->Clone(
"simSig"));
1413 simSig->SetTitle(
"Signal (simulated)");
1414 simSig->SetDirectory(0);
1416 simSig->Add(simBg, -1);
1419 TH2* alpha =
static_cast<TH2*
>(trueGen->Clone(
"alpha"));
1420 alpha->SetTitle(
"#alpha");
1421 alpha->SetDirectory(0);
1422 alpha->Divide(simSig);
1424 TH2* fiducial =
static_cast<TH2*
>(alpha->Clone(
"fiducial"));
1425 fiducial->SetTitle(
"Fiducial cut");
1426 fiducial->SetDirectory(0);
1427 for (
Int_t i = 1; i <= fiducial->GetNbinsX(); i++) {
1428 for (
Int_t j = 1; j <= fiducial->GetNbinsY(); j++) {
1429 Double_t a = fiducial->GetBinContent(i,j);
1430 fiducial->SetBinError(i,j,0);
1431 fiducial->SetBinContent(i,j, (a > 0 && a <= 2.5));
1434 alpha ->Multiply(fiducial);
1437 TH2* result =
static_cast<TH2*
>(realSig->Clone(
"result"));
1438 result->SetTitle(
"Result");
1439 result->SetDirectory(0);
1440 result->Multiply(alpha);
1445 TH1* truth =
Avg(trueGen, trueIPz, mode,
"truth");
1447 SetAttr(truth,
cc[b%11], 24, 1.5, 0, 3, 2);
1450 TH1* dndeta =
Avg(result, realIPz, mode,
"dndeta");
1452 SetAttr(dndeta,
cc[b%11], 20, 1.2, 0, 1, 1);
1454 for (
Int_t i = 1; i <= dndeta->GetNbinsX(); i++) {
1455 dndeta->SetBinError(i,1./3*dndeta->GetBinError(i));
1459 TF1* f =
new TF1(
"ff",
"pol0", -flim, flim);
1460 TFitResultPtr r = dndeta->Fit(f,
"QN0RS",
"", -flim, flim);
1461 Printf(
"dNch/deta %20s: %6.1f +/- %6.1f (%5.2f - %5.3f+/-%5.3f)",
1463 r->Parameter(0), r->ParError(0), r->Chi2()/r->Ndf(),
1468 realMeas->Multiply(fiducial);
1469 simMeas ->Multiply(fiducial);
1470 realSig ->Multiply(fiducial);
1471 simSig ->Multiply(fiducial);
1472 realBg ->Multiply(fiducial);
1473 simBg ->Multiply(fiducial);
1474 realSig ->Multiply(fiducial);
1475 simSig ->Multiply(fiducial);
1476 TH1* realAvgMeas =
Avg(realMeas,realIPz, mode,
"realAvgMeas",result);
1477 TH1* simAvgMeas =
Avg(simMeas, simIPz, mode,
"simAvgMeas", result);
1478 TH1* realAvgBg =
Avg(realBg, realIPz, mode,
"realAvgBg", result);
1479 TH1* simAvgBg =
Avg(simBg, simIPz, mode,
"simAvgBg", result);
1480 TH1* realAvgSig =
Avg(realSig, realIPz, mode,
"realAvgSig", result);
1481 TH1* simAvgSig =
Avg(simSig, simIPz, mode,
"simAvgSig", result);
1482 SetAttr(realAvgMeas, kGreen+2, 21);
1483 SetAttr(realAvgBg, kGreen+2, 22);
1484 SetAttr(realAvgSig, kGreen+2, 23);
1485 SetAttr(simAvgMeas, kRed+2, 22);
1486 SetAttr(simAvgBg, kBlue+2, 21);
1487 SetAttr(simAvgSig, kBlue+2, 22);
1490 THStack* summary =
new THStack(
"summary", dndeta->GetYaxis()->GetTitle());
1491 summary->Add(dndeta,
"e2");
1492 summary->Add(truth,
"E");
1493 summary->Add(realAvgMeas);
1494 summary->Add(simAvgMeas);
1495 summary->Add(realAvgBg);
1496 summary->Add(simAvgBg);
1497 summary->Add(realAvgSig);
1498 summary->Add(simAvgSig);
1502 summary->Add(histK);
1507 realIPz->Scale(1/realNev);
1508 simIPz ->Scale(1/simNev);
1509 trueIPz->Scale(1/trueNev);
1512 TCanvas*
c =
new TCanvas(Form(
"c%02d",b),
1513 Form(
"Centrality bin %d",b),
1516 c->cd(1); realMeas->DrawCopy(
"lego2 e");
1517 c->cd(2); realBg ->DrawCopy(
"lego2 e");
1518 c->cd(3); realSig ->DrawCopy(
"lego2 e");
1519 c->cd(4); realIPz ->DrawCopy();
1520 c->cd(4); simIPz ->DrawCopy(
"same");
1521 c->cd(4); trueIPz ->DrawCopy(
"same");
1522 c->cd(5); simMeas ->DrawCopy(
"lego2 e");
1523 c->cd(6); simBg ->DrawCopy(
"lego2 e");
1524 c->cd(7); simSig ->DrawCopy(
"lego2 e");
1525 c->cd(8); alpha ->DrawCopy(
"lego2 e");
1526 c->cd(9); fiducial->DrawCopy(
"lego2 e");
1527 c->cd(10); trueGen ->DrawCopy(
"lego2 e");
1528 c->cd(11); result ->DrawCopy(
"lego2 e");
1529 c->cd(12); summary ->DrawClone(
"nostack");
1530 c->GetPad(4)->Modified();
1531 c->GetPad(12)->Modified();
1532 c->GetPad(4)->BuildLegend();
1533 c->GetPad(12)->BuildLegend();
1534 TLatex* ltx =
new TLatex(.52, .3,
1535 Form(
"d#it{N}_{ch}/d#eta|_{|#eta|<%3.1f}"
1536 "=%6.1f#pm%6.1f (#chi^{2}/#nu=%5.2f)",
1537 flim, r->Parameter(0),
1538 r->ParError(0), r->Chi2()/r->Ndf()));
1539 ltx->SetTextAlign(23);
1541 ltx->SetTextFont(42);
1542 ltx->SetTextSize(0.035);
1543 c->cd(12); ltx->Draw();
1544 ltx->DrawLatex(.52,.25,
1545 Form(
"%sk%s=%5.3f#pm%5.3f",
1552 stack->Add(dndeta,
"e2");
1560 if (histK) histK ->Write();
1564 TDirectory* det = dir->mkdir(
"details");
1582 if (histK) histK->Write();
1583 THStack* deltas =
new THStack(
"deltas",
"#Delta distributions");
1584 deltas->Add(
GetDelta(realList,
"Data",b,realNev,0, 0, lim));
1585 deltas->Add(
GetDelta(simList,
"Data",b,simNev, scale,scaleE,lim));
1586 deltas->Add(
GetDelta(realList,
"Inj", b,realNev,0, 0, lim));
1587 deltas->Add(
GetDelta(simList,
"Inj", b,simNev, scale,scaleE,lim));
1588 deltas->Add(
GetDelta(simList,
"Comb",b,simNev, scale,scaleE,lim));
1589 deltas->Add(
GetDelta(simList,
"Prim",b,trueNev,scale,scaleE,lim));
1590 deltas->Add(
GetDelta(simList,
"Sec", b,trueNev,scale,scaleE,lim));
1591 TH1* sum =
GetDelta(simList,
"Prim", b,trueNev,0, 0, lim);
1592 sum->Add(
GetDelta(simList,
"Sec", b,trueNev,0, 0, lim));
1593 sum->Add(
GetDelta(simList,
"Comb", b,simNev, 0, 0, lim));
1594 sum->SetTitle(Form(
"%5.3f#pm%5.3f #times (Prim+Sec+Comb)",scale,scaleE));
1595 Scale(sum,scale,scaleE);
1600 TDirectory* avs = dir->mkdir(
"averages");
1602 realAvgMeas->Write();
1603 realAvgSig ->Write();
1604 realAvgBg ->Write();
1605 simAvgMeas ->Write();
1606 simAvgSig ->Write();
1630 Printf(
"=== %-15s - %s ===", title, path);
1633 Warning(
"",
"No stats");
1636 TH1* stats =
static_cast<TH1*
>(tmp->Clone());
1637 stats->SetDirectory(0);
1638 TAxis* axis = stats->GetXaxis();
1639 for (
Int_t i = 1; i <= 4; i++) {
1640 Printf(
"%31s: %f", axis->GetBinLabel(i), stats->GetBinContent(i));
1642 Double_t scale = stats->GetBinContent(4);
1643 stats->Scale(1/scale);
1644 for (
Int_t i = 6; i <= 23; i++) {
1645 Printf(
"%31s: %f", axis->GetBinLabel(i), stats->GetBinContent(i));
1680 const char* realFileName=
"trdt.root",
1681 const char* simFileName=
"trmc.root",
1683 const char* otherFileName=
"")
1685 TFile* realFile = TFile::Open(realFileName,
"READ");
1686 TFile* simFile = TFile::Open(simFileName,
"READ");
1687 TFile* otherFile = 0;
1688 if (otherFileName && otherFileName[0] !=
'\0')
1689 otherFile = TFile::Open(otherFileName,
"READ");
1690 TList* realList =
static_cast<TList*
>(realFile->Get(
"clist"));
1691 TList* simList =
static_cast<TList*
>(simFile ->Get(
"clist"));
1694 otherList =
static_cast<TObjArray*
>(otherFile->Get(
"TObjArray"));
1696 ShowStats(realList,
"Real data", realFile->GetName());
1697 ShowStats(simList,
"Simulated data",simFile ->GetName());
1699 TH1* realCent =
GetH1(realList,
"EvCentr");
1700 TH1* simCent =
GetH1(simList,
"EvCentr");
1702 Warning(
"SimpleCorrection",
"Inconsistent centralities");
1705 TH2* realIpz =
GetH2(realList,
"b0_TrData_ZvEtaCutT");
1706 TH2* simIpz =
GetH2(simList,
"b0_TrData_ZvEtaCutT");
1708 Warning(
"SimpleCorrection",
"Inconsistent IPz or eta axis");
1714 Bool_t realComb = flags & 0x1;
1715 Bool_t simComb = flags & 0x2;
1716 Bool_t preScale = flags & 0x4;
1717 Bool_t full = flags & 0x8;
1718 Bool_t showOne = flags & 0x20;
1720 Printf(
"Real background: %s\n"
1721 "Simulated background: %s\n"
1723 "Full pre-scale: %s\n"
1724 "Show each bin: %s\n"
1725 "Delta upper limit: %g\n"
1727 (realComb ?
"comb" :
"inj"),
1728 (simComb ?
"comb" :
"inj"),
1729 (preScale ?
"yes" :
"no"),
1730 (full ?
"yes" :
"no"),
1731 (showOne ?
"yes" :
"no"),
1733 (scaleMode==0 ?
"fixed" :
1734 scaleMode==1 ?
"c" :
1735 scaleMode==2 ?
"eta-c" :
1739 if (flags & 0x10) realList = simList;
1741 tit.Form(
"Real BG: %s Simulated BG: %s%s",
1742 realComb ?
"MC-labels" :
"Injection",
1743 simComb ?
"MC-labels" :
"Injection",
1744 preScale ?
" (pre-scaled)" :
"");
1746 TFile* output = TFile::Open(Form(
"result_0x%x.root", flags&0x3),
1748 realCent->Write(
"cent");
1749 THStack* res =
new THStack(
"result", tit);
1750 TH1* cent =
GetH1(realList,
"EvCentr");
1751 for (
Int_t b = 0; b < nBins; b++) {
1752 Double_t c1 = cent->GetXaxis()->GetBinLowEdge(b+1);
1753 Double_t c2 = cent->GetXaxis()->GetBinUpEdge (b+1);
1755 name.Form(
"cent%03dd%02d_%03dd%02d",
1758 TDirectory*
dir = output->mkdir(name);
1759 CorrectOneBin(res,dir,b,realList,simList,realComb,simComb,preScale,full,
1760 showOne, scaleMode, deltaLim);
1763 TH1* od =
GetH1(otherList, Form(
"bin%d_DataCorrSignal_Bla",b));
1764 TH1* ot =
GetH1(otherList, Form(
"bin%d_MCTruth",b));
1765 TH1* td =
static_cast<TH1*
>(od->Clone(
"otherdNdeta"));
1766 TH1* tt =
static_cast<TH1*
>(ot->Clone(
"otherTruth"));
1767 td->SetDirectory(0);
1768 tt->SetDirectory(0);
1769 SetAttr(td,
cc[b%11], 21, 1.2, 3004, 1, 1);
1770 SetAttr(tt,
cc[b%11], 25, 1.4, 0, 7, 3);
1771 if (td->GetMinimum() < 1e-6) { td->SetMinimum(); td->SetMaximum(); }
1772 if (tt->GetMinimum() < 1e-6) { tt->SetMinimum(); tt->SetMaximum(); }
1782 res->SetMaximum(1.2*res->GetMaximum(
"nostack"));
Bool_t CheckAxisNBins(const char *which, const TAxis *a1, const TAxis *a2)
Int_t color[]
print message on plot with ok/not ok
Bool_t CheckAxisLabels(const char *which, const TAxis *a1, const TAxis *a2)
const Bool_t kSimpleCorrectLoaded
Bool_t CheckAxisBins(const char *which, const TAxis *a1, const TAxis *a2)
TH1 * Avg2(TH2 *h, TH1 *ipz, UShort_t mode, const char *name)
TH1 * GetH1(TSeqCollection *l, const char *name)
TH1 * SetAttr(TH1 *h, Color_t color, Style_t marker=20, Double_t size=1., Style_t fill=0, Style_t line=1, Width_t width=1)
void ShowStats(TList *list, const char *title, const char *path)
Double_t Integrate(TH1 *h, Double_t min, Double_t max, Double_t &err)
TH2 * GetH2(TSeqCollection *l, const char *name)
TH2 * CorrectSignal(Bool_t comb, TH2 *meas, TH2 *bg, TH2 *simMeas)
Bool_t CheckConsistency(const TH1 *h1, const TH1 *h2)
TH1 * GetDeltaDist(TH2 *deltaReal, TH2 *deltaSim, Double_t nEvReal, Double_t nEvSim, Double_t lim=1e9)
TH1 * Avg(TH2 *h, TH1 *ipz, UShort_t mode, const char *name, TH2 *other=0)
TH2 * CopyH2(TSeqCollection *l, const char *name, const char *newName=0)
Double_t RatioE(Double_t n, Double_t en, Double_t d, Double_t ed, Double_t &er)
TH2 * GetBg(Bool_t comb, TSeqCollection *list, Int_t b, Double_t nEv=1, Double_t lim=1e9)
TObject * CopyO(TSeqCollection *l, const char *name, const char *newName, TClass *cls)
Bool_t CheckAxisLimits(const char *which, const TAxis *a1, const TAxis *a2)
void PrintH(TH2 *h, Int_t prec=2)
void SimpleCorrect(UShort_t flags=0x3, UShort_t scaleMode=1, const char *realFileName="trdt.root", const char *simFileName="trmc.root", Int_t nBins=9, const char *otherFileName="")
TH2 * Scale(TH2 *h, TH1 *g)
TObject * GetO(TSeqCollection *l, const char *name, TClass *cls)
TH1 * CopyH1(TSeqCollection *l, const char *name, const char *newName=0)
Bool_t CheckAxis(const char *which, const TAxis *a1, const TAxis *a2, Bool_t alsoLbls)
Double_t GetDeltaScale(TH1 *deltaRec, TH1 *deltaBg, Double_t &eR, Double_t nEvRec=1, Double_t nEvBg=1, Double_t lim=1e9)
TH1 * Coerce(const TH1 *templ, TH1 *src)
TH2 * CorrectForIPz(TH2 *in, TH1 *ipz, Bool_t full=true)
void CorrectOneBin(THStack *stack, TDirectory *dir, Int_t b, TSeqCollection *realList, TSeqCollection *simList, Bool_t realComb, Bool_t simComb, Bool_t preScale=false, Bool_t full=false, Bool_t showOne=false, UShort_t scaleMode=2, Double_t lim=1e9)
TH1 * GetDelta(TSeqCollection *l, const char *which, Int_t b, Int_t nEv, Double_t scale=0, Double_t scaleE=0, Double_t lim=1e9)