REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestSensitivity.h
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
23#ifndef REST_TRestSensitivity
24#define REST_TRestSensitivity
25
26#include "TRestExperiment.h"
27
30 private:
32 std::vector<TRestExperiment*> fExperiments;
33
35 std::vector<Double_t> fParameterizationNodes; //<
36
38 std::vector<std::vector<Double_t>> fCurves; //<
39
41 Bool_t fFrozen = false; //< Only needed if we add experiments by other means than RML
42
44 TH1D* fSignalTest = nullptr; //<
45
47 TCanvas* fCanvas = nullptr;
48
49 protected:
50 void InitFromConfigFile() override;
51
52 Double_t UnbinnedLogLikelihood(const TRestExperiment* experiment, Double_t node, Double_t g4 = 0);
53 Double_t ApproachByFactor(Double_t node, Double_t g4, Double_t chi0, Double_t target, Double_t factor);
54
55 public:
56 void Initialize() override;
57
58 void ExtractExperimentParameterizationNodes(Bool_t rescan = false);
59 std::vector<Double_t> GetParameterizationNodes() { return fParameterizationNodes; }
60 void PrintParameterizationNodes();
61
62 Double_t GetCoupling(Double_t node, Double_t sigma = 2, Double_t precision = 0.01);
63 void AddCurve(const std::vector<Double_t>& curve) { fCurves.push_back(curve); }
64 void ImportCurve(const std::vector<Double_t>& curve) { AddCurve(curve); }
65 void GenerateCurve();
66 void GenerateCurves(Int_t N);
67
68 std::vector<Double_t> GetCurve(size_t n = 0);
69 std::vector<Double_t> GetAveragedCurve();
70 std::vector<std::vector<Double_t>> GetLevelCurves(const std::vector<Double_t>& levels);
71
72 void ExportCurve(std::string fname, Double_t factor = 1.e-10, Double_t power = 0.25, int n = 0);
73 void ExportAveragedCurve(std::string fname, Double_t factor = 1.e-10, Double_t power = 0.25);
74
75 TH1D* SignalStatisticalTest(Double_t node, Int_t N);
76
77 void Freeze() { fFrozen = true; }
78
79 std::vector<TRestExperiment*> GetExperiments() { return fExperiments; }
80 TRestExperiment* GetExperiment(const size_t& n) {
81 if (n >= GetNumberOfExperiments())
82 return nullptr;
83 else
84 return fExperiments[n];
85 }
86
87 size_t GetNumberOfExperiments() { return fExperiments.size(); }
88 size_t GetNumberOfCurves() { return fCurves.size(); }
89 size_t GetNumberOfNodes() { return fParameterizationNodes.size(); }
90
91 void PrintMetadata() override;
92
93 TCanvas* DrawCurves();
94 TCanvas* DrawLevelCurves();
95
96 TRestSensitivity(const char* cfgFileName, const std::string& name = "");
99
100 ClassDefOverride(TRestSensitivity, 2);
101};
102#endif
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
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...
Bool_t fFrozen
A flag that will frozen adding more experiments in the future.
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...