FairRoot/PandaRoot
Public Member Functions | Private Attributes | List of all members
PndSimpleNtuple Class Reference

#include <PndSimpleNtuple.h>

Public Member Functions

 PndSimpleNtuple (TString name, TString title, TString precut="")
 
 ~PndSimpleNtuple ()
 
void Column (const TString &name, double value)
 
void Column (const TString &name, float value)
 
void Column (const TString &name, int value)
 
void Column (const TString &name, bool value)
 
void ColumnD (const TString &name, double value)
 
void ColumnF (const TString &name, float value)
 
void ColumnI (const TString &name, int value)
 
void ColumnB (const TString &name, bool value)
 
bool DumpData ()
 
bool Accept ()
 
bool AcceptData ()
 
TTree * GetTree ()
 
double GetCurrentValue (TString name)
 
void SetPrecut (TString precut)
 
bool BranchExists (TString name)
 
int ShowBranches ()
 

Private Attributes

TTree * fTree
 
TTree * fTmpTree
 
std::map< TString, double > fDValues
 
std::map< TString, float > fFValues
 
std::map< TString, int > fIValues
 
std::map< TString, char > fBValues
 
std::map< TString, int > fBrTypes
 
TString fPrecut
 
TTreeFormula * fFml
 

Detailed Description

Definition at line 40 of file PndSimpleNtuple.h.

Constructor & Destructor Documentation

PndSimpleNtuple::PndSimpleNtuple ( TString  name,
TString  title,
TString  precut = "" 
)

Definition at line 44 of file PndSimpleNtuple.cxx.

References fTmpTree, and fTree.

44  :
45  fTree(new TTree(name, title)),fTmpTree(new TTree("tmp", "")), fPrecut(precut), fFml(0)
46 {
47  fTree->SetDirectory(0);
48  fTmpTree->SetDirectory(0);
49 };
TTreeFormula * fFml
TString name
PndSimpleNtuple::~PndSimpleNtuple ( )
inline

Definition at line 44 of file PndSimpleNtuple.h.

References fFml.

44 {if (fFml!=0) delete fFml;}
TTreeFormula * fFml

Member Function Documentation

bool PndSimpleNtuple::Accept ( )

Definition at line 161 of file PndSimpleNtuple.cxx.

References fFml, fPrecut, and fTmpTree.

Referenced by AcceptData().

162 {
163  if (fPrecut=="") return true;
164 
165  if (fFml==0) fFml = new TTreeFormula("fFml", fPrecut, fTmpTree);
166 
167  fTmpTree->Fill();
168  fTmpTree->GetEntry(fTmpTree->GetEntriesFast()-1);
169  return fFml->EvalInstance();
170 }
TTreeFormula * fFml
bool PndSimpleNtuple::AcceptData ( )
inline

Definition at line 63 of file PndSimpleNtuple.h.

References Accept(), and DumpData().

63 {if (Accept()) {DumpData(); return true;} return false; } // writes current event if accepted
bool PndSimpleNtuple::BranchExists ( TString  name)
inline

Definition at line 68 of file PndSimpleNtuple.h.

References fBrTypes.

68 { return (fBrTypes.find(name)!=fBrTypes.end()); }
TString name
std::map< TString, int > fBrTypes
void PndSimpleNtuple::Column ( const TString name,
double  value 
)

Definition at line 54 of file PndSimpleNtuple.cxx.

References fBrTypes, fDValues, fTmpTree, fTree, name, and typearr.

Referenced by ColumnB(), ColumnD(), ColumnF(), and ColumnI().

55 {
56  // new branch?
57  if (fBrTypes.find(name) == fBrTypes.end())
58  {
59  fBrTypes[name] = 0; // type double
60  fDValues[name] = value;
61  fTree->Branch(name.Data(), &(fDValues[name]), (name+"/D").Data());
62  fTmpTree->Branch(name.Data(), &(fDValues[name]), (name+"/D").Data());
63  }
64  else if ( fBrTypes[name] != 0)
65  {
66  std::cout<<" - WARNING - Type mismatch of branch '"<<name<<"/"<<typearr[fBrTypes[name]]<<"' with type '"<<typearr[0]<<"'"<<std::endl;
67  return;
68  }
69  else
70  fDValues[name] = value;
71 }
std::map< TString, double > fDValues
TString typearr[8]
TString name
std::map< TString, int > fBrTypes
void PndSimpleNtuple::Column ( const TString name,
float  value 
)

Definition at line 76 of file PndSimpleNtuple.cxx.

References fBrTypes, fFValues, fTmpTree, fTree, name, and typearr.

77 {
78  // new branch?
79  if (fBrTypes.find(name) == fBrTypes.end())
80  {
81  fBrTypes[name] = 1; // type double
82  fFValues[name] = value;
83  fTree->Branch(name.Data(), &(fFValues[name]), (name+"/F").Data());
84  fTmpTree->Branch(name.Data(), &(fFValues[name]), (name+"/F").Data());
85  }
86  else if ( fBrTypes[name] != 1)
87  {
88  std::cout<<" - WARNING - Type mismatch of branch '"<<name<<"/"<<typearr[fBrTypes[name]]<<"' with type '"<<typearr[1]<<"'"<<std::endl;
89  return;
90  }
91  else
92  fFValues[name] = value;
93 }
TString typearr[8]
std::map< TString, float > fFValues
TString name
std::map< TString, int > fBrTypes
void PndSimpleNtuple::Column ( const TString name,
int  value 
)

Definition at line 98 of file PndSimpleNtuple.cxx.

References fBrTypes, fIValues, fTmpTree, fTree, name, and typearr.

99 {
100  // new branch?
101  if (fBrTypes.find(name) == fBrTypes.end())
102  {
103  fBrTypes[name] = 2; // type double
104  fIValues[name] = value;
105  fTree->Branch(name.Data(), &(fIValues[name]), (name+"/I").Data());
106  fTmpTree->Branch(name.Data(), &(fIValues[name]), (name+"/I").Data());
107  }
108  else if ( fBrTypes[name] != 2)
109  {
110  std::cout<<" - WARNING - Type mismatch of branch '"<<name<<"/"<<typearr[fBrTypes[name]]<<"' with type '"<<typearr[2]<<"'"<<std::endl;
111  return;
112  }
113  else
114  fIValues[name] = value;
115 }
TString typearr[8]
std::map< TString, int > fIValues
TString name
std::map< TString, int > fBrTypes
void PndSimpleNtuple::Column ( const TString name,
bool  value 
)

Definition at line 120 of file PndSimpleNtuple.cxx.

References fBrTypes, fBValues, fTmpTree, fTree, name, and typearr.

121 {
122  // new branch?
123  if (fBrTypes.find(name) == fBrTypes.end())
124  {
125  fBrTypes[name] = 3; // type double
126  fBValues[name] = value;
127  fTree->Branch(name.Data(), &(fBValues[name]), (name+"/B").Data());
128  fTmpTree->Branch(name.Data(), &(fBValues[name]), (name+"/B").Data());
129  }
130  else if ( fBrTypes[name] != 3)
131  {
132  std::cout<<" - WARNING - Type mismatch of branch '"<<name<<"/"<<typearr[fBrTypes[name]]<<"' with type '"<<typearr[3]<<"'"<<std::endl;
133  return;
134  }
135  else
136  fBValues[name] = value;
137 }
TString typearr[8]
TString name
std::map< TString, char > fBValues
std::map< TString, int > fBrTypes
void PndSimpleNtuple::ColumnB ( const TString name,
bool  value 
)
inline

Definition at line 54 of file PndSimpleNtuple.h.

References Column().

54 {Column(name, value);} // sets bool value of a variable 'name'
TString name
void Column(const TString &name, double value)
void PndSimpleNtuple::ColumnD ( const TString name,
double  value 
)
inline

Definition at line 51 of file PndSimpleNtuple.h.

References Column().

51 {Column(name, value);} // sets double value of a variable 'name'
TString name
void Column(const TString &name, double value)
void PndSimpleNtuple::ColumnF ( const TString name,
float  value 
)
inline

Definition at line 52 of file PndSimpleNtuple.h.

References Column().

52 {Column(name, value);} // sets float value of a variable 'name'
TString name
void Column(const TString &name, double value)
void PndSimpleNtuple::ColumnI ( const TString name,
int  value 
)
inline

Definition at line 53 of file PndSimpleNtuple.h.

References Column().

53 {Column(name, value);} // sets integer value of a variable 'name'
TString name
void Column(const TString &name, double value)
bool PndSimpleNtuple::DumpData ( )
inline

Definition at line 61 of file PndSimpleNtuple.h.

References fTmpTree, and fTree.

Referenced by AcceptData().

61 { fTree->Fill(); fTmpTree->Reset(); return true; } // writes the current event
double PndSimpleNtuple::GetCurrentValue ( TString  name)

Definition at line 141 of file PndSimpleNtuple.cxx.

References fBrTypes, fBValues, fDValues, fFValues, fIValues, name, and sqrt().

142 {
143  // branch does not exist
144  if (fBrTypes.find(name)==fBrTypes.end()) return sqrt(-1.);
145 
146  int brtype = fBrTypes[name];
147 
148  switch (brtype)
149  {
150  case 0: return fDValues[name];
151  case 1: return (double) fFValues[name];
152  case 2: return (double) fIValues[name];
153  case 3: return (double) fBValues[name];
154  }
155 
156  return sqrt(-1);
157 }
std::map< TString, double > fDValues
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
std::map< TString, float > fFValues
std::map< TString, int > fIValues
TString name
std::map< TString, char > fBValues
std::map< TString, int > fBrTypes
TTree* PndSimpleNtuple::GetTree ( )
inline

Definition at line 65 of file PndSimpleNtuple.h.

References fTree.

65 {return fTree;}
void PndSimpleNtuple::SetPrecut ( TString  precut)

Definition at line 187 of file PndSimpleNtuple.cxx.

References fFml, and fPrecut.

188 {
189  if (precut!=fPrecut)
190  {
191  if (fFml!=0) {delete fFml; fFml=0;}
192  fPrecut=precut;
193  }
194 }
TTreeFormula * fFml
int PndSimpleNtuple::ShowBranches ( )

Definition at line 174 of file PndSimpleNtuple.cxx.

References cnt, fBrTypes, typearr, and x.

175 {
176  int cnt=0;
177  for (auto x: fBrTypes)
178  {
179  cout <<cnt++<<" : "<<x.first<<"/"<<typearr[x.second]<<endl;
180  }
181 
182  return fBrTypes.size();
183 }
TString typearr[8]
Double_t x
Int_t cnt
Definition: hist-t7.C:106
std::map< TString, int > fBrTypes

Member Data Documentation

std::map<TString, int> PndSimpleNtuple::fBrTypes
private

Definition at line 80 of file PndSimpleNtuple.h.

Referenced by BranchExists(), Column(), GetCurrentValue(), and ShowBranches().

std::map<TString, char> PndSimpleNtuple::fBValues
private

Definition at line 78 of file PndSimpleNtuple.h.

Referenced by Column(), and GetCurrentValue().

std::map<TString, double> PndSimpleNtuple::fDValues
private

Definition at line 75 of file PndSimpleNtuple.h.

Referenced by Column(), and GetCurrentValue().

TTreeFormula* PndSimpleNtuple::fFml
private

Definition at line 83 of file PndSimpleNtuple.h.

Referenced by Accept(), SetPrecut(), and ~PndSimpleNtuple().

std::map<TString, float> PndSimpleNtuple::fFValues
private

Definition at line 76 of file PndSimpleNtuple.h.

Referenced by Column(), and GetCurrentValue().

std::map<TString, int> PndSimpleNtuple::fIValues
private

Definition at line 77 of file PndSimpleNtuple.h.

Referenced by Column(), and GetCurrentValue().

TString PndSimpleNtuple::fPrecut
private

Definition at line 82 of file PndSimpleNtuple.h.

Referenced by Accept(), and SetPrecut().

TTree* PndSimpleNtuple::fTmpTree
private

Definition at line 73 of file PndSimpleNtuple.h.

Referenced by Accept(), Column(), DumpData(), and PndSimpleNtuple().

TTree* PndSimpleNtuple::fTree
private

Definition at line 72 of file PndSimpleNtuple.h.

Referenced by Column(), DumpData(), GetTree(), and PndSimpleNtuple().


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