laura is hosted by Hepforge, IPPP Durham
Laura++  v3r0p1
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauFitDataTree.hh
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2013.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
21 #ifndef LAU_FIT_DATA_TREE
22 #define LAU_FIT_DATA_TREE
23 
24 #include <vector>
25 #include <map>
26 
27 #include "TEventList.h"
28 #include "TTree.h"
29 
30 class TFile;
31 class TLeaf;
32 
34 typedef std::map<TString,Double_t> LauFitData;
35 
36 
38 
39  public:
41 
45  LauFitDataTree(const TString& rootFileName, const TString& rootTreeName);
46 
48  virtual ~LauFitDataTree();
49 
51 
54  const TString& fileName() const {return rootFileName_;}
55 
57 
60  const TString& treeName() const {return rootTreeName_;}
61 
63 
66  Bool_t findBranches();
67 
69  void readAllData();
70 
72 
75  void readExperimentData( UInt_t iExpt );
76 
78 
84  void appendFakePoints(const std::vector<Double_t>& xCoords, const std::vector<Double_t>& yCoords);
85 
87 
90  UInt_t nBranches() const;
91 
93 
96  UInt_t nTreeEvents() const {return rootTree_ ? static_cast<UInt_t>(rootTree_->GetEntries()) : 0;}
97 
99 
102  UInt_t nEvents() const {return eventList_ ? static_cast<UInt_t>(eventList_->GetN()) : this->nTreeEvents();}
103 
105 
108  UInt_t nFakeEvents() const {return fakeEvents_.size();}
109 
111 
114  UInt_t nTotalEvents() const {return this->nTreeEvents()+this->nFakeEvents();}
115 
117 
121  Bool_t haveBranch(const TString& name) const;
122 
124 
128  const LauFitData& getData(UInt_t iEvt) const;
129 
131  void disableAllBranches() const;
132 
134  void enableAllBranches() const;
135 
137 
140  void enableBranch(const TString& name) const;
141 
143 
146  void disableBranch(const TString& name) const;
147 
148  protected:
150  void openFileAndTree();
151 
153 
156  void loadData();
157 
158  private:
160  LauFitDataTree(const LauFitDataTree& rhs);
161 
164 
166  typedef std::map<TString,UInt_t> LauNameIndexMap;
167 
169  typedef std::vector<Double_t> LauEventData;
170 
172  typedef std::vector<TLeaf*> LauLeafList;
173 
175  TString rootFileName_;
176 
178  TString rootTreeName_;
179 
181  TFile* rootFile_;
182 
184  TTree* rootTree_;
185 
187  TEventList* eventList_;
188 
191 
194 
197 
200 
202  std::vector<LauEventData> treeEvents_;
203 
205  std::vector<LauEventData> fakeEvents_;
206 
207  ClassDef(LauFitDataTree, 0)
208 };
209 
210 #endif
Bool_t findBranches()
Find all of the branches in the tree.
void loadData()
Load events from the tree.
LauEventData eventData_
Stores the current event.
std::map< TString, UInt_t > LauNameIndexMap
The type used to map the leaf names to the vector indices.
void openFileAndTree()
Open the file and tree.
std::vector< LauEventData > fakeEvents_
The fake events, which are not from the tree.
void readAllData()
Read all events from the tree.
TFile * rootFile_
The file containing the data.
TString rootFileName_
The name of the file containing the data.
LauFitData eventDataOut_
Stores the current event (for external use)
void enableAllBranches() const
Enable all branches.
std::map< TString, Double_t > LauFitData
Type for holding event data.
std::vector< LauEventData > treeEvents_
The events read from the tree.
UInt_t nBranches() const
Obtain the number of branches in the tree.
TString rootTreeName_
The name of the tree ocntaining the data.
void readExperimentData(UInt_t iExpt)
Read events only for the given experiment.
LauNameIndexMap leafNames_
Stores the mapping from the leaf names to the vector indices.
UInt_t nTotalEvents() const
Retrieve the total number of events.
TEventList * eventList_
A list of the events in the current experiment.
std::vector< TLeaf * > LauLeafList
The type used to hold the leaves.
const TString & treeName() const
Retrieve the tree name.
void disableAllBranches() const
Disable all branches.
TTree * rootTree_
The tree containing the data.
LauFitDataTree & operator=(const LauFitDataTree &rhs)
Copy assignment operator (not implemented)
virtual ~LauFitDataTree()
Destructor.
const TString & fileName() const
Retrieve the file name.
LauFitDataTree(const TString &rootFileName, const TString &rootTreeName)
Constructor.
void appendFakePoints(const std::vector< Double_t > &xCoords, const std::vector< Double_t > &yCoords)
Add fake events to the data.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
UInt_t nTreeEvents() const
Retrieve the number of events in the tree.
UInt_t nFakeEvents() const
Retrieve the number of fake events.
std::vector< Double_t > LauEventData
The type used to contain the data for each event.
Bool_t haveBranch(const TString &name) const
Check if the named branch is stored.
UInt_t nEvents() const
Retrieve the number of events.
void enableBranch(const TString &name) const
Enable the named branch.
void disableBranch(const TString &name) const
Disable the named branch.
LauLeafList leaves_
The leaf objects.
Class to store the input fit variables.