laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauFitDataTree.hh
Go to the documentation of this file.
1 
2 /*
3 Copyright 2004 University of Warwick
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 */
17 
18 /*
19 Laura++ package authors:
20 John Back
21 Paul Harrison
22 Thomas Latham
23 */
24 
35 #ifndef LAU_FIT_DATA_TREE
36 #define LAU_FIT_DATA_TREE
37 
38 #include "TEventList.h"
39 #include "TTree.h"
40 
41 #include <map>
42 #include <vector>
43 
44 class TFile;
45 class TLeaf;
46 
48 typedef std::map<TString, Double_t> LauFitData;
49 
51 
52  public:
54 
58  LauFitDataTree( const TString& rootFileName, const TString& rootTreeName );
59 
61  virtual ~LauFitDataTree();
62 
64 
67  const TString& fileName() const { return rootFileName_; }
68 
70 
73  const TString& treeName() const { return rootTreeName_; }
74 
76 
79  Bool_t findBranches();
80 
82  void readAllData();
83 
85 
88  void readExperimentData( UInt_t iExpt );
89 
91 
97  void appendFakePoints( const std::vector<Double_t>& xCoords,
98  const std::vector<Double_t>& yCoords );
99 
101 
104  UInt_t nBranches() const;
105 
107 
110  UInt_t nTreeEvents() const
111  {
112  return rootTree_ ? static_cast<UInt_t>( rootTree_->GetEntries() ) : 0;
113  }
114 
116 
119  UInt_t nEvents() const
120  {
121  return eventList_ ? static_cast<UInt_t>( eventList_->GetN() ) : this->nTreeEvents();
122  }
123 
125 
128  UInt_t nFakeEvents() const { return fakeEvents_.size(); }
129 
131 
134  UInt_t nTotalEvents() const { return this->nTreeEvents() + this->nFakeEvents(); }
135 
137 
141  Bool_t haveBranch( const TString& name ) const;
142 
144 
148  const LauFitData& getData( UInt_t iEvt ) const;
149 
151  void disableAllBranches() const;
152 
154  void enableAllBranches() const;
155 
157 
160  void enableBranch( const TString& name ) const;
161 
163 
166  void disableBranch( const TString& name ) const;
167 
168  protected:
170  void openFileAndTree();
171 
173 
176  void loadData();
177 
178  private:
181 
184 
186  typedef std::map<TString, UInt_t> LauNameIndexMap;
187 
189  typedef std::vector<Double_t> LauEventData;
190 
192  typedef std::vector<TLeaf*> LauLeafList;
193 
195  TString rootFileName_;
196 
198  TString rootTreeName_;
199 
201  TFile* rootFile_;
202 
204  TTree* rootTree_;
205 
207  TEventList* eventList_;
208 
211 
214 
217 
220 
222  std::vector<LauEventData> treeEvents_;
223 
225  std::vector<LauEventData> fakeEvents_;
226 
227  ClassDef( LauFitDataTree, 0 )
228 };
229 
230 #endif
LauFitDataTree(const TString &rootFileName, const TString &rootTreeName)
Constructor.
TString rootTreeName_
The name of the tree ocntaining the data.
void readExperimentData(UInt_t iExpt)
Read events only for the given experiment.
UInt_t nBranches() const
Obtain the number of branches in the tree.
const LauFitData & getData(UInt_t iEvt) const
Retrieve the data for a given event.
LauLeafList leaves_
The leaf objects.
void appendFakePoints(const std::vector< Double_t > &xCoords, const std::vector< Double_t > &yCoords)
Add fake events to the data.
Bool_t findBranches()
Find all of the branches in the tree.
LauEventData eventData_
Stores the current event.
const TString & treeName() const
Retrieve the tree name.
void disableBranch(const TString &name) const
Disable the named branch.
void enableBranch(const TString &name) const
Enable the named branch.
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.
std::map< TString, UInt_t > LauNameIndexMap
The type used to map the leaf names to the vector indices.
TString rootFileName_
The name of the file containing the data.
Class to store the input fit variables.
LauFitDataTree(const LauFitDataTree &rhs)
Copy constructor (not implemented)
UInt_t nFakeEvents() const
Retrieve the number of fake events.
void readAllData()
Read all events from the tree.
UInt_t nTreeEvents() const
Retrieve the number of events in the tree.
LauFitData eventDataOut_
Stores the current event (for external use)
void loadData()
Load events from the tree.
TEventList * eventList_
A list of the events in the current experiment.
void enableAllBranches() const
Enable all branches.
std::vector< TLeaf * > LauLeafList
The type used to hold the leaves.
const TString & fileName() const
Retrieve the file name.
void openFileAndTree()
Open the file and tree.
UInt_t nTotalEvents() const
Retrieve the total number of events.
void disableAllBranches() const
Disable all branches.
std::vector< LauEventData > fakeEvents_
The fake events, which are not from the tree.
std::map< TString, Double_t > LauFitData
Type for holding event data.
LauNameIndexMap leafNames_
Stores the mapping from the leaf names to the vector indices.
std::vector< LauEventData > treeEvents_
The events read from the tree.
TTree * rootTree_
The tree containing the data.
TFile * rootFile_
The file containing the data.
virtual ~LauFitDataTree()
Destructor.
LauFitDataTree & operator=(const LauFitDataTree &rhs)
Copy assignment operator (not implemented)
UInt_t nEvents() const
Retrieve the number of events.