ROOT  6.06/09
Reference Guide
DataSet.h
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Peter Speckmayer, Joerg Stelzer, Helge Voss
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : DataSet *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Contains all the data information *
12  * *
13  * Authors (alphabetical): *
14  * Andreas Hoecker <Andreas.Hocker@cern.ch> - CERN, Switzerland *
15  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - CERN, Switzerland *
16  * Peter Speckmayer <Peter.Speckmayer@cern.ch> - CERN, Switzerland *
17  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
18  * *
19  * Copyright (c) 2006: *
20  * CERN, Switzerland *
21  * U. of Victoria, Canada *
22  * MPI-K Heidelberg, Germany *
23  * *
24  * Redistribution and use in source and binary forms, with or without *
25  * modification, are permitted according to the terms listed in LICENSE *
26  * (http://tmva.sourceforge.net/LICENSE) *
27  **********************************************************************************/
28 
29 #ifndef ROOT_TMVA_DataSet
30 #define ROOT_TMVA_DataSet
31 
32 //////////////////////////////////////////////////////////////////////////
33 // //
34 // DataSet //
35 // //
36 // Class that contains all the data information //
37 // //
38 //////////////////////////////////////////////////////////////////////////
39 
40 #include <vector>
41 #include <map>
42 #include <string>
43 
44 #ifndef ROOT_TObject
45 #include "TObject.h"
46 #endif
47 #ifndef ROOT_TString
48 #include "TString.h"
49 #endif
50 #ifndef ROOT_TTree
51 #include "TTree.h"
52 #endif
53 //#ifndef ROOT_TCut
54 //#include "TCut.h"
55 //#endif
56 //#ifndef ROOT_TMatrixDfwd
57 //#include "TMatrixDfwd.h"
58 //#endif
59 //#ifndef ROOT_TPrincipal
60 //#include "TPrincipal.h"
61 //#endif
62 #ifndef ROOT_TRandom3
63 #include "TRandom3.h"
64 #endif
65 
66 #ifndef ROOT_TMVA_Types
67 #include "TMVA/Types.h"
68 #endif
69 #ifndef ROOT_TMVA_VariableInfo
70 #include "TMVA/VariableInfo.h"
71 #endif
72 
73 namespace TMVA {
74 
75  class Event;
76  class DataSetInfo;
77  class MsgLogger;
78  class Results;
79 
80  class DataSet {
81 
82  public:
83 
84  DataSet(const DataSetInfo&);
85  virtual ~DataSet();
86 
87  void AddEvent( Event *, Types::ETreeType );
88 
92 
93  // const getters
94  const Event* GetEvent() const; // returns event without transformations
95  const Event* GetEvent ( Long64_t ievt ) const { fCurrentEventIdx = ievt; return GetEvent(); } // returns event without transformations
96  const Event* GetTrainingEvent( Long64_t ievt ) const { return GetEvent(ievt, Types::kTraining); }
97  const Event* GetTestEvent ( Long64_t ievt ) const { return GetEvent(ievt, Types::kTesting); }
98  const Event* GetEvent ( Long64_t ievt, Types::ETreeType type ) const
99  {
100  fCurrentTreeIdx = TreeIndex(type); fCurrentEventIdx = ievt; return GetEvent();
101  }
102 
103 
104 
105 
106  UInt_t GetNVariables() const;
107  UInt_t GetNTargets() const;
108  UInt_t GetNSpectators() const;
109 
110  void SetCurrentEvent( Long64_t ievt ) const { fCurrentEventIdx = ievt; }
113 
114  void SetEventCollection( std::vector<Event*>*, Types::ETreeType );
115  const std::vector<Event*>& GetEventCollection( Types::ETreeType type = Types::kMaxTreeType ) const;
117 
122 
124 
125  Results* GetResults ( const TString &,
127  Types::EAnalysisType analysistype );
128  void DeleteResults ( const TString &,
130  Types::EAnalysisType analysistype );
131 
132  void SetVerbose( Bool_t ) {}
133 
134  // sets the number of blocks to which the training set is divided,
135  // some of which are given to the Validation sample. As default they belong all to Training set.
136  void DivideTrainingSet( UInt_t blockNum );
137 
138  // sets a certrain block from the origin training set to belong to either Training or Validation set
139  void MoveTrainingBlock( Int_t blockInd,Types::ETreeType dest, Bool_t applyChanges = kTRUE );
140 
141  void IncrementNClassEvents( Int_t type, UInt_t classNumber );
142  Long64_t GetNClassEvents ( Int_t type, UInt_t classNumber );
143  void ClearNClassEvents ( Int_t type );
144 
146 
147  // accessors for random and importance sampling
148  void InitSampling( Float_t fraction, Float_t weight, UInt_t seed = 0 );
149  void EventResult( Bool_t successful, Long64_t evtNumber = -1 );
150  void CreateSampling() const;
151 
153 
154  private:
155 
156  // data members
157  DataSet();
158  void DestroyCollection( Types::ETreeType type, Bool_t deleteEvents );
159 
160  const DataSetInfo& fdsi; //! datasetinfo that created this dataset
161 
162  std::vector<Event*>::iterator fEvtCollIt;
163  std::vector< std::vector<Event*>* > fEventCollection; //! list of events for training/testing/...
164 
165  std::vector< std::map< TString, Results* > > fResults; //! [train/test/...][method-identifier]
166 
169 
170  // event sampling
171  std::vector<Char_t> fSampling; // random or importance sampling (not all events are taken) !! Bool_t are stored ( no std::vector<bool> taken for speed (performance) issues )
172  std::vector<Int_t> fSamplingNEvents; // number of events which should be sampled
173  std::vector<Float_t> fSamplingWeight; // weight change factor [weight is indicating if sampling is random (1.0) or importance (<1.0)]
174  mutable std::vector< std::vector< std::pair< Float_t, Long64_t >* > > fSamplingEventList; // weights and indices for sampling
175  mutable std::vector< std::vector< std::pair< Float_t, Long64_t >* > > fSamplingSelected; // selected events
176  TRandom3 *fSamplingRandom; // random generator for sampling
177 
178 
179  // further things
180  std::vector< std::vector<Long64_t> > fClassEvents; //! number of events of class 0,1,2,... in training[0]
181  // and testing[1] (+validation, trainingoriginal)
182 
183  Bool_t fHasNegativeEventWeights; // true if at least one signal or bkg event has negative weight
184 
185  mutable MsgLogger* fLogger; // message logger
186  MsgLogger& Log() const { return *fLogger; }
187  std::vector<Char_t> fBlockBelongToTraining; // when dividing the dataset to blocks, sets whether
188  // the certain block is in the Training set or else
189  // in the validation set
190  // boolean are stored, taken std::vector<Char_t> for performance reasons (instead of std::vector<Bool_t>)
191  Long64_t fTrainingBlockSize; // block size into which the training dataset is divided
192 
195  };
196 }
197 
198 
199 //_______________________________________________________________________
201 {
202  switch (type) {
203  case Types::kMaxTreeType : return fCurrentTreeIdx;
204  case Types::kTraining : return 0;
205  case Types::kTesting : return 1;
206  case Types::kValidation : return 2;
207  case Types::kTrainingOriginal : return 3;
208  default : return fCurrentTreeIdx;
209  }
210 }
211 
212 //_______________________________________________________________________
214 {
215  switch (fCurrentTreeIdx) {
216  case 0: return Types::kTraining;
217  case 1: return Types::kTesting;
218  case 2: return Types::kValidation;
219  case 3: return Types::kTrainingOriginal;
220  }
221  return Types::kMaxTreeType;
222 }
223 
224 //_______________________________________________________________________
226 {
227  Int_t treeIdx = TreeIndex(type);
228  if (fSampling.size() > UInt_t(treeIdx) && fSampling.at(treeIdx)) {
229  return fSamplingSelected.at(treeIdx).size();
230  }
231  return GetEventCollection(type).size();
232 }
233 
234 //_______________________________________________________________________
235 inline const std::vector<TMVA::Event*>& TMVA::DataSet::GetEventCollection( TMVA::Types::ETreeType type ) const
236 {
237  return *(fEventCollection.at(TreeIndex(type)));
238 }
239 
240 
241 #endif
UInt_t GetNSpectators() const
access the number of targets through the datasetinfo
Definition: DataSet.cxx:210
void SetEventCollection(std::vector< Event * > *, Types::ETreeType)
Sets the event collection (by DataSetFactory)
Definition: DataSet.cxx:229
Long64_t GetNTestEvents() const
Definition: DataSet.h:91
Random number generator class based on M.
Definition: TRandom3.h:29
long long Long64_t
Definition: RtypesCore.h:69
Long64_t fTrainingBlockSize
Definition: DataSet.h:191
void AddEvent(Event *, Types::ETreeType)
add event to event list after which the event is owned by the dataset
Definition: DataSet.cxx:219
float Float_t
Definition: RtypesCore.h:53
MsgLogger & Log() const
Definition: DataSet.h:186
const Event * GetTrainingEvent(Long64_t ievt) const
Definition: DataSet.h:96
TRandom3 * fSamplingRandom
Definition: DataSet.h:176
EAnalysisType
Definition: Types.h:124
UInt_t TreeIndex(Types::ETreeType type) const
Definition: DataSet.h:200
const Event * GetEvent(Long64_t ievt) const
Definition: DataSet.h:95
Basic string class.
Definition: TString.h:137
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
std::vector< Char_t > fBlockBelongToTraining
Definition: DataSet.h:187
void ClearNClassEvents(Int_t type)
Definition: DataSet.cxx:138
std::vector< std::vector< std::pair< Float_t, Long64_t > * > > fSamplingEventList
Definition: DataSet.h:174
Long64_t GetNEvtBkgdTrain()
return number of background training events in dataset
Definition: DataSet.cxx:414
const Event * GetTestEvent(Long64_t ievt) const
Definition: DataSet.h:97
Bool_t HasNegativeEventWeights() const
Definition: DataSet.h:123
virtual ~DataSet()
destructor
Definition: DataSet.cxx:94
TTree * GetTree(Types::ETreeType type)
create the test/trainings tree with all the variables, the weights, the classes, the targets...
Definition: DataSet.cxx:573
void SetCurrentEvent(Long64_t ievt) const
Definition: DataSet.h:110
Bool_t fHasNegativeEventWeights
number of events of class 0,1,2,... in training[0]
Definition: DataSet.h:183
void MoveTrainingBlock(Int_t blockInd, Types::ETreeType dest, Bool_t applyChanges=kTRUE)
move training block
Definition: DataSet.cxx:378
void ApplyTrainingSetDivision()
apply division of data set
Definition: DataSet.cxx:358
const DataSetInfo & fdsi
Definition: DataSet.h:160
Results * GetResults(const TString &, Types::ETreeType type, Types::EAnalysisType analysistype)
TString info(resultsName+"/"); switch(type) { case Types::kTraining: info += "kTraining/"; break; cas...
Definition: DataSet.cxx:257
std::vector< std::vector< Long64_t > > fClassEvents
Definition: DataSet.h:180
Long64_t GetNEvtSigTest()
return number of signal test events in dataset
Definition: DataSet.cxx:390
std::vector< std::vector< std::pair< Float_t, Long64_t > * > > fSamplingSelected
Definition: DataSet.h:175
void DeleteResults(const TString &, Types::ETreeType type, Types::EAnalysisType analysistype)
delete the results stored for this particulary Method instance (here appareantly called resultsName i...
Definition: DataSet.cxx:308
unsigned int UInt_t
Definition: RtypesCore.h:42
void DivideTrainingSet(UInt_t blockNum)
divide training set
Definition: DataSet.cxx:334
void DestroyCollection(Types::ETreeType type, Bool_t deleteEvents)
destroys the event collection (events + vector)
Definition: DataSet.cxx:167
Types::ETreeType GetCurrentType() const
Definition: DataSet.h:213
Long64_t GetNEvtBkgdTest()
return number of background test events in dataset
Definition: DataSet.cxx:398
void CreateSampling() const
create an event sampling (random or importance sampling)
Definition: DataSet.cxx:472
std::vector< Event * >::iterator fEvtCollIt
datasetinfo that created this dataset
Definition: DataSet.h:162
void IncrementNClassEvents(Int_t type, UInt_t classNumber)
Definition: DataSet.cxx:129
const Event * GetEvent() const
Definition: DataSet.cxx:180
void SetCurrentType(Types::ETreeType type) const
Definition: DataSet.h:111
std::vector< Char_t > fSampling
Definition: DataSet.h:171
void EventResult(Bool_t successful, Long64_t evtNumber=-1)
increase the importance sampling weight of the event when not successful and decrease it when success...
Definition: DataSet.cxx:536
std::vector< Float_t > fSamplingWeight
Definition: DataSet.h:173
Long64_t fCurrentEventIdx
Definition: DataSet.h:168
UInt_t fCurrentTreeIdx
[train/test/...][method-identifier]
Definition: DataSet.h:167
Long64_t GetNEvtSigTrain()
return number of signal training events in dataset
Definition: DataSet.cxx:406
int type
Definition: TGX11.cxx:120
void SetVerbose(Bool_t)
Definition: DataSet.h:132
Long64_t GetNEvents(Types::ETreeType type=Types::kMaxTreeType) const
Definition: DataSet.h:225
UInt_t GetNTargets() const
access the number of targets through the datasetinfo
Definition: DataSet.cxx:202
UInt_t GetNVariables() const
access the number of variables through the datasetinfo
Definition: DataSet.cxx:194
void ApplyTrainingBlockDivision()
std::vector< std::vector< Event * > * > fEventCollection
Definition: DataSet.h:163
std::vector< Int_t > fSamplingNEvents
Definition: DataSet.h:172
const Event * GetEvent(Long64_t ievt, Types::ETreeType type) const
Definition: DataSet.h:98
Long64_t GetNClassEvents(Int_t type, UInt_t classNumber)
Definition: DataSet.cxx:146
MsgLogger * fLogger
Definition: DataSet.h:185
Abstract ClassifierFactory template that handles arbitrary types.
#define dest(otri, vertexptr)
Definition: triangle.c:1040
const std::vector< Event * > & GetEventCollection(Types::ETreeType type=Types::kMaxTreeType) const
Definition: DataSet.h:235
Long64_t GetNTrainingEvents() const
Definition: DataSet.h:90
A TTree object has a header with a name and a title.
Definition: TTree.h:94
const TTree * GetEventCollectionAsTree()
const Bool_t kTRUE
Definition: Rtypes.h:91
void InitSampling(Float_t fraction, Float_t weight, UInt_t seed=0)
initialize random or importance sampling
Definition: DataSet.cxx:422
std::vector< std::map< TString, Results * > > fResults
list of events for training/testing/...
Definition: DataSet.h:165