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