17 #include "TMatrixDSym.h"
23 #include "TObjString.h"
39 std::vector <double>& Prior,
40 std::vector <bool>&
isFlat,
41 std::vector<std::string>& ParamNames,
42 std::vector<std::vector<double>>& invCovMatrix,
46 TFile *TempFile =
M3::Open(inputFile,
"open", __FILE__, __LINE__);
49 TDirectory* CovarianceFolder = TempFile->Get<TDirectory>(
"CovarianceFolder");
53 TMacro *Config = TempFile->Get<TMacro>(
"MaCh3_Config");
54 if (Config ==
nullptr) {
55 MACH3LOG_ERROR(
"Didn't find MaCh3_Config tree in MCMC file! {}", inputFile);
63 std::vector<std::string> CovPos = GetFromManager<std::vector<std::string>>(Settings[
"General"][
"Systematics"][
"XsecCovFile"], {
"none"});
64 if(CovPos.back() ==
"none")
72 if (std::getenv(
"MACH3") !=
nullptr) {
73 MACH3LOG_INFO(
"Found MACH3 environment variable: {}", std::getenv(
"MACH3"));
74 for(
unsigned int i = 0; i < CovPos.size(); i++)
75 CovPos[i].insert(0, std::string(std::getenv(
"MACH3"))+
"/");
79 CovFile[
"Systematics"] = YAML::Node(YAML::NodeType::Sequence);
80 for(
unsigned int i = 0; i < CovPos.size(); i++)
83 for (
const auto& item : YAMLDocTemp[
"Systematics"]) {
84 CovFile[
"Systematics"].push_back(item);
88 nParams = CovMatrix->GetNrows();
90 auto systematics = CovFile[
"Systematics"];
91 for (
auto it = systematics.begin(); it != systematics.end(); ++it)
93 auto const ¶m = *it;
95 ParamNames.push_back(param[
"Systematic"][
"Names"][
"FancyName"].as<std::string>());
96 Prior.push_back( param[
"Systematic"][
"ParameterValues"][
"PreFitValue"].as<double>() );
99 if (param[
"Systematic"][
"FlatPrior"]) { flat = param[
"Systematic"][
"FlatPrior"].as<
bool>(); }
105 invCovMatrix.resize(nParams, std::vector<double>(nParams, -999));
108 #pragma omp parallel for collapse(2)
110 for (
int i = 0; i < nParams; i++)
112 for (
int j = 0; j < nParams; ++j)
114 invCovMatrix[i][j] = (*CovMatrix)(i,j);
123 std::vector<std::string>& SetsNames,
124 std::vector<std::string>& FancyTitle,
125 std::vector<std::vector<bool>>& isRelevantParam,
126 const std::vector<std::string>& ParamNames,
129 std::vector<std::string> node = Settings[
"GetPenaltyTerm"][
"PenaltySets"].as<std::vector<std::string>>();
130 std::vector<std::vector<std::string>> RemoveNames;
131 std::vector<bool> Exclude;
133 for (
unsigned int i = 0; i < node.size(); i++)
135 std::string ParName = node[i];
136 SetsNames.push_back(ParName);
138 const auto& Set = Settings[
"GetPenaltyTerm"][ParName];
140 RemoveNames.push_back(Set[0].as<std::vector<std::string>>());
141 Exclude.push_back(Set[1].as<bool>());
142 FancyTitle.push_back(Set[2].as<std::string>());
145 const int NSets = int(SetsNames.size());
147 isRelevantParam.resize(NSets);
149 for(
int i = 0; i < NSets; i++)
151 isRelevantParam[i].resize(nParams);
154 for (
int j = 0; j < nParams; j++)
156 isRelevantParam[i][j] =
false;
162 for (
unsigned int k = 0; k < RemoveNames[i].size(); k++)
164 if (ParamNames[j].rfind(RemoveNames[i][k], 0) == 0)
171 isRelevantParam[i][j] =
true;
178 for (
unsigned int k = 0; k < RemoveNames[i].size(); k++)
180 if (ParamNames[j].rfind(RemoveNames[i][k], 0) == 0)
182 isRelevantParam[i][j] =
true;
189 MACH3LOG_INFO(
" Found {} params for set {}", counter, SetsNames[i]);
195 auto canvas = std::make_unique<TCanvas>(
"canvas",
"canvas", 0, 0, 1024, 1024);
200 canvas->SetBottomMargin(0.1f);
201 canvas->SetTopMargin(0.02f);
202 canvas->SetRightMargin(0.08f);
203 canvas->SetLeftMargin(0.15f);
205 gStyle->SetOptTitle(0);
206 gStyle->SetOptStat(0);
207 gStyle->SetPalette(51);
209 std::vector <double> Prior;
210 std::vector <bool>
isFlat;
211 std::vector<std::string> ParamNames;
212 std::vector<std::vector<double>> invCovMatrix;
219 TChain* Chain =
new TChain(
"posteriors",
"");
220 Chain->Add(inputFile.c_str());
223 TObjArray* brlis = Chain->GetListOfBranches();
226 int nBranches = brlis->GetEntries();
227 int RelevantBranches = 0;
228 for (
int i = 0; i < nBranches; i++)
231 TBranch* br =
static_cast<TBranch*
>(brlis->At(i));
236 TString bname = br->GetName();
239 if(bname.BeginsWith(
"param_"))
247 Chain->SetBranchStatus(
"*",
false);
249 std::vector<double> fParProp(RelevantBranches);
251 for (
int i = 0; i < RelevantBranches; ++i)
253 Chain->SetBranchStatus(
BranchNames[i].Data(),
true);
254 Chain->SetBranchAddress(
BranchNames[i].Data(), &fParProp[i]);
258 std::vector<std::string> SetsNames;
259 std::vector<std::string> FancyTitle;
260 std::vector<std::vector<bool>> isRelevantParam;
262 LoadSettings(Settings, SetsNames, FancyTitle, isRelevantParam, ParamNames, nParams);
264 const int NSets = int(SetsNames.size());
265 int AllEvents = int(Chain->GetEntries());
266 std::vector<std::unique_ptr<TH1D>> hLogL(NSets);
267 for (
int i = 0; i < NSets; i++) {
268 std::string NameTemp =
"LogL_" + SetsNames[i];
269 hLogL[i] = std::make_unique<TH1D>(NameTemp.c_str(), NameTemp.c_str(), AllEvents, 0, AllEvents);
270 hLogL[i]->SetLineColor(kBlue);
272 std::vector<double> logL(NSets, 0.0);
273 for(
int n = 0; n < AllEvents; ++n)
279 for(
int k = 0; k < NSets; ++k) logL[k] = 0.;
282 double *logL_private =
nullptr;
286 #pragma omp parallel private(logL_private)
288 logL_private =
new double[NSets];
289 for(
int k = 0; k < NSets; ++k) logL_private[k] = 0.;
292 for (
int i = 0; i < nParams; i++)
294 for (
int j = 0; j <= i; ++j)
299 for(
int k = 0; k < NSets; ++k)
302 if (isRelevantParam[k][i] && isRelevantParam[k][j])
306 if(i != j) scale = 2;
307 logL_private[k] += scale * 0.5*(fParProp[i] - Prior[i])*(fParProp[j] - Prior[j])*invCovMatrix[i][j];
314 for(
int k = 0; k < NSets; ++k)
317 logL[k] += logL_private[k];
320 delete[] logL_private;
324 for (
int i = 0; i < nParams; i++)
326 for (
int j = 0; j <= i; ++j)
331 for(
int k = 0; k < NSets; ++k)
334 if (isRelevantParam[k][i] && isRelevantParam[k][j])
338 if(i != j) scale = 2;
339 logL[k] += scale * 0.5*(fParProp[i] - Prior[i])*(fParProp[j] - Prior[j])*invCovMatrix[i][j];
346 for(
int k = 0; k < NSets; ++k)
348 hLogL[k]->SetBinContent(n, logL[k]);
353 std::string OutputName = inputFile +
"_PenaltyTerm" +
".root";
354 TFile *OutputFile =
M3::Open(OutputName,
"recreate", __FILE__, __LINE__);
355 TDirectory *PenaltyTermDir = OutputFile->mkdir(
"PenaltyTerm");
357 canvas->Print(Form(
"%s_PenaltyTerm.pdf[",inputFile.c_str()),
"pdf");
358 for(
int i = 0; i < NSets; i++)
360 const double Maximum = hLogL[i]->GetMaximum();
361 hLogL[i]->GetYaxis()->SetRangeUser(0., Maximum*1.2);
362 hLogL[i]->SetTitle(FancyTitle[i].c_str());
363 hLogL[i]->GetXaxis()->SetTitle(
"Step");
364 hLogL[i]->GetYaxis()->SetTitle(FancyTitle[i].c_str());
365 hLogL[i]->GetYaxis()->SetTitleOffset(1.4f);
369 PenaltyTermDir->cd();
372 canvas->Print(Form(
"%s_PenaltyTerm.pdf",inputFile.c_str()),
"pdf");
374 canvas->Print(Form(
"%s_PenaltyTerm.pdf]",inputFile.c_str()),
"pdf");
381 int main(
int argc,
char *argv[])
391 std::string filename = argv[1];
392 std::string
config = argv[2];
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
int main(int argc, char *argv[])
void ReadCovFile(const std::string &inputFile, std::vector< double > &Prior, std::vector< bool > &isFlat, std::vector< std::string > &ParamNames, std::vector< std::vector< double >> &invCovMatrix, int &nParams)
void GetPenaltyTerm(const std::string &inputFile, const std::string &configFile)
void LoadSettings(YAML::Node &Settings, std::vector< std::string > &SetsNames, std::vector< std::string > &FancyTitle, std::vector< std::vector< bool >> &isRelevantParam, const std::vector< std::string > &ParamNames, const int nParams)
void SetMaCh3LoggerFormat()
Set messaging format of the logger.
std::vector< TString > BranchNames
bool isFlat(TSpline3_red *&spl)
CW: Helper function used in the constructor, tests to see if the spline is flat.
YAML::Node TMacroToYAML(const TMacro ¯o)
KS: Convert a ROOT TMacro object to a YAML node.
#define M3OpenConfig(filename)
Macro to simplify calling LoadYaml with file and line info.
Custom exception class used throughout MaCh3.
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.
TMatrixDSym * GetCovMatrixFromChain(TDirectory *TempFile)
KS: Retrieve the cross-section covariance matrix from the given TDirectory. Historically,...
void PrintProgressBar(const Long64_t Done, const Long64_t All)
KS: Simply print progress bar.
void PrintConfig(const YAML::Node &node)
KS: Print Yaml config using logger.
void MaCh3Welcome()
KS: Prints welcome message with MaCh3 logo.