Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
DataLoader.cxx
Go to the documentation of this file.
1// @(#)root/tmva $Id$
2// Author: Omar Zapata
3// Mentors: Lorenzo Moneta, Sergei Gleyzer
4//NOTE: Based on TMVA::Factory
5
6/**********************************************************************************
7 * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
8 * Package: TMVA *
9 * Class : DataLoader *
10 * *
11 * *
12 * Description: *
13 * This is a class to load datasets into every booked method *
14 * *
15 * Authors (alphabetical): *
16 * Lorenzo Moneta <Lorenzo.Moneta@cern.ch> - CERN, Switzerland *
17 * Omar Zapata <Omar.Zapata@cern.ch> - ITM/UdeA, Colombia *
18 * Sergei Gleyzer<sergei.gleyzer@cern.ch> - CERN, Switzerland *
19 * *
20 * Copyright (c) 2005-2015: *
21 * CERN, Switzerland *
22 * ITM/UdeA, Colombia *
23 * *
24 * Redistribution and use in source and binary forms, with or without *
25 * modification, are permitted according to the terms listed in LICENSE *
26 * (see tmva/doc/LICENSE) *
27 **********************************************************************************/
28
29
30/*! \class TMVA::DataLoader
31\ingroup TMVA
32
33*/
34
35#include "TTree.h"
36#include "TH2.h"
37#include "TMatrixD.h"
38
39#include "TMVA/DataLoader.h"
40#include "TMVA/Config.h"
41#include "TMVA/CvSplit.h"
42#include "TMVA/Tools.h"
43#include "TMVA/IMethod.h"
44#include "TMVA/MethodBase.h"
46#include "TMVA/DataSetManager.h"
47#include "TMVA/DataSetInfo.h"
48#include "TMVA/MethodBoost.h"
49#include "TMVA/MethodCategory.h"
50
51#include "TMVA/VariableInfo.h"
58
60
61////////////////////////////////////////////////////////////////////////////////
62/*** Create a data loader
63 \param[in] thedlName name of DataLoader object. This name will be used as the
64 top directory name where the training results
65 (weights, i.e .XML and .C files) will be stored.
66 The results will be stored by default in the `theDlName/weights`
67 directory and relative to the current directory. If the directory is not existing,
68 a new one will be created automatically.
69 For using a different location (i.e. a different path to the current directory) one
70 can set an absolute path location in `TMVA::gConfig()::GetIONames().fWeightFileDirPrefix`
71 For example, by setting
72~~~~~~~~~~~~~~~{.cpp}
73 TMVA::gConfig()::GetIONames().fWeightFileDirPrefix = "/tmp";
74 TMVA::gConfig()::GetIONames().fWeightFileDir = "myTrainingResults";
75~~~~~~~~~~~~~~~
76 The training results will be stored in the `/tmp/thedlName/myTrainingResults`
77 directory.
78**/
79
81: Configurable( ),
82 fDataSetManager ( NULL ), //DSMTEST
83 fDataInputHandler ( new DataInputHandler ),
84 fTransformations ( "I" ),
85 fVerbose ( kFALSE ),
86 fDataAssignType ( kAssignEvents ),
87 fATreeEvent (0)
88{
90 SetName(thedlName.Data());
91 fLogger->SetSource("DataLoader");
92}
93
94////////////////////////////////////////////////////////////////////////////////
95
97{
98 // destructor
99
100 std::vector<TMVA::VariableTransformBase*>::iterator trfIt = fDefaultTrfs.begin();
101 for (;trfIt != fDefaultTrfs.end(); ++trfIt) delete (*trfIt);
102
103 delete fDataInputHandler;
104
105 // destroy singletons
106 // DataSetManager::DestroyInstance(); // DSMTEST replaced by following line
107 delete fDataSetManager; // DSMTEST
108
109 // problem with call of REGISTER_METHOD macro ...
110 // ClassifierDataLoader::DestroyInstance();
111 // Types::DestroyInstance();
112 //Tools::DestroyInstance();
113 //Config::DestroyInstance();
114}
115
116
117////////////////////////////////////////////////////////////////////////////////
118
120{
121 return fDataSetManager->AddDataSetInfo(dsi); // DSMTEST
122}
123
124////////////////////////////////////////////////////////////////////////////////
125
127{
128 DataSetInfo* dsi = fDataSetManager->GetDataSetInfo(dsiName); // DSMTEST
129
130 if (dsi!=0) return *dsi;
131
132 return fDataSetManager->AddDataSetInfo(*(new DataSetInfo(dsiName))); // DSMTEST
133}
134
135////////////////////////////////////////////////////////////////////////////////
136
138{
139 return DefaultDataSetInfo(); // DSMTEST
140}
141
142////////////////////////////////////////////////////////////////////////////////
143/// Transforms the variables and return a new DataLoader with the transformed
144/// variables
145
147{
148 TString trOptions = "0";
149 TString trName = "None";
150 if (trafoDefinition.Contains("(")) {
151
152 // contains transformation parameters
153 Ssiz_t parStart = trafoDefinition.Index( "(" );
155
158 trOptions.Remove(parLen-1,1);
159 trOptions.Remove(0,1);
160 }
161 else
163
164 VarTransformHandler* handler = new VarTransformHandler(this);
165 // variance threshold variable transformation
166 if (trName == "VT") {
167
168 // find threshold value from given input
169 Double_t threshold = 0.0;
170 if (!trOptions.IsFloat()){
171 Log() << kFATAL << " VT transformation must be passed a floating threshold value" << Endl;
172 delete handler;
173 return this;
174 }
175 else
176 threshold = trOptions.Atof();
178 delete handler;
179 return transformedLoader;
180 }
181 else {
182 delete handler;
183 Log() << kFATAL << "Incorrect transformation string provided, please check" << Endl;
184 }
185 Log() << kINFO << "No transformation applied, returning original loader" << Endl;
186 return this;
187}
188
189////////////////////////////////////////////////////////////////////////////////
190// the next functions are to assign events directly
191
192////////////////////////////////////////////////////////////////////////////////
193/// create the data assignment tree (for event-wise data assignment by user)
194
196{
197 TTree * assignTree = new TTree( name, name );
198 assignTree->SetDirectory(nullptr);
199 assignTree->Branch( "type", &fATreeType, "ATreeType/I" );
200 assignTree->Branch( "weight", &fATreeWeight, "ATreeWeight/F" );
201
202 std::vector<VariableInfo>& vars = DefaultDataSetInfo().GetVariableInfos();
203 std::vector<VariableInfo>& tgts = DefaultDataSetInfo().GetTargetInfos();
204 std::vector<VariableInfo>& spec = DefaultDataSetInfo().GetSpectatorInfos();
205
206 if (fATreeEvent.size()==0) fATreeEvent.resize(vars.size()+tgts.size()+spec.size());
207 // add variables
208 for (UInt_t ivar=0; ivar<vars.size(); ivar++) {
209 TString vname = vars[ivar].GetExpression();
210 assignTree->Branch( vname, &fATreeEvent[ivar], vname + "/F" );
211 }
212 // add targets
213 for (UInt_t itgt=0; itgt<tgts.size(); itgt++) {
214 TString vname = tgts[itgt].GetExpression();
215 assignTree->Branch( vname, &fATreeEvent[vars.size()+itgt], vname + "/F" );
216 }
217 // add spectators
218 for (UInt_t ispc=0; ispc<spec.size(); ispc++) {
219 TString vname = spec[ispc].GetExpression();
220 assignTree->Branch( vname, &fATreeEvent[vars.size()+tgts.size()+ispc], vname + "/F" );
221 }
222 return assignTree;
223}
224
225////////////////////////////////////////////////////////////////////////////////
226/// add signal training event
227
228void TMVA::DataLoader::AddSignalTrainingEvent( const std::vector<Double_t>& event, Double_t weight )
229{
230 AddEvent( "Signal", Types::kTraining, event, weight );
231}
232
233////////////////////////////////////////////////////////////////////////////////
234/// add signal testing event
235
236void TMVA::DataLoader::AddSignalTestEvent( const std::vector<Double_t>& event, Double_t weight )
237{
238 AddEvent( "Signal", Types::kTesting, event, weight );
239}
240
241////////////////////////////////////////////////////////////////////////////////
242/// add signal training event
243
244void TMVA::DataLoader::AddBackgroundTrainingEvent( const std::vector<Double_t>& event, Double_t weight )
245{
246 AddEvent( "Background", Types::kTraining, event, weight );
247}
248
249////////////////////////////////////////////////////////////////////////////////
250/// add signal training event
251
252void TMVA::DataLoader::AddBackgroundTestEvent( const std::vector<Double_t>& event, Double_t weight )
253{
254 AddEvent( "Background", Types::kTesting, event, weight );
255}
256
257////////////////////////////////////////////////////////////////////////////////
258/// add signal training event
259
260void TMVA::DataLoader::AddTrainingEvent( const TString& className, const std::vector<Double_t>& event, Double_t weight )
261{
262 AddEvent( className, Types::kTraining, event, weight );
263}
264
265////////////////////////////////////////////////////////////////////////////////
266/// add signal test event
267
268void TMVA::DataLoader::AddTestEvent( const TString& className, const std::vector<Double_t>& event, Double_t weight )
269{
270 AddEvent( className, Types::kTesting, event, weight );
271}
272
273////////////////////////////////////////////////////////////////////////////////
274/// add event
275/// vector event : the order of values is: variables + targets + spectators
276
278 const std::vector<Double_t>& event, Double_t weight )
279{
280 ClassInfo* theClass = DefaultDataSetInfo().AddClass(className); // returns class (creates it if necessary)
281 UInt_t clIndex = theClass->GetNumber();
282
283
284 // set analysistype to "kMulticlass" if more than two classes and analysistype == kNoAnalysisType
285 if( fAnalysisType == Types::kNoAnalysisType && DefaultDataSetInfo().GetNClasses() > 2 )
286 fAnalysisType = Types::kMulticlass;
287
288
289 if (clIndex>=fTrainAssignTree.size()) {
290 fTrainAssignTree.resize(clIndex+1, 0);
291 fTestAssignTree.resize(clIndex+1, 0);
292 }
293
294 if (fTrainAssignTree[clIndex]==0) { // does not exist yet
295 fTrainAssignTree[clIndex] = CreateEventAssignTrees( TString::Format("TrainAssignTree_%s", className.Data()).Data() );
296 fTestAssignTree[clIndex] = CreateEventAssignTrees( TString::Format("TestAssignTree_%s", className.Data()).Data() );
297 }
298
299 fATreeType = clIndex;
300 fATreeWeight = weight;
301 if (event.size() > fATreeEvent.size()) {
302 Error("AddEvent",
303 "Number of variables defined through DataLoader::AddVariable (%zu) is inconsistent"
304 " with number of variables given to DataLoader::Add*Event (%zu)."
305 " Please check your variable definitions and statement ordering."
306 " This event will not be added.", fATreeEvent.size(), event.size());
307 return;
308 }
309 for (UInt_t ivar=0; ivar<event.size(); ivar++) fATreeEvent[ivar] = event[ivar];
310
311 if(tt==Types::kTraining) fTrainAssignTree[clIndex]->Fill();
312 else fTestAssignTree[clIndex]->Fill();
313
314}
315
316////////////////////////////////////////////////////////////////////////////////
317///
318
320{
321 return fTrainAssignTree[clIndex]!=0;
322}
323
324////////////////////////////////////////////////////////////////////////////////
325/// assign event-wise local trees to data set
326
328{
329 UInt_t size = fTrainAssignTree.size();
330 for(UInt_t i=0; i<size; i++) {
331 if(!UserAssignEvents(i)) continue;
332 const TString& className = DefaultDataSetInfo().GetClassInfo(i)->GetName();
333 SetWeightExpression( "weight", className );
334 AddTree(fTrainAssignTree[i], className, 1.0, TCut(""), Types::kTraining );
335 AddTree(fTestAssignTree[i], className, 1.0, TCut(""), Types::kTesting );
336 }
337}
338
339////////////////////////////////////////////////////////////////////////////////
340/// number of signal events (used to compute significance)
341
342void TMVA::DataLoader::AddTree( TTree* tree, const TString& className, Double_t weight,
343 const TCut& cut, const TString& treetype )
344{
347 if (tmpTreeType.Contains( "train" ) && tmpTreeType.Contains( "test" )) tt = Types::kMaxTreeType;
348 else if (tmpTreeType.Contains( "train" )) tt = Types::kTraining;
349 else if (tmpTreeType.Contains( "test" )) tt = Types::kTesting;
350 else {
351 Log() << kFATAL << "<AddTree> cannot interpret tree type: \"" << treetype
352 << "\" should be \"Training\" or \"Test\" or \"Training and Testing\"" << Endl;
353 }
354 AddTree( tree, className, weight, cut, tt );
355}
356
357////////////////////////////////////////////////////////////////////////////////
358
359void TMVA::DataLoader::AddTree( TTree* tree, const TString& className, Double_t weight,
360 const TCut& cut, Types::ETreeType tt )
361{
362 if(!tree)
363 Log() << kFATAL << "Tree does not exist (empty pointer)." << Endl;
364
365 DefaultDataSetInfo().AddClass( className );
366
367 // set analysistype to "kMulticlass" if more than two classes and analysistype == kNoAnalysisType
368 if( fAnalysisType == Types::kNoAnalysisType && DefaultDataSetInfo().GetNClasses() > 2 )
369 fAnalysisType = Types::kMulticlass;
370
371 Log() << kINFO<< "Add Tree " << tree->GetName() << " of type " << className
372 << " with " << tree->GetEntries() << " events" << Endl;
373 DataInput().AddTree( tree, className, weight, cut, tt );
374}
375
376////////////////////////////////////////////////////////////////////////////////
377/// number of signal events (used to compute significance)
378
380{
381 AddTree( signal, "Signal", weight, TCut(""), treetype );
382}
383
384////////////////////////////////////////////////////////////////////////////////
385/// add signal tree from text file
386
388{
389 // create trees from these ascii files
390 TTree* signalTree = new TTree( "TreeS", "Tree (S)" );
391 signalTree->ReadFile( datFileS );
392
393 Log() << kINFO << "Create TTree objects from ASCII input files ... \n- Signal file : \""
394 << datFileS << Endl;
395
396 // number of signal events (used to compute significance)
397 AddTree( signalTree, "Signal", weight, TCut(""), treetype );
398}
399
400////////////////////////////////////////////////////////////////////////////////
401
403{
404 AddTree( signal, "Signal", weight, TCut(""), treetype );
405}
406
407////////////////////////////////////////////////////////////////////////////////
408/// number of signal events (used to compute significance)
409
411{
412 AddTree( signal, "Background", weight, TCut(""), treetype );
413}
414
415////////////////////////////////////////////////////////////////////////////////
416/// add background tree from text file
417
419{
420 // create trees from these ascii files
421 TTree* bkgTree = new TTree( "TreeB", "Tree (B)" );
422 bkgTree->ReadFile( datFileB );
423
424 Log() << kINFO << "Create TTree objects from ASCII input files ... \n- Background file : \""
425 << datFileB << Endl;
426
427 // number of signal events (used to compute significance)
428 AddTree( bkgTree, "Background", weight, TCut(""), treetype );
429}
430
431////////////////////////////////////////////////////////////////////////////////
432
434{
435 AddTree( signal, "Background", weight, TCut(""), treetype );
436}
437
438////////////////////////////////////////////////////////////////////////////////
439
441{
442 AddTree( tree, "Signal", weight );
443}
444
445////////////////////////////////////////////////////////////////////////////////
446
448{
449 AddTree( tree, "Background", weight );
450}
451
452////////////////////////////////////////////////////////////////////////////////
453/// set background tree
454
455void TMVA::DataLoader::SetTree( TTree* tree, const TString& className, Double_t weight )
456{
457 AddTree( tree, className, weight, TCut(""), Types::kMaxTreeType );
458}
459
460////////////////////////////////////////////////////////////////////////////////
461/// define the input trees for signal and background; no cuts are applied
462
469
470////////////////////////////////////////////////////////////////////////////////
471
474{
475 DataInput().AddTree( datFileS, "Signal", signalWeight );
476 DataInput().AddTree( datFileB, "Background", backgroundWeight );
477}
478
479////////////////////////////////////////////////////////////////////////////////
480/// define the input trees for signal and background from single input tree,
481/// containing both signal and background events distinguished by the type
482/// identifiers: SigCut and BgCut
483
485{
486 AddTree( inputTree, "Signal", 1.0, SigCut, Types::kMaxTreeType );
487 AddTree( inputTree, "Background", 1.0, BgCut , Types::kMaxTreeType );
488}
489
490////////////////////////////////////////////////////////////////////////////////
491/// user inserts discriminating variable in data set info
492
493void TMVA::DataLoader::AddVariable( const TString& expression, const TString& title, const TString& unit,
494 char type, Double_t min, Double_t max )
495{
496 DefaultDataSetInfo().AddVariable( expression, title, unit, min, max, type );
497}
498
499////////////////////////////////////////////////////////////////////////////////
500/// user inserts discriminating variable in data set info
501
502void TMVA::DataLoader::AddVariable( const TString& expression, char type,
503 Double_t min, Double_t max )
504{
505 DefaultDataSetInfo().AddVariable( expression, "", "", min, max, type );
506}
507
508////////////////////////////////////////////////////////////////////////////////
509/// user inserts discriminating array of variables in data set info
510/// in case input tree provides an array of values
511
512void TMVA::DataLoader::AddVariablesArray(const TString &expression, int size, char type,
513 Double_t min, Double_t max)
514{
515 DefaultDataSetInfo().AddVariablesArray(expression, size, "", "", min, max, type);
516}
517////////////////////////////////////////////////////////////////////////////////
518/// user inserts target in data set info
519
520void TMVA::DataLoader::AddTarget( const TString& expression, const TString& title, const TString& unit,
521 Double_t min, Double_t max )
522{
523 if( fAnalysisType == Types::kNoAnalysisType )
524 fAnalysisType = Types::kRegression;
525
526 DefaultDataSetInfo().AddTarget( expression, title, unit, min, max );
527}
528
529////////////////////////////////////////////////////////////////////////////////
530/// user inserts target in data set info
531
532void TMVA::DataLoader::AddSpectator( const TString& expression, const TString& title, const TString& unit,
533 Double_t min, Double_t max )
534{
535 DefaultDataSetInfo().AddSpectator( expression, title, unit, min, max );
536}
537
538////////////////////////////////////////////////////////////////////////////////
539/// default creation
540
542{
543 return AddDataSet( fName );
544}
545
546////////////////////////////////////////////////////////////////////////////////
547/// fill input variables in data set
548
550{
551 for (std::vector<TString>::iterator it=theVariables->begin();
552 it!=theVariables->end(); ++it) AddVariable(*it);
553}
554
555////////////////////////////////////////////////////////////////////////////////
556
558{
559 DefaultDataSetInfo().SetWeightExpression(variable, "Signal");
560}
561
562////////////////////////////////////////////////////////////////////////////////
563
565{
566 DefaultDataSetInfo().SetWeightExpression(variable, "Background");
567}
568
569////////////////////////////////////////////////////////////////////////////////
570
571void TMVA::DataLoader::SetWeightExpression( const TString& variable, const TString& className )
572{
573 //Log() << kWarning << DefaultDataSetInfo().GetNClasses() /*fClasses.size()*/ << Endl;
574 if (className=="") {
575 SetSignalWeightExpression(variable);
576 SetBackgroundWeightExpression(variable);
577 }
578 else DefaultDataSetInfo().SetWeightExpression( variable, className );
579}
580
581////////////////////////////////////////////////////////////////////////////////
582
583void TMVA::DataLoader::SetCut( const TString& cut, const TString& className ) {
584 SetCut( TCut(cut), className );
585}
586
587////////////////////////////////////////////////////////////////////////////////
588
589void TMVA::DataLoader::SetCut( const TCut& cut, const TString& className )
590{
591 DefaultDataSetInfo().SetCut( cut, className );
592}
593
594////////////////////////////////////////////////////////////////////////////////
595
596void TMVA::DataLoader::AddCut( const TString& cut, const TString& className )
597{
598 AddCut( TCut(cut), className );
599}
600
601////////////////////////////////////////////////////////////////////////////////
602void TMVA::DataLoader::AddCut( const TCut& cut, const TString& className )
603{
604 DefaultDataSetInfo().AddCut( cut, className );
605}
606
607////////////////////////////////////////////////////////////////////////////////
608/// prepare the training and test trees
609
612 const TString& otherOpt )
613{
614 SetInputTreesFromEventAssignTrees();
615
616 AddCut( cut );
617
618 DefaultDataSetInfo().SetSplitOptions( TString::Format("nTrain_Signal=%i:nTrain_Background=%i:nTest_Signal=%i:nTest_Background=%i:%s",
619 NsigTrain, NbkgTrain, NsigTest, NbkgTest, otherOpt.Data()).Data() );
620}
621
622////////////////////////////////////////////////////////////////////////////////
623/// prepare the training and test trees
624/// kept for backward compatibility
625
627{
628 SetInputTreesFromEventAssignTrees();
629
630 AddCut( cut );
631
632 DefaultDataSetInfo().SetSplitOptions( TString::Format("nTrain_Signal=%i:nTrain_Background=%i:nTest_Signal=%i:nTest_Background=%i:SplitMode=Random:EqualTrainSample:!V",
633 Ntrain, Ntrain, Ntest, Ntest).Data() );
634}
635
636////////////////////////////////////////////////////////////////////////////////
637/// prepare the training and test trees
638/// -> same cuts for signal and background
639
641{
642 SetInputTreesFromEventAssignTrees();
643
644 DefaultDataSetInfo().PrintClasses();
645 AddCut( cut );
646 DefaultDataSetInfo().SetSplitOptions( opt );
647}
648
649////////////////////////////////////////////////////////////////////////////////
650/// prepare the training and test trees
651
653{
654 // if event-wise data assignment, add local trees to dataset first
655 SetInputTreesFromEventAssignTrees();
656
657 //Log() << kINFO <<"Preparing trees for training and testing..."<< Endl;
658 AddCut( sigcut, "Signal" );
659 AddCut( bkgcut, "Background" );
660
661 DefaultDataSetInfo().SetSplitOptions( splitOpt );
662}
663
664////////////////////////////////////////////////////////////////////////////////
665/// Function required to split the training and testing datasets into a
666/// number of folds. Required by the CrossValidation and HyperParameterOptimisation
667/// classes. The option to split the training dataset into a training set and
668/// a validation set is implemented but not currently used.
669
671{
672 s.MakeKFoldDataSet( DefaultDataSetInfo() );
673}
674
675////////////////////////////////////////////////////////////////////////////////
676/// Function for assigning the correct folds to the testing or training set.
677
682
683
684////////////////////////////////////////////////////////////////////////////////
685/// Recombines the dataset. The precise semantics depend on the actual split.
686///
687/// Similar to the inverse operation of `MakeKFoldDataSet` but _will_ differ.
688/// See documentation for each particular split for more information.
689///
690
695
696////////////////////////////////////////////////////////////////////////////////
697/// Copy method use in VI and CV
698
705
706////////////////////////////////////////////////////////////////////////////////
707///Loading Dataset from DataInputHandler for subseed
708
710{
711 for( std::vector<TreeInfo>::const_iterator treeinfo=src->DataInput().Sbegin();treeinfo!=src->DataInput().Send();++treeinfo)
712 {
713 des->AddSignalTree( (*treeinfo).GetTree(), (*treeinfo).GetWeight(),(*treeinfo).GetTreeType());
714 }
715
716 for( std::vector<TreeInfo>::const_iterator treeinfo=src->DataInput().Bbegin();treeinfo!=src->DataInput().Bend();++treeinfo)
717 {
718 des->AddBackgroundTree( (*treeinfo).GetTree(), (*treeinfo).GetWeight(),(*treeinfo).GetTreeType());
719 }
720}
721
722////////////////////////////////////////////////////////////////////////////////
723/// returns the correlation matrix of datasets
724
726{
727 const TMatrixD * m = DefaultDataSetInfo().CorrelationMatrix(className);
728 return DefaultDataSetInfo().CreateCorrelationMatrixHist(m,
729 "CorrelationMatrix"+className, "Correlation Matrix ("+className+")");
730}
bool fVerbose
The verbosity flag.
Definition HLFactory.h:70
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
#define ClassImp(name)
Definition Rtypes.h:374
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 src
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
TDataSetManager * fDataSetManager
Definition TProofLite.h:23
const_iterator begin() const
const_iterator end() const
A specialized string object used for TTree selections.
Definition TCut.h:25
Service class for 2-D histogram classes.
Definition TH2.h:30
Class that contains all the information of a class.
Definition ClassInfo.h:49
MsgLogger * fLogger
! message logger
virtual void RecombineKFoldDataSet(DataSetInfo &dsi, Types::ETreeType tt=Types::kTraining)
Definition CvSplit.cxx:114
virtual void MakeKFoldDataSet(DataSetInfo &dsi)=0
virtual void PrepareFoldDataSet(DataSetInfo &dsi, UInt_t foldNumber, Types::ETreeType tt)
Set training and test set vectors of dataset described by dsi.
Definition CvSplit.cxx:57
Class that contains all the data information.
DataInputHandler * fDataInputHandler
->
Definition DataLoader.h:189
TTree * CreateEventAssignTrees(const TString &name)
create the data assignment tree (for event-wise data assignment by user)
void AddVariablesArray(const TString &expression, int size, char type='F', Double_t min=0, Double_t max=0)
user inserts discriminating array of variables in data set info in case input tree provides an array ...
void SetBackgroundTree(TTree *background, Double_t weight=1.0)
void AddSignalTree(TTree *signal, Double_t weight=1.0, Types::ETreeType treetype=Types::kMaxTreeType)
number of signal events (used to compute significance)
DataSetInfo & AddDataSet(DataSetInfo &)
void AddSpectator(const TString &expression, const TString &title="", const TString &unit="", Double_t min=0, Double_t max=0)
user inserts target in data set info
void SetInputTreesFromEventAssignTrees()
assign event-wise local trees to data set
void AddTrainingEvent(const TString &className, const std::vector< Double_t > &event, Double_t weight)
add signal training event
void SetTree(TTree *tree, const TString &className, Double_t weight)
deprecated
void AddSignalTestEvent(const std::vector< Double_t > &event, Double_t weight=1.0)
add signal testing event
DataSetInfo & DefaultDataSetInfo()
default creation
void AddBackgroundTestEvent(const std::vector< Double_t > &event, Double_t weight=1.0)
add signal training event
DataSetManager * fDataSetManager
Definition DataLoader.h:186
DataLoader * MakeCopy(TString name)
Copy method use in VI and CV.
void SetSignalWeightExpression(const TString &variable)
void MakeKFoldDataSet(CvSplit &s)
Function required to split the training and testing datasets into a number of folds.
void SetWeightExpression(const TString &variable, const TString &className="")
void AddBackgroundTrainingEvent(const std::vector< Double_t > &event, Double_t weight=1.0)
add signal training event
void RecombineKFoldDataSet(CvSplit &s, Types::ETreeType tt=Types::kTraining)
Recombines the dataset.
DataLoader * VarTransform(TString trafoDefinition)
Transforms the variables and return a new DataLoader with the transformed variables.
void SetBackgroundWeightExpression(const TString &variable)
void AddCut(const TString &cut, const TString &className="")
void AddEvent(const TString &className, Types::ETreeType tt, const std::vector< Double_t > &event, Double_t weight)
add event vector event : the order of values is: variables + targets + spectators
DataLoader(TString thedlName="default")
void PrepareTrainingAndTestTree(const TCut &cut, const TString &splitOpt)
prepare the training and test trees -> same cuts for signal and background
void AddBackgroundTree(TTree *background, Double_t weight=1.0, Types::ETreeType treetype=Types::kMaxTreeType)
number of signal events (used to compute significance)
DataSetInfo & GetDataSetInfo()
void AddTarget(const TString &expression, const TString &title="", const TString &unit="", Double_t min=0, Double_t max=0)
user inserts target in data set info
TH2 * GetCorrelationMatrix(const TString &className)
returns the correlation matrix of datasets
Bool_t UserAssignEvents(UInt_t clIndex)
void AddSignalTrainingEvent(const std::vector< Double_t > &event, Double_t weight=1.0)
add signal training event
void AddTestEvent(const TString &className, const std::vector< Double_t > &event, Double_t weight)
add signal test event
void SetSignalTree(TTree *signal, Double_t weight=1.0)
void SetInputTrees(const TString &signalFileName, const TString &backgroundFileName, Double_t signalWeight=1.0, Double_t backgroundWeight=1.0)
virtual ~DataLoader()
void AddTree(TTree *tree, const TString &className, Double_t weight=1.0, const TCut &cut="", Types::ETreeType tt=Types::kMaxTreeType)
void SetInputVariables(std::vector< TString > *theVariables)
deprecated
void SetCut(const TString &cut, const TString &className="")
void AddVariable(const TString &expression, const TString &title, const TString &unit, char type='F', Double_t min=0, Double_t max=0)
user inserts discriminating variable in data set info
void PrepareFoldDataSet(CvSplit &s, UInt_t foldNumber, Types::ETreeType tt=Types::kTraining)
Function for assigning the correct folds to the testing or training set.
Class that contains all the data information.
Definition DataSetInfo.h:62
Class that contains all the data information.
void SetSource(const std::string &source)
Definition MsgLogger.h:68
@ kMulticlass
Definition Types.h:129
@ kNoAnalysisType
Definition Types.h:130
@ kRegression
Definition Types.h:128
@ kMaxTreeType
also used as temporary storage for trees not yet assigned for testing;training...
Definition Types.h:145
@ kTraining
Definition Types.h:143
TMVA::DataLoader * VarianceThreshold(Double_t threshold)
Computes variance of all the variables and returns a new DataLoader with the selected variables whose...
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:150
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
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
A TTree represents a columnar dataset.
Definition TTree.h:84
void DataLoaderCopy(TMVA::DataLoader *des, TMVA::DataLoader *src)
MsgLogger & Endl(MsgLogger &ml)
Definition MsgLogger.h:148
TMarker m
Definition textangle.C:8
auto * tt
Definition textangle.C:16