REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestExperimentList.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 
43 ClassImp(TRestExperimentList);
44 
49 
54 
68 TRestExperimentList::TRestExperimentList(const char* cfgFileName, const std::string& name)
69  : TRestMetadata(cfgFileName) {
71 }
72 
77 void TRestExperimentList::Initialize() { SetSectionName(this->ClassName()); }
78 
84 
85  if (!fExperimentsFile.empty() && fExperiments.empty()) {
87 
88  for (auto& row : fExperimentsTable)
89  for (auto& el : row) el = REST_StringHelper::ReplaceMathematicalExpressions(el);
90 
91  if (fExperimentsTable.empty()) {
92  RESTError << "TRestExperimentList::InitFromConfigFile. The experiments table is empty!"
93  << RESTendl;
94  return;
95  }
96 
97  Int_t nTableColumns = fExperimentsTable[0].size();
98 
99  int cont = 0;
100  TRestComponent* comp = (TRestComponent*)this->InstantiateChildMetadata(cont, "Component");
101  while (comp != nullptr) {
102  if (ToLower(comp->GetNature()) == "background")
103  fBackground = comp;
104  else if (ToLower(comp->GetNature()) == "signal")
105  fSignal = comp;
106  else
107  RESTWarning << "TRestExperimentList::InitFromConfigFile. Unknown component!" << RESTendl;
108 
109  cont++;
110  comp = (TRestComponent*)this->InstantiateChildMetadata(cont, "Component");
111  }
112 
113  Int_t nExpectedColumns = 3;
114  if (fSignal) nExpectedColumns--;
115  if (fBackground) nExpectedColumns--;
116  if (fExposureTime > 0) nExpectedColumns--;
117 
118  if (nExpectedColumns == 0) {
119  RESTError << "TRestExperimentList::InitFromConfigFile. At least one free parameter required! "
120  "(Exposure/Background/Signal)"
121  << RESTendl;
122  return;
123  }
124 
125  if (nExpectedColumns != nTableColumns) {
126  std::string expectedColumns = "";
127  if (fExposureTime == 0) expectedColumns += "exposure";
128  if (!fSignal) {
129  if (fExposureTime == 0) expectedColumns += "/";
130  expectedColumns += "signal";
131  }
132  if (!fBackground) {
133  if (fExposureTime == 0 || !fSignal) expectedColumns += "/";
134  expectedColumns += "background";
135  }
136 
137  RESTError << "TRestExperimentList::InitFromConfigFile. Number of expected columns does not match "
138  "the number of table columns"
139  << RESTendl;
140  RESTError << "Number of table columns : " << nTableColumns << RESTendl;
141  RESTError << "Number of expected columns : " << nExpectedColumns << RESTendl;
142  RESTError << "Expected columns : " << expectedColumns << RESTendl;
143  return;
144  }
145 
147 
148  Bool_t generateMockData = false;
149  for (const auto& experimentRow : fExperimentsTable) {
150  TRestExperiment* experiment = new TRestExperiment();
151 
152  std::string rowStr = "";
153  for (const auto& el : experimentRow) {
154  rowStr += el + " ";
155  }
156 
157  RESTInfo << "TRestExperimentList. Loading experiment: " << rowStr << RESTendl;
158 
159  int column = 0;
160  if (fExposureTime == 0) {
161  if (REST_StringHelper::isANumber(experimentRow[column])) {
162  experiment->SetExposureInSeconds(
163  REST_StringHelper::StringToDouble(experimentRow[column]));
164  // We will generate mock data once we load the background component
165  generateMockData = true;
166  } else if (TRestTools::isRootFile(experimentRow[column])) {
167  // We load the file with the dataset into the experimental data
168  std::string fname = SearchFile(experimentRow[column]);
169  experiment->SetExperimentalDataSet(fname);
170  RESTWarning << "Loading experimental data havent been tested yet!" << RESTendl;
171  RESTWarning
172  << "It might require further development. Remove these lines once it works smooth!"
173  << RESTendl;
174  } else {
175  RESTError << experimentRow[column] << " is not a exposure time or an experimental dataset"
176  << RESTendl;
177  }
178  column++;
179  } else {
180  if (ToLower(fExposureStrategy) == "unique") {
181  experiment->SetExposureInSeconds(fExposureTime * units("s"));
182  // We will generate mock data once we load the background component
183  generateMockData = true;
184  }
185  }
186 
187  if (!fSignal) {
188  if (GetComponent(experimentRow[column])) {
189  TRestComponent* sgnl = (TRestComponent*)GetComponent(experimentRow[column])->Clone();
190  sgnl->SetName((TString)experimentRow[column]);
191  experiment->SetSignal(sgnl);
192  } else {
193  RESTError << "TRestExperimentList. Signal component : " << experimentRow[column]
194  << " not found!" << RESTendl;
195  }
196  column++;
197  } else {
198  experiment->SetSignal(fSignal);
199  }
200 
201  if (!fBackground) {
202  if (GetComponent(experimentRow[column])) {
203  TRestComponent* bck = (TRestComponent*)GetComponent(experimentRow[column])->Clone();
204  experiment->SetBackground(bck);
205  } else {
206  RESTError << "TRestExperimentList. Background component : " << experimentRow[column]
207  << " not found!" << RESTendl;
208  }
209  } else {
210  experiment->SetBackground(fBackground);
211  }
212 
213  if (generateMockData) {
214  RESTInfo << "TRestExperimentList. Generating mock dataset" << RESTendl;
215  experiment->GenerateMockDataSet(fUseAverage);
216  }
217 
218  if (experiment->GetSignal() && experiment->GetBackground()) {
219  experiment->SetName(experiment->GetSignal()->GetName());
220  fExperiments.push_back(experiment);
221  }
222  }
223 
224  if (fExposureTime > 0 && ToLower(fExposureStrategy) == "exp") {
226 
227  Double_t sum = 0;
228  for (size_t n = 0; n < fExperiments.size(); n++) sum += TMath::Exp((double)n * fExposureFactor);
229 
230  Double_t A = fExposureTime * units("s") / sum;
231  for (size_t n = 0; n < fExperiments.size(); n++) {
232  fExperiments[n]->SetExposureInSeconds(A * TMath::Exp((double)n * fExposureFactor));
233  fExperiments[n]->GenerateMockDataSet(fUseAverage);
234  }
235  }
236 
237  if (fExposureTime > 0 && ToLower(fExposureStrategy) == "power") {
239 
240  Double_t sum = 0;
241  for (size_t n = 0; n < fExperiments.size(); n++) sum += TMath::Power((double)n, fExposureFactor);
242 
243  Double_t A = fExposureTime * units("s") / sum;
244  for (size_t n = 0; n < fExperiments.size(); n++) {
245  fExperiments[n]->SetExposureInSeconds(A * TMath::Power((double)n, fExposureFactor));
246  fExperiments[n]->GenerateMockDataSet(fUseAverage);
247  }
248  }
249 
250  if (fExposureTime > 0 && ToLower(fExposureStrategy) == "equal") {
252 
253  for (size_t n = 0; n < fExperiments.size(); n++) {
254  fExperiments[n]->SetExposureInSeconds(fExposureTime * units("s") / fExperiments.size());
255  fExperiments[n]->GenerateMockDataSet(fUseAverage);
256  }
257  }
258  }
259 }
260 
269  fParameterizationNodes.clear();
270 
271  for (const auto& experiment : fExperiments) {
272  std::vector<Double_t> nodes = experiment->GetSignal()->GetParameterizationNodes();
273  fParameterizationNodes.insert(fParameterizationNodes.end(), nodes.begin(), nodes.end());
274 
275  std::sort(fParameterizationNodes.begin(), fParameterizationNodes.end());
276  auto last = std::unique(fParameterizationNodes.begin(), fParameterizationNodes.end());
278  }
279 }
280 
281 void TRestExperimentList::PrintParameterizationNodes() {
282  std::cout << "Experiment list nodes: ";
283  for (const auto& node : fParameterizationNodes) std::cout << node << "\t";
284  std::cout << std::endl;
285 }
286 
287 TRestComponent* TRestExperimentList::GetComponent(std::string compName) {
288  TRestComponent* component = nullptr;
289  for (const auto& c : fComponentFiles) {
290  TFile* f = TFile::Open(c.c_str(), "READ");
291  TObject* obj = f->Get((TString)compName);
292 
293  if (!obj) continue;
294 
295  if (obj->InheritsFrom("TRestComponent")) {
296  return (TRestComponent*)obj;
297  } else {
298  RESTError << "An object named : " << compName
299  << " exists inside the file, but it does not inherit from TRestComponent" << RESTendl;
300  }
301  }
302 
303  return component;
304 }
305 
311 
312  RESTMetadata << "Number of experiments loaded: " << fExperiments.size() << RESTendl;
313 
314  if (fUseAverage) RESTMetadata << "Average MonteCarlo counts generation was enabled" << RESTendl;
315 
316  RESTMetadata << "----" << RESTendl;
317 }
It defines a background/signal model distribution in a given parameter space (tipically x,...
A helper metadata class to create a list of TRestExperiment instances.
std::string fExperimentsFile
A file where we define experiment components, exposureTime, and tracking data of each experiment.
std::vector< TRestExperiment * > fExperiments
A vector with a list of experiments includes the background components in this model.
std::vector< Double_t > fParameterizationNodes
The fusioned list of parameterization nodes found at each experiment signal.
std::vector< std::string > fComponentFiles
A vector with filenames containing the components.
void InitFromConfigFile() override
It customizes the retrieval of XML data values of this class.
TRestComponent * fBackground
If not null this will be the common signal used in each experiment.
Bool_t fUseAverage
The mock dataset will be generated using the mean counts instead of a real MonteCarlo.
void Initialize() override
It will initialize the data frame with the filelist and column names (or observables) that have been ...
std::string fExposureStrategy
In case an exposure time is given it defines how to assign time to each experiment (equal/ksvz).
void ExtractExperimentParameterizationNodes()
It scans all the experiment signals parametric nodes to build a complete list of nodes used to build ...
Double_t fExposureFactor
The factor used on the exponential exposure time as a function of the experiment number.
std::string fComponentPattern
A fullpath filename pattern helping to initialize the component files vector.
TRestExperimentList()
Default constructor.
TRestComponent * fSignal
If not null this will be the common signal used in each experiment.
std::vector< std::vector< std::string > > fExperimentsTable
A table with the experiment file information.
~TRestExperimentList()
Default destructor.
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestAxionSolarFlux.
Double_t fExposureTime
If not zero this will be the common exposure time in micro-seconds (standard REST units)
It includes a model definition and experimental data used to obtain a final experimental sensitivity.
void SetExperimentalDataSet(const std::string &filename)
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 SearchFile(std::string filename)
Search files in current directory and directories specified in "searchPath" section.
std::string fConfigFileName
Full name of the rml file.
static std::vector< std::string > GetFilesMatchingPattern(std::string pattern, bool unlimited=false)
Returns a list of files whose name match the pattern string. Key word is "*". e.g....
Definition: TRestTools.cxx:976
static bool isRootFile(const std::string &filename)
Returns true if the filename has *.root* extension.
Definition: TRestTools.cxx:733
static int ReadASCIITable(std::string fName, std::vector< std::vector< Double_t >> &data, Int_t skipLines=0, std::string separator="\t")
Reads an ASCII file containing a table with values.
Definition: TRestTools.cxx:577
Double_t StringToDouble(std::string in)
Gets a double from a string.
std::string ToLower(std::string in)
Convert string to its lower case. Alternative of TString::ToLower.
Int_t isANumber(std::string in)
Returns 1 only if a valid number is found in the string in. If not it returns 0.
std::string ReplaceMathematicalExpressions(std::string buffer, Int_t precision=0, std::string errorMessage="")
Evaluates and replaces valid mathematical expressions found in the input string buffer.