Logo ROOT   6.14/05
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 /*! \class TMVA::DataSetInfo
28 \ingroup TMVA
29 
30 Class that contains all the data information.
31 
32 */
33 
34 #include <vector>
35 
36 #include "TEventList.h"
37 #include "TFile.h"
38 #include "TH1.h"
39 #include "TH2.h"
40 #include "TProfile.h"
41 #include "TRandom3.h"
42 #include "TMatrixF.h"
43 #include "TVectorF.h"
44 #include "TMath.h"
45 #include "TROOT.h"
46 #include "TObjString.h"
47 
48 #include "TMVA/MsgLogger.h"
49 #include "TMVA/Tools.h"
50 #include "TMVA/DataSet.h"
51 #include "TMVA/DataSetInfo.h"
52 #include "TMVA/DataSetManager.h"
53 #include "TMVA/Event.h"
54 
55 #include "TMVA/Types.h"
56 #include "TMVA/VariableInfo.h"
57 
58 ////////////////////////////////////////////////////////////////////////////////
59 /// constructor
60 
62  : TObject(),
63  fDataSetManager(NULL),
64  fName(name),
65  fDataSet( 0 ),
66  fNeedsRebuilding( kTRUE ),
67  fVariables(),
68  fTargets(),
69  fSpectators(),
70  fClasses( 0 ),
71  fNormalization( "NONE" ),
72  fSplitOptions(""),
73  fTrainingSumSignalWeights(-1),
74  fTrainingSumBackgrWeights(-1),
75  fTestingSumSignalWeights (-1),
76  fTestingSumBackgrWeights (-1),
77  fOwnRootDir(0),
78  fVerbose( kFALSE ),
79  fSignalClass(0),
80  fTargetsForMulticlass(0),
81  fLogger( new MsgLogger("DataSetInfo", kINFO) )
82 {
83 }
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 /// destructor
87 
89 {
90  ClearDataSet();
91 
92  for(UInt_t i=0, iEnd = fClasses.size(); i<iEnd; ++i) {
93  delete fClasses[i];
94  }
95 
96  delete fTargetsForMulticlass;
97 
98  delete fLogger;
99 }
100 
101 ////////////////////////////////////////////////////////////////////////////////
102 
104 {
105  if(fDataSet!=0) { delete fDataSet; fDataSet=0; }
106 }
107 
108 ////////////////////////////////////////////////////////////////////////////////
109 
110 void
111 TMVA::DataSetInfo::SetMsgType( EMsgType t ) const
112 {
113  fLogger->SetMinType(t);
114 }
115 
116 ////////////////////////////////////////////////////////////////////////////////
117 
119 {
120  ClassInfo* theClass = GetClassInfo(className);
121  if (theClass) return theClass;
122 
123 
124  fClasses.push_back( new ClassInfo(className) );
125  fClasses.back()->SetNumber(fClasses.size()-1);
126 
127  //Log() << kHEADER << Endl;
128 
129  Log() << kHEADER << Form("[%s] : ",fName.Data()) << "Added class \"" << className << "\""<< Endl;
130 
131  Log() << kDEBUG <<"\t with internal class number " << fClasses.back()->GetNumber() << Endl;
132 
133 
134  if (className == "Signal") fSignalClass = fClasses.size()-1; // store the signal class index ( for comparison reasons )
135 
136  return fClasses.back();
137 }
138 
139 ////////////////////////////////////////////////////////////////////////////////
140 
142 {
143  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
144  if ((*it)->GetName() == name) return (*it);
145  }
146  return 0;
147 }
148 
149 ////////////////////////////////////////////////////////////////////////////////
150 
152 {
153  try {
154  return fClasses.at(cls);
155  }
156  catch(...) {
157  return 0;
158  }
159 }
160 
161 ////////////////////////////////////////////////////////////////////////////////
162 
164 {
165  for (UInt_t cls = 0; cls < GetNClasses() ; cls++) {
166  Log() << kINFO << Form("Dataset[%s] : ",fName.Data()) << "Class index : " << cls << " name : " << GetClassInfo(cls)->GetName() << Endl;
167  }
168 }
169 
170 ////////////////////////////////////////////////////////////////////////////////
171 
173 {
174  return (ev->GetClass() == fSignalClass);
175 }
176 
177 ////////////////////////////////////////////////////////////////////////////////
178 
179 std::vector<Float_t>* TMVA::DataSetInfo::GetTargetsForMulticlass( const TMVA::Event* ev )
180 {
181  if( !fTargetsForMulticlass ) fTargetsForMulticlass = new std::vector<Float_t>( GetNClasses() );
182  // fTargetsForMulticlass->resize( GetNClasses() );
183  fTargetsForMulticlass->assign( GetNClasses(), 0.0 );
184  fTargetsForMulticlass->at( ev->GetClass() ) = 1.0;
185  return fTargetsForMulticlass;
186 }
187 
188 
189 ////////////////////////////////////////////////////////////////////////////////
190 
192 {
193  Bool_t hasCuts = kFALSE;
194  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
195  if( TString((*it)->GetCut()) != TString("") ) hasCuts = kTRUE;
196  }
197  return hasCuts;
198 }
199 
200 ////////////////////////////////////////////////////////////////////////////////
201 
202 const TMatrixD* TMVA::DataSetInfo::CorrelationMatrix( const TString& className ) const
203 {
204  ClassInfo* ptr = GetClassInfo(className);
205  return ptr?ptr->GetCorrelationMatrix():0;
206 }
207 
208 ////////////////////////////////////////////////////////////////////////////////
209 /// add a variable (can be a complex expression) to the set of
210 /// variables used in the MV analysis
211 
213  const TString& title,
214  const TString& unit,
215  Double_t min, Double_t max,
216  char varType,
217  Bool_t normalized,
218  void* external )
219 {
220  TString regexpr = expression; // remove possible blanks
221  regexpr.ReplaceAll(" ", "" );
222  fVariables.push_back(VariableInfo( regexpr, title, unit,
223  fVariables.size()+1, varType, external, min, max, normalized ));
225  return fVariables.back();
226 }
227 
228 ////////////////////////////////////////////////////////////////////////////////
229 /// add variable with given VariableInfo
230 
232  fVariables.push_back(VariableInfo( varInfo ));
234  return fVariables.back();
235 }
236 
237 ////////////////////////////////////////////////////////////////////////////////
238 /// add a variable (can be a complex expression) to the set of
239 /// variables used in the MV analysis
240 
242  const TString& title,
243  const TString& unit,
244  Double_t min, Double_t max,
245  Bool_t normalized,
246  void* external )
247 {
248  TString regexpr = expression; // remove possible blanks
249  regexpr.ReplaceAll(" ", "" );
250  char type='F';
251  fTargets.push_back(VariableInfo( regexpr, title, unit,
252  fTargets.size()+1, type, external, min,
253  max, normalized ));
255  return fTargets.back();
256 }
257 
258 ////////////////////////////////////////////////////////////////////////////////
259 /// add target with given VariableInfo
260 
262  fTargets.push_back(VariableInfo( varInfo ));
264  return fTargets.back();
265 }
266 
267 ////////////////////////////////////////////////////////////////////////////////
268 /// add a spectator (can be a complex expression) to the set of spectator variables used in
269 /// the MV analysis
270 
272  const TString& title,
273  const TString& unit,
274  Double_t min, Double_t max, char type,
275  Bool_t normalized, void* external )
276 {
277  TString regexpr = expression; // remove possible blanks
278  regexpr.ReplaceAll(" ", "" );
279  fSpectators.push_back(VariableInfo( regexpr, title, unit,
280  fSpectators.size()+1, type, external, min, max, normalized ));
282  return fSpectators.back();
283 }
284 
285 ////////////////////////////////////////////////////////////////////////////////
286 /// add spectator with given VariableInfo
287 
289  fSpectators.push_back(VariableInfo( varInfo ));
291  return fSpectators.back();
292 }
293 
294 ////////////////////////////////////////////////////////////////////////////////
295 /// find variable by name
296 
298 {
299  for (UInt_t ivar=0; ivar<GetNVariables(); ivar++)
300  if (var == GetVariableInfo(ivar).GetInternalName()) return ivar;
301 
302  for (UInt_t ivar=0; ivar<GetNVariables(); ivar++)
303  Log() << kINFO << Form("Dataset[%s] : ",fName.Data()) << GetVariableInfo(ivar).GetInternalName() << Endl;
304 
305  Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << "<FindVarIndex> Variable \'" << var << "\' not found." << Endl;
306 
307  return -1;
308 }
309 
310 ////////////////////////////////////////////////////////////////////////////////
311 /// set the weight expressions for the classes
312 /// if class name is specified, set only for this class
313 /// if class name is unknown, register new class with this name
314 
315 void TMVA::DataSetInfo::SetWeightExpression( const TString& expr, const TString& className )
316 {
317  if (className != "") {
318  TMVA::ClassInfo* ci = AddClass(className);
319  ci->SetWeight( expr );
320  }
321  else {
322  // no class name specified, set weight for all classes
323  if (fClasses.empty()) {
324  Log() << kWARNING << Form("Dataset[%s] : ",fName.Data()) << "No classes registered yet, cannot specify weight expression!" << Endl;
325  }
326  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
327  (*it)->SetWeight( expr );
328  }
329  }
330 }
331 
332 ////////////////////////////////////////////////////////////////////////////////
333 
335 {
336  GetClassInfo(className)->SetCorrelationMatrix(matrix);
337 }
338 
339 ////////////////////////////////////////////////////////////////////////////////
340 /// set the cut for the classes
341 
342 void TMVA::DataSetInfo::SetCut( const TCut& cut, const TString& className )
343 {
344  if (className == "") { // if no className has been given set the cut for all the classes
345  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
346  (*it)->SetCut( cut );
347  }
348  }
349  else {
350  TMVA::ClassInfo* ci = AddClass(className);
351  ci->SetCut( cut );
352  }
353 }
354 
355 ////////////////////////////////////////////////////////////////////////////////
356 /// set the cut for the classes
357 
358 void TMVA::DataSetInfo::AddCut( const TCut& cut, const TString& className )
359 {
360  if (className == "") { // if no className has been given set the cut for all the classes
361  for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
362  const TCut& oldCut = (*it)->GetCut();
363  (*it)->SetCut( oldCut+cut );
364  }
365  }
366  else {
367  TMVA::ClassInfo* ci = AddClass(className);
368  ci->SetCut( ci->GetCut()+cut );
369  }
370 }
371 
372 ////////////////////////////////////////////////////////////////////////////////
373 /// returns list of variables
374 
375 std::vector<TString> TMVA::DataSetInfo::GetListOfVariables() const
376 {
377  std::vector<TString> vNames;
378  std::vector<TMVA::VariableInfo>::const_iterator viIt = GetVariableInfos().begin();
379  for(;viIt != GetVariableInfos().end(); ++viIt) vNames.push_back( (*viIt).GetExpression() );
380 
381  return vNames;
382 }
383 
384 ////////////////////////////////////////////////////////////////////////////////
385 /// calculates the correlation matrices for signal and background,
386 /// prints them to standard output, and fills 2D histograms
387 
389 {
390 
391  Log() << kHEADER //<< Form("Dataset[%s] : ",fName.Data())
392  << "Correlation matrix (" << className << "):" << Endl;
394 }
395 
396 ////////////////////////////////////////////////////////////////////////////////
397 
399  const TString& hName,
400  const TString& hTitle ) const
401 {
402  if (m==0) return 0;
403 
404  const UInt_t nvar = GetNVariables();
405 
406  // workaround till the TMatrix templates are commonly used
407  // this keeps backward compatibility
408  TMatrixF* tm = new TMatrixF( nvar, nvar );
409  for (UInt_t ivar=0; ivar<nvar; ivar++) {
410  for (UInt_t jvar=0; jvar<nvar; jvar++) {
411  (*tm)(ivar, jvar) = (*m)(ivar,jvar);
412  }
413  }
414 
415  TH2F* h2 = new TH2F( *tm );
416  h2->SetNameTitle( hName, hTitle );
417 
418  for (UInt_t ivar=0; ivar<nvar; ivar++) {
419  h2->GetXaxis()->SetBinLabel( ivar+1, GetVariableInfo(ivar).GetTitle() );
420  h2->GetYaxis()->SetBinLabel( ivar+1, GetVariableInfo(ivar).GetTitle() );
421  }
422 
423  // present in percent, and round off digits
424  // also, use absolute value of correlation coefficient (ignore sign)
425  h2->Scale( 100.0 );
426  for (UInt_t ibin=1; ibin<=nvar; ibin++) {
427  for (UInt_t jbin=1; jbin<=nvar; jbin++) {
428  h2->SetBinContent( ibin, jbin, Int_t(h2->GetBinContent( ibin, jbin )) );
429  }
430  }
431 
432  // style settings
433  const Float_t labelSize = 0.055;
434  h2->SetStats( 0 );
435  h2->GetXaxis()->SetLabelSize( labelSize );
436  h2->GetYaxis()->SetLabelSize( labelSize );
437  h2->SetMarkerSize( 1.5 );
438  h2->SetMarkerColor( 0 );
439  h2->LabelsOption( "d" ); // diagonal labels on x axis
440  h2->SetLabelOffset( 0.011 );// label offset on x axis
441  h2->SetMinimum( -100.0 );
442  h2->SetMaximum( +100.0 );
443 
444  // -------------------------------------------------------------------------------------
445  // just in case one wants to change the position of the color palette axis
446  // -------------------------------------------------------------------------------------
447  // gROOT->SetStyle("Plain");
448  // TStyle* gStyle = gROOT->GetStyle( "Plain" );
449  // gStyle->SetPalette( 1, 0 );
450  // TPaletteAxis* paletteAxis
451  // = (TPaletteAxis*)h2->GetListOfFunctions()->FindObject( "palette" );
452  // -------------------------------------------------------------------------------------
453 
454  Log() << kDEBUG << Form("Dataset[%s] : ",fName.Data()) << "Created correlation matrix as 2D histogram: " << h2->GetName() << Endl;
455 
456  return h2;
457 }
458 
459 ////////////////////////////////////////////////////////////////////////////////
460 /// returns data set
461 
463 {
464  if (fDataSet==0 || fNeedsRebuilding) {
465  if(fDataSet!=0) ClearDataSet();
466  // fDataSet = DataSetManager::Instance().CreateDataSet(GetName()); //DSMTEST replaced by following lines
467  if( !fDataSetManager )
468  Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << "DataSetManager has not been set in DataSetInfo (GetDataSet() )." << Endl;
470 
472  }
473  return fDataSet;
474 }
475 
476 ////////////////////////////////////////////////////////////////////////////////
477 
479 {
480  if(all)
481  return fSpectators.size();
482  UInt_t nsp(0);
483  for(std::vector<VariableInfo>::const_iterator spit=fSpectators.begin(); spit!=fSpectators.end(); ++spit) {
484  if(spit->GetVarType()!='C') nsp++;
485  }
486  return nsp;
487 }
488 
489 ////////////////////////////////////////////////////////////////////////////////
490 
492 {
493  Int_t maxL = 0;
494  for (UInt_t cl = 0; cl < GetNClasses(); cl++) {
495  if (TString(GetClassInfo(cl)->GetName()).Length() > maxL) maxL = TString(GetClassInfo(cl)->GetName()).Length();
496  }
497 
498  return maxL;
499 }
500 
501 ////////////////////////////////////////////////////////////////////////////////
502 
504 {
505  Int_t maxL = 0;
506  for (UInt_t i = 0; i < GetNVariables(); i++) {
507  if (TString(GetVariableInfo(i).GetExpression()).Length() > maxL) maxL = TString(GetVariableInfo(i).GetExpression()).Length();
508  }
509 
510  return maxL;
511 }
512 
513 ////////////////////////////////////////////////////////////////////////////////
514 
516 {
517  Int_t maxL = 0;
518  for (UInt_t i = 0; i < GetNTargets(); i++) {
519  if (TString(GetTargetInfo(i).GetExpression()).Length() > maxL) maxL = TString(GetTargetInfo(i).GetExpression()).Length();
520  }
521 
522  return maxL;
523 }
524 
525 ////////////////////////////////////////////////////////////////////////////////
526 
528  if (fTrainingSumSignalWeights<0) Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << " asking for the sum of training signal event weights which is not initialized yet" << Endl;
530 }
531 
532 ////////////////////////////////////////////////////////////////////////////////
533 
535  if (fTrainingSumBackgrWeights<0) Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << " asking for the sum of training backgr event weights which is not initialized yet" << Endl;
537 }
538 
539 ////////////////////////////////////////////////////////////////////////////////
540 
542  if (fTestingSumSignalWeights<0) Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << " asking for the sum of testing signal event weights which is not initialized yet" << Endl;
543  return fTestingSumSignalWeights ;
544 }
545 
546 ////////////////////////////////////////////////////////////////////////////////
547 
549  if (fTestingSumBackgrWeights<0) Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << " asking for the sum of testing backgr event weights which is not initialized yet" << Endl;
550  return fTestingSumBackgrWeights ;
551 }
552 
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual void SetNameTitle(const char *name, const char *title)
Change the name and title of this histogram.
Definition: TH1.cxx:8268
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6101
UInt_t GetNVariables() const
Definition: DataSetInfo.h:110
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:158
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:390
auto * m
Definition: textangle.C:8
virtual void LabelsOption(Option_t *option="h", Option_t *axis="X")
Set option(s) to draw axis with labels.
Definition: TH1.cxx:5085
const TString & GetInternalName() const
Definition: VariableInfo.h:58
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:22
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
std::vector< TString > GetListOfVariables() const
returns list of variables
std::vector< VariableInfo > fTargets
Definition: DataSetInfo.h:193
DataSet * CreateDataSet(const TString &dsiName)
Creates the singleton dataset.
Int_t GetTargetNameMaxLength() const
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:391
Basic string class.
Definition: TString.h:131
std::vector< ClassInfo * > fClasses
Definition: DataSetInfo.h:197
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
void ClearDataSet() const
UInt_t GetNClasses() const
Definition: DataSetInfo.h:136
std::vector< VariableInfo > fSpectators
Definition: DataSetInfo.h:194
void AddCut(const TCut &cut, const TString &className)
set the cut for the classes
std::vector< VariableInfo > fVariables
Definition: DataSetInfo.h:192
Double_t GetTrainingSumSignalWeights()
const TString & GetExpression() const
Definition: VariableInfo.h:57
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...
Class that contains all the information of a class.
Definition: ClassInfo.h:49
DataSet * fDataSet
Definition: DataSetInfo.h:188
TH2 * CreateCorrelationMatrixHist(const TMatrixD *m, const TString &hName, const TString &hTitle) const
void SetWeight(const TString &weight)
Definition: ClassInfo.h:57
virtual ~DataSetInfo()
destructor
Definition: DataSetInfo.cxx:88
Double_t fTrainingSumSignalWeights
Definition: DataSetInfo.h:202
UInt_t GetClass() const
Definition: Event.h:81
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
MsgLogger & Log() const
message logger
Definition: DataSetInfo.h:217
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:72
Double_t fTestingSumSignalWeights
Definition: DataSetInfo.h:204
DataSetInfo(const TString &name="Default")
constructor
Definition: DataSetInfo.cxx:61
A specialized string object used for TTree selections.
Definition: TCut.h:25
void SetCorrelationMatrix(const TString &className, TMatrixD *matrix)
void SetCorrelationMatrix(TMatrixD *matrix)
Definition: ClassInfo.h:60
Class that contains all the data information.
Definition: DataSet.h:69
MsgLogger * fLogger
Definition: DataSetInfo.h:216
Double_t fTrainingSumBackgrWeights
Definition: DataSetInfo.h:203
UInt_t GetNTargets() const
Definition: DataSetInfo.h:111
Double_t GetTestingSumSignalWeights()
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
ClassInfo * GetClassInfo(Int_t clNum) const
const TMatrixD * CorrelationMatrix(const TString &className) const
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:250
VariableInfo & GetTargetInfo(Int_t i)
Definition: DataSetInfo.h:101
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
char * Form(const char *fmt,...)
Ssiz_t Length() const
Definition: TString.h:405
UInt_t GetNSpectators(bool all=kTRUE) const
TAxis * GetYaxis()
Definition: TH1.h:316
void SetMsgType(EMsgType t) const
Tools & gTools()
void PrintClasses() const
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:41
Int_t FindVarIndex(const TString &) const
find variable by name
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
Int_t GetVariableNameMaxLength() const
Double_t GetTrainingSumBackgrWeights()
int type
Definition: TGX11.cxx:120
const TMatrixD * GetCorrelationMatrix() const
Definition: ClassInfo.h:66
void SetCut(const TCut &cut)
Definition: ClassInfo.h:58
VariableInfo & GetVariableInfo(Int_t i)
Definition: DataSetInfo.h:96
virtual void SetBinLabel(Int_t bin, const char *label)
Set label for bin.
Definition: TAxis.cxx:809
ClassInfo * AddClass(const TString &className)
Int_t GetClassNameMaxLength() const
virtual const char * GetName() const
Returns name of object.
Definition: DataSetInfo.h:67
ostringstream derivative to redirect and format output
Definition: MsgLogger.h:59
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:899
Mother of all ROOT objects.
Definition: TObject.h:37
virtual const char * GetTitle() const
Returns title of object.
Definition: TObject.cxx:401
std::vector< Float_t > * fTargetsForMulticlass
Definition: DataSetInfo.h:214
const TCut & GetCut() const
Definition: ClassInfo.h:64
TMVA::DataSetManager * fDataSetManager
Definition: DataSetInfo.h:178
Bool_t HasCuts() const
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:84
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 ...
Bool_t IsSignal(const Event *ev) const
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2500
Bool_t fNeedsRebuilding
Definition: DataSetInfo.h:189
Class for type info of MVA input variable.
Definition: VariableInfo.h:47
Double_t fTestingSumBackgrWeights
Definition: DataSetInfo.h:205
THist< 2, float, THistStatContent, THistStatUncertainty > TH2F
Definition: THist.hxx:291
const Bool_t kTRUE
Definition: RtypesCore.h:87
DataSet * GetDataSet() const
returns data set
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8284
std::vector< VariableInfo > & GetVariableInfos()
Definition: DataSetInfo.h:94
Double_t GetTestingSumBackgrWeights()
char name[80]
Definition: TGX11.cxx:109
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:315
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
const char * Data() const
Definition: TString.h:364