REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestAnalysisTree.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 http://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 http://www.gnu.org/licenses/. *
20  * For the list of contributors see $REST_PATH/CREDITS. *
21  *************************************************************************/
22 
23 #ifndef RestCore_TRestAnalysisTree
24 #define RestCore_TRestAnalysisTree
25 
26 #include <TChain.h>
27 #include <TTree.h>
28 
29 #include <limits>
30 
31 #include "TRestEvent.h"
32 #include "TRestReflector.h"
33 
35 class TRestAnalysisTree : public TTree {
36  private:
37  Int_t fEventID;
38  Int_t fSubEventID;
39  Double_t fTimeStamp;
40  TString* fSubEventTag;
41  Int_t fRunOrigin;
42  Int_t fSubRunOrigin;
43 
44  //
45  Int_t fStatus = 0;
46  Int_t fSetObservableCalls = 0;
47  Int_t fSetObservableIndex = 0;
48  Bool_t fQuickSetObservableValue = true;
49  std::vector<RESTValue> fObservables;
50  std::map<std::string, int> fObservableIdMap;
51  std::map<std::string, int> fObservableIdSearchMap;
52  TChain* fChain = nullptr;
53 
54  // for storage
56  std::vector<TString> fObservableNames;
57  std::vector<TString> fObservableDescriptions;
58  std::vector<TString> fObservableTypes;
59 
60  void Initialize();
61  int EvaluateStatus();
62  void UpdateObservables();
63  void UpdateBranches();
64  void InitObservables();
65  void MakeObservableIdMap();
66  void ReadLeafValueToObservable(TLeaf* lf, RESTValue& obs);
67  bool BranchesExist() { return GetListOfBranches()->GetEntriesFast() > 0; }
68 
71  Error = -1,
73  None = 0,
78  Created = 1,
84  Retrieved = 2,
88  EmptyCloned = 3,
91  Connected = 4,
94  Filled = 5,
95  };
96 
97  TRestAnalysisTree(TTree* tree);
98 
99  protected:
100  public:
102 
103  // Get the status of this tree. This call will not evaluate the status.
104  inline int GetStatus() const { return fStatus; }
105  Int_t GetObservableID(const std::string& obsName);
106  Int_t GetMatchedObservableID(const std::string& obsName);
107  Bool_t ObservableExists(const std::string& obsName);
108  // six basic event parameters
109  Int_t GetEventID() { return fChain ? ((TRestAnalysisTree*)fChain->GetTree())->GetEventID() : fEventID; }
110  Int_t GetSubEventID() {
111  return fChain ? ((TRestAnalysisTree*)fChain->GetTree())->GetSubEventID() : fSubEventID;
112  }
113  Double_t GetTimeStamp() {
114  return fChain ? ((TRestAnalysisTree*)fChain->GetTree())->GetTimeStamp() : fTimeStamp;
115  }
116  TString GetSubEventTag() {
117  return fChain ? ((TRestAnalysisTree*)fChain->GetTree())->GetSubEventTag() : *fSubEventTag;
118  }
119  // we suppose all the chained trees have same run and sub run id.
120  // so there is no need to call fChain->GetTree()
121  Int_t GetRunOrigin() { return fRunOrigin; }
122  Int_t GetSubRunOrigin() { return fSubRunOrigin; }
123  Int_t GetNumberOfObservables() { return fNObservables; }
124 
125  // observable method
126  RESTValue GetObservable(const std::string& obsName);
127  RESTValue GetObservable(Int_t n);
128  TString GetObservableName(Int_t n);
129  TString GetObservableDescription(Int_t n);
130  TString GetObservableType(Int_t n);
131  TString GetObservableType(const std::string& obsName);
132  Double_t GetDblObservableValue(const std::string& obsName);
133  Double_t GetDblObservableValue(Int_t n);
134 
137  template <class T>
138  T GetObservableValue(Int_t n) {
139  // id check
140  if (n >= fNObservables) {
141  std::cout << "Error! TRestAnalysisTree::GetObservableValue(): index outside limits!" << std::endl;
142  return std::numeric_limits<T>::quiet_NaN();
143  }
144  if (fChain != nullptr) {
145  return ((TRestAnalysisTree*)fChain->GetTree())->GetObservableValue<T>(n);
146  }
147  return fObservables[n].GetValue<T>();
148  }
157  template <class T>
158  T GetObservableValue(std::string obsName) {
159  Int_t id = GetObservableID(obsName);
160  if (id == -1) {
161  // try to find matched observables
162  id = GetMatchedObservableID(obsName);
163  if (id == -1) {
164  return std::numeric_limits<T>::quiet_NaN();
165  }
166  }
167  return GetObservableValue<T>(id);
168  }
169 
172  template <class T>
173  void SetObservableValue(const Int_t& id, const T& value) {
174  // id check
175  if (id >= fNObservables) {
176  std::cout << "Error! TRestAnalysisTree::SetObservableValue(): index outside limits!" << std::endl;
177  return;
178  }
179  if (fChain != nullptr) {
180  std::cout << "Error! cannot set observable! AnalysisTree is in chain state" << std::endl;
181  return;
182  }
183  fObservables[id].SetValue(value);
184  }
223  template <class T>
224  void SetObservableValue(const std::string& name, const T& value) {
225  if (fQuickSetObservableValue && fStatus == Filled && fSetObservableCalls == fNObservables) {
226  SetObservableValue(fSetObservableIndex, value);
227  fSetObservableIndex++;
228  return;
229  }
230 
231  int id = GetObservableID(name);
232  if (id != -1) {
233  SetObservableValue(id, value);
234  } else {
235  AddObservable<T>(name);
236  id = GetObservableID(name);
237  if (id != -1) {
238  SetObservableValue(id, value);
239  fSetObservableCalls++;
240  }
241  }
242  }
243 
244  void SetObservable(Int_t id, RESTValue obs);
245  void SetObservable(std::string name, RESTValue value);
246 
247  void PrintObservables();
248  void PrintObservable(int N);
249 
250  void SetRunOrigin(Int_t run_origin) { fRunOrigin = run_origin; }
251  void SetSubRunOrigin(Int_t sub_run_origin) { fSubRunOrigin = sub_run_origin; }
252 
253  void SetEventInfo(TRestAnalysisTree* tree);
254  void SetEventInfo(TRestEvent* evt);
255  Int_t Fill();
256 
257  RESTValue AddObservable(const TString& observableName, const TString& observableType = "double",
258  const TString& description = "");
259  template <typename T>
260  T& AddObservable(TString observableName, TString description = "") {
261  return *(T*)AddObservable(observableName, REST_Reflection::GetTypeName<T>(), description);
262  }
263 
264  Int_t GetEntry(Long64_t entry = 0, Int_t getall = 0);
265 
266  Bool_t EvaluateCuts(const std::string& expression);
267  Bool_t EvaluateCut(const std::string& expression);
268 
269  std::vector<std::string> GetObservableNames();
270 
271  std::vector<std::string> GetCutObservables(const std::string& cut_str);
272 
273  void EnableBranches(std::vector<std::string> obsNames);
274  void DisableBranches(std::vector<std::string> obsNames);
275 
276  void EnableAllBranches();
277  void DisableAllBranches();
278 
281 
282  Double_t GetIntegral(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1) {
283  return GetObservableIntegral(obsName, xLow, xHigh);
284  }
285 
286  Double_t GetAverage(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1) {
287  return GetObservableAverage(obsName, xLow, xHigh);
288  }
289 
290  Double_t GetRMS(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1) {
291  return GetObservableRMS(obsName, xLow, xHigh);
292  }
293 
294  Double_t GetMinimum(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1) {
295  return GetObservableMinimum(obsName, xLow, xHigh);
296  }
297 
298  Double_t GetMaximum(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1) {
299  return GetObservableMaximum(obsName, xLow, xHigh);
300  }
301 
302  Double_t GetContour(const TString& obsName, const TString& obsIndexer, Double_t level = 0.5,
303  Int_t nBins = -1, Double_t xLow = -1, Double_t xHigh = -1) {
304  return GetObservableContour(obsName, obsIndexer, level, nBins, xLow, xHigh);
305  }
306 
307  Double_t GetObservableIntegral(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1);
308 
309  Double_t GetObservableAverage(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1);
310 
311  Double_t GetObservableRMS(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1);
312 
313  Double_t GetObservableMinimum(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1);
314  Double_t GetObservableMaximum(const TString& obsName, Double_t xLow = -1, Double_t xHigh = -1);
315 
316  Double_t GetObservableContour(const TString& obsName, const TString& obsIndexer, Double_t level = 0.5,
317  Int_t nBins = -1, Double_t xLow = -1, Double_t xHigh = -1);
318 
319  Int_t WriteAsTTree(const char* name = 0, Int_t option = 0, Int_t bufsize = 0);
320 
321  Bool_t AddChainFile(const std::string& file);
322 
323  TTree* GetTree() const;
324 
325  TChain* GetChain() { return fChain; }
326 
327  Long64_t LoadTree(Long64_t entry);
328 
329  Long64_t GetEntries() const;
330 
331  Long64_t GetEntries(const char* sel);
332 
333  void Browse(TBrowser* b) { fChain ? fChain->Browse(b) : TTree::Browse(b); }
334  Long64_t Draw(const char* varexp, const TCut& selection, Option_t* option = "",
335  Long64_t nentries = kMaxEntries, Long64_t firstentry = 0) {
336  return fChain ? fChain->Draw(varexp, selection, option, nentries, firstentry)
337  : TTree::Draw(varexp, selection, option, nentries, firstentry);
338  }
339  Long64_t Draw(const char* varexp, const char* selection, Option_t* option = "",
340  Long64_t nentries = kMaxEntries, Long64_t firstentry = 0) {
341  return fChain ? fChain->Draw(varexp, selection, option, nentries, firstentry)
342  : TTree::Draw(varexp, selection, option, nentries, firstentry);
343  }
344  void Draw(Option_t* opt) { fChain ? fChain->Draw(opt) : TTree::Draw(opt); }
345  TBranch* FindBranch(const char* name) {
346  return fChain ? fChain->FindBranch(name) : TTree::FindBranch(name);
347  }
348  TLeaf* FindLeaf(const char* name) { return fChain ? fChain->FindLeaf(name) : TTree::FindLeaf(name); }
349  TBranch* GetBranch(const char* name) { return fChain ? fChain->GetBranch(name) : TTree::GetBranch(name); }
350  Bool_t GetBranchStatus(const char* branchname) const {
351  return fChain ? fChain->GetBranchStatus(branchname) : TTree::GetBranchStatus(branchname);
352  }
353  Long64_t GetCacheSize() const { return fChain ? fChain->GetCacheSize() : TTree::GetCacheSize(); }
354  Long64_t GetChainEntryNumber(Long64_t entry) const {
355  return fChain ? fChain->GetChainEntryNumber(entry) : TTree::GetChainEntryNumber(entry);
356  }
357  Long64_t GetEntryNumber(Long64_t entry) const {
358  return fChain ? fChain->GetEntryNumber(entry) : TTree::GetEntryNumber(entry);
359  }
360  Long64_t GetReadEntry() const { return fChain ? fChain->GetReadEntry() : TTree::GetReadEntry(); }
361 
362  TH1* GetHistogram() { return fChain ? fChain->GetHistogram() : TTree::GetHistogram(); }
363  TLeaf* GetLeaf(const char* branchname, const char* leafname) {
364  return fChain ? fChain->GetLeaf(branchname, leafname) : TTree::GetLeaf(branchname, leafname);
365  }
366  TLeaf* GetLeaf(const char* name) { return fChain ? fChain->GetLeaf(name) : TTree::GetLeaf(name); }
367  TObjArray* GetListOfBranches() {
368  return fChain ? fChain->GetListOfBranches() : TTree::GetListOfBranches();
369  }
370  TObjArray* GetListOfLeaves() { return fChain ? fChain->GetListOfLeaves() : TTree::GetListOfLeaves(); }
371  Long64_t Process(const char* filename, Option_t* option = "", Long64_t nentries = kMaxEntries,
372  Long64_t firstentry = 0) {
373  return fChain ? fChain->Process(filename, option, nentries, firstentry)
374  : TTree::Process(filename, option, nentries, firstentry);
375  }
376  Long64_t Process(TSelector* selector, Option_t* option = "", Long64_t nentries = kMaxEntries,
377  Long64_t firstentry = 0) {
378  return fChain ? fChain->Process(selector, option, nentries, firstentry)
379  : TTree::Process(selector, option, nentries, firstentry);
380  }
381  Long64_t Scan(const char* varexp = "", const char* selection = "", Option_t* option = "",
382  Long64_t nentries = kMaxEntries, Long64_t firstentry = 0) {
383  return fChain ? fChain->Scan(varexp, selection, option, nentries, firstentry)
384  : TTree::Scan(varexp, selection, option, nentries, firstentry);
385  }
386  Int_t SetBranchAddress(const char* bname, void* add, TBranch** ptr = 0) {
387  return fChain ? fChain->SetBranchAddress(bname, add, ptr) : TTree::SetBranchAddress(bname, add, ptr);
388  }
389  Int_t SetBranchAddress(const char* bname, void* add, TBranch** ptr, TClass* realClass, EDataType datatype,
390  Bool_t isptr) {
391  return fChain ? fChain->SetBranchAddress(bname, add, ptr, realClass, datatype, isptr)
392  : TTree::SetBranchAddress(bname, add, ptr, realClass, datatype, isptr);
393  }
394  Int_t SetBranchAddress(const char* bname, void* add, TClass* realClass, EDataType datatype,
395  Bool_t isptr) {
396  return fChain ? fChain->SetBranchAddress(bname, add, realClass, datatype, isptr)
397  : TTree::SetBranchAddress(bname, add, realClass, datatype, isptr);
398  }
399  void SetBranchStatus(const char* bname, Bool_t status = 1, UInt_t* found = 0) {
400  fChain ? fChain->SetBranchStatus(bname, status, found) : TTree::SetBranchStatus(bname, status, found);
401  }
402  void SetDirectory(TDirectory* dir) { fChain ? fChain->SetDirectory(dir) : TTree::SetDirectory(dir); }
403 
404  void ResetBranchAddress(TBranch* br) {
405  fChain ? fChain->ResetBranchAddress(br) : TTree::ResetBranchAddress(br);
406  }
407  void ResetBranchAddresses() { fChain ? fChain->ResetBranchAddresses() : TTree::ResetBranchAddresses(); }
408 
409  // Constructor
411  TRestAnalysisTree(const TString& name, const TString& title);
412  static TRestAnalysisTree* ConvertFromTTree(TTree* tree);
413  // Destructor
414  virtual ~TRestAnalysisTree();
415 
416  ClassDef(TRestAnalysisTree, 3);
417 };
418 
419 #endif
REST core data-saving helper based on TTree.
void SetObservableValue(const std::string &name, const T &value)
Set the value of observable. May not check the name.
TString GetObservableDescription(Int_t n)
Get the observable description according to id.
T GetObservableValue(Int_t n)
Get observable in a given type, according to its id.
void DisableQuickObservableValueSetting()
It will disable quick observable value setting.
std::vector< std::string > GetCutObservables(const std::string &cut_str)
It will return a list with the names found in a string with conditions, as given in methods as Evalua...
void UpdateObservables()
Update observables stored in the tree.
Double_t GetObservableContour(const TString &obsName, const TString &obsIndexer, Double_t level=0.5, Int_t nBins=-1, Double_t xLow=-1, Double_t xHigh=-1)
This method generates a histogram of the the observable obsName given in the argument weighting it wi...
void EnableAllBranches()
It will enable all branches in the tree.
void EnableQuickObservableValueSetting()
It will enable quick observable value setting.
Long64_t LoadTree(Long64_t entry)
Overrides TTree::LoadTree(), to set current tree according to the given entry number,...
Bool_t AddChainFile(const std::string &file)
Add a series output file like TChain.
TString GetObservableName(Int_t n)
Get the observable name according to id.
void UpdateBranches()
Update branches in the tree.
Double_t GetObservableIntegral(const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
It returns the integral of the observable considering the given range. If no range is given the full ...
TChain * fChain
used for quick search of certain observables
Bool_t EvaluateCuts(const std::string &expression)
This method will receive a string with several analysis observables/cuts in the same fashion as used ...
TTree * GetTree() const
Overrides TTree::GetTree(), to get the actual tree used in case of chain operation(fCurrentTree !...
Double_t GetObservableMaximum(const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
It returns the maximum value of obsName considering the given range. If no range is given the full hi...
Int_t GetObservableID(const std::string &obsName)
Get the index of the specified observable.
void MakeObservableIdMap()
Update the map of observable name to observable id.
Int_t GetMatchedObservableID(const std::string &obsName)
Get the index of the observable, erasing the prefix if exists.
Double_t GetObservableAverage(const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
It returns the average of the observable considering the given range. If no range is given the full h...
void DisableBranches(std::vector< std::string > obsNames)
It will disable the branches given by argument.
Double_t GetObservableRMS(const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
It returns the RMS of the observable considering the given range. If no range is given the full histo...
Double_t GetDblObservableValue(const std::string &obsName)
Get double value of the observable, according to the name.
TString GetObservableType(Int_t n)
Get the observable type according to id.
void SetObservableValue(const Int_t &id, const T &value)
Set the value of observable whose index is id.
T GetObservableValue(std::string obsName)
Get observable in a given type, according to its name.
std::vector< std::string > GetObservableNames()
It returns a vector with strings containing all the observables that exist in the analysis tree.
@ Error
Error state.
void DisableAllBranches()
It will disable all branches in the tree.
Double_t GetObservableMinimum(const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
It returns the minimum value of obsName considering the given range. If no range is given the full hi...
int EvaluateStatus()
Evaluate the Status of this tree.
Bool_t ObservableExists(const std::string &obsName)
Get if the specified observable exists.
Int_t fNObservables
in case multiple files for reading
void EnableBranches(std::vector< std::string > obsNames)
It will enable the branches given by argument.
Bool_t EvaluateCut(const std::string &expression)
This method will evaluate a condition on one analysis observable constructed as "observable>value"....
void SetObservable(Int_t id, RESTValue obs)
Set the value of observable whose id is as specified.
Int_t WriteAsTTree(const char *name=0, Int_t option=0, Int_t bufsize=0)
It writes TRestAnalysisTree as TTree, in order to migrate files to pure ROOT users.
A base class for any REST event.
Definition: TRestEvent.h:38