7 #include "TObjString.h"
16 int FileROOTVersion = file->GetVersion();
17 int MainFileROOTVersion = FileROOTVersion;
21 while (MainFileROOTVersion >= 10)
22 MainFileROOTVersion /= 10;
24 std::string SystemROOTVersion = std::string(ROOT_RELEASE);
25 int MainSystemROOTVersion = SystemROOTVersion.at(0) -
'0';
27 if(MainFileROOTVersion != MainSystemROOTVersion)
29 MACH3LOG_ERROR(
"File was produced with: {} ROOT version", FileROOTVersion);
30 MACH3LOG_ERROR(
"Found: {} ROOT version in the system", SystemROOTVersion);
42 for(
int iOverflow = -1; iOverflow > -10; iOverflow--)
44 overflow += poly->GetBinContent(iOverflow);
56 for(
int i = 1; i < poly->GetNumberOfBins()+1; i++)
58 integral += poly->GetBinContent(i);
66 TH1D*
PolyProjectionX(TObject* poly,
const std::string& TempName,
const std::vector<double>& xbins,
const bool computeErrors) {
68 TH1D* hProjX =
new TH1D((TempName+
"_x").c_str(),(TempName+
"_x").c_str(),
int(xbins.size()-1), &xbins[0]);
71 std::vector<double> hProjX_Error(hProjX->GetXaxis()->GetNbins() + 1, 0.0);
72 double xlow, xup, frac = 0.;
75 for (
int i = 0; i < static_cast<TH2Poly*>(poly)->GetNumberOfBins(); i++)
78 TH2PolyBin* bin =
static_cast<TH2PolyBin*
>(
static_cast<TH2Poly*
>(poly)->GetBins()->At(i));
79 xlow = bin->GetXMin();
83 for(
int dx = 0; dx < int(xbins.size()); dx++)
85 if(xbins[dx+1] <= xlow || xbins[dx] >= xup)
89 else if(xbins[dx] <= xlow && xbins[dx+1] >= xup)
93 else if(xbins[dx] <= xlow && xbins[dx+1] <= xup)
95 frac = (xbins[dx+1]-xlow)/(xup-xlow);
97 else if(xbins[dx] >= xlow && xbins[dx+1] >= xup)
99 frac = (xup-xbins[dx])/(xup-xlow);
101 else if(xbins[dx] >= xlow && xbins[dx+1] <= xup)
103 frac = (xbins[dx+1]-xbins[dx])/(xup-xlow);
109 hProjX->SetBinContent(dx+1,hProjX->GetBinContent(dx+1)+frac*bin->GetContent());
115 const double Temp_Err = frac*
static_cast<TH2Poly*
>(poly)->GetBinError(i+1)*frac*
static_cast<TH2Poly*
>(poly)->GetBinError(i+1);
116 hProjX_Error[dx+1] += Temp_Err;
123 for (
int i = 1; i <= hProjX->GetXaxis()->GetNbins(); ++i)
125 const double Error = std::sqrt(hProjX_Error[i]);
126 hProjX->SetBinError(i, Error);
134 TH1D*
PolyProjectionY(TObject* poly,
const std::string& TempName,
const std::vector<double>& ybins,
const bool computeErrors) {
136 TH1D* hProjY =
new TH1D((TempName+
"_y").c_str(),(TempName+
"_y").c_str(),
int(ybins.size()-1),&ybins[0]);
138 std::vector<double> hProjY_Error(hProjY->GetXaxis()->GetNbins() + 1, 0.0);
139 double ylow, yup, frac = 0.;
142 for (
int i = 0; i < static_cast<TH2Poly*>(poly)->GetNumberOfBins(); i++)
145 TH2PolyBin* bin =
static_cast<TH2PolyBin*
>(
static_cast<TH2Poly*
>(poly)->GetBins()->At(i));
146 ylow = bin->GetYMin();
147 yup = bin->GetYMax();
150 for(
int dy = 0; dy < int(ybins.size()); dy++)
152 if(ybins[dy+1]<=ylow || ybins[dy] >= yup)
156 else if(ybins[dy] <= ylow && ybins[dy+1] >= yup)
160 else if(ybins[dy] <= ylow && ybins[dy+1] <= yup)
162 frac = (ybins[dy+1]-ylow)/(yup-ylow);
164 else if(ybins[dy] >= ylow && ybins[dy+1] >= yup)
166 frac = (yup-ybins[dy])/(yup-ylow);
168 else if(ybins[dy] >= ylow && ybins[dy+1] <= yup)
170 frac = (ybins[dy+1]-ybins[dy])/(yup-ylow);
176 hProjY->SetBinContent(dy+1,hProjY->GetBinContent(dy+1)+frac*bin->GetContent());
182 const double Temp_Err = frac*
static_cast<TH2Poly*
>(poly)->GetBinError(i+1)*frac*
static_cast<TH2Poly*
>(poly)->GetBinError(i+1);
183 hProjY_Error[dy+1] += Temp_Err;
190 for (
int i = 1; i <= hProjY->GetXaxis()->GetNbins(); ++i)
192 const double Error = std::sqrt(hProjY_Error[i]);
193 hProjY->SetBinError(i, Error);
203 TH2Poly* HistCopy =
static_cast<TH2Poly*
>(Histogram->Clone());
206 std::string title = std::string(HistCopy->GetName())+
"_norm";
207 HistCopy->SetNameTitle(title.c_str(), title.c_str());
217 for(
int j = 1; j < Histogram->GetNumberOfBins()+1; j++)
219 Histogram->SetBinContent(j, Histogram->GetBinContent(j)/Integral);
225 template<
class HistType>
226 HistType*
RatioHists(HistType *NumHist, HistType *DenomHist) {
228 HistType *NumCopy =
static_cast<HistType*
>(NumHist->Clone());
229 std::string title = std::string(DenomHist->GetName()) +
"_ratio";
230 NumCopy->SetNameTitle(title.c_str(), title.c_str());
231 NumCopy->Divide(DenomHist);
240 TH2Poly *NumCopy =
static_cast<TH2Poly*
>(NumHist->Clone());
241 std::string title = std::string(DenomHist->GetName()) +
"_ratio";
242 NumCopy->SetNameTitle(title.c_str(), title.c_str());
244 for(
int i = 1; i < NumCopy->GetNumberOfBins()+1; ++i) {
245 NumCopy->SetBinContent(i,NumHist->GetBinContent(i)/DenomHist->GetBinContent(i));
253 double xlow, xup, ylow, yup;
254 std::string HistTempName = poly->GetName();
258 TH2D *hist =
static_cast<TH2D*
>(h2dhist->Clone());
259 hist->SetNameTitle(HistTempName.c_str(), HistTempName.c_str());
261 for(
int ix = 0; ix < hist->GetNbinsX() + 2; ix++) {
262 for(
int iy = 0; iy < hist->GetNbinsY() + 2; iy++) {
263 hist->SetBinContent(ix, iy, 0);
267 for(
int i = 0; i< poly->GetNumberOfBins(); i++){
268 TH2PolyBin & polybin =
static_cast<TH2PolyBin &
>(*poly->GetBins()->At(i));
269 xlow = polybin.GetXMin();
270 xup = polybin.GetXMax();
271 ylow = polybin.GetYMin();
272 yup = polybin.GetYMax();
275 xbin = hist->GetXaxis()->FindBin(xlow+(xup-xlow)/2);
276 ybin = hist->GetYaxis()->FindBin(ylow+(yup-ylow)/2);
278 MACH3LOG_TRACE(
"Poly bin {}, xlow: {}, xup: {}, ylow: {}, yup: {}. Finding bin for ({}, {}). Found Bin ({}, {}) with content {}. But Poly content: {}",
279 i, xlow, xup, ylow, yup, (xlow + (xup - xlow) / 2), (ylow + (yup - ylow) / 2), xbin, ybin, polybin.GetContent(), poly->GetBinContent(i));
280 hist->SetBinContent(xbin, ybin, polybin.GetContent());
288 TAxis* xaxis = hist->GetXaxis();
290 TAxis* yaxis = hist->GetYaxis();
292 TString histname = hist->GetName();
294 TH2Poly* poly =
new TH2Poly();
295 poly->SetName(histname);
296 poly->SetTitle(histname);
299 poly->GetXaxis()->SetTitle(xaxis->GetTitle());
300 poly->GetYaxis()->SetTitle(yaxis->GetTitle());
302 double xmax, xmin, ymax, ymin;
303 for (
int iy = 1; iy <= yaxis->GetNbins(); iy++) {
304 ymax = yaxis->GetBinUpEdge(iy);
305 ymin = yaxis->GetBinLowEdge(iy);
306 for (
int ix = 1; ix <= xaxis->GetNbins(); ix++) {
307 xmax = xaxis->GetBinUpEdge(ix);
308 xmin = xaxis->GetBinLowEdge(ix);
309 double binofx[] = {xmin, xmax, xmax, xmin};
310 double binofy[] = {ymin, ymin, ymax, ymax};
311 poly->AddBin(4, binofx, binofy);
315 for (
int iy = 1; iy <= yaxis->GetNbins(); iy++) {
316 ymax = yaxis->GetBinUpEdge(iy);
317 ymin = yaxis->GetBinLowEdge(iy);
318 for (
int ix = 1; ix <= xaxis->GetNbins(); ix++) {
319 xmax = xaxis->GetBinUpEdge(ix);
320 xmin = xaxis->GetBinLowEdge(ix);
323 int bin = hist->GetBin(ix, iy);
324 double content = hist->GetBinContent(bin);
325 poly->SetBinContent(poly->FindBin((xmin + xmax) / 2, (ymin + ymax) / 2), content);
336 TH2Poly* HistCopy =
static_cast<TH2Poly*
>(Histogram->Clone());
337 double xlow, xup, ylow, yup, area;
339 for(
int i = 1; i < HistCopy->GetNumberOfBins()+1; i++)
341 TH2PolyBin* bin =
static_cast<TH2PolyBin*
>(Histogram->GetBins()->At(i - 1));
342 xlow = bin->GetXMin();
343 xup = bin->GetXMax();
344 ylow = bin->GetYMin();
345 yup = bin->GetYMax();
346 area = (xup-xlow)*(yup-ylow);
347 HistCopy->SetBinContent(i, Histogram->GetBinContent(i)/(area*scale));
357 double xlow, xup, ylow, yup, area;
359 for(
int i = 1; i < Histogram->GetNumberOfBins()+1; i++)
361 TH2PolyBin* bin =
static_cast<TH2PolyBin*
>(Histogram->GetBins()->At(i - 1));
362 xlow = bin->GetXMin();
363 xup = bin->GetXMax();
364 ylow = bin->GetYMin();
365 yup = bin->GetYMax();
366 area = (xup-xlow)*(yup-ylow);
367 integral += Histogram->GetBinContent(i)*area;
373 TH2Poly*
MakePolyHist(
const std::string& name,
const std::vector<double>& BinArray_x,
const std::vector<double>& BinArray_y) {
375 TH2Poly* poly =
new TH2Poly();
376 poly->SetName(name.c_str());
377 poly->SetTitle(name.c_str());
378 double xmax, xmin, ymax, ymin;
379 for(
unsigned int iy = 0; iy < BinArray_y.size()-1; iy++)
381 ymax = BinArray_y[iy+1];
382 ymin = BinArray_y[iy];
383 for(
unsigned int ix = 0; ix < BinArray_x.size()-1; ix++)
385 xmax = BinArray_x[ix+1];
386 xmin = BinArray_x[ix];
387 double binofx[] = {xmin, xmax, xmax, xmin};
388 double binofy[] = {ymin, ymin, ymax, ymax};
389 poly->AddBin(4,binofx,binofy);
399 TList *listOfFunctions = hist->GetListOfFunctions();
400 TF1 *fitter =
dynamic_cast<TF1*
>(listOfFunctions->FindObject(name.c_str()));
402 listOfFunctions->Remove(fitter);
411 FluctHist->Reset(
"");
412 FluctHist->Fill(0.0, 0.0);
414 for (
int i = 1; i <= PolyHist->GetXaxis()->GetNbins(); ++i)
417 const double MeanContent = PolyHist->GetBinContent(i);
419 const double Random = rand->PoissonD(MeanContent);
421 FluctHist->SetBinContent(i,Random);
430 FluctHist->Reset(
"");
431 FluctHist->Fill(0.0, 0.0, 0.0);
433 for (
int i = 1; i < FluctHist->GetNumberOfBins()+1; ++i)
436 const double MeanContent = PolyHist->GetBinContent(i);
438 const double Random = rand->PoissonD(MeanContent);
440 FluctHist->SetBinContent(i,Random);
449 FluctHist->Reset(
"");
450 FluctHist->Fill(0.0, 0.0, 0.0);
453 const int nBinsX = Hist->GetXaxis()->GetNbins();
454 const int nBinsY = Hist->GetYaxis()->GetNbins();
455 for (
int ix = 1; ix <= nBinsX; ++ix) {
456 for (
int iy = 1; iy <= nBinsY; ++iy) {
458 const double MeanContent = Hist->GetBinContent(ix, iy);
460 const double Random = rand->PoissonD(MeanContent);
462 FluctHist->SetBinContent(ix, iy, Random);
472 FluctHist->Reset(
"");
473 FluctHist->Fill(0.0, 0.0);
475 const double evrate = PolyHist->Integral();
476 #pragma GCC diagnostic push
477 #pragma GCC diagnostic ignored "-Wconversion"
478 const int num = rand->Poisson(evrate);
479 #pragma GCC diagnostic pop
483 const double candidate = PolyHist->GetRandom();
484 FluctHist->Fill(candidate);
494 FluctHist->Fill(0.0, 0.0, 0.0);
496 const double evrate = PolyHist->Integral();
497 #pragma GCC diagnostic push
498 #pragma GCC diagnostic ignored "-Wconversion"
499 const int num = rand->Poisson(evrate);
500 #pragma GCC diagnostic pop
503 for (
int count = 0; count < num; ++count) {
504 PolyHist->GetRandom2(x, y);
505 FluctHist->Fill(x, y);
513 const int nbins = PolyHist->GetNumberOfBins();
514 const double r1 = rand->Rndm();
516 double* fIntegral =
new double[nbins+2];
520 for (
int i = 1; i < nbins+1; ++i)
523 const double content = PolyHist->GetBinContent(i);
524 fIntegral[i] += fIntegral[i - 1] + content;
526 for (Int_t bin = 1; bin < nbins+1; ++bin) fIntegral[bin] /= fIntegral[nbins];
527 fIntegral[nbins+1] = PolyHist->GetEntries();
530 int iBin = int(TMath::BinarySearch(nbins, fIntegral, r1));
543 FluctHist->Reset(
"");
544 FluctHist->Fill(0.0, 0.0, 0.0);
547 #pragma GCC diagnostic push
548 #pragma GCC diagnostic ignored "-Wconversion"
549 const int num = rand->Poisson(evrate);
550 #pragma GCC diagnostic pop
555 FluctHist->SetBinContent(iBin, FluctHist->GetBinContent(iBin) + 1);
564 for (
int x = 0; x < violin->GetXaxis()->GetNbins(); ++x)
566 const int y = violin->GetYaxis()->FindBin(hist_1d->GetBinContent(x+1));
567 violin->SetBinContent(x+1, y, violin->GetBinContent(x+1, y)+1);
576 constexpr
double refractiveIndex = 1.334;
578 double momentumThreshold = mass/sqrt(refractiveIndex*refractiveIndex-1.);
579 return momentumThreshold;
584 double CalculateQ2(
double PLep,
double PUpd,
double EnuTrue,
double InitialQ2){
586 constexpr
double MLep = 0.10565837;
587 constexpr
double MLep2 = MLep * MLep;
590 double ELep = sqrt((MLep2)+(PLep*PLep));
592 double CosTh = (2*EnuTrue*ELep - MLep2 - InitialQ2)/(2*EnuTrue*PLep);
594 ELep = sqrt((MLep2)+(PUpd*PUpd));
597 double Q2Upd = -(MLep2) + 2.0*EnuTrue*(ELep - PUpd*CosTh);
599 return Q2Upd - InitialQ2;
604 double CalculateEnu(
double PLep,
double costh,
double Eb,
bool neutrino){
606 double mNeff = 0.93956536 - Eb / 1000.;
607 double mNoth = 0.93827203;
610 mNeff = 0.93827203 - Eb / 1000.;
614 constexpr
double mLep = 0.10565837;
615 constexpr
double mLep2 = mLep * mLep;
617 const double eLep = sqrt(PLep * PLep + mLep2);
618 double Enu = (2 * mNeff * eLep - mLep2 + mNoth * mNoth - mNeff * mNeff) /(2 * (mNeff - eLep + PLep * costh));
625 const std::string& name) {
627 const int nBinsX = Spectra->GetNbinsX();
629 std::vector<double> edges(nBinsX + 1);
630 for (
int i = 0; i < nBinsX; ++i) {
631 edges[i] = Spectra->GetXaxis()->GetBinLowEdge(i+1);
633 edges[nBinsX] = Spectra->GetXaxis()->GetBinUpEdge(nBinsX);
635 auto h1 = std::make_unique<TH1D>(name.c_str(), name.c_str(), nBinsX, edges.data());
636 h1->SetDirectory(
nullptr);
640 for (
int ix = 1; ix <= nBinsX; ++ix) {
642 std::unique_ptr<TH1D> slice(Spectra->ProjectionY(Form(
"_py_%d", ix), ix, ix));
643 slice->SetDirectory(
nullptr);
644 if (slice->GetEntries() == 0) {
645 h1->SetBinContent(ix, 0.0);
646 h1->SetBinError(ix, 0.0);
650 const double mean = slice->GetMean();
651 const double rms = slice->GetRMS();
653 h1->SetBinContent(ix, mean);
654 h1->SetBinError(ix, rms);
657 h1->GetXaxis()->SetTitle(Spectra->GetXaxis()->GetTitle());
658 h1->GetYaxis()->SetTitle(
"Events");
667 YAML::Node bin_range_specifier;
668 if (bin_range[
"linspace"]) {
669 bin_range_specifier = bin_range[
"linspace"];
670 }
else if (bin_range[
"logspace"]) {
672 bin_range_specifier = bin_range[
"logspace"];
674 std::stringstream ss;
676 MACH3LOG_ERROR(
"When parsing binning, expected bin range specifier with "
677 "key linspace or logspace, but found,\n{}",
682 auto nb = Get<int>(bin_range_specifier[
"nb"], __FILE__, __LINE__);
683 auto low = Get<double>(bin_range_specifier[
"low"], __FILE__, __LINE__);
684 auto up = Get<double>(bin_range_specifier[
"up"], __FILE__, __LINE__);
686 std::vector<double> edges(nb + 1, low);
694 double bw = (up - low) / nb;
695 for (
int i = 0; i < (nb - 1); ++i) {
696 edges[i + 1] = edges[i] + bw;
699 double llow = std::log10(low);
700 double lup = std::log10(up);
701 double lbw = (lup - llow) / nb;
702 for (
int i = 0; i < (nb - 1); ++i) {
703 edges[i + 1] = std::pow(10, llow + (i + 1) * lbw);
719 bool &found_range_specifier) {
721 if (bin_edges_node.IsMap()) {
722 found_range_specifier =
true;
725 std::vector<double> edges_builder;
726 if (!bin_edges_node.IsSequence()) {
727 std::stringstream ss;
728 ss << bin_edges_node;
730 "When parsing binning, expected to find a YAML map or sequence, "
735 for (
auto const &it : bin_edges_node) {
737 edges_builder.push_back(it.as<
double>());
738 }
else if (it.IsMap()) {
739 found_range_specifier =
true;
741 std::copy(range_edges.begin(), range_edges.end(),
742 std::back_inserter(edges_builder));
744 std::stringstream ss;
745 ss << bin_edges_node;
747 "When parsing binning, expected elements in outer sequence to all be "
748 "either scalars or maps, but found:\n{}",
755 std::vector<double> edges;
756 for (
size_t eb_it = 0; eb_it < edges_builder.size(); ++eb_it) {
758 if (edges_builder[eb_it] == edges.back()) {
760 }
else if (edges_builder[eb_it] < edges.back()) {
761 std::stringstream ss;
763 for(
auto const & e : edges_builder){
764 ss << fmt::format(
"{:.3g} ", e);
768 "When parsing binning, found edges that were not monotonically "
769 "increasing, problem bin at index: {}:\n{}",
774 edges.push_back(edges_builder[eb_it]);
782 TFile*
Open(
const std::string& Name,
const std::string& Type,
const std::string& File,
const int Line) {
784 TFile* OutFile =
new TFile(Name.c_str(), Type.c_str());
787 if (OutFile->IsZombie()) {
789 std::string lowerType = Type;
790 std::transform(lowerType.begin(), lowerType.end(), lowerType.begin(), ::tolower);
791 if (lowerType ==
"recreate") {
803 for (
int j = 0; j <= Sample_Hist->GetNbinsX(); ++j)
805 double num = Sample_Hist->GetBinContent(j);
806 double numErr = Sample_Hist->GetBinError(j);
807 double den = Sample_Hist->GetBinWidth(j);
809 double valueErr = 0.;
812 value = num/(den/scale);
813 valueErr = numErr/(den/scale);
814 Sample_Hist->SetBinContent(j,value);
815 Sample_Hist->SetBinError(j,valueErr);
821 void CheckBinningMatch(
const TH1D* Hist1,
const TH1D* Hist2,
const std::string& File,
const int Line) {
823 if (Hist1->GetNbinsX() != Hist2->GetNbinsX()) {
824 MACH3LOG_ERROR(
"Number of bins does not match for TH1D: {} vs {}", Hist1->GetNbinsX(), Hist2->GetNbinsX());
827 for (
int i = 1; i <= Hist1->GetNbinsX(); ++i) {
828 if (std::fabs(Hist1->GetXaxis()->GetBinLowEdge(i) - Hist2->GetXaxis()->GetBinLowEdge(i)) > 0.001 ||
829 std::fabs(Hist1->GetXaxis()->GetBinUpEdge(i) - Hist2->GetXaxis()->GetBinUpEdge(i)) > 0.001) {
837 void CheckBinningMatch(
const TH2D* Hist1,
const TH2D* Hist2,
const std::string& File,
const int Line) {
839 if (Hist1->GetNbinsX() != Hist2->GetNbinsX() || Hist1->GetNbinsY() != Hist2->GetNbinsY()) {
844 for (
int i = 1; i <= Hist1->GetNbinsX(); ++i) {
845 if (std::fabs(Hist1->GetXaxis()->GetBinLowEdge(i) - Hist2->GetXaxis()->GetBinLowEdge(i)) > 0.001 ||
846 std::fabs(Hist1->GetXaxis()->GetBinUpEdge(i) - Hist2->GetXaxis()->GetBinUpEdge(i)) > 0.001) {
851 for (
int j = 1; j <= Hist1->GetNbinsY(); ++j) {
852 if (std::fabs(Hist1->GetYaxis()->GetBinLowEdge(j) - Hist2->GetYaxis()->GetBinLowEdge(j)) > 0.001 ||
853 std::fabs(Hist1->GetYaxis()->GetBinUpEdge(j) - Hist2->GetYaxis()->GetBinUpEdge(j)) > 0.001) {
862 void CheckBinningMatch(TH2Poly* Hist1, TH2Poly* Hist2,
const std::string& File,
const int Line) {
864 int NBins1 = Hist1->GetNumberOfBins();
865 int NBins2 = Hist2->GetNumberOfBins();
866 if (NBins1 != NBins2) {
867 MACH3LOG_ERROR(
"Number of bins does not match for TH2Poly: {} vs {}", NBins1, NBins2);
870 for (
int j = 1; j <= NBins1; j++)
873 TH2PolyBin* polybin1 =
static_cast<TH2PolyBin*
>(Hist1->GetBins()->At(j - 1));
874 TH2PolyBin* polybin2 =
static_cast<TH2PolyBin*
>(Hist2->GetBins()->At(j - 1));
876 if( std::fabs(polybin2->GetXMin() - polybin1->GetXMin()) > 0.001 ||
877 std::fabs(polybin2->GetXMax() - polybin1->GetXMax()) > 0.001 ||
878 std::fabs(polybin2->GetYMin() - polybin1->GetYMin()) > 0.001 ||
879 std::fabs(polybin2->GetYMax() - polybin1->GetYMax()) > 0.001 )
882 MACH3LOG_ERROR(
"data x min {} x max {} y min {} y max {}", polybin2->GetXMin(), polybin2->GetXMax(), polybin2->GetYMin(), polybin2->GetYMax());
883 MACH3LOG_ERROR(
"mc x min {} x max {} y min {} y max {}", polybin1->GetXMin(), polybin1->GetXMax(), polybin1->GetYMin(), polybin1->GetYMax());
891 YAML::Node
PolyToYaml(TH2Poly* Hist,
const std::string& YamlName,
const std::string& File,
const int Line) {
898 YAML::Node bins(YAML::NodeType::Sequence);
899 bins.SetStyle(YAML::EmitterStyle::Flow);
900 const int NBins = Hist->GetNumberOfBins();
902 for (
int j = 1; j <= NBins; j++)
904 TH2PolyBin* polybin =
static_cast<TH2PolyBin*
>(Hist->GetBins()->At(j - 1));
906 double xmin = polybin->GetXMin();
907 double xmax = polybin->GetXMax();
908 double ymin = polybin->GetYMin();
909 double ymax = polybin->GetYMax();
911 YAML::Node xNode(YAML::NodeType::Sequence);
912 xNode.SetStyle(YAML::EmitterStyle::Flow);
913 xNode.push_back(xmin);
914 xNode.push_back(xmax);
916 YAML::Node yNode(YAML::NodeType::Sequence);
917 yNode.SetStyle(YAML::EmitterStyle::Flow);
918 yNode.push_back(ymin);
919 yNode.push_back(ymax);
921 YAML::Node bin(YAML::NodeType::Sequence);
922 bin.SetStyle(YAML::EmitterStyle::Flow);
923 bin.push_back(xNode);
924 bin.push_back(yNode);
930 result[YamlName] = bins;
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
#define _MaCh3_Safe_Include_End_
void MakeFluctuatedHistogramAlternative(TH1D *FluctHist, TH1D *PolyHist, TRandom3 *rand)
Make Poisson fluctuation of TH1D hist using slow method which is only for cross-check.
TH2D * ConvertTH2PolyToTH2D(TH2Poly *poly, TH2D *h2dhist)
KS: Convert TH2D to TH2Poly.
void NormaliseTH2Poly(TH2Poly *Histogram)
Helper to Normalise histograms.
double returnCherenkovThresholdMomentum(const int PDG)
DB Get the Cherenkov momentum threshold in MeV.
double OverflowIntegral(TH2Poly *poly)
WP: Helper function for calculating unbinned Integral of TH2Poly i.e including overflow.
void MakeFluctuatedHistogramStandard(TH1D *FluctHist, TH1D *PolyHist, TRandom3 *rand)
Make Poisson fluctuation of TH1D hist using default fast method.
TH1D * PolyProjectionX(TObject *poly, const std::string &TempName, const std::vector< double > &xbins, const bool computeErrors)
WP: Poly Projectors.
int GetRandomPoly2(const TH2Poly *PolyHist, TRandom3 *rand)
KS: ROOT developers were too lazy do develop getRanom2 for TH2Poly, this implementation is based on l...
std::vector< double > BuildBinEdgesFromNode(YAML::Node const &bin_edges_node, bool &found_range_specifier)
Builds a single dimension's bin edges from YAML::Node.
void FastViolinFill(TH2D *violin, TH1D *hist_1d)
KS: Fill Violin histogram with entry from a toy.
TH2Poly * NormalisePoly(TH2Poly *Histogram)
WP: Helper to Normalise histograms.
HistType * RatioHists(HistType *NumHist, HistType *DenomHist)
Helper to make ratio histograms.
std::vector< double > BinRangeToBinEdges(YAML::Node const &bin_range)
Converts a range (linspace/logspace) to a std::vector<double>
TH2Poly * RatioPolys(TH2Poly *NumHist, TH2Poly *DenomHist)
Helper to make ratio of TH2Polys.
TH2Poly * PolyScaleWidth(TH2Poly *Histogram, double scale)
WP: Helper to scale th2poly analogous to th2d scale with option "width".
void RemoveFitter(TH1D *hist, const std::string &name)
KS: Remove fitted TF1 from hist to make comparison easier.
double CalculateEnu(double PLep, double costh, double Eb, bool neutrino)
Recalculate Enu after Eb shift. Takes in shifted lepton momentum, lepton angle, and binding energy ch...
double PolyIntegralWidth(TH2Poly *Histogram)
WP: Helper to calc integral of th2poly analogous to th2d integra; with option "width".
_MaCh3_Safe_Include_Start_ _MaCh3_Safe_Include_End_ void CheckTH2PolyFileVersion(TFile *file)
KS: ROOT changes something with binning when moving from ROOT 5 to ROOT 6. If you open ROOT5 produced...
double NoOverflowIntegral(TH2Poly *poly)
WP: Helper function for calculating binned Integral of TH2Poly i.e not including overflow.
TH2Poly * ConvertTH2DtoTH2Poly(TH2D *hist)
KS: Convert TH2Poly to TH2D.
TH2Poly * MakePolyHist(const std::string &name, const std::vector< double > &BinArray_x, const std::vector< double > &BinArray_y)
WP: Helper function to create TH2Poly histogram with uniform binning.
TH1D * PolyProjectionY(TObject *poly, const std::string &TempName, const std::vector< double > &ybins, const bool computeErrors)
WP: Poly Projectors.
std::unique_ptr< TH1D > MakeSummaryFromSpectra(const TH2D *Spectra, const std::string &name)
Build a 1D posterior-predictive summary from a violin spectrum.
double CalculateQ2(double PLep, double PUpd, double EnuTrue, double InitialQ2)
Recalculate Q^2 after Eb shift. Takes in shifted lepton momentum, lepton angle, and true neutrino ene...
Custom exception class used throughout MaCh3.
double GetMassFromPDG(const int PDG)
Return mass for given PDG.
Main namespace for MaCh3 software.
void CheckBinningMatch(const TH1D *Hist1, const TH1D *Hist2, const std::string &File, const int Line)
KS: Helper function check if data and MC binning matches.
void ScaleHistogram(TH1 *Sample_Hist, const double scale)
Scale histogram to get divided by bin width.
TFile * Open(const std::string &Name, const std::string &Type, const std::string &File, const int Line)
Opens a ROOT file with the given name and mode.
YAML::Node PolyToYaml(TH2Poly *Hist, const std::string &YamlName, const std::string &File, const int Line)
KS: Convert TH2Poly into yaml config accepted by MaCh3.