Logo ROOT   6.07/09
Reference Guide
DataSetInfo.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Joerg Stelzer, Peter Speckmeier
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : DataSetInfo *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementation (see header for description) *
12  * *
13  * Authors (alphabetical): *
14  * Peter Speckmayer <speckmay@mail.cern.ch> - CERN, Switzerland *
15  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - DESY, Germany *
16  * *
17  * Copyright (c) 2008: *
18  * CERN, Switzerland *
19  * MPI-K Heidelberg, Germany *
20  * DESY Hamburg, Germany *
21  * *
22  * Redistribution and use in source and binary forms, with or without *
23  * modification, are permitted according to the terms listed in LICENSE *
24  * (http://tmva.sourceforge.net/LICENSE) *
25  **********************************************************************************/
26 
27 #include <vector>
28 
29 #include "TEventList.h"
30 #include "TFile.h"
31 #include "TH1.h"
32 #include "TH2.h"
33 #include "TProfile.h"
34 #include "TRandom3.h"
35 #include "TMatrixF.h"
36 #include "TVectorF.h"
37 #include "TMath.h"
38 #include "TROOT.h"
39 #include "TObjString.h"
40 
41 #ifndef ROOT_TMVA_MsgLogger
42 #include "TMVA/MsgLogger.h"
43 #endif
44 #ifndef ROOT_TMVA_Tools
45 #include "TMVA/Tools.h"
46 #endif
47 #ifndef ROOT_TMVA_DataSet
48 #include "TMVA/DataSet.h"
49 #endif
50 #ifndef ROOT_TMVA_DataSetInfo
51 #include "TMVA/DataSetInfo.h"
52 #endif
53 #ifndef ROOT_TMVA_DataSetManager
54 #include "TMVA/DataSetManager.h"
55 #endif
56 #ifndef ROOT_TMVA_Event
57 #include "TMVA/Event.h"
58 #endif
59 
60 #include "TMVA/Types.h"
61 #include "TMVA/VariableInfo.h"
62 
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 /// constructor
66 
68  : TObject(),
69  fDataSetManager(NULL),
70  fName(name),
71  fDataSet( 0 ),
72  fNeedsRebuilding( kTRUE ),
73  fVariables(),
74  fTargets(),
75  fSpectators(),
76  fClasses( 0 ),
77  fNormalization( "NONE" ),
78  fSplitOptions(""),
79  fTrainingSumSignalWeights(-1),
80  fTrainingSumBackgrWeights(-1),
81  fTestingSumSignalWeights (-1),
82  fTestingSumBackgrWeights (-1),
83  fOwnRootDir(0),
84  fVerbose( kFALSE ),
85  fSignalClass(0),
86  fTargetsForMulticlass(0),
87  fLogger( new MsgLogger("DataSetInfo", kINFO) )
88 {
89 }
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// destructor
93 
95 {
96  ClearDataSet();
97 
98  for(UInt_t i=0, iEnd = fClasses.size(); i<iEnd; ++i) {
99  delete fClasses[i];
100  }
101 
102  delete fTargetsForMulticlass;
103 
104  delete fLogger;
105 }
106 
107 ////////////////////////////////////////////////////////////////////////////////
108 
110 {
111  if(fDataSet!=0) { delete fDataSet; fDataSet=0; }
112 }
113 
114 void
116 {
117  fLogger->SetMinType(t);
118 }
119 
120 ////////////////////////////////////////////////////////////////////////////////
121 
123 {
124  ClassInfo* theClass = GetClassInfo(className);
125  if (theClass) return theClass;
126 
127 
128  fClasses.push_back( new ClassInfo(className) );
129  fClasses.back()->SetNumber(fClasses.size()-1);
130 
131  //Log() << kHEADER << Endl;
132 
133  Log() << kHEADER << Form("[%s] : ",fName.Data()) << "Added class \"" << className << "\""<< Endl;
134 
135  Log() << kDEBUG <<"\t with internal class number " << fClasses.back()->GetNumber() << Endl;
136 
137 
138  if (className == "Signal") fSignalClass = fClasses.size()-1; // store the signal class index ( for comparison reasons )
139 
140  return fClasses.back();
141 }
142 
143 ////////////////////////////////////////////////////////////////////////////////
144 
146 {
147  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); it++) {
148  if ((*it)->GetName() == name) return (*it);
149  }
150  return 0;
151 }
152 
153 ////////////////////////////////////////////////////////////////////////////////
154 
156 {
157  try {
158  return fClasses.at(cls);
159  }
160  catch(...) {
161  return 0;
162  }
163 }
164 
165 ////////////////////////////////////////////////////////////////////////////////
166 
168 {
169  for (UInt_t cls = 0; cls < GetNClasses() ; cls++) {
170  Log() << kINFO << Form("Dataset[%s] : ",fName.Data()) << "Class index : " << cls << " name : " << GetClassInfo(cls)->GetName() << Endl;
171  }
172 }
173 
174 ////////////////////////////////////////////////////////////////////////////////
175 
177 {
178  return (ev->GetClass() == fSignalClass);
179 }
180 
181 ////////////////////////////////////////////////////////////////////////////////
182 
183 std::vector<Float_t>* TMVA::DataSetInfo::GetTargetsForMulticlass( const TMVA::Event* ev )
184 {
185  if( !fTargetsForMulticlass ) fTargetsForMulticlass = new std::vector<Float_t>( GetNClasses() );
186  // fTargetsForMulticlass->resize( GetNClasses() );
187  fTargetsForMulticlass->assign( GetNClasses(), 0.0 );
188  fTargetsForMulticlass->at( ev->GetClass() ) = 1.0;
189  return fTargetsForMulticlass;
190 }
191 
192 
193 ////////////////////////////////////////////////////////////////////////////////
194 
196 {
197  Bool_t hasCuts = kFALSE;
198  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); it++) {
199  if( TString((*it)->GetCut()) != TString("") ) hasCuts = kTRUE;
200  }
201  return hasCuts;
202 }
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 
206 const TMatrixD* TMVA::DataSetInfo::CorrelationMatrix( const TString& className ) const
207 {
208  ClassInfo* ptr = GetClassInfo(className);
209  return ptr?ptr->GetCorrelationMatrix():0;
210 }
211 
212 ////////////////////////////////////////////////////////////////////////////////
213 /// add a variable (can be a complex expression) to the set of
214 /// variables used in the MV analysis
215 
217  const TString& title,
218  const TString& unit,
219  Double_t min, Double_t max,
220  char varType,
221  Bool_t normalized,
222  void* external )
223 {
224  TString regexpr = expression; // remove possible blanks
225  regexpr.ReplaceAll(" ", "" );
226  fVariables.push_back(VariableInfo( regexpr, title, unit,
227  fVariables.size()+1, varType, external, min, max, normalized ));
229  return fVariables.back();
230 }
231 
232 ////////////////////////////////////////////////////////////////////////////////
233 /// add variable with given VariableInfo
234 
236  fVariables.push_back(VariableInfo( varInfo ));
238  return fVariables.back();
239 }
240 
241 ////////////////////////////////////////////////////////////////////////////////
242 /// add a variable (can be a complex expression) to the set of
243 /// variables used in the MV analysis
244 
246  const TString& title,
247  const TString& unit,
248  Double_t min, Double_t max,
249  Bool_t normalized,
250  void* external )
251 {
252  TString regexpr = expression; // remove possible blanks
253  regexpr.ReplaceAll(" ", "" );
254  char type='F';
255  fTargets.push_back(VariableInfo( regexpr, title, unit,
256  fTargets.size()+1, type, external, min,
257  max, normalized ));
259  return fTargets.back();
260 }
261 
262 ////////////////////////////////////////////////////////////////////////////////
263 /// add target with given VariableInfo
264 
266  fTargets.push_back(VariableInfo( varInfo ));
268  return fTargets.back();
269 }
270 
271 ////////////////////////////////////////////////////////////////////////////////
272 /// add a spectator (can be a complex expression) to the set of spectator variables used in
273 /// the MV analysis
274 
276  const TString& title,
277  const TString& unit,
278  Double_t min, Double_t max, char type,
279  Bool_t normalized, void* external )
280 {
281  TString regexpr = expression; // remove possible blanks
282  regexpr.ReplaceAll(" ", "" );
283  fSpectators.push_back(VariableInfo( regexpr, title, unit,
284  fSpectators.size()+1, type, external, min, max, normalized ));
286  return fSpectators.back();
287 }
288 
289 ////////////////////////////////////////////////////////////////////////////////
290 /// add spectator with given VariableInfo
291 
293  fSpectators.push_back(VariableInfo( varInfo ));
295  return fSpectators.back();
296 }
297 
298 ////////////////////////////////////////////////////////////////////////////////
299 /// find variable by name
300 
302 {
303  for (UInt_t ivar=0; ivar<GetNVariables(); ivar++)
304  if (var == GetVariableInfo(ivar).GetInternalName()) return ivar;
305 
306  for (UInt_t ivar=0; ivar<GetNVariables(); ivar++)
307  Log() << kINFO << Form("Dataset[%s] : ",fName.Data()) << GetVariableInfo(ivar).GetInternalName() << Endl;
308 
309  Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << "<FindVarIndex> Variable \'" << var << "\' not found." << Endl;
310 
311  return -1;
312 }
313 
314 ////////////////////////////////////////////////////////////////////////////////
315 /// set the weight expressions for the classes
316 /// if class name is specified, set only for this class
317 /// if class name is unknown, register new class with this name
318 
319 void TMVA::DataSetInfo::SetWeightExpression( const TString& expr, const TString& className )
320 {
321  if (className != "") {
322  TMVA::ClassInfo* ci = AddClass(className);
323  ci->SetWeight( expr );
324  }
325  else {
326  // no class name specified, set weight for all classes
327  if (fClasses.empty()) {
328  Log() << kWARNING << Form("Dataset[%s] : ",fName.Data()) << "No classes registered yet, cannot specify weight expression!" << Endl;
329  }
330  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); it++) {
331  (*it)->SetWeight( expr );
332  }
333  }
334 }
335 
336 ////////////////////////////////////////////////////////////////////////////////
337 
338 void TMVA::DataSetInfo::SetCorrelationMatrix( const TString& className, TMatrixD* matrix )
339 {
340  GetClassInfo(className)->SetCorrelationMatrix(matrix);
341 }
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 /// set the cut for the classes
345 
346 void TMVA::DataSetInfo::SetCut( const TCut& cut, const TString& className )
347 {
348  if (className == "") { // if no className has been given set the cut for all the classes
349  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); it++) {
350  (*it)->SetCut( cut );
351  }
352  }
353  else {
354  TMVA::ClassInfo* ci = AddClass(className);
355  ci->SetCut( cut );
356  }
357 }
358 
359 ////////////////////////////////////////////////////////////////////////////////
360 /// set the cut for the classes
361 
362 void TMVA::DataSetInfo::AddCut( const TCut& cut, const TString& className )
363 {
364  if (className == "") { // if no className has been given set the cut for all the classes
365  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); it++) {
366  const TCut& oldCut = (*it)->GetCut();
367  (*it)->SetCut( oldCut+cut );
368  }
369  }
370  else {
371  TMVA::ClassInfo* ci = AddClass(className);
372  ci->SetCut( ci->GetCut()+cut );
373  }
374 }
375 
376 ////////////////////////////////////////////////////////////////////////////////
377 /// returns list of variables
378 
379 std::vector<TString> TMVA::DataSetInfo::GetListOfVariables() const
380 {
381  std::vector<TString> vNames;
382  std::vector<TMVA::VariableInfo>::const_iterator viIt = GetVariableInfos().begin();
383  for(;viIt != GetVariableInfos().end(); viIt++) vNames.push_back( (*viIt).GetExpression() );
384 
385  return vNames;
386 }
387 
388 ////////////////////////////////////////////////////////////////////////////////
389 /// calculates the correlation matrices for signal and background,
390 /// prints them to standard output, and fills 2D histograms
391 
393 {
394 
395  Log() << kHEADER //<< Form("Dataset[%s] : ",fName.Data())
396  << "Correlation matrix (" << className << "):" << Endl;
398 }
399 
400 ////////////////////////////////////////////////////////////////////////////////
401 
403  const TString& hName,
404  const TString& hTitle ) const
405 {
406  if (m==0) return 0;
407 
408  const UInt_t nvar = GetNVariables();
409 
410  // workaround till the TMatrix templates are comonly used
411  // this keeps backward compatibility
412  TMatrixF* tm = new TMatrixF( nvar, nvar );
413  for (UInt_t ivar=0; ivar<nvar; ivar++) {
414  for (UInt_t jvar=0; jvar<nvar; jvar++) {
415  (*tm)(ivar, jvar) = (*m)(ivar,jvar);
416  }
417  }
418 
419  TH2F* h2 = new TH2F( *tm );
420  h2->SetNameTitle( hName, hTitle );
421 
422  for (UInt_t ivar=0; ivar<nvar; ivar++) {
423  h2->GetXaxis()->SetBinLabel( ivar+1, GetVariableInfo(ivar).GetTitle() );
424  h2->GetYaxis()->SetBinLabel( ivar+1, GetVariableInfo(ivar).GetTitle() );
425  }
426 
427  // present in percent, and round off digits
428  // also, use absolute value of correlation coefficient (ignore sign)
429  h2->Scale( 100.0 );
430  for (UInt_t ibin=1; ibin<=nvar; ibin++) {
431  for (UInt_t jbin=1; jbin<=nvar; jbin++) {
432  h2->SetBinContent( ibin, jbin, Int_t(h2->GetBinContent( ibin, jbin )) );
433  }
434  }
435 
436  // style settings
437  const Float_t labelSize = 0.055;
438  h2->SetStats( 0 );
439  h2->GetXaxis()->SetLabelSize( labelSize );
440  h2->GetYaxis()->SetLabelSize( labelSize );
441  h2->SetMarkerSize( 1.5 );
442  h2->SetMarkerColor( 0 );
443  h2->LabelsOption( "d" ); // diagonal labels on x axis
444  h2->SetLabelOffset( 0.011 );// label offset on x axis
445  h2->SetMinimum( -100.0 );
446  h2->SetMaximum( +100.0 );
447 
448  // -------------------------------------------------------------------------------------
449  // just in case one wants to change the position of the color palette axis
450  // -------------------------------------------------------------------------------------
451  // gROOT->SetStyle("Plain");
452  // TStyle* gStyle = gROOT->GetStyle( "Plain" );
453  // gStyle->SetPalette( 1, 0 );
454  // TPaletteAxis* paletteAxis
455  // = (TPaletteAxis*)h2->GetListOfFunctions()->FindObject( "palette" );
456  // -------------------------------------------------------------------------------------
457 
458  Log() << kDEBUG << Form("Dataset[%s] : ",fName.Data()) << "Created correlation matrix as 2D histogram: " << h2->GetName() << Endl;
459 
460  return h2;
461 }
462 
463 ////////////////////////////////////////////////////////////////////////////////
464 /// returns data set
465 
467 {
468  if (fDataSet==0 || fNeedsRebuilding) {
469  if(fDataSet!=0) ClearDataSet();
470  // fDataSet = DataSetManager::Instance().CreateDataSet(GetName()); //DSMTEST replaced by following lines
471  if( !fDataSetManager )
472  Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << "DataSetManager has not been set in DataSetInfo (GetDataSet() )." << Endl;
474 
476  }
477  return fDataSet;
478 }
479 
480 ////////////////////////////////////////////////////////////////////////////////
481 
483 {
484  if(all)
485  return fSpectators.size();
486  UInt_t nsp(0);
487  for(std::vector<VariableInfo>::const_iterator spit=fSpectators.begin(); spit!=fSpectators.end(); ++spit) {
488  if(spit->GetVarType()!='C') nsp++;
489  }
490  return nsp;
491 }
492 
493 ////////////////////////////////////////////////////////////////////////////////
494 
496 {
497  Int_t maxL = 0;
498  for (UInt_t cl = 0; cl < GetNClasses(); cl++) {
499  if (TString(GetClassInfo(cl)->GetName()).Length() > maxL) maxL = TString(GetClassInfo(cl)->GetName()).Length();
500  }
501 
502  return maxL;
503 }
504 
505 ////////////////////////////////////////////////////////////////////////////////
506 
508 {
509  Int_t maxL = 0;
510  for (UInt_t i = 0; i < GetNVariables(); i++) {
511  if (TString(GetVariableInfo(i).GetExpression()).Length() > maxL) maxL = TString(GetVariableInfo(i).GetExpression()).Length();
512  }
513 
514  return maxL;
515 }
516 
517 ////////////////////////////////////////////////////////////////////////////////
518 
520 {
521  Int_t maxL = 0;
522  for (UInt_t i = 0; i < GetNTargets(); i++) {
523  if (TString(GetTargetInfo(i).GetExpression()).Length() > maxL) maxL = TString(GetTargetInfo(i).GetExpression()).Length();
524  }
525 
526  return maxL;
527 }
528 
530  if (fTrainingSumSignalWeights<0) Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << " asking for the sum of training signal event weights which is not initicalised yet" << Endl;
532 }
534  if (fTrainingSumBackgrWeights<0) Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << " asking for the sum of training backgr event weights which is not initicalised yet" << Endl;
536 }
538  if (fTestingSumSignalWeights<0) Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << " asking for the sum of testing signal event weights which is not initicalised yet" << Endl;
539  return fTestingSumSignalWeights ;
540 }
542  if (fTestingSumBackgrWeights<0) Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << " asking for the sum of testing backgr event weights which is not initicalised yet" << Endl;
543  return fTestingSumBackgrWeights ;
544 }
545 
virtual void SetNameTitle(const char *name, const char *title)
Change the name and title of this histogram.
Definition: TH1.cxx:8042
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:5893
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
VariableInfo & AddTarget(const TString &expression, const TString &title, const TString &unit, Double_t min, Double_t max, Bool_t normalized=kTRUE, void *external=0)
add a variable (can be a complex expression) to the set of variables used in the MV analysis ...
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:399
virtual void LabelsOption(Option_t *option="h", Option_t *axis="X")
Set option(s) to draw axis with labels.
Definition: TH1.cxx:4908
Ssiz_t Length() const
Definition: TString.h:390
float Float_t
Definition: RtypesCore.h:53
void SetCut(const TCut &cut, const TString &className)
set the cut for the classes
TMatrixT< Float_t > TMatrixF
Definition: TMatrixFfwd.h:24
const TString & GetExpression() const
Definition: VariableInfo.h:65
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
Int_t GetClassNameMaxLength() const
UInt_t GetNClasses() const
Definition: DataSetInfo.h:154
std::vector< VariableInfo > fTargets
Definition: DataSetInfo.h:211
DataSet * CreateDataSet(const TString &dsiName)
Creates the singleton dataset.
UInt_t GetNTargets() const
Definition: DataSetInfo.h:129
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:400
Basic string class.
Definition: TString.h:137
std::vector< ClassInfo * > fClasses
Definition: DataSetInfo.h:215
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
std::vector< VariableInfo > fSpectators
Definition: DataSetInfo.h:212
void AddCut(const TCut &cut, const TString &className)
set the cut for the classes
std::vector< VariableInfo > fVariables
Definition: DataSetInfo.h:210
UInt_t GetNVariables() const
Definition: DataSetInfo.h:128
Double_t GetTrainingSumSignalWeights()
DataSet * GetDataSet() const
returns data set
const TString & GetInternalName() const
Definition: VariableInfo.h:66
virtual const char * GetName() const
Returns name of object.
Definition: DataSetInfo.h:85
const char * Data() const
Definition: TString.h:349
Tools & gTools()
Definition: Tools.cxx:79
DataSet * fDataSet
Definition: DataSetInfo.h:206
void SetWeight(const TString &weight)
Definition: ClassInfo.h:65
Bool_t IsSignal(const Event *ev) const
virtual ~DataSetInfo()
destructor
Definition: DataSetInfo.cxx:94
Double_t fTrainingSumSignalWeights
Definition: DataSetInfo.h:220
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:43
void PrintCorrelationMatrix(const TString &className)
calculates the correlation matrices for signal and background, prints them to standard output...
void SetMinType(EMsgType minType)
Definition: MsgLogger.h:76
Double_t fTestingSumSignalWeights
Definition: DataSetInfo.h:222
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:90
DataSetInfo(const TString &name="Default")
constructor
Definition: DataSetInfo.cxx:67
A specialized string object used for TTree selections.
Definition: TCut.h:27
void SetMsgType(EMsgType t) const
void SetCorrelationMatrix(const TString &className, TMatrixD *matrix)
void SetCorrelationMatrix(TMatrixD *matrix)
Definition: ClassInfo.h:68
MsgLogger * fLogger
Definition: DataSetInfo.h:234
Double_t fTrainingSumBackgrWeights
Definition: DataSetInfo.h:221
Double_t GetTestingSumSignalWeights()
Service class for 2-Dim histogram classes.
Definition: TH2.h:36
tomato 2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:255
VariableInfo & GetTargetInfo(Int_t i)
Definition: DataSetInfo.h:119
ClassInfo * GetClassInfo(Int_t clNum) const
EMsgType
Definition: Types.h:61
const TMatrixD * GetCorrelationMatrix() const
Definition: ClassInfo.h:74
const TMatrixD * CorrelationMatrix(const TString &className) const
void SetWeightExpression(const TString &exp, const TString &className="")
set the weight expressions for the classes if class name is specified, set only for this class if cla...
unsigned int UInt_t
Definition: RtypesCore.h:42
TMarker * m
Definition: textangle.C:8
char * Form(const char *fmt,...)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
TAxis * GetYaxis()
Definition: TH1.h:325
const TCut & GetCut() const
Definition: ClassInfo.h:72
virtual void SetLabelSize(Float_t size=0.04)
Set size of axis labels The size is expressed in per cent of the pad width.
Definition: TAttAxis.cxx:204
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition: TAttMarker.h:46
void PrintClasses() const
TH2 * CreateCorrelationMatrixHist(const TMatrixD *m, const TString &hName, const TString &hTitle) const
MsgLogger & Log() const
Definition: DataSetInfo.h:235
Int_t GetTargetNameMaxLength() const
double Double_t
Definition: RtypesCore.h:55
Double_t GetTrainingSumBackgrWeights()
int type
Definition: TGX11.cxx:120
void SetCut(const TCut &cut)
Definition: ClassInfo.h:66
VariableInfo & GetVariableInfo(Int_t i)
Definition: DataSetInfo.h:114
virtual void SetBinLabel(Int_t bin, const char *label)
Set label for bin.
Definition: TAxis.cxx:808
void ClearDataSet() const
ClassInfo * AddClass(const TString &className)
UInt_t GetClass() const
Definition: Event.h:89
void FormattedOutput(const std::vector< Double_t > &, const std::vector< TString > &, const TString titleVars, const TString titleValues, MsgLogger &logger, TString format="%+1.3f")
formatted output of simple table
Definition: Tools.cxx:896
Mother of all ROOT objects.
Definition: TObject.h:44
std::vector< Float_t > * fTargetsForMulticlass
Definition: DataSetInfo.h:232
VariableInfo & AddSpectator(const TString &expression, const TString &title, const TString &unit, Double_t min, Double_t max, char type= 'F', Bool_t normalized=kTRUE, void *external=0)
add a spectator (can be a complex expression) to the set of spectator variables used in the MV analys...
TMVA::DataSetManager * fDataSetManager
Definition: DataSetInfo.h:196
UInt_t GetNSpectators(bool all=kTRUE) const
Int_t FindVarIndex(const TString &) const
find variable by name
#define NULL
Definition: Rtypes.h:82
Int_t GetVariableNameMaxLength() const
VariableInfo & AddVariable(const TString &expression, const TString &title="", const TString &unit="", Double_t min=0, Double_t max=0, char varType='F', Bool_t normalized=kTRUE, void *external=0)
add a variable (can be a complex expression) to the set of variables used in the MV analysis ...
std::vector< TString > GetListOfVariables() const
returns list of variables
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2475
Bool_t fNeedsRebuilding
Definition: DataSetInfo.h:207
Double_t fTestingSumBackgrWeights
Definition: DataSetInfo.h:223
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual const char * GetTitle() const
Returns title of object.
Definition: TObject.cxx:460
THist< 2, float, THistStatContent, THistStatUncertainty > TH2F
Definition: THist.hxx:308
Bool_t HasCuts() const
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8058
std::vector< VariableInfo > & GetVariableInfos()
Definition: DataSetInfo.h:112
Double_t GetTestingSumBackgrWeights()
char name[80]
Definition: TGX11.cxx:109
TAxis * GetXaxis()
Definition: TH1.h:324
std::vector< Float_t > * GetTargetsForMulticlass(const Event *ev)
virtual void SetLabelOffset(Float_t offset=0.005, Option_t *axis="X")
Set offset between axis and axis&#39; labels.
Definition: Haxis.cxx:267