REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes
TRestAnalysisTree Class Reference

Detailed Description

REST core data-saving helper based on TTree.

TRestAnalysisTree is a TTree but with custom objects for the branches that will be filled. The user will decide in each event data processing chain which branches/observables/variables will be finally added to the analysis tree. Inside a TRestAnalysisTree we find always the following six branches containing event information: runOrigin, subRunOrigin, eventID, subEventID, subEventTag and timeStamp. Those branches point to the corresponding class members inside TRestAnalysisTree, we name those branches the event branches. Additional branches can be added by the user in a processing chain by any class inheriting by TRestEventProcess. Those process generated branches will also point to some objects whose addresses are also stored in this class. Those objects are called observables.

In the traditional TTree case, the user defines multiple global variables, and adds branches with the address of these variables to the tree. Then the user changes the value of those variables somewhere in the code, and calls TTree::Fill() to create and save a new entry inside the tree.

In TRestAnalysisTree, the concept of "Branch" is weakened. We update the variables by invoking the SetObservableValue() method and then TRestAnalysisTree::Fill() to generate a new entry inside the tree. As soon as TRestEventProcess::SetObservable method is invoked, a new branch will be generated inside this tree. The code inside REST processes will be simplified while sacrificing a little performance. We can use temporary variable to set observable value directly. We can the focus on the analysis code inside each process, without caring about variable initialization before that.

As soon as TRestEventProcess::SetObservable method is invoked, a new branch will be generated inside this tree.

The following is a summary of speed of filling 1000000 entries for TTree and TRestAnalysisTree. Four observables and six event branches are added. We take the average of 3 tests as the result. See the file pipeline/analysistree/testspeed.cpp. for more details.

Condition time(us)
A. Do not use observable 846,522
B. Use quick observable (default) 1,188,232
C. Do not use quick observable 2,014,646
D. Use reflected observable 8,425,772
TTree 841,744

RESTsoft - Software for Rare Event Searches with TPCs

History of developments:

2016-Mar: First implementation 2019-May: Updated to support any type of observables 2020-Oct: Updated to be free from "Branch" concept

Author
Javier Galan, Ni Kaixiang

Definition at line 35 of file TRestAnalysisTree.h.

#include <TRestAnalysisTree.h>

Inheritance diagram for TRestAnalysisTree:

Public Member Functions

Bool_t AddChainFile (const std::string &file)
 Add a series output file like TChain. More...
 
RESTValue AddObservable (const TString &observableName, const TString &observableType="double", const TString &description="")
 
template<typename T >
T & AddObservable (TString observableName, TString description="")
 
void Browse (TBrowser *b)
 
void DisableAllBranches ()
 It will disable all branches in the tree. More...
 
void DisableBranches (std::vector< std::string > obsNames)
 It will disable the branches given by argument. More...
 
void DisableQuickObservableValueSetting ()
 It will disable quick observable value setting. More...
 
Long64_t Draw (const char *varexp, const char *selection, Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
 
Long64_t Draw (const char *varexp, const TCut &selection, Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
 
void Draw (Option_t *opt)
 
void EnableAllBranches ()
 It will enable all branches in the tree. More...
 
void EnableBranches (std::vector< std::string > obsNames)
 It will enable the branches given by argument. More...
 
void EnableQuickObservableValueSetting ()
 It will enable quick observable value setting. More...
 
Bool_t EvaluateCut (const std::string &expression)
 This method will evaluate a condition on one analysis observable constructed as "observable>value". For example, "rawAna_NumberOfSignals>10". More...
 
Bool_t EvaluateCuts (const std::string &expression)
 This method will receive a string with several analysis observables/cuts in the same fashion as used by ROOT. The string must be constructed as follows "obsName1>value1 && obsName2>value2 && ...". More...
 
Int_t Fill ()
 
TBranch * FindBranch (const char *name)
 
TLeaf * FindLeaf (const char *name)
 
Double_t GetAverage (const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
 
TBranch * GetBranch (const char *name)
 
Bool_t GetBranchStatus (const char *branchname) const
 
Long64_t GetCacheSize () const
 
TChain * GetChain ()
 
Long64_t GetChainEntryNumber (Long64_t entry) const
 
Double_t GetContour (const TString &obsName, const TString &obsIndexer, Double_t level=0.5, Int_t nBins=-1, Double_t xLow=-1, Double_t xHigh=-1)
 
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 EvaluateCuts. I.e. a construction as "obsName1==value1&&obsName2<=value2" will return {obsName1,obsName2}. More...
 
Double_t GetDblObservableValue (const std::string &obsName)
 Get double value of the observable, according to the name. More...
 
Double_t GetDblObservableValue (Int_t n)
 Get double value of the observable, according to the id. More...
 
Long64_t GetEntries () const
 
Long64_t GetEntries (const char *sel)
 
Int_t GetEntry (Long64_t entry=0, Int_t getall=0)
 
Long64_t GetEntryNumber (Long64_t entry) const
 
Int_t GetEventID ()
 
TH1 * GetHistogram ()
 
Double_t GetIntegral (const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
 
TLeaf * GetLeaf (const char *branchname, const char *leafname)
 
TLeaf * GetLeaf (const char *name)
 
TObjArray * GetListOfBranches ()
 
TObjArray * GetListOfLeaves ()
 
Int_t GetMatchedObservableID (const std::string &obsName)
 Get the index of the observable, erasing the prefix if exists. More...
 
Double_t GetMaximum (const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
 
Double_t GetMinimum (const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
 
Int_t GetNumberOfObservables ()
 
RESTValue GetObservable (const std::string &obsName)
 
RESTValue GetObservable (Int_t n)
 Get the observable object wrapped with "RESTValue" class, according to the id. More...
 
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 histogram range will be considered. More...
 
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 with a second observable also given by argument as obsWeight. More...
 
TString GetObservableDescription (Int_t n)
 Get the observable description according to id. More...
 
Int_t GetObservableID (const std::string &obsName)
 Get the index of the specified observable. More...
 
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 histogram range will be considered. More...
 
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 histogram range will be considered. More...
 
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 histogram range will be considered. More...
 
TString GetObservableName (Int_t n)
 Get the observable name according to id. More...
 
std::vector< std::string > GetObservableNames ()
 It returns a vector with strings containing all the observables that exist in the analysis tree. More...
 
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 histogram range will be considered. More...
 
TString GetObservableType (const std::string &obsName)
 Get the observable type according to its name. More...
 
TString GetObservableType (Int_t n)
 Get the observable type according to id. More...
 
template<class T >
GetObservableValue (Int_t n)
 Get observable in a given type, according to its id. More...
 
template<class T >
GetObservableValue (std::string obsName)
 Get observable in a given type, according to its name. More...
 
Long64_t GetReadEntry () const
 
Double_t GetRMS (const TString &obsName, Double_t xLow=-1, Double_t xHigh=-1)
 
Int_t GetRunOrigin ()
 
int GetStatus () const
 
Int_t GetSubEventID ()
 
TString GetSubEventTag ()
 
Int_t GetSubRunOrigin ()
 
Double_t GetTimeStamp ()
 
TTree * GetTree () const
 Overrides TTree::GetTree(), to get the actual tree used in case of chain operation(fCurrentTree != nullptr) More...
 
Long64_t LoadTree (Long64_t entry)
 Overrides TTree::LoadTree(), to set current tree according to the given entry number, in case of chain operation More...
 
Bool_t ObservableExists (const std::string &obsName)
 Get if the specified observable exists. More...
 
void PrintObservable (int N)
 
void PrintObservables ()
 
Long64_t Process (const char *filename, Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
 
Long64_t Process (TSelector *selector, Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
 
void ResetBranchAddress (TBranch *br)
 
void ResetBranchAddresses ()
 
Long64_t Scan (const char *varexp="", const char *selection="", Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
 
Int_t SetBranchAddress (const char *bname, void *add, TBranch **ptr, TClass *realClass, EDataType datatype, Bool_t isptr)
 
Int_t SetBranchAddress (const char *bname, void *add, TBranch **ptr=0)
 
Int_t SetBranchAddress (const char *bname, void *add, TClass *realClass, EDataType datatype, Bool_t isptr)
 
void SetBranchStatus (const char *bname, Bool_t status=1, UInt_t *found=0)
 
void SetDirectory (TDirectory *dir)
 
void SetEventInfo (TRestAnalysisTree *tree)
 
void SetEventInfo (TRestEvent *evt)
 
void SetObservable (Int_t id, RESTValue obs)
 Set the value of observable whose id is as specified. More...
 
void SetObservable (std::string name, RESTValue value)
 Set the value of observable whose name is as specified. More...
 
template<class T >
void SetObservableValue (const Int_t &id, const T &value)
 Set the value of observable whose index is id. More...
 
template<class T >
void SetObservableValue (const std::string &name, const T &value)
 Set the value of observable. May not check the name. More...
 
void SetRunOrigin (Int_t run_origin)
 
void SetSubRunOrigin (Int_t sub_run_origin)
 
 TRestAnalysisTree (const TString &name, const TString &title)
 
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. More...
 

Static Public Member Functions

static TRestAnalysisTreeConvertFromTTree (TTree *tree)
 

Private Types

enum  TRestAnalysisTree_Status {
  Error = -1 , None = 0 , Created = 1 , Retrieved = 2 ,
  EmptyCloned = 3 , Connected = 4 , Filled = 5
}
 

Private Member Functions

bool BranchesExist ()
 
int EvaluateStatus ()
 Evaluate the Status of this tree. More...
 
void Initialize ()
 
void InitObservables ()
 
void MakeObservableIdMap ()
 Update the map of observable name to observable id. More...
 
void ReadLeafValueToObservable (TLeaf *lf, RESTValue &obs)
 
 TRestAnalysisTree (TTree *tree)
 
void UpdateBranches ()
 Update branches in the tree. More...
 
void UpdateObservables ()
 Update observables stored in the tree. More...
 

Private Attributes

TChain * fChain = nullptr
 used for quick search of certain observables More...
 
Int_t fEventID
 
Int_t fNObservables
 in case multiple files for reading More...
 
std::vector< TString > fObservableDescriptions
 
std::map< std::string, int > fObservableIdMap
 
std::map< std::string, int > fObservableIdSearchMap
 
std::vector< TString > fObservableNames
 
std::vector< RESTValuefObservables
 
std::vector< TString > fObservableTypes
 
Bool_t fQuickSetObservableValue = true
 
Int_t fRunOrigin
 
Int_t fSetObservableCalls = 0
 
Int_t fSetObservableIndex = 0
 
Int_t fStatus = 0
 
Int_t fSubEventID
 
TString * fSubEventTag
 
Int_t fSubRunOrigin
 
Double_t fTimeStamp
 

Member Enumeration Documentation

◆ TRestAnalysisTree_Status

Enumerator
Error 

Error state.

First status when constructed as empty tree

None 

Status during the first loop of SetObservableValue() before Fill() and after construction. No branches are created and filled. Some observables may be added. We will create branches when we call Fill(). We will update branches when we call PrintObservables(). The status will become "Filled" and "Connected" respectivelly

Created 

There are branches in the tree with data, but no observables are added. This means the tree was just retrieved from root file. In this case, when calling GetEntry() or Fill(), we will first recover the observables and connect them to the branches. We forbid to add new observable to the tree. After reconnection, the status will become "Filled"

Retrieved 

There are branches in the tree but with no data(i.e. entry number is 0). No observables are added either. This means the tree could be an empty clone of another TRestAnalysisTree. Adding new observable is allowed in this case.

EmptyCloned 

There are branches in the tree with no data. There are observables in the list. same as "Created"

Connected 

There are branches in the tree with data. There are observables in the list. Once filled, it is forbidden to add new observable to the tree.

Definition at line 69 of file TRestAnalysisTree.h.

Constructor & Destructor Documentation

◆ TRestAnalysisTree() [1/3]

TRestAnalysisTree::TRestAnalysisTree ( TTree *  tree)
private

Definition at line 107 of file TRestAnalysisTree.cxx.

◆ TRestAnalysisTree() [2/3]

TRestAnalysisTree::TRestAnalysisTree ( )

Definition at line 100 of file TRestAnalysisTree.cxx.

◆ TRestAnalysisTree() [3/3]

TRestAnalysisTree::TRestAnalysisTree ( const TString &  name,
const TString &  title 
)

Definition at line 133 of file TRestAnalysisTree.cxx.

◆ ~TRestAnalysisTree()

TRestAnalysisTree::~TRestAnalysisTree ( )
virtual

Definition at line 1300 of file TRestAnalysisTree.cxx.

Member Function Documentation

◆ AddChainFile()

Bool_t TRestAnalysisTree::AddChainFile ( const std::string &  file)

Add a series output file like TChain.

Parameters
fileThe input file that contains another AnalysisTree with same run id
Returns

Definition at line 1219 of file TRestAnalysisTree.cxx.

◆ AddObservable() [1/2]

RESTValue TRestAnalysisTree::AddObservable ( const TString &  observableName,
const TString &  observableType = "double",
const TString &  description = "" 
)

Definition at line 634 of file TRestAnalysisTree.cxx.

◆ AddObservable() [2/2]

template<typename T >
T & TRestAnalysisTree::AddObservable ( TString  observableName,
TString  description = "" 
)
inline

Definition at line 260 of file TRestAnalysisTree.h.

◆ BranchesExist()

bool TRestAnalysisTree::BranchesExist ( )
inlineprivate

Definition at line 67 of file TRestAnalysisTree.h.

◆ Browse()

void TRestAnalysisTree::Browse ( TBrowser *  b)
inline

Definition at line 333 of file TRestAnalysisTree.h.

◆ ConvertFromTTree()

TRestAnalysisTree * TRestAnalysisTree::ConvertFromTTree ( TTree *  tree)
static

Definition at line 125 of file TRestAnalysisTree.cxx.

◆ DisableAllBranches()

void TRestAnalysisTree::DisableAllBranches ( )

It will disable all branches in the tree.

Definition at line 975 of file TRestAnalysisTree.cxx.

◆ DisableBranches()

void TRestAnalysisTree::DisableBranches ( std::vector< std::string >  obsNames)

It will disable the branches given by argument.

Definition at line 963 of file TRestAnalysisTree.cxx.

◆ DisableQuickObservableValueSetting()

void TRestAnalysisTree::DisableQuickObservableValueSetting ( )

It will disable quick observable value setting.

Definition at line 988 of file TRestAnalysisTree.cxx.

◆ Draw() [1/3]

Long64_t TRestAnalysisTree::Draw ( const char *  varexp,
const char *  selection,
Option_t *  option = "",
Long64_t  nentries = kMaxEntries,
Long64_t  firstentry = 0 
)
inline

Definition at line 339 of file TRestAnalysisTree.h.

◆ Draw() [2/3]

Long64_t TRestAnalysisTree::Draw ( const char *  varexp,
const TCut &  selection,
Option_t *  option = "",
Long64_t  nentries = kMaxEntries,
Long64_t  firstentry = 0 
)
inline

Definition at line 334 of file TRestAnalysisTree.h.

◆ Draw() [3/3]

void TRestAnalysisTree::Draw ( Option_t *  opt)
inline

Definition at line 344 of file TRestAnalysisTree.h.

◆ EnableAllBranches()

void TRestAnalysisTree::EnableAllBranches ( )

It will enable all branches in the tree.

Definition at line 970 of file TRestAnalysisTree.cxx.

◆ EnableBranches()

void TRestAnalysisTree::EnableBranches ( std::vector< std::string >  obsNames)

It will enable the branches given by argument.

Definition at line 956 of file TRestAnalysisTree.cxx.

◆ EnableQuickObservableValueSetting()

void TRestAnalysisTree::EnableQuickObservableValueSetting ( )

It will enable quick observable value setting.

When enabled, when calling SetObservableValue(string, value), we skip the index searching and directly set the observable at n-th call of the method. This is possible when multiple SetObservableValue() is called in sequence in a loop. This will significantally improve the speed of that method.

Definition at line 984 of file TRestAnalysisTree.cxx.

◆ EvaluateCut()

Bool_t TRestAnalysisTree::EvaluateCut ( const std::string &  expression)

This method will evaluate a condition on one analysis observable constructed as "observable>value". For example, "rawAna_NumberOfSignals>10".

It will evaluate the given conditions and return the result. Valid operators are "==", "<=", ">=", "!=", "=", ">" and "<".

Definition at line 900 of file TRestAnalysisTree.cxx.

◆ EvaluateCuts()

Bool_t TRestAnalysisTree::EvaluateCuts ( const std::string &  expression)

This method will receive a string with several analysis observables/cuts in the same fashion as used by ROOT. The string must be constructed as follows "obsName1>value1 && obsName2>value2 && ...".

It will evaluate the given conditions and return the result. Valid operators are "==", "<=", ">=", "!=", "=", ">" and "<".

Definition at line 882 of file TRestAnalysisTree.cxx.

◆ EvaluateStatus()

int TRestAnalysisTree::EvaluateStatus ( )
private

Evaluate the Status of this tree.

Definition at line 270 of file TRestAnalysisTree.cxx.

◆ Fill()

Int_t TRestAnalysisTree::Fill ( )

Definition at line 768 of file TRestAnalysisTree.cxx.

◆ FindBranch()

TBranch * TRestAnalysisTree::FindBranch ( const char *  name)
inline

Definition at line 345 of file TRestAnalysisTree.h.

◆ FindLeaf()

TLeaf * TRestAnalysisTree::FindLeaf ( const char *  name)
inline

Definition at line 348 of file TRestAnalysisTree.h.

◆ GetAverage()

Double_t TRestAnalysisTree::GetAverage ( const TString &  obsName,
Double_t  xLow = -1,
Double_t  xHigh = -1 
)
inline

Definition at line 286 of file TRestAnalysisTree.h.

◆ GetBranch()

TBranch * TRestAnalysisTree::GetBranch ( const char *  name)
inline

Definition at line 349 of file TRestAnalysisTree.h.

◆ GetBranchStatus()

Bool_t TRestAnalysisTree::GetBranchStatus ( const char *  branchname) const
inline

Definition at line 350 of file TRestAnalysisTree.h.

◆ GetCacheSize()

Long64_t TRestAnalysisTree::GetCacheSize ( ) const
inline

Definition at line 353 of file TRestAnalysisTree.h.

◆ GetChain()

TChain * TRestAnalysisTree::GetChain ( )
inline

Definition at line 325 of file TRestAnalysisTree.h.

◆ GetChainEntryNumber()

Long64_t TRestAnalysisTree::GetChainEntryNumber ( Long64_t  entry) const
inline

Definition at line 354 of file TRestAnalysisTree.h.

◆ GetContour()

Double_t TRestAnalysisTree::GetContour ( const TString &  obsName,
const TString &  obsIndexer,
Double_t  level = 0.5,
Int_t  nBins = -1,
Double_t  xLow = -1,
Double_t  xHigh = -1 
)
inline

Definition at line 302 of file TRestAnalysisTree.h.

◆ GetCutObservables()

vector< string > TRestAnalysisTree::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 EvaluateCuts. I.e. a construction as "obsName1==value1&&obsName2<=value2" will return {obsName1,obsName2}.

Definition at line 935 of file TRestAnalysisTree.cxx.

◆ GetDblObservableValue() [1/2]

Double_t TRestAnalysisTree::GetDblObservableValue ( const std::string &  obsName)

Get double value of the observable, according to the name.

It assumes the observable is in double/int type. If not it will print error and return 0

Definition at line 612 of file TRestAnalysisTree.cxx.

◆ GetDblObservableValue() [2/2]

Double_t TRestAnalysisTree::GetDblObservableValue ( Int_t  n)

Get double value of the observable, according to the id.

It assumes the observable is in double/int type. If not it will print error and return 0

Definition at line 619 of file TRestAnalysisTree.cxx.

◆ GetEntries() [1/2]

Long64_t TRestAnalysisTree::GetEntries ( ) const

Definition at line 1286 of file TRestAnalysisTree.cxx.

◆ GetEntries() [2/2]

Long64_t TRestAnalysisTree::GetEntries ( const char *  sel)

Definition at line 1293 of file TRestAnalysisTree.cxx.

◆ GetEntry()

Int_t TRestAnalysisTree::GetEntry ( Long64_t  entry = 0,
Int_t  getall = 0 
)

Definition at line 716 of file TRestAnalysisTree.cxx.

◆ GetEntryNumber()

Long64_t TRestAnalysisTree::GetEntryNumber ( Long64_t  entry) const
inline

Definition at line 357 of file TRestAnalysisTree.h.

◆ GetEventID()

Int_t TRestAnalysisTree::GetEventID ( )
inline

Definition at line 109 of file TRestAnalysisTree.h.

◆ GetHistogram()

TH1 * TRestAnalysisTree::GetHistogram ( )
inline

Definition at line 362 of file TRestAnalysisTree.h.

◆ GetIntegral()

Double_t TRestAnalysisTree::GetIntegral ( const TString &  obsName,
Double_t  xLow = -1,
Double_t  xHigh = -1 
)
inline

Definition at line 282 of file TRestAnalysisTree.h.

◆ GetLeaf() [1/2]

TLeaf * TRestAnalysisTree::GetLeaf ( const char *  branchname,
const char *  leafname 
)
inline

Definition at line 363 of file TRestAnalysisTree.h.

◆ GetLeaf() [2/2]

TLeaf * TRestAnalysisTree::GetLeaf ( const char *  name)
inline

Definition at line 366 of file TRestAnalysisTree.h.

◆ GetListOfBranches()

TObjArray * TRestAnalysisTree::GetListOfBranches ( )
inline

Definition at line 367 of file TRestAnalysisTree.h.

◆ GetListOfLeaves()

TObjArray * TRestAnalysisTree::GetListOfLeaves ( )
inline

Definition at line 370 of file TRestAnalysisTree.h.

◆ GetMatchedObservableID()

Int_t TRestAnalysisTree::GetMatchedObservableID ( const std::string &  obsName)

Get the index of the observable, erasing the prefix if exists.

Ignores prefix like "sAna_". Case sensitive, misspelling prompted. If not exist, it will return -1.

Definition at line 181 of file TRestAnalysisTree.cxx.

◆ GetMaximum()

Double_t TRestAnalysisTree::GetMaximum ( const TString &  obsName,
Double_t  xLow = -1,
Double_t  xHigh = -1 
)
inline

Definition at line 298 of file TRestAnalysisTree.h.

◆ GetMinimum()

Double_t TRestAnalysisTree::GetMinimum ( const TString &  obsName,
Double_t  xLow = -1,
Double_t  xHigh = -1 
)
inline

Definition at line 294 of file TRestAnalysisTree.h.

◆ GetNumberOfObservables()

Int_t TRestAnalysisTree::GetNumberOfObservables ( )
inline

Definition at line 123 of file TRestAnalysisTree.h.

◆ GetObservable() [1/2]

RESTValue TRestAnalysisTree::GetObservable ( const std::string &  obsName)

Definition at line 554 of file TRestAnalysisTree.cxx.

◆ GetObservable() [2/2]

RESTValue TRestAnalysisTree::GetObservable ( Int_t  n)

Get the observable object wrapped with "RESTValue" class, according to the id.

Definition at line 562 of file TRestAnalysisTree.cxx.

◆ GetObservableAverage()

Double_t TRestAnalysisTree::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 histogram range will be considered.

Definition at line 1019 of file TRestAnalysisTree.cxx.

◆ GetObservableContour()

Double_t TRestAnalysisTree::GetObservableContour ( const TString &  obsName,
const TString &  obsWeight,
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 with a second observable also given by argument as obsWeight.

This method will return the value at which obsName integral reaches a fraction of the total integral defined by the argument level. E.g. if level=0.5, then the value of obsName at which the histogram obsName reaches half the integral is returned.

If not given the default level value is 0.5.

Optionally we may define the parameters of the histogram. If not, ROOT will use the default values defined by the user.

For example, we could bin the variable final_R between 0 and 1cm in 1000 bins, and get the value of final_R where optics_efficiency integrated events is 80% of the total.

analysisTree->GetObservableContour("final_R", "optics_efficiency", 0.8, 1000, 0, 1)
Returns
Returns the value at which obsName contains a level fraction of the integral of the obsWeight observable.

Definition at line 1144 of file TRestAnalysisTree.cxx.

◆ GetObservableDescription()

TString TRestAnalysisTree::GetObservableDescription ( Int_t  n)

Get the observable description according to id.

Definition at line 585 of file TRestAnalysisTree.cxx.

◆ GetObservableID()

Int_t TRestAnalysisTree::GetObservableID ( const std::string &  obsName)

Get the index of the specified observable.

If not exist, it will return -1. It will call MakeObservableIdMap() to update observable id map before searching

Definition at line 165 of file TRestAnalysisTree.cxx.

◆ GetObservableIntegral()

Double_t TRestAnalysisTree::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 histogram range will be considered.

Definition at line 994 of file TRestAnalysisTree.cxx.

◆ GetObservableMaximum()

Double_t TRestAnalysisTree::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 histogram range will be considered.

Definition at line 1074 of file TRestAnalysisTree.cxx.

◆ GetObservableMinimum()

Double_t TRestAnalysisTree::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 histogram range will be considered.

Definition at line 1099 of file TRestAnalysisTree.cxx.

◆ GetObservableName()

TString TRestAnalysisTree::GetObservableName ( Int_t  n)

Get the observable name according to id.

Definition at line 576 of file TRestAnalysisTree.cxx.

◆ GetObservableNames()

std::vector< std::string > TRestAnalysisTree::GetObservableNames ( )

It returns a vector with strings containing all the observables that exist in the analysis tree.

Definition at line 1176 of file TRestAnalysisTree.cxx.

◆ GetObservableRMS()

Double_t TRestAnalysisTree::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 histogram range will be considered.

Definition at line 1048 of file TRestAnalysisTree.cxx.

◆ GetObservableType() [1/2]

TString TRestAnalysisTree::GetObservableType ( const std::string &  obsName)

Get the observable type according to its name.

Definition at line 603 of file TRestAnalysisTree.cxx.

◆ GetObservableType() [2/2]

TString TRestAnalysisTree::GetObservableType ( Int_t  n)

Get the observable type according to id.

Definition at line 594 of file TRestAnalysisTree.cxx.

◆ GetObservableValue() [1/2]

template<class T >
T TRestAnalysisTree::GetObservableValue ( Int_t  n)
inline

Get observable in a given type, according to its id.

Definition at line 138 of file TRestAnalysisTree.h.

◆ GetObservableValue() [2/2]

template<class T >
T TRestAnalysisTree::GetObservableValue ( std::string  obsName)
inline

Get observable in a given type, according to its name.

The returned value is directly the type.

Example: std::vector<int> v = AnalysisTree->GetObservableValue<std::vector<int>>("myvec1"); double a = AnalysisTree->GetObservableValue<double>("myval");

Definition at line 158 of file TRestAnalysisTree.h.

◆ GetReadEntry()

Long64_t TRestAnalysisTree::GetReadEntry ( ) const
inline

Definition at line 360 of file TRestAnalysisTree.h.

◆ GetRMS()

Double_t TRestAnalysisTree::GetRMS ( const TString &  obsName,
Double_t  xLow = -1,
Double_t  xHigh = -1 
)
inline

Definition at line 290 of file TRestAnalysisTree.h.

◆ GetRunOrigin()

Int_t TRestAnalysisTree::GetRunOrigin ( )
inline

Definition at line 121 of file TRestAnalysisTree.h.

◆ GetStatus()

int TRestAnalysisTree::GetStatus ( ) const
inline

Definition at line 104 of file TRestAnalysisTree.h.

◆ GetSubEventID()

Int_t TRestAnalysisTree::GetSubEventID ( )
inline

Definition at line 110 of file TRestAnalysisTree.h.

◆ GetSubEventTag()

TString TRestAnalysisTree::GetSubEventTag ( )
inline

Definition at line 116 of file TRestAnalysisTree.h.

◆ GetSubRunOrigin()

Int_t TRestAnalysisTree::GetSubRunOrigin ( )
inline

Definition at line 122 of file TRestAnalysisTree.h.

◆ GetTimeStamp()

Double_t TRestAnalysisTree::GetTimeStamp ( )
inline

Definition at line 113 of file TRestAnalysisTree.h.

◆ GetTree()

TTree * TRestAnalysisTree::GetTree ( ) const

Overrides TTree::GetTree(), to get the actual tree used in case of chain operation(fCurrentTree != nullptr)

Returns

Definition at line 1263 of file TRestAnalysisTree.cxx.

◆ Initialize()

void TRestAnalysisTree::Initialize ( )
private

Definition at line 140 of file TRestAnalysisTree.cxx.

◆ InitObservables()

void TRestAnalysisTree::InitObservables ( )
private

Definition at line 445 of file TRestAnalysisTree.cxx.

◆ LoadTree()

Long64_t TRestAnalysisTree::LoadTree ( Long64_t  entry)

Overrides TTree::LoadTree(), to set current tree according to the given entry number, in case of chain operation

Parameters
entry
Returns

Definition at line 1276 of file TRestAnalysisTree.cxx.

◆ MakeObservableIdMap()

void TRestAnalysisTree::MakeObservableIdMap ( )
private

Update the map of observable name to observable id.

Using map will improve the speed of "SetObservableValue"

Definition at line 459 of file TRestAnalysisTree.cxx.

◆ ObservableExists()

Bool_t TRestAnalysisTree::ObservableExists ( const std::string &  obsName)

Get if the specified observable exists.

It will call MakeObservableIdMap() to update observable id map before searching

Definition at line 243 of file TRestAnalysisTree.cxx.

◆ PrintObservable()

void TRestAnalysisTree::PrintObservable ( int  N)

Definition at line 705 of file TRestAnalysisTree.cxx.

◆ PrintObservables()

void TRestAnalysisTree::PrintObservables ( )

Definition at line 676 of file TRestAnalysisTree.cxx.

◆ Process() [1/2]

Long64_t TRestAnalysisTree::Process ( const char *  filename,
Option_t *  option = "",
Long64_t  nentries = kMaxEntries,
Long64_t  firstentry = 0 
)
inline

Definition at line 371 of file TRestAnalysisTree.h.

◆ Process() [2/2]

Long64_t TRestAnalysisTree::Process ( TSelector *  selector,
Option_t *  option = "",
Long64_t  nentries = kMaxEntries,
Long64_t  firstentry = 0 
)
inline

Definition at line 376 of file TRestAnalysisTree.h.

◆ ReadLeafValueToObservable()

void TRestAnalysisTree::ReadLeafValueToObservable ( TLeaf *  lf,
RESTValue obs 
)
private

Definition at line 471 of file TRestAnalysisTree.cxx.

◆ ResetBranchAddress()

void TRestAnalysisTree::ResetBranchAddress ( TBranch *  br)
inline

Definition at line 404 of file TRestAnalysisTree.h.

◆ ResetBranchAddresses()

void TRestAnalysisTree::ResetBranchAddresses ( )
inline

Definition at line 407 of file TRestAnalysisTree.h.

◆ Scan()

Long64_t TRestAnalysisTree::Scan ( const char *  varexp = "",
const char *  selection = "",
Option_t *  option = "",
Long64_t  nentries = kMaxEntries,
Long64_t  firstentry = 0 
)
inline

Definition at line 381 of file TRestAnalysisTree.h.

◆ SetBranchAddress() [1/3]

Int_t TRestAnalysisTree::SetBranchAddress ( const char *  bname,
void *  add,
TBranch **  ptr,
TClass *  realClass,
EDataType  datatype,
Bool_t  isptr 
)
inline

Definition at line 389 of file TRestAnalysisTree.h.

◆ SetBranchAddress() [2/3]

Int_t TRestAnalysisTree::SetBranchAddress ( const char *  bname,
void *  add,
TBranch **  ptr = 0 
)
inline

Definition at line 386 of file TRestAnalysisTree.h.

◆ SetBranchAddress() [3/3]

Int_t TRestAnalysisTree::SetBranchAddress ( const char *  bname,
void *  add,
TClass *  realClass,
EDataType  datatype,
Bool_t  isptr 
)
inline

Definition at line 394 of file TRestAnalysisTree.h.

◆ SetBranchStatus()

void TRestAnalysisTree::SetBranchStatus ( const char *  bname,
Bool_t  status = 1,
UInt_t *  found = 0 
)
inline

Definition at line 399 of file TRestAnalysisTree.h.

◆ SetDirectory()

void TRestAnalysisTree::SetDirectory ( TDirectory *  dir)
inline

Definition at line 402 of file TRestAnalysisTree.h.

◆ SetEventInfo() [1/2]

void TRestAnalysisTree::SetEventInfo ( TRestAnalysisTree tree)

Definition at line 736 of file TRestAnalysisTree.cxx.

◆ SetEventInfo() [2/2]

void TRestAnalysisTree::SetEventInfo ( TRestEvent evt)

Definition at line 752 of file TRestAnalysisTree.cxx.

◆ SetObservable() [1/2]

void TRestAnalysisTree::SetObservable ( Int_t  id,
RESTValue  obs 
)

Set the value of observable whose id is as specified.

Definition at line 807 of file TRestAnalysisTree.cxx.

◆ SetObservable() [2/2]

void TRestAnalysisTree::SetObservable ( std::string  name,
RESTValue  value 
)

Set the value of observable whose name is as specified.

The input type is "RESTValue". This class is able to be constructed from any object, Therefore the input could be equal to SetObservableValue(). But in this case it would be much slower than SetObservableValue(), since converting to "RESTValue" contains a type name reflection which takes time. This method is better used to copy directly "RESTValue" objects, e.g., observables from another TRestAnalysisTree.

Example:

TFile* fin = new TFile("simpleTree.root", "open");
TRestAnalysisTree* treein = (TRestAnalysisTree*)fin->Get("AnalysisTree");
TRestAnalysisTree* treeout = (TRestAnalysisTree*)tmp_tree->CloneTree(0);
RESTValue c = treein->GetObservable("vec");
for (int i = 0; i < 10; i++) {
treein->GetEntry(i);
treeout->SetObservable("interitedvec", c);
treeout->Fill();
}
REST core data-saving helper based on TTree.
void SetObservable(Int_t id, RESTValue obs)
Set the value of observable whose id is as specified.

Definition at line 870 of file TRestAnalysisTree.cxx.

◆ SetObservableValue() [1/2]

template<class T >
void TRestAnalysisTree::SetObservableValue ( const Int_t &  id,
const T &  value 
)
inline

Set the value of observable whose index is id.

Definition at line 173 of file TRestAnalysisTree.h.

◆ SetObservableValue() [2/2]

template<class T >
void TRestAnalysisTree::SetObservableValue ( const std::string &  name,
const T &  value 
)
inline

Set the value of observable. May not check the name.

Any type of input value is accepted. We can directly set value from a std::vector or std::map object. If the observable does not exist, it will create a new one if the tree is not filled yet. If fQuickSetObservableValue == true we directly set the observable whose index is equal to fSetObservableIndex, and increase fSetObservableIndex at each call of this method. Otherwise we first find the index of the observable with given name, and then set its value.

If the method is not called in the loop with equal times(e.g. when manually typing the code in ROOT prompt), one needs to disable quick observable value setting before calling the method, or consider to use SetObservable()

Example:

for (int i = 0; i < 1000000; i++) {
tree->SetObservableValue("gaus", gRandom->Gaus(100, 20));
tree->SetObservableValue("poisson", gRandom->Poisson(36));
tree->SetObservableValue("rndm", gRandom->Rndm());
tree->SetObservableValue("landau", gRandom->Landau(10, 2));
tree->Fill();
}
tree->SetObservableValue("myval", 20);
tree->SetObservableValue("myvec", std::vector<int>{11,23,37,41});
tree->Fill();
tree->SetObservableValue("myvec", std::vector<int>{2,3,5,7});
tree->SetObservableValue("myval", 30);
tree->Fill();
tree->SetObservableValue("myval", 40);
tree->Fill();
void DisableQuickObservableValueSetting()
It will disable quick observable value setting.
void SetObservableValue(const Int_t &id, const T &value)
Set the value of observable whose index is id.

Definition at line 224 of file TRestAnalysisTree.h.

◆ SetRunOrigin()

void TRestAnalysisTree::SetRunOrigin ( Int_t  run_origin)
inline

Definition at line 250 of file TRestAnalysisTree.h.

◆ SetSubRunOrigin()

void TRestAnalysisTree::SetSubRunOrigin ( Int_t  sub_run_origin)
inline

Definition at line 251 of file TRestAnalysisTree.h.

◆ UpdateBranches()

void TRestAnalysisTree::UpdateBranches ( )
private

Update branches in the tree.

It will create branches if not exist, both for event data members and observables. Note that this method can be called multiple times during the first loop of observable setting. After branch creation, this method will change status 1->4, or stay 4.

Definition at line 391 of file TRestAnalysisTree.cxx.

◆ UpdateObservables()

void TRestAnalysisTree::UpdateObservables ( )
private

Update observables stored in the tree.

It will first connect the six event data members(i.e. runid, eventid, etc.) to the existing TTree branches. Then it will create new observable objects by reflection, and connect them also to the existing TTree branches. After re-connection, this method will change status 2->5, 3->4

Definition at line 327 of file TRestAnalysisTree.cxx.

◆ WriteAsTTree()

Int_t TRestAnalysisTree::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.

Definition at line 1190 of file TRestAnalysisTree.cxx.

Field Documentation

◆ fChain

TChain* TRestAnalysisTree::fChain = nullptr
private

used for quick search of certain observables

Definition at line 52 of file TRestAnalysisTree.h.

◆ fEventID

Int_t TRestAnalysisTree::fEventID
private

Definition at line 37 of file TRestAnalysisTree.h.

◆ fNObservables

Int_t TRestAnalysisTree::fNObservables
private

in case multiple files for reading

Definition at line 55 of file TRestAnalysisTree.h.

◆ fObservableDescriptions

std::vector<TString> TRestAnalysisTree::fObservableDescriptions
private

Definition at line 57 of file TRestAnalysisTree.h.

◆ fObservableIdMap

std::map<std::string, int> TRestAnalysisTree::fObservableIdMap
private

Definition at line 50 of file TRestAnalysisTree.h.

◆ fObservableIdSearchMap

std::map<std::string, int> TRestAnalysisTree::fObservableIdSearchMap
private

Definition at line 51 of file TRestAnalysisTree.h.

◆ fObservableNames

std::vector<TString> TRestAnalysisTree::fObservableNames
private

Definition at line 56 of file TRestAnalysisTree.h.

◆ fObservables

std::vector<RESTValue> TRestAnalysisTree::fObservables
private

Definition at line 49 of file TRestAnalysisTree.h.

◆ fObservableTypes

std::vector<TString> TRestAnalysisTree::fObservableTypes
private

Definition at line 58 of file TRestAnalysisTree.h.

◆ fQuickSetObservableValue

Bool_t TRestAnalysisTree::fQuickSetObservableValue = true
private

Definition at line 48 of file TRestAnalysisTree.h.

◆ fRunOrigin

Int_t TRestAnalysisTree::fRunOrigin
private

Definition at line 41 of file TRestAnalysisTree.h.

◆ fSetObservableCalls

Int_t TRestAnalysisTree::fSetObservableCalls = 0
private

Definition at line 46 of file TRestAnalysisTree.h.

◆ fSetObservableIndex

Int_t TRestAnalysisTree::fSetObservableIndex = 0
private

Definition at line 47 of file TRestAnalysisTree.h.

◆ fStatus

Int_t TRestAnalysisTree::fStatus = 0
private

Definition at line 45 of file TRestAnalysisTree.h.

◆ fSubEventID

Int_t TRestAnalysisTree::fSubEventID
private

Definition at line 38 of file TRestAnalysisTree.h.

◆ fSubEventTag

TString* TRestAnalysisTree::fSubEventTag
private

Definition at line 40 of file TRestAnalysisTree.h.

◆ fSubRunOrigin

Int_t TRestAnalysisTree::fSubRunOrigin
private

Definition at line 42 of file TRestAnalysisTree.h.

◆ fTimeStamp

Double_t TRestAnalysisTree::fTimeStamp
private

Definition at line 39 of file TRestAnalysisTree.h.


The documentation for this class was generated from the following files: