Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 * *
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 * (see tmva/doc/LICENSE) *
25 **********************************************************************************/
26
27/*! \class TMVA::DataSetInfo
28\ingroup TMVA
29
30Class that contains all the data information.
31
32*/
33
34#include <vector>
35
36#include "TEventList.h"
37#include "TH2.h"
38#include "TRandom3.h"
39#include "TMatrixF.h"
40#include "TVectorF.h"
41#include "TROOT.h"
42
43#include "TMVA/MsgLogger.h"
44#include "TMVA/Tools.h"
45#include "TMVA/DataSet.h"
46#include "TMVA/DataSetInfo.h"
47#include "TMVA/DataSetManager.h"
48#include "TMVA/Event.h"
49
50#include "TMVA/Types.h"
51#include "TMVA/VariableInfo.h"
52
53////////////////////////////////////////////////////////////////////////////////
54/// constructor
55
57 : TObject(),
58 fDataSetManager(NULL),
59 fName(name),
60 fDataSet( 0 ),
61 fNeedsRebuilding( kTRUE ),
62 fVariables(),
63 fTargets(),
64 fSpectators(),
65 fClasses( 0 ),
66 fNormalization( "NONE" ),
67 fSplitOptions(""),
68 fTrainingSumSignalWeights(-1),
69 fTrainingSumBackgrWeights(-1),
70 fTestingSumSignalWeights (-1),
71 fTestingSumBackgrWeights (-1),
72 fOwnRootDir(0),
73 fVerbose( kFALSE ),
74 fSignalClass(0),
75 fTargetsForMulticlass(0),
76 fLogger( new MsgLogger("DataSetInfo", kINFO) )
77{
78}
79
80////////////////////////////////////////////////////////////////////////////////
81/// destructor
82
84{
85 ClearDataSet();
86
87 for(UInt_t i=0, iEnd = fClasses.size(); i<iEnd; ++i) {
88 if (fClasses[i]) delete fClasses[i];
89 }
90
91 if (fTargetsForMulticlass) delete fTargetsForMulticlass;
92
93 delete fLogger;
94}
95
96////////////////////////////////////////////////////////////////////////////////
97
99{
100 if(fDataSet) { delete fDataSet; fDataSet=nullptr; }
101}
102
103////////////////////////////////////////////////////////////////////////////////
104
105void
107{
108 fLogger->SetMinType(t);
109}
110
111////////////////////////////////////////////////////////////////////////////////
112
114{
115 ClassInfo* theClass = GetClassInfo(className);
116 if (theClass) return theClass;
117
118
119 fClasses.push_back( new ClassInfo(className) );
120 fClasses.back()->SetNumber(fClasses.size()-1);
121
122 //Log() << kHEADER << Endl;
123
124 Log() << kHEADER << Form("[%s] : ",fName.Data()) << "Added class \"" << className << "\""<< Endl;
125
126 Log() << kDEBUG <<"\t with internal class number " << fClasses.back()->GetNumber() << Endl;
127
128
129 if (className == "Signal") fSignalClass = fClasses.size()-1; // store the signal class index ( for comparison reasons )
130
131 return fClasses.back();
132}
133
134////////////////////////////////////////////////////////////////////////////////
135
137{
138 for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
139 if ((*it)->GetName() == name) return (*it);
140 }
141 return 0;
142}
143
144////////////////////////////////////////////////////////////////////////////////
145
147{
148 try {
149 return fClasses.at(cls);
150 }
151 catch(...) {
152 return 0;
153 }
154}
155
156////////////////////////////////////////////////////////////////////////////////
157
159{
160 for (UInt_t cls = 0; cls < GetNClasses() ; cls++) {
161 Log() << kINFO << Form("Dataset[%s] : ",fName.Data()) << "Class index : " << cls << " name : " << GetClassInfo(cls)->GetName() << Endl;
162 }
163}
164
165////////////////////////////////////////////////////////////////////////////////
166
168{
169 return (ev->GetClass() == fSignalClass);
170}
171
172////////////////////////////////////////////////////////////////////////////////
173
175{
176 if( !fTargetsForMulticlass ) fTargetsForMulticlass = new std::vector<Float_t>( GetNClasses() );
177 // fTargetsForMulticlass->resize( GetNClasses() );
178 fTargetsForMulticlass->assign( GetNClasses(), 0.0 );
179 fTargetsForMulticlass->at( ev->GetClass() ) = 1.0;
180 return fTargetsForMulticlass;
181}
182
183
184////////////////////////////////////////////////////////////////////////////////
185
187{
188 Bool_t hasCuts = kFALSE;
189 for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
190 if( TString((*it)->GetCut()) != TString("") ) hasCuts = kTRUE;
191 }
192 return hasCuts;
193}
194
195////////////////////////////////////////////////////////////////////////////////
196
198{
199 ClassInfo* ptr = GetClassInfo(className);
200 return ptr?ptr->GetCorrelationMatrix():0;
201}
202
203////////////////////////////////////////////////////////////////////////////////
204/// add a variable (can be a complex expression) to the set of
205/// variables used in the MV analysis
206
208 const TString& title,
209 const TString& unit,
210 Double_t min, Double_t max,
211 char varType,
212 Bool_t normalized,
213 void* external )
214{
215 TString regexpr = expression; // remove possible blanks
216 regexpr.ReplaceAll(" ", "" );
217 fVariables.push_back(VariableInfo( regexpr, title, unit,
218 fVariables.size()+1, varType, external, min, max, normalized ));
219 fNeedsRebuilding = kTRUE;
220 return fVariables.back();
221}
222
223////////////////////////////////////////////////////////////////////////////////
224/// add variable with given VariableInfo
225
227 fVariables.push_back(VariableInfo( varInfo ));
228 fNeedsRebuilding = kTRUE;
229 return fVariables.back();
230}
231
232////////////////////////////////////////////////////////////////////////////////
233/// add an array of variables identified by an expression corresponding to an array entry in the tree
234
235void TMVA::DataSetInfo::AddVariablesArray(const TString &expression, Int_t size, const TString &title, const TString &unit,
236 Double_t min, Double_t max, char varType, Bool_t normalized,
237 void *external)
238{
239 TString regexpr = expression; // remove possible blanks
240 regexpr.ReplaceAll(" ", "");
241 fVariables.reserve(fVariables.size() + size);
242 for (int i = 0; i < size; ++i) {
243 TString newTitle = title + TString::Format("[%d]", i);
244
245 fVariables.emplace_back(regexpr, newTitle, unit, fVariables.size() + 1, varType, external, min, max, normalized);
246 // set corresponding bit indicating is a variable from an array
247 fVariables.back().SetBit(kIsArrayVariable);
248 TString newVarName = fVariables.back().GetInternalName() + TString::Format("[%d]", i);
249 fVariables.back().SetInternalName(newVarName);
250
251 // move "external" pointer to the next variable in the array
252 if (varType == 'F') {
253 float *ptr = (float *)external;
254 ++ptr;
255 external = (void *)ptr;
256 } else if (varType == 'I') {
257 int *ptr = (int *)external;
258 ++ptr;
259 external = (void *)ptr;
260 } else {
261 Error("TMVA::DataSetInfo::AddVariablesArray", "'%c' variable type is not supported", varType);
262 }
263 }
264 fVarArrays[regexpr] = size;
265 fNeedsRebuilding = kTRUE;
266}
267
268////////////////////////////////////////////////////////////////////////////////
269/// add a variable (can be a complex expression) to the set of
270/// variables used in the MV analysis
271
273 const TString& title,
274 const TString& unit,
275 Double_t min, Double_t max,
276 Bool_t normalized,
277 void* external )
278{
279 TString regexpr = expression; // remove possible blanks
280 regexpr.ReplaceAll(" ", "" );
281 char type='F';
282 fTargets.push_back(VariableInfo( regexpr, title, unit,
283 fTargets.size()+1, type, external, min,
284 max, normalized ));
285 fNeedsRebuilding = kTRUE;
286 return fTargets.back();
287}
288
289////////////////////////////////////////////////////////////////////////////////
290/// add target with given VariableInfo
291
293 fTargets.push_back(VariableInfo( varInfo ));
294 fNeedsRebuilding = kTRUE;
295 return fTargets.back();
296}
297
298////////////////////////////////////////////////////////////////////////////////
299/// add a spectator (can be a complex expression) to the set of spectator variables used in
300/// the MV analysis
301
303 const TString& title,
304 const TString& unit,
305 Double_t min, Double_t max, char type,
306 Bool_t normalized, void* external )
307{
308 TString regexpr = expression; // remove possible blanks
309 regexpr.ReplaceAll(" ", "" );
310 fSpectators.push_back(VariableInfo( regexpr, title, unit,
311 fSpectators.size()+1, type, external, min, max, normalized ));
312 fNeedsRebuilding = kTRUE;
313 return fSpectators.back();
314}
315
316////////////////////////////////////////////////////////////////////////////////
317/// add spectator with given VariableInfo
318
320 fSpectators.push_back(VariableInfo( varInfo ));
321 fNeedsRebuilding = kTRUE;
322 return fSpectators.back();
323}
324
325////////////////////////////////////////////////////////////////////////////////
326/// find variable by name
327
329{
330 for (UInt_t ivar=0; ivar<GetNVariables(); ivar++)
331 if (var == GetVariableInfo(ivar).GetInternalName()) return ivar;
332
333 for (UInt_t ivar=0; ivar<GetNVariables(); ivar++)
334 Log() << kINFO << Form("Dataset[%s] : ",fName.Data()) << GetVariableInfo(ivar).GetInternalName() << Endl;
335
336 Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << "<FindVarIndex> Variable \'" << var << "\' not found." << Endl;
337
338 return -1;
339}
340
341////////////////////////////////////////////////////////////////////////////////
342/// set the weight expressions for the classes
343/// if class name is specified, set only for this class
344/// if class name is unknown, register new class with this name
345
346void TMVA::DataSetInfo::SetWeightExpression( const TString& expr, const TString& className )
347{
348 if (className != "") {
349 TMVA::ClassInfo* ci = AddClass(className);
350 ci->SetWeight( expr );
351 }
352 else {
353 // no class name specified, set weight for all classes
354 if (fClasses.empty()) {
355 Log() << kWARNING << Form("Dataset[%s] : ",fName.Data()) << "No classes registered yet, cannot specify weight expression!" << Endl;
356 }
357 for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
358 (*it)->SetWeight( expr );
359 }
360 }
361}
362
363////////////////////////////////////////////////////////////////////////////////
364
366{
367 GetClassInfo(className)->SetCorrelationMatrix(matrix);
368}
369
370////////////////////////////////////////////////////////////////////////////////
371/// set the cut for the classes
372
373void TMVA::DataSetInfo::SetCut( const TCut& cut, const TString& className )
374{
375 if (className == "") { // if no className has been given set the cut for all the classes
376 for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
377 (*it)->SetCut( cut );
378 }
379 }
380 else {
381 TMVA::ClassInfo* ci = AddClass(className);
382 ci->SetCut( cut );
383 }
384}
385
386////////////////////////////////////////////////////////////////////////////////
387/// set the cut for the classes
388
389void TMVA::DataSetInfo::AddCut( const TCut& cut, const TString& className )
390{
391 if (className == "") { // if no className has been given set the cut for all the classes
392 for (std::vector<ClassInfo*>::iterator it = fClasses.begin(); it < fClasses.end(); ++it) {
393 const TCut& oldCut = (*it)->GetCut();
394 (*it)->SetCut( oldCut+cut );
395 }
396 }
397 else {
398 TMVA::ClassInfo* ci = AddClass(className);
399 ci->SetCut( ci->GetCut()+cut );
400 }
401}
402
403////////////////////////////////////////////////////////////////////////////////
404/// returns list of variables
405
406std::vector<TString> TMVA::DataSetInfo::GetListOfVariables() const
407{
408 std::vector<TString> vNames;
409 std::vector<TMVA::VariableInfo>::const_iterator viIt = GetVariableInfos().begin();
410 for(;viIt != GetVariableInfos().end(); ++viIt) vNames.push_back( (*viIt).GetInternalName() );
411
412 return vNames;
413}
414
415////////////////////////////////////////////////////////////////////////////////
416/// calculates the correlation matrices for signal and background,
417/// prints them to standard output, and fills 2D histograms
418
420{
421
422 Log() << kHEADER //<< Form("Dataset[%s] : ",fName.Data())
423 << "Correlation matrix (" << className << "):" << Endl;
424 gTools().FormattedOutput( *CorrelationMatrix( className ), GetListOfVariables(), Log() );
425}
426
427////////////////////////////////////////////////////////////////////////////////
428
430 const TString& hName,
431 const TString& hTitle ) const
432{
433 if (m==0) return 0;
434
435 const UInt_t nvar = GetNVariables();
436
437 // workaround till the TMatrix templates are commonly used
438 // this keeps backward compatibility
439 TMatrixF* tm = new TMatrixF( nvar, nvar );
440 for (UInt_t ivar=0; ivar<nvar; ivar++) {
441 for (UInt_t jvar=0; jvar<nvar; jvar++) {
442 (*tm)(ivar, jvar) = (*m)(ivar,jvar);
443 }
444 }
445
446 TH2F* h2 = new TH2F( *tm );
447 h2->SetNameTitle( hName, hTitle );
448
449 for (UInt_t ivar=0; ivar<nvar; ivar++) {
450 h2->GetXaxis()->SetBinLabel( ivar+1, GetVariableInfo(ivar).GetTitle() );
451 h2->GetYaxis()->SetBinLabel( ivar+1, GetVariableInfo(ivar).GetTitle() );
452 }
453
454 // present in percent, and round off digits
455 // also, use absolute value of correlation coefficient (ignore sign)
456 h2->Scale( 100.0 );
457 for (UInt_t ibin=1; ibin<=nvar; ibin++) {
458 for (UInt_t jbin=1; jbin<=nvar; jbin++) {
459 h2->SetBinContent( ibin, jbin, Int_t(h2->GetBinContent( ibin, jbin )) );
460 }
461 }
462
463 // style settings
464 const Float_t labelSize = 0.055;
465 h2->SetStats( 0 );
466 h2->GetXaxis()->SetLabelSize( labelSize );
467 h2->GetYaxis()->SetLabelSize( labelSize );
468 h2->SetMarkerSize( 1.5 );
469 h2->SetMarkerColor( 0 );
470 h2->LabelsOption( "d" ); // diagonal labels on x axis
471 h2->SetLabelOffset( 0.011 );// label offset on x axis
472 h2->SetMinimum( -100.0 );
473 h2->SetMaximum( +100.0 );
474
475 // -------------------------------------------------------------------------------------
476 // just in case one wants to change the position of the color palette axis
477 // -------------------------------------------------------------------------------------
478 // gROOT->SetStyle("Plain");
479 // TStyle* gStyle = gROOT->GetStyle( "Plain" );
480 // gStyle->SetPalette( 1, 0 );
481 // TPaletteAxis* paletteAxis
482 // = (TPaletteAxis*)h2->GetListOfFunctions()->FindObject( "palette" );
483 // -------------------------------------------------------------------------------------
484
485 Log() << kDEBUG << Form("Dataset[%s] : ",fName.Data()) << "Created correlation matrix as 2D histogram: " << h2->GetName() << Endl;
486
487 return h2;
488}
489
490////////////////////////////////////////////////////////////////////////////////
491/// returns data set
492
494{
495 if (fDataSet==0 || fNeedsRebuilding) {
496 if (fNeedsRebuilding) Log() << kINFO << "Rebuilding Dataset " << fName << Endl;
497 if (fDataSet != 0)
498 ClearDataSet();
499 // fDataSet = DataSetManager::Instance().CreateDataSet(GetName()); //DSMTEST replaced by following lines
500 if( !fDataSetManager )
501 Log() << kFATAL << Form("Dataset[%s] : ",fName.Data()) << "DataSetManager has not been set in DataSetInfo (GetDataSet() )." << Endl;
502 fDataSet = fDataSetManager->CreateDataSet(GetName());
503
504 fNeedsRebuilding = kFALSE;
505 }
506 return fDataSet;
507}
508
509////////////////////////////////////////////////////////////////////////////////
510
512{
513 if(all)
514 return fSpectators.size();
515 UInt_t nsp(0);
516 for(std::vector<VariableInfo>::const_iterator spit=fSpectators.begin(); spit!=fSpectators.end(); ++spit) {
517 if(spit->GetVarType()!='C') nsp++;
518 }
519 return nsp;
520}
521
522////////////////////////////////////////////////////////////////////////////////
523
525{
526 Int_t maxL = 0;
527 for (UInt_t cl = 0; cl < GetNClasses(); cl++) {
528 if (TString(GetClassInfo(cl)->GetName()).Length() > maxL) maxL = TString(GetClassInfo(cl)->GetName()).Length();
529 }
530
531 return maxL;
532}
533
534////////////////////////////////////////////////////////////////////////////////
535
537{
538 Int_t maxL = 0;
539 for (UInt_t i = 0; i < GetNVariables(); i++) {
540 if (TString(GetVariableInfo(i).GetExpression()).Length() > maxL) maxL = TString(GetVariableInfo(i).GetExpression()).Length();
541 }
542
543 return maxL;
544}
545
546////////////////////////////////////////////////////////////////////////////////
547
549{
550 Int_t maxL = 0;
551 for (UInt_t i = 0; i < GetNTargets(); i++) {
552 if (TString(GetTargetInfo(i).GetExpression()).Length() > maxL) maxL = TString(GetTargetInfo(i).GetExpression()).Length();
553 }
554
555 return maxL;
556}
557
558////////////////////////////////////////////////////////////////////////////////
559
561 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;
562 return fTrainingSumSignalWeights;
563}
564
565////////////////////////////////////////////////////////////////////////////////
566
568 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;
569 return fTrainingSumBackgrWeights;
570}
571
572////////////////////////////////////////////////////////////////////////////////
573
575 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;
576 return fTestingSumSignalWeights ;
577}
578
579////////////////////////////////////////////////////////////////////////////////
580
582 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;
583 return fTestingSumBackgrWeights ;
584}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
TMatrixT< Float_t > TMatrixF
Definition TMatrixFfwd.h:23
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
A specialized string object used for TTree selections.
Definition TCut.h:25
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:308
Service class for 2-D histogram classes.
Definition TH2.h:30
Class that contains all the information of a class.
Definition ClassInfo.h:49
const TMatrixD * GetCorrelationMatrix() const
Definition ClassInfo.h:66
const TCut & GetCut() const
Definition ClassInfo.h:64
void SetCut(const TCut &cut)
Definition ClassInfo.h:58
void SetWeight(const TString &weight)
Definition ClassInfo.h:57
void SetNumber(const UInt_t index)
Definition ClassInfo.h:59
Bool_t HasCuts() const
UInt_t GetNSpectators(bool all=kTRUE) 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=nullptr)
add a variable (can be a complex expression) to the set of variables used in the MV analysis
ClassInfo * AddClass(const TString &className)
const TMatrixD * CorrelationMatrix(const TString &className) const
Int_t GetTargetNameMaxLength() const
virtual ~DataSetInfo()
destructor
Double_t GetTestingSumBackgrWeights()
void SetMsgType(EMsgType t) const
DataSet * GetDataSet() const
returns data set
DataSetInfo(const TString &name="Default")
constructor
TH2 * CreateCorrelationMatrixHist(const TMatrixD *m, const TString &hName, const TString &hTitle) const
std::vector< TString > GetListOfVariables() const
returns list of variables
ClassInfo * GetClassInfo(Int_t clNum) const
Double_t GetTrainingSumSignalWeights()
void PrintClasses() const
Int_t GetClassNameMaxLength() const
Double_t GetTrainingSumBackgrWeights()
VariableInfo & AddTarget(const TString &expression, const TString &title, const TString &unit, Double_t min, Double_t max, Bool_t normalized=kTRUE, void *external=nullptr)
add a variable (can be a complex expression) to the set of variables used in the MV analysis
void PrintCorrelationMatrix(const TString &className)
calculates the correlation matrices for signal and background, prints them to standard output,...
void SetCut(const TCut &cut, const TString &className)
set the cut for the classes
Double_t GetTestingSumSignalWeights()
Int_t FindVarIndex(const TString &) const
find variable by name
Int_t GetVariableNameMaxLength() const
Bool_t IsSignal(const Event *ev) const
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=nullptr)
add a spectator (can be a complex expression) to the set of spectator variables used in the MV analys...
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...
void AddCut(const TCut &cut, const TString &className)
set the cut for the classes
void AddVariablesArray(const TString &expression, Int_t size, const TString &title="", const TString &unit="", Double_t min=0, Double_t max=0, char type='F', Bool_t normalized=kTRUE, void *external=nullptr)
add an array of variables identified by an expression corresponding to an array entry in the tree
std::vector< Float_t > * GetTargetsForMulticlass(const Event *ev)
void SetCorrelationMatrix(const TString &className, TMatrixD *matrix)
void ClearDataSet() const
Class that contains all the data information.
Definition DataSet.h:58
UInt_t GetClass() const
Definition Event.h:86
ostringstream derivative to redirect and format output
Definition MsgLogger.h:57
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:887
Class for type info of MVA input variable.
Mother of all ROOT objects.
Definition TObject.h:41
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:417
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
Tools & gTools()
MsgLogger & Endl(MsgLogger &ml)
Definition MsgLogger.h:148
TMarker m
Definition textangle.C:8