hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
MutableTTree.h
Go to the documentation of this file.
1#ifndef MUTABLE_TTREE_H
2#define MUTABLE_TTREE_H
3
4#include <iostream>
5#include <fstream>
6#include <string>
7#include <cstdlib>
8#include <TTree.h>
9#include <TFile.h>
10#include <TBranch.h>
11#include <functional>
12
17
18 public:
19
20 MutableTTree(TFile* infile, std::string tree_name);
21
25 int GetEntries(){return tree_->GetEntries();}
26
31 void GetEntry(int entry){tree_->GetEntry(entry);}
32
36 void Fill();
37
43 void shiftVariable(std::string variable, double shift);
44
50 double getValue(std::string branch_name);
51
55 void printEvent();
56
62 void setBranchValue(std::string branch_name, double value){*tuple_[branch_name] = value;}
63
68 void addNewBranch(std::string branch);
69
75 void defineMassWindow(double lowMass, double highMass);
76
80 std::vector<std::string> getAllVariables();
81
88 bool variableExists(std::string variable);
89
90 virtual void addVariable(std::string variableName, double param)=0;
91
92 void addVariableToTBranch(const std::string& variableName);
93
95
96 protected:
97 TTree* tree_{nullptr};
98 TTree* newtree_{nullptr};
99 std::map<std::string,double*> tuple_;
100 std::map<std::string,TBranch*> new_branches;
101 std::map<std::string, double*> new_variables_;
102 std::map<std::string,std::function<double()>> functions_;
103 std::map<std::string,std::function<double()>> variable_shifts_;
104
105 private:
111 void initializeFlatTuple(TTree* tree, std::map<std::string, double*> &tuple_map);
112
116 void copyTTree();
117
118 double lowMass_{-999.9};
119 double highMass_{-999.9};
120};
121
122#endif // __MUTABLE_TTREE_H
Reads flat TTree and allows user to create new variables in the TTree.
void initializeFlatTuple(TTree *tree, std::map< std::string, double * > &tuple_map)
read in the initial flat TTree
void GetEntry(int entry)
get tree entry
void defineMassWindow(double lowMass, double highMass)
Set the mass window within which to read the input ttree.
double getValue(std::string branch_name)
Get the value of a flat tuple variable.
void shiftVariable(std::string variable, double shift)
Apply any corrections to specified variable.
int GetEntries()
return number of entries in tree
void setBranchValue(std::string branch_name, double value)
Set branch value.
TTree * newtree_
temporary ttree used to create and fill new branches
std::vector< std::string > getAllVariables()
Get list of all variables defined in ttree.
void addNewBranch(std::string branch)
Add new branch to ttree.
TTree * tree_
flat ttree
void copyTTree()
copy the TTree
void printEvent()
Print TTree Event.
virtual void addVariable(std::string variableName, double param)=0
std::map< std::string, double * > tuple_
holds all variables and values
void addVariableToTBranch(const std::string &variableName)
std::map< std::string, TBranch * > new_branches
list of new branches added to ttree
bool variableExists(std::string variable)
Check if a variable exists in the ttree.
double lowMass_
mass window low
double highMass_
mass window high
std::map< std::string, std::function< double()> > variable_shifts_
variable corrections
std::map< std::string, std::function< double()> > functions_
functions that calculate new variables
void Fill()
Fill ttree with new variables included.
std::map< std::string, double * > new_variables_
list of new variables