41#include <TRestExperimentList.h>
42#include <TRestSensitivity.h>
99 }
while (Chi2 - chi0 < target);
108 }
while (Chi2 - chi0 > target);
113void TRestSensitivity::GenerateCurve() {
116 if (GetNumberOfCurves() > 0)
118 exp->GenerateMockDataSet();
121 RESTInfo <<
"Generating sensitivity curve" <<
RESTendl;
122 std::vector<Double_t> curve;
124 RESTInfo <<
"Generating node : " << node <<
RESTendl;
129 RESTInfo <<
"Curve has been generated. You may use now TRestSensitivity::ExportCurve( fname.txt )."
133void TRestSensitivity::GenerateCurves(Int_t N) {
141 for (
int n = 0; n < N; n++) GenerateCurve();
144std::vector<Double_t> TRestSensitivity::GetCurve(
size_t n) {
145 if (n >= GetNumberOfCurves()) {
146 RESTWarning <<
"Requested curve number : " << n <<
" but only " << GetNumberOfCurves() <<
" generated"
148 return std::vector<Double_t>();
153std::vector<Double_t> TRestSensitivity::GetAveragedCurve() {
154 if (GetNumberOfCurves() <= 0)
return std::vector<Double_t>();
156 std::vector<double> averagedCurve(
fCurves[0].size(), 0.0);
158 for (
const auto& row :
fCurves) {
159 for (
size_t i = 0; i < row.size(); ++i) {
160 averagedCurve[i] += row[i];
164 for (
double& avg : averagedCurve) {
165 avg /=
static_cast<double>(
fCurves.size());
168 return averagedCurve;
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;
177 std::ofstream outputFile(fname);
181 RESTError <<
"TRestSensitivity::ExportCurve. Error opening file for writing!" <<
RESTendl;
186 RESTError <<
"TRestSensitivity::ExportCurve. Curve has not been properly generated" <<
RESTendl;
188 RESTError <<
"Try invoking TRestSensitivity::GenerateCurve" <<
RESTendl;
194 outputFile << node <<
" " << factor * TMath::Power(curve[m], power) << std::endl;
200 RESTInfo <<
"TRestSensitivity::ExportCurve. File has been written successfully!" <<
RESTendl;
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;
208 std::ofstream outputFile(fname);
212 RESTError <<
"TRestSensitivity::ExportCurve. Error opening file for writing!" <<
RESTendl;
217 RESTError <<
"TRestSensitivity::ExportCurve. Curve has not been properly generated" <<
RESTendl;
218 RESTError <<
"Try invoking TRestSensitivity::GenerateCurve" <<
RESTendl;
224 outputFile << node <<
" " << factor * TMath::Power(curve[m], power) << std::endl;
230 RESTInfo <<
"TRestSensitivity::ExportCurve. File has been written successfully!" <<
RESTendl;
242 Double_t target = sigma * sigma;
260 if (!experiment->IsDataReady()) {
261 RESTError <<
"TRestSensitivity::UnbinnedLogLikelihood. Experiment " << experiment->GetName()
266 if (!experiment->GetSignal()->
HasNodes()) {
267 RESTError <<
"Experiment signal : " << experiment->GetSignal()->GetName() <<
" has no nodes!"
278 RESTWarning <<
"Node : " << node <<
" not found in signal : " << experiment->GetSignal()->GetName()
288 if (experiment->GetBackground()->
HasNodes()) {
290 <<
"TRestSensitivity::UnbinnedLogLikelihood is not ready to have background parameter nodes!"
295 Double_t signal = g4 * experiment->GetSignal()->
GetTotalRate() * experiment->GetExposureInSeconds();
299 if (experiment->GetExperimentalCounts() == 0)
return lhood;
301 if (ROOT::IsImplicitMTEnabled()) ROOT::DisableImplicitMT();
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);
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]);
314 Double_t bckRate = experiment->GetBackground()->
GetRate(point);
315 Double_t sgnlRate = experiment->GetSignal()->
GetRate(point);
317 Double_t expectedRate = bckRate + g4 * sgnlRate;
318 lhood += TMath::Log(expectedRate);
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();
332 Double_t coupling = TMath::Sqrt(TMath::Sqrt(
GetCoupling(node)));
333 couplings.push_back(coupling);
337 double min_value = *std::min_element(couplings.begin(), couplings.end());
338 double max_value = *std::max_element(couplings.begin(), couplings.end());
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);
355 while (metadata !=
nullptr) {
357 if (metadata->InheritsFrom(
"TRestExperimentList")) {
359 std::vector<TRestExperiment*> exList = experimentsList->GetExperiments();
361 }
else if (metadata->InheritsFrom(
"TRestExperiment")) {
383 std::vector<std::vector<Double_t>> curves(levels.size());
385 for (
const auto& l : levels) {
386 if (l >= 1 || l <= 0) {
387 RESTError <<
"The level value should be between 0 and 1" <<
RESTendl;
392 std::vector<int> intLevels;
393 for (
const auto& l : levels) {
394 int val = (int)round(l *
fCurves.size());
396 if (val < 0) val = 0;
398 intLevels.push_back(val);
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]);
405 std::sort(v.begin(), v.begin() + v.size());
407 for (
size_t n = 0; n < intLevels.size(); n++) curves[n].push_back(v[intLevels[n]]);
413TCanvas* TRestSensitivity::DrawCurves() {
418 fCanvas =
new TCanvas(
"canv",
"This is the canvas title", 600, 450);
421 TPad* pad1 =
new TPad(
"pad1",
"This is pad1", 0.01, 0.02, 0.99, 0.97);
427 pad1->cd()->SetRightMargin(0.09);
428 pad1->cd()->SetLeftMargin(0.15);
429 pad1->cd()->SetBottomMargin(0.15);
431 std::vector<TGraph*> graphs;
433 for (
size_t n = 0; n < 20; n++) {
435 TGraph* gr =
new TGraph();
436 gr->SetName(grname.c_str());
437 for (
size_t m = 0; m < this->GetCurve(n).size(); m++)
439 TMath::Sqrt(TMath::Sqrt(this->GetCurve(n)[m])));
441 gr->SetLineColorAlpha(kBlue + n, 0.3);
443 graphs.push_back(gr);
446 TGraph* avGr =
new TGraph();
447 std::vector<Double_t> avCurve = GetAveragedCurve();
448 for (
size_t m = 0; m < avCurve.size(); m++)
450 avGr->SetLineColor(kBlack);
451 avGr->SetLineWidth(2);
453 graphs[0]->GetXaxis()->SetLimits(0, 0.25);
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);
466 graphs[0]->Draw(
"AL");
467 for (
unsigned int n = 1; n < graphs.size(); n++) graphs[n]->Draw(
"L");
494TCanvas* TRestSensitivity::DrawLevelCurves() {
499 fCanvas =
new TCanvas(
"canv",
"This is the canvas title", 500, 400);
505 std::vector<std::vector<Double_t>> levelCurves =
GetLevelCurves({0.025, 0.16, 0.375, 0.625, 0.84, 0.975});
507 std::vector<TGraph*> graphs;
508 for (
size_t n = 0; n < levelCurves.size(); n++) {
510 TGraph* gr =
new TGraph();
511 gr->SetName(grname.c_str());
512 for (
size_t m = 0; m < levelCurves[n].size(); m++)
515 gr->SetLineColor(kBlack);
517 graphs.push_back(gr);
520 TGraph* centralGr =
new TGraph();
522 for (
size_t m = 0; m < centralCurve.size(); m++)
524 TMath::Sqrt(TMath::Sqrt(centralCurve[m])));
525 centralGr->SetLineColor(kBlack);
526 centralGr->SetLineWidth(2);
527 centralGr->SetMarkerSize(0.1);
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);
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);
544 graphs[0]->Draw(
"AL");
546 TGraph* randomGr =
new TGraph();
547 std::vector<Double_t> randomCurve =
fCurves[13];
548 for (
size_t m = 0; m < randomCurve.size(); 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);
556 std::vector<TGraph*> shadeGraphs;
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++)
564 TMath::Sqrt(TMath::Sqrt(levelCurves[x][m])));
565 for (
int m = N - 1; m >= 0; --m)
567 TMath::Sqrt(TMath::Sqrt(levelCurves[M - 1 - x][m])));
568 shade->SetFillColorAlpha(kRed, 0.25);
570 shadeGraphs.push_back(shade);
573 for (
unsigned int n = 1; n < graphs.size(); n++) graphs[n]->Draw(
"Lsame");
574 randomGr->Draw(
"LPsame");
592 std::vector<Double_t> nodes = experiment->GetSignal()->GetParameterizationNodes();
602void TRestSensitivity::PrintParameterizationNodes() {
603 std::cout <<
"Curve sensitivity nodes: ";
605 std::cout << std::endl;
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;
618 RESTMetadata <<
" You may access experiment info using TRestSensitivity::GetExperiment(n)" <<
RESTendl;
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.
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.