REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestComponent.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_TRestComponent
24#define REST_TRestComponent
25
26#include <TCanvas.h>
27#include <THn.h>
28#include <TRandom3.h>
29
30#include <ROOT/RDataFrame.hxx>
31#include <ROOT/RVec.hxx>
32
33#include "TRestDataSet.h"
34#include "TRestMetadata.h"
35#include "TRestResponse.h"
36
39 protected:
41 std::string fNature = "unknown"; //<
42
44 std::vector<std::string> fVariables; //<
45
47 std::vector<TVector2> fRanges; //<
48
50 std::vector<Int_t> fNbins; //<
51
53 std::string fParameter = ""; //<
54
56 std::vector<Double_t> fParameterizationNodes; //<
57
59 Double_t fFirstParameterValue = 0; //<
60
62 Double_t fLastParameterValue = 0; //<
63
65 Double_t fStepParameterValue = 0; //<
66
68 Bool_t fExponential = false; //<
69
71 Int_t fActiveNode = -1; //<
72
74 std::vector<THnD*> fNodeDensity; //<
75
77 Int_t fSamples = 0; //<
78
80 Bool_t fInterpolation = true; //<
81
83 TRestResponse* fResponse = nullptr; //<
84
86 Float_t fPrecision = 0.01; //<
87
89 TRandom3* fRandom = nullptr;
90
92 UInt_t fSeed = 0; //<
93
95 TCanvas* fCanvas = nullptr;
96
97 Bool_t HasDensity() { return !fNodeDensity.empty(); }
98
100 Bool_t ValidNode(Double_t node) {
101 SetActiveNode(node);
102 if (GetActiveNode() >= 0) return true;
103 return false;
104 }
105
106 Int_t GetVariableIndex(std::string varName);
107
108 void InitFromConfigFile() override;
109
110 virtual void FillHistograms() = 0;
111
112 public:
113 void Initialize() override;
114 void RegenerateHistograms(UInt_t seed = 0);
115
116 void RegenerateParametricNodes(Double_t from, Double_t to, Double_t step, Bool_t expIncrease = false);
117
119 Bool_t HasNodes() { return !fParameterizationNodes.empty(); }
120
121 virtual void RegenerateActiveNodeDensity() {}
122
123 std::string GetNature() const { return fNature; }
124 TRestResponse* GetResponse() const { return fResponse; }
125 Float_t GetPrecision() { return fPrecision; }
126 size_t GetDimensions() { return fVariables.size(); }
127 Int_t GetSamples() { return fSamples; }
128 Int_t GetActiveNode() { return fActiveNode; }
129 Double_t GetActiveNodeValue() {
130 if (fActiveNode >= 0 && fActiveNode < (Int_t)fParameterizationNodes.size())
132 return 0;
133 }
134 std::vector<Double_t> GetParameterizationNodes() { return fParameterizationNodes; }
135
136 std::vector<std::string> GetVariables() const { return fVariables; }
137 std::vector<TVector2> GetRanges() const { return fRanges; }
138 std::vector<Int_t> GetNbins() const { return fNbins; }
139
140 Double_t GetRawRate(std::vector<Double_t> point);
141 Double_t GetTotalRate();
142 Double_t GetMaxRate();
143 Double_t GetAllNodesIntegratedRate();
144 Double_t GetNormalizedRate(std::vector<Double_t> point);
145 Double_t GetRate(std::vector<Double_t> point);
146
147 Double_t GetBinCenter(Int_t nDim, const Int_t bin);
148
149 void SetPrecision(const Float_t& pr) { fPrecision = pr; }
150
151 Int_t FindActiveNode(Double_t node);
152 Int_t SetActiveNode(Double_t node);
153 Int_t SetActiveNode(Int_t n) {
154 fActiveNode = n;
155 return fActiveNode;
156 }
157
158 void SetSamples(Int_t samples) { fSamples = samples; }
159
160 Bool_t Interpolation() { return fInterpolation; }
161 void EnableInterpolation() { fInterpolation = true; }
162 void DisableInterpolation() { fInterpolation = false; }
163
164 THnD* GetDensityForNode(Double_t value);
165 THnD* GetDensityForActiveNode();
166 THnD* GetDensity() { return GetDensityForActiveNode(); }
167
168 TH1D* GetHistogram(Double_t node, std::string varName);
169 TH2D* GetHistogram(Double_t node, std::string varName1, std::string varName2);
170 TH3D* GetHistogram(Double_t node, std::string varName1, std::string varName2, std::string varName3);
171
172 TH1D* GetHistogram(std::string varName);
173 TH2D* GetHistogram(std::string varName1, std::string varName2);
174 TH3D* GetHistogram(std::string varName1, std::string varName2, std::string varName3);
175
176 ROOT::RVecD GetRandom();
177
178 ROOT::RDF::RNode GetMonteCarloDataFrame(Int_t N = 100);
179
180 TCanvas* DrawComponent(std::vector<std::string> drawVariables, std::vector<std::string> scanVariables,
181 Int_t binScanSize = 1, TString drawOption = "");
182
183 void LoadResponse(const TRestResponse& resp);
184
185 void PrintMetadata() override;
186
187 void PrintStatistics();
188 void PrintNodes();
189
190 TRestComponent(const char* cfgFileName, const std::string& name = "");
193
194 ClassDefOverride(TRestComponent, 6);
195};
196#endif
It defines a background/signal model distribution in a given parameter space (tipically x,...
Double_t fFirstParameterValue
It defines the first parametric node value in case of automatic parameter generation.
Double_t fLastParameterValue
It defines the upper limit for the automatic parametric node values generation.
UInt_t fSeed
Seed used in random generator.
Double_t fStepParameterValue
It defines the increasing step for automatic parameter list generation.
void InitFromConfigFile() override
It customizes the retrieval of XML data values of this class.
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestAxionSolarFlux.
Double_t GetAllNodesIntegratedRate()
This method returns the integrated total rate for all the nodes The result will be returned in s-1.
Int_t fActiveNode
It is used to define the node that will be accessed for rate retrieval.
Int_t FindActiveNode(Double_t node)
It returns the position of the fParameterizationNodes element for the variable name given by argument...
void RegenerateParametricNodes(Double_t from, Double_t to, Double_t step, Bool_t expIncrease=false)
It allows to produce a parameter nodes list providing the initial value, the final value and the step...
Int_t fSamples
It introduces a fixed number of samples (if 0 it will take all available samples)
Int_t GetVariableIndex(std::string varName)
It returns the position of the fVariable element for the variable name given by argument.
Double_t GetBinCenter(Int_t nDim, const Int_t bin)
It returns the bin center of the given component dimension.
void PrintNodes()
It prints out on screen the values of the parametric node.
Double_t GetRawRate(std::vector< Double_t > point)
It returns the intensity/rate (in seconds) corresponding to the generated distribution or formula eva...
Bool_t fExponential
It true the parametric values automatically generated will grow exponentially.
Bool_t fInterpolation
Enables or disables the interpolation at TRestComponentDataSet::GetRawRate.
std::string fParameter
It is used to parameterize a set of distribution densities (e.g. WIMP or axion mass)
Double_t GetNormalizedRate(std::vector< Double_t > point)
It returns the intensity/rate (in seconds) corresponding to the generated distribution or formula eva...
TRestComponent()
Default constructor.
TH1D * GetHistogram(Double_t node, std::string varName)
It returns a 1-dimensional projected histogram for the variable names provided in the argument.
Float_t fPrecision
A precision used to select the node value with a given range defined as a fraction of the value.
std::vector< Int_t > fNbins
The number of bins in which we should divide each variable.
Double_t GetMaxRate()
This method returns the total rate for the node that has the highest contribution The result will be ...
TCanvas * DrawComponent(std::vector< std::string > drawVariables, std::vector< std::string > scanVariables, Int_t binScanSize=1, TString drawOption="")
A method allowing to draw a series of plots representing the density distributions.
void Initialize() override
It initializes the random number. We avoid to define the section name here since we will never define...
TRestResponse * fResponse
A pointer to the detector response.
Double_t GetTotalRate()
This method integrates the rate to all the parameter space defined in the density function....
std::vector< TVector2 > fRanges
The range of each of the variables used to create the PDF distribution.
std::vector< Double_t > fParameterizationNodes
It defines the nodes of the parameterization (Initialized by the dataset)
Bool_t HasNodes()
It returns true if any nodes have been defined.
TCanvas * fCanvas
A canvas for drawing the active node component.
Double_t GetRate(std::vector< Double_t > point)
It returns the intensity/rate (in seconds) corresponding to the generated distribution or formula eva...
TRandom3 * fRandom
Internal process random generator.
std::vector< std::string > fVariables
A list with the branches that will be used to create the distribution space.
~TRestComponent()
Default destructor.
void RegenerateHistograms(UInt_t seed=0)
It will produce a histogram with the distribution defined using the variables and the weights for eac...
Bool_t ValidNode(Double_t node)
It returns true if the node has been properly identified.
Int_t SetActiveNode(Double_t node)
It returns the position of the fParameterizationNodes element for the variable name given by argument...
std::vector< THnD * > fNodeDensity
The generated N-dimensional variable space density for a given node.
std::string fNature
It defines the component type (unknown/signal/background)
A base class for any REST metadata class.
Definition: TRestMetadata.h:74
A response matrix that might be applied to a given component inside a TRestComponent.
Definition: TRestResponse.h:29