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 Int_t fActiveNode = -1; //<
60
62 std::vector<THnD*> fNodeDensity; //<
63
65 Int_t fSamples = 0; //<
66
68 Bool_t fInterpolation = true; //<
69
71 TRestResponse* fResponse = nullptr; //<
72
74 Float_t fPrecision = 0.01; //<
75
77 TRandom3* fRandom = nullptr;
78
80 UInt_t fSeed = 0; //<
81
83 TCanvas* fCanvas = nullptr;
84
85 Bool_t HasDensity() { return !fNodeDensity.empty(); }
86
88 Bool_t ValidNode(Double_t node) {
89 SetActiveNode(node);
90 if (GetActiveNode() >= 0) return true;
91 return false;
92 }
93
94 Int_t GetVariableIndex(std::string varName);
95
96 void InitFromConfigFile() override;
97
98 virtual void FillHistograms() = 0;
99
100 public:
101 void Initialize() override;
102 void RegenerateHistograms(UInt_t seed = 0);
103
105 Bool_t HasNodes() { return !fParameterizationNodes.empty(); }
106
107 virtual void RegenerateActiveNodeDensity() {}
108
109 std::string GetNature() const { return fNature; }
110 TRestResponse* GetResponse() const { return fResponse; }
111 Float_t GetPrecision() { return fPrecision; }
112 size_t GetDimensions() { return fVariables.size(); }
113 Int_t GetSamples() { return fSamples; }
114 Int_t GetActiveNode() { return fActiveNode; }
115 Double_t GetActiveNodeValue() { return fParameterizationNodes[fActiveNode]; }
116 std::vector<Double_t> GetParameterizationNodes() { return fParameterizationNodes; }
117
118 std::vector<std::string> GetVariables() const { return fVariables; }
119 std::vector<TVector2> GetRanges() const { return fRanges; }
120 std::vector<Int_t> GetNbins() const { return fNbins; }
121
122 Double_t GetRawRate(std::vector<Double_t> point);
123 Double_t GetTotalRate();
124 Double_t GetNormalizedRate(std::vector<Double_t> point);
125 Double_t GetRate(std::vector<Double_t> point);
126
127 Double_t GetBinCenter(Int_t nDim, const Int_t bin);
128
129 void SetPrecision(const Float_t& pr) { fPrecision = pr; }
130
131 Int_t FindActiveNode(Double_t node);
132 Int_t SetActiveNode(Double_t node);
133 Int_t SetActiveNode(Int_t n) {
134 fActiveNode = n;
135 return fActiveNode;
136 }
137
138 void SetSamples(Int_t samples) { fSamples = samples; }
139
140 Bool_t Interpolation() { return fInterpolation; }
141 void EnableInterpolation() { fInterpolation = true; }
142 void DisableInterpolation() { fInterpolation = false; }
143
144 THnD* GetDensityForNode(Double_t value);
145 THnD* GetDensityForActiveNode();
146 THnD* GetDensity() { return GetDensityForActiveNode(); }
147
148 TH1D* GetHistogram(Double_t node, std::string varName);
149 TH2D* GetHistogram(Double_t node, std::string varName1, std::string varName2);
150 TH3D* GetHistogram(Double_t node, std::string varName1, std::string varName2, std::string varName3);
151
152 TH1D* GetHistogram(std::string varName);
153 TH2D* GetHistogram(std::string varName1, std::string varName2);
154 TH3D* GetHistogram(std::string varName1, std::string varName2, std::string varName3);
155
156 ROOT::RVecD GetRandom();
157
158 ROOT::RDF::RNode GetMonteCarloDataFrame(Int_t N = 100);
159
160 TCanvas* DrawComponent(std::vector<std::string> drawVariables, std::vector<std::string> scanVariables,
161 Int_t binScanSize = 1, TString drawOption = "");
162
163 void LoadResponse(const TRestResponse& resp);
164
165 void PrintMetadata() override;
166
167 void PrintStatistics();
168 void PrintNodes();
169
170 TRestComponent(const char* cfgFileName, const std::string& name = "");
173
174 ClassDefOverride(TRestComponent, 5);
175};
176#endif
It defines a background/signal model distribution in a given parameter space (tipically x,...
UInt_t fSeed
Seed used in random generator.
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.
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...
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 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.
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:70
A response matrix that might be applied to a given component inside a TRestComponent.
Definition: TRestResponse.h:29