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.cxx
Go to the documentation of this file.
1#include <MutableTTree.h>
2
3MutableTTree::MutableTTree(TFile* infile, std::string tree_name){
4 std::cout << "[MutableTTree]::Reading in tree: " << tree_name << std::endl;
5 tree_ = (TTree*)infile->Get(tree_name.c_str());
6 if(tree_ == nullptr)
7 std::cout << "[MutableTTree]::ERROR READING TREE " << tree_name << " from file " << std::endl;
9 newtree_ = new TTree();
10 copyTTree();
11}
12
13double MutableTTree::getValue(std::string branch_name){
14 if(tuple_.find(branch_name) == tuple_.end()){
15 return -9999.9;
16 }
17 else
18 return *tuple_[branch_name];
19}
20
22 //Make new TTree. Copy branches from original. Add new branch. Fill new TTree with old TTree
23 int nBr = tree_->GetListOfBranches()->GetEntries();
24 for(int iBr = 0; iBr < nBr; iBr++){
25 TBranch *br = dynamic_cast<TBranch*>(tree_->GetListOfBranches()->At(iBr));
26 std::string varname = (std::string)br->GetFullName();
27 newtree_->Branch(varname.c_str(),tuple_[varname],(varname+"/D").c_str());
28 }
29}
30
31void MutableTTree::defineMassWindow(double lowMass, double highMass){
32 lowMass_ = lowMass;
33 highMass_ = highMass;
34}
35
36bool MutableTTree::variableExists(std::string variable){
37 if(tuple_.find(variable) != tuple_.end())
38 return true;
39 else
40 return false;
41}
42
43std::vector<std::string> MutableTTree::getAllVariables(){
44 std::vector<std::string> variables;
45 for(std::map<std::string,double*>::iterator it = tuple_.begin(); it != tuple_.end(); it++){
46 variables.push_back(it->first);
47 }
48
49 return variables;
50}
51
53
54 for(int e=0; e < tree_->GetEntries(); e++){
55 tree_->GetEntry(e);
56
57 //Mass Window (if set)
58 if(lowMass_ != -999.9 && highMass_ != -999.9){
59 if(getValue("unc_vtx_mass")*1000.0 > highMass_) continue;
60 if(getValue("unc_vtx_mass")*1000.0 < lowMass_) continue;
61 }
62
63 //Apply varible shifts here
64 for(std::map<std::string,std::function<double()>>::iterator it = variable_shifts_.begin();
65 it != variable_shifts_.end(); it++){
66 variable_shifts_[it->first]();
67 }
68
69 //Add new variables here
70 for(std::map<std::string,double*>::iterator it = new_variables_.begin(); it != new_variables_.end(); it++){
71 *new_variables_[it->first] = functions_[it->first]();
72 }
73
74 newtree_->Fill();
75 }
76
77 delete tree_;
79}
80
81void MutableTTree::shiftVariable(std::string variable, double shift){
82 std::cout << "[MutableTTree]::Shifting Variable " << variable << " by " << shift << std::endl;
83 std::function<double()> shiftVariableFunc = [&,shift]()->double{
84 return *tuple_[variable] = *tuple_[variable] + shift;
85 };
86 variable_shifts_[variable] = shiftVariableFunc;
87}
88
89void MutableTTree::addNewBranch(std::string branch){
90 double* value = new double{999.9};
91 tuple_[branch] = value;
92 tree_->Branch(branch.c_str(),tuple_[branch],(branch+"/D").c_str());
93}
94
96 for(std::map<std::string,double*>::iterator it = tuple_.begin(); it != tuple_.end(); it ++){
97 std::cout << it->first << ": " << *it->second << std::endl;
98 }
99 std::cout << "[MutableTree}::End print" << std::endl;
100}
101
102void MutableTTree::initializeFlatTuple(TTree* tree, std::map<std::string, double*> &tuple_map){
103 int nBr = tree->GetListOfBranches()->GetEntries();
104 for(int iBr = 0; iBr < nBr; iBr++){
105 TBranch *br = dynamic_cast<TBranch*>(tree->GetListOfBranches()->At(iBr));
106 double* value = new double;
107 std::string varname = (std::string)br->GetFullName();
108 tuple_map[varname] = value;
109 tree->SetBranchAddress(varname.c_str(),tuple_map[varname]);
110 }
111}
112
113void MutableTTree::addVariableToTBranch(const std::string& variableName){
114 double* variable = new double{999.9};
115 tuple_[variableName] = variable;
116 newtree_->Branch(variableName.c_str(), tuple_[variableName], (variableName+"/D").c_str());
117 new_variables_[variableName] = variable;
118}
119
void initializeFlatTuple(TTree *tree, std::map< std::string, double * > &tuple_map)
read in the initial flat TTree
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.
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
MutableTTree(TFile *infile, std::string tree_name)
void copyTTree()
copy the TTree
void printEvent()
Print TTree Event.
std::map< std::string, double * > tuple_
holds all variables and values
void addVariableToTBranch(const std::string &variableName)
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