REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestSensitivity.cxx
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see https://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see https://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
41#include <TRestExperimentList.h>
42#include <TRestSensitivity.h>
43
44ClassImp(TRestSensitivity);
45
50
55
70TRestSensitivity::TRestSensitivity(const char* cfgFileName, const std::string& name)
71 : TRestMetadata(cfgFileName) {
73}
74
79void TRestSensitivity::Initialize() { SetSectionName(this->ClassName()); }
80
86Double_t TRestSensitivity::ApproachByFactor(Double_t node, Double_t g4, Double_t chi0, Double_t target,
87 Double_t factor) {
88 if (factor == 1) {
89 return 0;
90 }
91
93 Double_t Chi2 = 0;
94 do {
95 Chi2 = 0;
96 for (const auto& exp : fExperiments) Chi2 += -2 * UnbinnedLogLikelihood(exp, node, g4);
97
98 g4 = factor * g4;
99 } while (Chi2 - chi0 < target);
100 g4 = g4 / factor;
101
103 do {
104 Chi2 = 0;
105 for (const auto& exp : fExperiments) Chi2 += -2 * UnbinnedLogLikelihood(exp, node, g4);
106
107 g4 = g4 / factor;
108 } while (Chi2 - chi0 > target);
109
110 return g4 * factor;
111}
112
113void TRestSensitivity::GenerateCurve() {
115
116 if (GetNumberOfCurves() > 0)
117 for (const auto& exp : fExperiments) {
118 exp->GenerateMockDataSet();
119 }
120
121 RESTInfo << "Generating sensitivity curve" << RESTendl;
122 std::vector<Double_t> curve;
123 for (const auto& node : fParameterizationNodes) {
124 RESTInfo << "Generating node : " << node << RESTendl;
125 curve.push_back(GetCoupling(node));
126 }
127 fCurves.push_back(curve);
128
129 RESTInfo << "Curve has been generated. You may use now TRestSensitivity::ExportCurve( fname.txt )."
130 << RESTendl;
131}
132
133void TRestSensitivity::GenerateCurves(Int_t N) {
134 /*
135 std::cout << "TRestSensitivity::GenerateCurves." << std::endl;
136 std::cout << "There is a potential memory leak when generating several curves." << std::endl;
137 std::cout << "This code needs to be reviewed" << std::endl;
138 return;
139 */
140
141 for (int n = 0; n < N; n++) GenerateCurve();
142}
143
144std::vector<Double_t> TRestSensitivity::GetCurve(size_t n) {
145 if (n >= GetNumberOfCurves()) {
146 RESTWarning << "Requested curve number : " << n << " but only " << GetNumberOfCurves() << " generated"
147 << RESTendl;
148 return std::vector<Double_t>();
149 }
150 return fCurves[n];
151}
152
153std::vector<Double_t> TRestSensitivity::GetAveragedCurve() {
154 if (GetNumberOfCurves() <= 0) return std::vector<Double_t>();
155
156 std::vector<double> averagedCurve(fCurves[0].size(), 0.0); // Initialize with zeros
157
158 for (const auto& row : fCurves) {
159 for (size_t i = 0; i < row.size(); ++i) {
160 averagedCurve[i] += row[i];
161 }
162 }
163
164 for (double& avg : averagedCurve) {
165 avg /= static_cast<double>(fCurves.size());
166 }
167
168 return averagedCurve;
169}
170
171void TRestSensitivity::ExportAveragedCurve(std::string fname, Double_t factor, Double_t power) {
172 std::vector<Double_t> curve = GetAveragedCurve();
173 if (curve.empty()) std::cout << "Curve is empty" << std::endl;
174 if (curve.empty()) return;
175
176 // Open a file for writing
177 std::ofstream outputFile(fname);
178
179 // Check if the file is opened successfully
180 if (!outputFile) {
181 RESTError << "TRestSensitivity::ExportCurve. Error opening file for writing!" << RESTendl;
182 return;
183 }
184
185 if (fParameterizationNodes.size() != curve.size()) {
186 RESTError << "TRestSensitivity::ExportCurve. Curve has not been properly generated" << RESTendl;
187 RESTError << "Parameterization nodes: " << fParameterizationNodes.size() << RESTendl;
188 RESTError << "Try invoking TRestSensitivity::GenerateCurve" << RESTendl;
189 return;
190 }
191
192 int m = 0;
193 for (const auto& node : fParameterizationNodes) {
194 outputFile << node << " " << factor * TMath::Power(curve[m], power) << std::endl;
195 m++;
196 }
197
198 outputFile.close();
199
200 RESTInfo << "TRestSensitivity::ExportCurve. File has been written successfully!" << RESTendl;
201}
202
203void TRestSensitivity::ExportCurve(std::string fname, Double_t factor, Double_t power, int n) {
204 std::vector<Double_t> curve = GetCurve(n);
205 if (curve.empty()) return;
206
207 // Open a file for writing
208 std::ofstream outputFile(fname);
209
210 // Check if the file is opened successfully
211 if (!outputFile) {
212 RESTError << "TRestSensitivity::ExportCurve. Error opening file for writing!" << RESTendl;
213 return;
214 }
215
216 if (fParameterizationNodes.size() != curve.size()) {
217 RESTError << "TRestSensitivity::ExportCurve. Curve has not been properly generated" << RESTendl;
218 RESTError << "Try invoking TRestSensitivity::GenerateCurve" << RESTendl;
219 return;
220 }
221
222 int m = 0;
223 for (const auto& node : fParameterizationNodes) {
224 outputFile << node << " " << factor * TMath::Power(curve[m], power) << std::endl;
225 m++;
226 }
227
228 outputFile.close();
229
230 RESTInfo << "TRestSensitivity::ExportCurve. File has been written successfully!" << RESTendl;
231}
232
236Double_t TRestSensitivity::GetCoupling(Double_t node, Double_t sigma, Double_t precision) {
237 Double_t Chi2_0 = 0;
238 for (const auto& exp : fExperiments) {
239 Chi2_0 += -2 * UnbinnedLogLikelihood(exp, node, 0);
240 }
241
242 Double_t target = sigma * sigma;
243
244 Double_t g4 = 0.5;
245
246 g4 = ApproachByFactor(node, g4, Chi2_0, target, 2);
247 g4 = ApproachByFactor(node, g4, Chi2_0, target, 1.2);
248 g4 = ApproachByFactor(node, g4, Chi2_0, target, 1.02);
249 g4 = ApproachByFactor(node, g4, Chi2_0, target, 1.0002);
250
251 return g4;
252}
253
257Double_t TRestSensitivity::UnbinnedLogLikelihood(const TRestExperiment* experiment, Double_t node,
258 Double_t g4) {
259 Double_t lhood = 0;
260 if (!experiment->IsDataReady()) {
261 RESTError << "TRestSensitivity::UnbinnedLogLikelihood. Experiment " << experiment->GetName()
262 << " is not ready!" << RESTendl;
263 return lhood;
264 }
265
266 if (!experiment->GetSignal()->HasNodes()) {
267 RESTError << "Experiment signal : " << experiment->GetSignal()->GetName() << " has no nodes!"
268 << RESTendl;
269 return lhood;
270 }
271
274 Int_t nd = experiment->GetSignal()->FindActiveNode(node);
275 if (nd >= 0)
276 experiment->GetSignal()->SetActiveNode(nd);
277 else {
278 RESTWarning << "Node : " << node << " not found in signal : " << experiment->GetSignal()->GetName()
279 << RESTendl;
280 return 0.0;
281 }
282
287
288 if (experiment->GetBackground()->HasNodes()) {
289 RESTWarning
290 << "TRestSensitivity::UnbinnedLogLikelihood is not ready to have background parameter nodes!"
291 << RESTendl;
292 return 0.0;
293 }
294
295 Double_t signal = g4 * experiment->GetSignal()->GetTotalRate() * experiment->GetExposureInSeconds();
296
297 lhood = -signal;
298
299 if (experiment->GetExperimentalCounts() == 0) return lhood;
300
301 if (ROOT::IsImplicitMTEnabled()) ROOT::DisableImplicitMT();
302
303 std::vector<std::vector<Double_t>> trackingData;
304 for (const auto& var : experiment->GetSignal()->GetVariables()) {
305 auto values = experiment->GetExperimentalDataFrame().Take<Double_t>(var);
306 std::vector<Double_t> vDbl = std::move(*values);
307 trackingData.push_back(vDbl);
308 }
309
310 for (size_t n = 0; n < trackingData[0].size(); n++) {
311 std::vector<Double_t> point;
312 for (size_t m = 0; m < trackingData.size(); m++) point.push_back(trackingData[m][n]);
313
314 Double_t bckRate = experiment->GetBackground()->GetRate(point);
315 Double_t sgnlRate = experiment->GetSignal()->GetRate(point);
316
317 Double_t expectedRate = bckRate + g4 * sgnlRate;
318 lhood += TMath::Log(expectedRate);
319 }
320
321 return lhood;
322}
323
327TH1D* TRestSensitivity::SignalStatisticalTest(Double_t node, Int_t N) {
328 std::vector<Double_t> couplings;
329 for (int n = 0; n < N; n++) {
330 for (const auto& exp : fExperiments) exp->GetSignal()->RegenerateActiveNodeDensity();
331
332 Double_t coupling = TMath::Sqrt(TMath::Sqrt(GetCoupling(node)));
333 couplings.push_back(coupling);
334 }
335
336 // Directly assign the minimum and maximum values
337 double min_value = *std::min_element(couplings.begin(), couplings.end());
338 double max_value = *std::max_element(couplings.begin(), couplings.end());
339
340 if (fSignalTest) delete fSignalTest;
341 fSignalTest = new TH1D("SignalTest", "A signal test", 100, 0.9 * min_value, 1.1 * max_value);
342 for (const auto& coup : couplings) fSignalTest->Fill(coup);
343
344 return fSignalTest;
345}
346
352
353 int cont = 0;
354 TRestMetadata* metadata = (TRestMetadata*)this->InstantiateChildMetadata(cont, "Experiment");
355 while (metadata != nullptr) {
356 cont++;
357 if (metadata->InheritsFrom("TRestExperimentList")) {
358 TRestExperimentList* experimentsList = (TRestExperimentList*)metadata;
359 std::vector<TRestExperiment*> exList = experimentsList->GetExperiments();
360 fExperiments.insert(fExperiments.end(), exList.begin(), exList.end());
361 } else if (metadata->InheritsFrom("TRestExperiment")) {
362 fExperiments.push_back((TRestExperiment*)metadata);
363 }
364
365 metadata = (TRestMetadata*)this->InstantiateChildMetadata(cont, "Experiment");
366 }
367
368 Initialize();
369}
370
375//
382std::vector<std::vector<Double_t>> TRestSensitivity::GetLevelCurves(const std::vector<Double_t>& levels) {
383 std::vector<std::vector<Double_t>> curves(levels.size());
384
385 for (const auto& l : levels) {
386 if (l >= 1 || l <= 0) {
387 RESTError << "The level value should be between 0 and 1" << RESTendl;
388 return curves;
389 }
390 }
391
392 std::vector<int> intLevels;
393 for (const auto& l : levels) {
394 int val = (int)round(l * fCurves.size());
395 if (val >= (int)fCurves.size()) val = fCurves.size() - 1;
396 if (val < 0) val = 0;
397
398 intLevels.push_back(val);
399 }
400
401 for (size_t m = 0; m < fCurves[0].size(); m++) {
402 std::vector<Double_t> v;
403 for (size_t n = 0; n < fCurves.size(); n++) v.push_back(fCurves[n][m]);
404
405 std::sort(v.begin(), v.begin() + v.size());
406
407 for (size_t n = 0; n < intLevels.size(); n++) curves[n].push_back(v[intLevels[n]]);
408 }
409
410 return curves;
411}
412
413TCanvas* TRestSensitivity::DrawCurves() {
414 if (fCanvas != NULL) {
415 delete fCanvas;
416 fCanvas = NULL;
417 }
418 fCanvas = new TCanvas("canv", "This is the canvas title", 600, 450);
419 fCanvas->Draw();
420
421 TPad* pad1 = new TPad("pad1", "This is pad1", 0.01, 0.02, 0.99, 0.97);
422 // pad1->Divide(2, 2);
423 pad1->Draw();
424
426 // pad1->cd()->SetLogx();
427 pad1->cd()->SetRightMargin(0.09);
428 pad1->cd()->SetLeftMargin(0.15);
429 pad1->cd()->SetBottomMargin(0.15);
430
431 std::vector<TGraph*> graphs;
432
433 for (size_t n = 0; n < 20; n++) {
434 std::string grname = "gr" + IntegerToString(n);
435 TGraph* gr = new TGraph();
436 gr->SetName(grname.c_str());
437 for (size_t m = 0; m < this->GetCurve(n).size(); m++)
438 gr->SetPoint(gr->GetN(), fParameterizationNodes[m],
439 TMath::Sqrt(TMath::Sqrt(this->GetCurve(n)[m])));
440
441 gr->SetLineColorAlpha(kBlue + n, 0.3);
442 gr->SetLineWidth(1);
443 graphs.push_back(gr);
444 }
445
446 TGraph* avGr = new TGraph();
447 std::vector<Double_t> avCurve = GetAveragedCurve();
448 for (size_t m = 0; m < avCurve.size(); m++)
449 avGr->SetPoint(avGr->GetN(), fParameterizationNodes[m], TMath::Sqrt(TMath::Sqrt(avCurve[m])));
450 avGr->SetLineColor(kBlack);
451 avGr->SetLineWidth(2);
452
453 graphs[0]->GetXaxis()->SetLimits(0, 0.25);
454 // graphs[0]->GetHistogram()->SetMaximum(1);
455 // graphs[0]->GetHistogram()->SetMinimum(lowRange);
456
457 graphs[0]->GetXaxis()->SetTitle("mass [eV]");
458 graphs[0]->GetXaxis()->SetTitleSize(0.05);
459 graphs[0]->GetXaxis()->SetLabelSize(0.05);
460 graphs[0]->GetYaxis()->SetTitle("g_{a#gamma} [10^{-10} GeV^{-1}]");
461 graphs[0]->GetYaxis()->SetTitleOffset(1.5);
462 graphs[0]->GetYaxis()->SetTitleSize(0.05);
463 // graphs[0]->GetYaxis()->SetLabelSize(0.05);
464 // graphs[0]->GetYaxis()->SetLabelOffset(0.0);
465 // pad1->cd()->SetLogy();
466 graphs[0]->Draw("AL");
467 for (unsigned int n = 1; n < graphs.size(); n++) graphs[n]->Draw("L");
468 avGr->Draw("L");
469
470 /*
471Double_t lx1 = 0.6, ly1 = 0.75, lx2 = 0.9, ly2 = 0.95;
472if (eLegendCoords.size() > 0) {
473 lx1 = eLegendCoords[0];
474 ly1 = eLegendCoords[1];
475 lx2 = eLegendCoords[2];
476 ly2 = eLegendCoords[3];
477}
478TLegend* legend = new TLegend(lx1, ly1, lx2, ly2);
479
480legend->SetTextSize(0.03);
481legend->SetHeader("Energies", "C"); // option "C" allows to center the header
482for (unsigned int n = 0; n < energies.size(); n++) {
483 std::string lname = "gr" + IntegerToString(n);
484 std::string ltitle = DoubleToString(energies[n]) + " keV";
485
486 legend->AddEntry(lname.c_str(), ltitle.c_str(), "l");
487}
488legend->Draw();
489 */
490
491 return fCanvas;
492}
493
494TCanvas* TRestSensitivity::DrawLevelCurves() {
495 if (fCanvas != NULL) {
496 delete fCanvas;
497 fCanvas = NULL;
498 }
499 fCanvas = new TCanvas("canv", "This is the canvas title", 500, 400);
500 fCanvas->Draw();
501 fCanvas->SetLeftMargin(0.15);
502 fCanvas->SetRightMargin(0.04);
503 fCanvas->SetLogx();
504
505 std::vector<std::vector<Double_t>> levelCurves = GetLevelCurves({0.025, 0.16, 0.375, 0.625, 0.84, 0.975});
506
507 std::vector<TGraph*> graphs;
508 for (size_t n = 0; n < levelCurves.size(); n++) {
509 std::string grname = "gr" + IntegerToString(n);
510 TGraph* gr = new TGraph();
511 gr->SetName(grname.c_str());
512 for (size_t m = 0; m < levelCurves[n].size(); m++)
513 gr->SetPoint(gr->GetN(), fParameterizationNodes[m], TMath::Sqrt(TMath::Sqrt(levelCurves[n][m])));
514
515 gr->SetLineColor(kBlack);
516 gr->SetLineWidth(1);
517 graphs.push_back(gr);
518 }
519
520 TGraph* centralGr = new TGraph();
521 std::vector<Double_t> centralCurve = GetLevelCurves({0.5})[0];
522 for (size_t m = 0; m < centralCurve.size(); m++)
523 centralGr->SetPoint(centralGr->GetN(), fParameterizationNodes[m],
524 TMath::Sqrt(TMath::Sqrt(centralCurve[m])));
525 centralGr->SetLineColor(kBlack);
526 centralGr->SetLineWidth(2);
527 centralGr->SetMarkerSize(0.1);
528
529 graphs[0]->GetYaxis()->SetRangeUser(0, 0.5);
530 graphs[0]->GetXaxis()->SetRangeUser(0.001, 0.25);
531 graphs[0]->GetXaxis()->SetLimits(0.0001, 0.25);
532 graphs[0]->GetXaxis()->SetTitle("mass [eV]");
533 graphs[0]->GetXaxis()->SetTitleSize(0.04);
534 graphs[0]->GetXaxis()->SetTitleOffset(1.15);
535 graphs[0]->GetXaxis()->SetLabelSize(0.04);
536
537 // graphs[0]->GetYaxis()->SetLabelFont(43);
538 graphs[0]->GetYaxis()->SetTitle("g_{a#gamma} [10^{-10} GeV^{-1}]");
539 graphs[0]->GetYaxis()->SetTitleOffset(1.5);
540 graphs[0]->GetYaxis()->SetTitleSize(0.04);
541 graphs[0]->GetYaxis()->SetLabelSize(0.04);
542 // graphs[0]->GetYaxis()->SetLabelOffset(0);
543 // graphs[0]->GetYaxis()->SetLabelFont(43);
544 graphs[0]->Draw("AL");
545
546 TGraph* randomGr = new TGraph();
547 std::vector<Double_t> randomCurve = fCurves[13];
548 for (size_t m = 0; m < randomCurve.size(); m++)
549 randomGr->SetPoint(randomGr->GetN(), fParameterizationNodes[m],
550 TMath::Sqrt(TMath::Sqrt(randomCurve[m])));
551 randomGr->SetLineColor(kBlack);
552 randomGr->SetLineWidth(1);
553 randomGr->SetMarkerSize(0.3);
554 randomGr->SetMarkerStyle(4);
555
556 std::vector<TGraph*> shadeGraphs;
557
558 int M = (int)levelCurves.size();
559 for (int x = 0; x < M / 2; x++) {
560 TGraph* shade = new TGraph();
561 int N = levelCurves[0].size();
562 for (size_t m = 0; m < levelCurves[0].size(); m++)
563 shade->SetPoint(shade->GetN(), fParameterizationNodes[m],
564 TMath::Sqrt(TMath::Sqrt(levelCurves[x][m])));
565 for (int m = N - 1; m >= 0; --m)
566 shade->SetPoint(shade->GetN(), fParameterizationNodes[m],
567 TMath::Sqrt(TMath::Sqrt(levelCurves[M - 1 - x][m])));
568 shade->SetFillColorAlpha(kRed, 0.25);
569 shade->Draw("f");
570 shadeGraphs.push_back(shade);
571 }
572
573 for (unsigned int n = 1; n < graphs.size(); n++) graphs[n]->Draw("Lsame");
574 randomGr->Draw("LPsame");
575 // centralGr->Draw("Lsame");
576
577 return fCanvas;
578}
579
588 if (fParameterizationNodes.empty() || rescan) {
590
591 for (const auto& experiment : fExperiments) {
592 std::vector<Double_t> nodes = experiment->GetSignal()->GetParameterizationNodes();
593 fParameterizationNodes.insert(fParameterizationNodes.end(), nodes.begin(), nodes.end());
594
595 std::sort(fParameterizationNodes.begin(), fParameterizationNodes.end());
596 auto last = std::unique(fParameterizationNodes.begin(), fParameterizationNodes.end());
598 }
599 }
600}
601
602void TRestSensitivity::PrintParameterizationNodes() {
603 std::cout << "Curve sensitivity nodes: ";
604 for (const auto& node : fParameterizationNodes) std::cout << node << "\t";
605 std::cout << std::endl;
606}
607
613
614 RESTMetadata << " - Number of parameterization nodes : " << GetNumberOfNodes() << RESTendl;
615 RESTMetadata << " - Number of experiments loaded : " << GetNumberOfExperiments() << RESTendl;
616 RESTMetadata << " - Number of sensitivity curves generated : " << GetNumberOfCurves() << RESTendl;
617 RESTMetadata << " " << RESTendl;
618 RESTMetadata << " You may access experiment info using TRestSensitivity::GetExperiment(n)" << RESTendl;
619
620 RESTMetadata << "----" << RESTendl;
621}
Int_t FindActiveNode(Double_t node)
It returns the position of the fParameterizationNodes element for the variable name given by argument...
Double_t GetTotalRate()
This method integrates the rate to all the parameter space defined in the density function....
Bool_t HasNodes()
It returns true if any nodes have been defined.
Double_t GetRate(std::vector< Double_t > point)
It returns the intensity/rate (in seconds) corresponding to the generated distribution or formula eva...
Int_t SetActiveNode(Double_t node)
It returns the position of the fParameterizationNodes element for the variable name given by argument...
A helper metadata class to create a list of TRestExperiment instances.
It includes a model definition and experimental data used to obtain a final experimental sensitivity.
A base class for any REST metadata class.
Definition: TRestMetadata.h:74
virtual void PrintMetadata()
Implemented it in the derived metadata class to print out specific metadata information.
endl_t RESTendl
Termination flag object for TRestStringOutput.
Int_t LoadConfigFromFile(const std::string &configFilename, const std::string &sectionName="")
Give the file name, find out the corresponding section. Then call the main starter.
TRestMetadata * InstantiateChildMetadata(int index, std::string pattern="")
This method will retrieve a new TRestMetadata instance of a child element of the present TRestMetadat...
virtual void InitFromConfigFile()
To make settings from rml file. This method must be implemented in the derived class.
void SetSectionName(std::string sName)
set the section name, clear the section content
std::string fConfigFileName
Full name of the rml file.
It combines a number of experimental conditions allowing to calculate a combined experimental sensiti...
void Initialize() override
It will initialize the data frame with the filelist and column names (or observables) that have been ...
~TRestSensitivity()
Default destructor.
std::vector< std::vector< Double_t > > fCurves
A vector of calculated sensitivity curves defined as a funtion of the parametric node.
std::vector< TRestExperiment * > fExperiments
A list of experimental conditions included to get a final sensitivity plot.
void InitFromConfigFile() override
It customizes the retrieval of XML data values of this class.
TCanvas * fCanvas
A canvas to draw.
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestAxionSolarFlux.
TH1D * fSignalTest
It is used to generate a histogram with the signal distribution produced with different signal sample...
void ExtractExperimentParameterizationNodes(Bool_t rescan=false)
It scans all the experiment signals parametric nodes to build a complete list of nodes used to build ...
TRestSensitivity()
Default constructor.
std::vector< std::vector< Double_t > > GetLevelCurves(const std::vector< Double_t > &levels)
This method is used to obtain the list of curves that satisfy that each value inside the curve is pla...
Double_t GetCoupling(Double_t node, Double_t sigma=2, Double_t precision=0.01)
It will return the coupling value for which Chi=sigma.
Double_t UnbinnedLogLikelihood(const TRestExperiment *experiment, Double_t node, Double_t g4=0)
It returns the Log(L) for the experiment and coupling given by argument.
std::vector< Double_t > fParameterizationNodes
The fusioned list of parameterization nodes found at each experiment signal.
Double_t ApproachByFactor(Double_t node, Double_t g4, Double_t chi0, Double_t target, Double_t factor)
It will return a value of the coupling, g4, such that (chi-chi0) gets closer to the target value give...
std::string IntegerToString(Int_t n, std::string format="%d")
Gets a string from an integer.