Logo ROOT   6.14/05
Reference Guide
DataSetFactory.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Peter Speckmayer, Joerg Stelzer, Eckhard von Toerne, Helge Voss
3 
4 /*****************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : DataSetFactory *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementation (see header for description) *
12  * *
13  * Authors (alphabetical): *
14  * Andreas Hoecker <Andreas.Hocker@cern.ch> - CERN, Switzerland *
15  * Peter Speckmayer <Peter.Speckmayer@cern.ch> - CERN, Switzerland *
16  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - MSU, USA *
17  * Eckhard von Toerne <evt@physik.uni-bonn.de> - U. of Bonn, Germany *
18  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
19  * *
20  * Copyright (c) 2009: *
21  * CERN, Switzerland *
22  * MPI-K Heidelberg, Germany *
23  * U. of Bonn, Germany *
24  * Redistribution and use in source and binary forms, with or without *
25  * modification, are permitted according to the terms listed in LICENSE *
26  * (http://tmva.sourceforge.net/LICENSE) *
27  *****************************************************************************/
28 
29 /*! \class TMVA::DataSetFactory
30 \ingroup TMVA
31 
32 Class that contains all the data information
33 
34 */
35 
36 #include <assert.h>
37 
38 #include <map>
39 #include <vector>
40 #include <iomanip>
41 #include <iostream>
42 
43 #include <algorithm>
44 #include <functional>
45 #include <numeric>
46 #include <random>
47 
48 #include "TMVA/DataSetFactory.h"
49 
50 #include "TEventList.h"
51 #include "TFile.h"
52 #include "TH1.h"
53 #include "TH2.h"
54 #include "TProfile.h"
55 #include "TRandom3.h"
56 #include "TMatrixF.h"
57 #include "TVectorF.h"
58 #include "TMath.h"
59 #include "TROOT.h"
60 
61 #include "TMVA/MsgLogger.h"
62 #include "TMVA/Configurable.h"
66 #include "TMVA/DataSet.h"
67 #include "TMVA/DataSetInfo.h"
68 #include "TMVA/DataInputHandler.h"
69 #include "TMVA/Event.h"
70 
71 #include "TMVA/Tools.h"
72 #include "TMVA/Types.h"
73 #include "TMVA/VariableInfo.h"
74 
75 using namespace std;
76 
77 //TMVA::DataSetFactory* TMVA::DataSetFactory::fgInstance = 0;
78 
79 namespace TMVA {
80  // calculate the largest common divider
81  // this function is not happy if numbers are negative!
83  {
84  if (a<b) {Int_t tmp = a; a=b; b=tmp; } // achieve a>=b
85  if (b==0) return a;
86  Int_t fullFits = a/b;
87  return LargestCommonDivider(b,a-b*fullFits);
88  }
89 }
90 
91 
92 ////////////////////////////////////////////////////////////////////////////////
93 /// constructor
94 
96  fVerbose(kFALSE),
97  fVerboseLevel(TString("Info")),
98  fScaleWithPreselEff(0),
99  fCurrentTree(0),
100  fCurrentEvtIdx(0),
101  fInputFormulas(0),
102  fLogger( new MsgLogger("DataSetFactory", kINFO) )
103 {
104 }
105 
106 ////////////////////////////////////////////////////////////////////////////////
107 /// destructor
108 
110 {
111  std::vector<TTreeFormula*>::const_iterator formIt;
112 
113  for (formIt = fInputFormulas.begin() ; formIt!=fInputFormulas.end() ; ++formIt) if (*formIt) delete *formIt;
114  for (formIt = fTargetFormulas.begin() ; formIt!=fTargetFormulas.end() ; ++formIt) if (*formIt) delete *formIt;
115  for (formIt = fCutFormulas.begin() ; formIt!=fCutFormulas.end() ; ++formIt) if (*formIt) delete *formIt;
116  for (formIt = fWeightFormula.begin() ; formIt!=fWeightFormula.end() ; ++formIt) if (*formIt) delete *formIt;
117  for (formIt = fSpectatorFormulas.begin(); formIt!=fSpectatorFormulas.end(); ++formIt) if (*formIt) delete *formIt;
118 
119  delete fLogger;
120 }
121 
122 ////////////////////////////////////////////////////////////////////////////////
123 /// steering the creation of a new dataset
124 
126  TMVA::DataInputHandler& dataInput )
127 {
128  // build the first dataset from the data input
129  DataSet * ds = BuildInitialDataSet( dsi, dataInput );
130 
131  if (ds->GetNEvents() > 1) {
132  CalcMinMax(ds,dsi);
133 
134  // from the the final dataset build the correlation matrix
135  for (UInt_t cl = 0; cl< dsi.GetNClasses(); cl++) {
136  const TString className = dsi.GetClassInfo(cl)->GetName();
137  dsi.SetCorrelationMatrix( className, CalcCorrelationMatrix( ds, cl ) );
138  if (fCorrelations) {
139  dsi.PrintCorrelationMatrix(className);
140  }
141  }
142  //Log() << kHEADER << Endl;
143  Log() << kHEADER << Form("[%s] : ",dsi.GetName()) << " " << Endl << Endl;
144  }
145 
146  return ds;
147 }
148 
149 ////////////////////////////////////////////////////////////////////////////////
150 
152 {
153  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "Build DataSet consisting of one Event with dynamically changing variables" << Endl;
154  DataSet* ds = new DataSet(dsi);
155 
156  // create a DataSet with one Event which uses dynamic variables
157  // (pointers to variables)
158  if(dsi.GetNClasses()==0){
159  dsi.AddClass( "data" );
160  dsi.GetClassInfo( "data" )->SetNumber(0);
161  }
162 
163  std::vector<Float_t*>* evdyn = new std::vector<Float_t*>(0);
164 
165  std::vector<VariableInfo>& varinfos = dsi.GetVariableInfos();
166 
167  if (varinfos.empty())
168  Log() << kFATAL << Form("Dataset[%s] : ",dsi.GetName()) << "Dynamic data set cannot be built, since no variable informations are present. Apparently no variables have been set. This should not happen, please contact the TMVA authors." << Endl;
169 
170  std::vector<VariableInfo>::iterator it = varinfos.begin(), itEnd=varinfos.end();
171  for (;it!=itEnd;++it) {
172  Float_t* external=(Float_t*)(*it).GetExternalLink();
173  if (external==0)
174  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "The link to the external variable is NULL while I am trying to build a dynamic data set. In this case fTmpEvent from MethodBase HAS TO BE USED in the method to get useful values in variables." << Endl;
175  else evdyn->push_back (external);
176  }
177 
178  std::vector<VariableInfo>& spectatorinfos = dsi.GetSpectatorInfos();
179  it = spectatorinfos.begin();
180  for (;it!=spectatorinfos.end();++it) evdyn->push_back( (Float_t*)(*it).GetExternalLink() );
181 
182  TMVA::Event * ev = new Event((const std::vector<Float_t*>*&)evdyn, varinfos.size());
183  std::vector<Event*>* newEventVector = new std::vector<Event*>;
184  newEventVector->push_back(ev);
185 
186  ds->SetEventCollection(newEventVector, Types::kTraining);
187  ds->SetCurrentType( Types::kTraining );
188  ds->SetCurrentEvent( 0 );
189 
190  delete newEventVector;
191  return ds;
192 }
193 
194 ////////////////////////////////////////////////////////////////////////////////
195 /// if no entries, than create a DataSet with one Event which uses
196 /// dynamic variables (pointers to variables)
197 
200  DataInputHandler& dataInput )
201 {
202  if (dataInput.GetEntries()==0) return BuildDynamicDataSet( dsi );
203  // -------------------------------------------------------------------------
204 
205  // register the classes in the datasetinfo-object
206  // information comes from the trees in the dataInputHandler-object
207  std::vector< TString >* classList = dataInput.GetClassList();
208  for (std::vector<TString>::iterator it = classList->begin(); it< classList->end(); ++it) {
209  dsi.AddClass( (*it) );
210  }
211  delete classList;
212 
213  EvtStatsPerClass eventCounts(dsi.GetNClasses());
214  TString normMode;
215  TString splitMode;
216  TString mixMode;
217  UInt_t splitSeed;
218 
219  InitOptions( dsi, eventCounts, normMode, splitSeed, splitMode , mixMode );
220  // ======= build event-vector from input, apply preselection ===============
221  EventVectorOfClassesOfTreeType tmpEventVector;
222  BuildEventVector( dsi, dataInput, tmpEventVector, eventCounts );
223 
224  DataSet* ds = MixEvents( dsi, tmpEventVector, eventCounts,
225  splitMode, mixMode, normMode, splitSeed );
226 
227  const Bool_t showCollectedOutput = kFALSE;
228  if (showCollectedOutput) {
229  Int_t maxL = dsi.GetClassNameMaxLength();
230  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "Collected:" << Endl;
231  for (UInt_t cl = 0; cl < dsi.GetNClasses(); cl++) {
232  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " "
233  << setiosflags(ios::left) << std::setw(maxL) << dsi.GetClassInfo(cl)->GetName()
234  << " training entries: " << ds->GetNClassEvents( 0, cl ) << Endl;
235  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " "
236  << setiosflags(ios::left) << std::setw(maxL) << dsi.GetClassInfo(cl)->GetName()
237  << " testing entries: " << ds->GetNClassEvents( 1, cl ) << Endl;
238  }
239  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " " << Endl;
240  }
241 
242  return ds;
243 }
244 
245 ////////////////////////////////////////////////////////////////////////////////
246 /// checks a TTreeFormula for problems
247 
249  const TString& expression,
250  Bool_t& hasDollar )
251 {
252  Bool_t worked = kTRUE;
253 
254  if( ttf->GetNdim() <= 0 )
255  Log() << kFATAL << "Expression " << expression.Data()
256  << " could not be resolved to a valid formula. " << Endl;
257  if( ttf->GetNdata() == 0 ){
258  Log() << kWARNING << "Expression: " << expression.Data()
259  << " does not provide data for this event. "
260  << "This event is not taken into account. --> please check if you use as a variable "
261  << "an entry of an array which is not filled for some events "
262  << "(e.g. arr[4] when arr has only 3 elements)." << Endl;
263  Log() << kWARNING << "If you want to take the event into account you can do something like: "
264  << "\"Alt$(arr[4],0)\" where in cases where arr doesn't have a 4th element, "
265  << " 0 is taken as an alternative." << Endl;
266  worked = kFALSE;
267  }
268  if( expression.Contains("$") )
269  hasDollar = kTRUE;
270  else
271  {
272  for (int i = 0, iEnd = ttf->GetNcodes (); i < iEnd; ++i)
273  {
274  TLeaf* leaf = ttf->GetLeaf (i);
275  if (!leaf->IsOnTerminalBranch())
276  hasDollar = kTRUE;
277  }
278  }
279  return worked;
280 }
281 
282 
283 ////////////////////////////////////////////////////////////////////////////////
284 /// While the data gets copied into the local training and testing
285 /// trees, the input tree can change (for instance when changing from
286 /// signal to background tree, or using TChains as input) The
287 /// TTreeFormulas, that hold the input expressions need to be
288 /// re-associated with the new tree, which is done here
289 
291 {
292  TTree *tr = tinfo.GetTree()->GetTree();
293 
294  tr->SetBranchStatus("*",1);
295  tr->ResetBranchAddresses();
296 
297  Bool_t hasDollar = kFALSE;
298 
299  // 1) the input variable formulas
300  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "transform input variables" << Endl;
301  std::vector<TTreeFormula*>::const_iterator formIt, formItEnd;
302  for (formIt = fInputFormulas.begin(), formItEnd=fInputFormulas.end(); formIt!=formItEnd; ++formIt) if (*formIt) delete *formIt;
303  fInputFormulas.clear();
304  TTreeFormula* ttf = 0;
305 
306  for (UInt_t i=0; i<dsi.GetNVariables(); i++) {
307  ttf = new TTreeFormula( Form( "Formula%s", dsi.GetVariableInfo(i).GetInternalName().Data() ),
308  dsi.GetVariableInfo(i).GetExpression().Data(), tr );
309  CheckTTreeFormula( ttf, dsi.GetVariableInfo(i).GetExpression(), hasDollar );
310  fInputFormulas.push_back( ttf );
311  }
312 
313  //
314  // targets
315  //
316  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "transform regression targets" << Endl;
317  for (formIt = fTargetFormulas.begin(), formItEnd = fTargetFormulas.end(); formIt!=formItEnd; ++formIt) if (*formIt) delete *formIt;
318  fTargetFormulas.clear();
319  for (UInt_t i=0; i<dsi.GetNTargets(); i++) {
320  ttf = new TTreeFormula( Form( "Formula%s", dsi.GetTargetInfo(i).GetInternalName().Data() ),
321  dsi.GetTargetInfo(i).GetExpression().Data(), tr );
322  CheckTTreeFormula( ttf, dsi.GetTargetInfo(i).GetExpression(), hasDollar );
323  fTargetFormulas.push_back( ttf );
324  }
325 
326  //
327  // spectators
328  //
329  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "transform spectator variables" << Endl;
330  for (formIt = fSpectatorFormulas.begin(), formItEnd = fSpectatorFormulas.end(); formIt!=formItEnd; ++formIt) if (*formIt) delete *formIt;
331  fSpectatorFormulas.clear();
332  for (UInt_t i=0; i<dsi.GetNSpectators(); i++) {
333  ttf = new TTreeFormula( Form( "Formula%s", dsi.GetSpectatorInfo(i).GetInternalName().Data() ),
334  dsi.GetSpectatorInfo(i).GetExpression().Data(), tr );
335  CheckTTreeFormula( ttf, dsi.GetSpectatorInfo(i).GetExpression(), hasDollar );
336  fSpectatorFormulas.push_back( ttf );
337  }
338 
339  //
340  // the cuts (one per class, if non-existent: formula pointer = 0)
341  //
342  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "transform cuts" << Endl;
343  for (formIt = fCutFormulas.begin(), formItEnd = fCutFormulas.end(); formIt!=formItEnd; ++formIt) if (*formIt) delete *formIt;
344  fCutFormulas.clear();
345  for (UInt_t clIdx=0; clIdx<dsi.GetNClasses(); clIdx++) {
346  const TCut& tmpCut = dsi.GetClassInfo(clIdx)->GetCut();
347  const TString tmpCutExp(tmpCut.GetTitle());
348  ttf = 0;
349  if (tmpCutExp!="") {
350  ttf = new TTreeFormula( Form("CutClass%i",clIdx), tmpCutExp, tr );
351  Bool_t worked = CheckTTreeFormula( ttf, tmpCutExp, hasDollar );
352  if( !worked ){
353  Log() << kWARNING << "Please check class \"" << dsi.GetClassInfo(clIdx)->GetName()
354  << "\" cut \"" << dsi.GetClassInfo(clIdx)->GetCut() << Endl;
355  }
356  }
357  fCutFormulas.push_back( ttf );
358  }
359 
360  //
361  // the weights (one per class, if non-existent: formula pointer = 0)
362  //
363  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "transform weights" << Endl;
364  for (formIt = fWeightFormula.begin(), formItEnd = fWeightFormula.end(); formIt!=formItEnd; ++formIt) if (*formIt) delete *formIt;
365  fWeightFormula.clear();
366  for (UInt_t clIdx=0; clIdx<dsi.GetNClasses(); clIdx++) {
367  const TString tmpWeight = dsi.GetClassInfo(clIdx)->GetWeight();
368 
369  if (dsi.GetClassInfo(clIdx)->GetName() != tinfo.GetClassName() ) { // if the tree is of another class
370  fWeightFormula.push_back( 0 );
371  continue;
372  }
373 
374  ttf = 0;
375  if (tmpWeight!="") {
376  ttf = new TTreeFormula( "FormulaWeight", tmpWeight, tr );
377  Bool_t worked = CheckTTreeFormula( ttf, tmpWeight, hasDollar );
378  if( !worked ){
379  Log() << kWARNING << Form("Dataset[%s] : ",dsi.GetName()) << "Please check class \"" << dsi.GetClassInfo(clIdx)->GetName()
380  << "\" weight \"" << dsi.GetClassInfo(clIdx)->GetWeight() << Endl;
381  }
382  }
383  else {
384  ttf = 0;
385  }
386  fWeightFormula.push_back( ttf );
387  }
388  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "enable branches" << Endl;
389  // now enable only branches that are needed in any input formula, target, cut, weight
390 
391  if (!hasDollar) {
392  tr->SetBranchStatus("*",0);
393  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "enable branches: input variables" << Endl;
394  // input vars
395  for (formIt = fInputFormulas.begin(); formIt!=fInputFormulas.end(); ++formIt) {
396  ttf = *formIt;
397  for (Int_t bi = 0; bi<ttf->GetNcodes(); bi++) {
398  tr->SetBranchStatus( ttf->GetLeaf(bi)->GetBranch()->GetName(), 1 );
399  }
400  }
401  // targets
402  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "enable branches: targets" << Endl;
403  for (formIt = fTargetFormulas.begin(); formIt!=fTargetFormulas.end(); ++formIt) {
404  ttf = *formIt;
405  for (Int_t bi = 0; bi<ttf->GetNcodes(); bi++)
406  tr->SetBranchStatus( ttf->GetLeaf(bi)->GetBranch()->GetName(), 1 );
407  }
408  // spectators
409  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "enable branches: spectators" << Endl;
410  for (formIt = fSpectatorFormulas.begin(); formIt!=fSpectatorFormulas.end(); ++formIt) {
411  ttf = *formIt;
412  for (Int_t bi = 0; bi<ttf->GetNcodes(); bi++)
413  tr->SetBranchStatus( ttf->GetLeaf(bi)->GetBranch()->GetName(), 1 );
414  }
415  // cuts
416  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "enable branches: cuts" << Endl;
417  for (formIt = fCutFormulas.begin(); formIt!=fCutFormulas.end(); ++formIt) {
418  ttf = *formIt;
419  if (!ttf) continue;
420  for (Int_t bi = 0; bi<ttf->GetNcodes(); bi++)
421  tr->SetBranchStatus( ttf->GetLeaf(bi)->GetBranch()->GetName(), 1 );
422  }
423  // weights
424  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "enable branches: weights" << Endl;
425  for (formIt = fWeightFormula.begin(); formIt!=fWeightFormula.end(); ++formIt) {
426  ttf = *formIt;
427  if (!ttf) continue;
428  for (Int_t bi = 0; bi<ttf->GetNcodes(); bi++)
429  tr->SetBranchStatus( ttf->GetLeaf(bi)->GetBranch()->GetName(), 1 );
430  }
431  }
432  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName()) << "tree initialized" << Endl;
433  return;
434 }
435 
436 ////////////////////////////////////////////////////////////////////////////////
437 /// compute covariance matrix
438 
440 {
441  const UInt_t nvar = ds->GetNVariables();
442  const UInt_t ntgts = ds->GetNTargets();
443  const UInt_t nvis = ds->GetNSpectators();
444 
445  Float_t *min = new Float_t[nvar];
446  Float_t *max = new Float_t[nvar];
447  Float_t *tgmin = new Float_t[ntgts];
448  Float_t *tgmax = new Float_t[ntgts];
449  Float_t *vmin = new Float_t[nvis];
450  Float_t *vmax = new Float_t[nvis];
451 
452  for (UInt_t ivar=0; ivar<nvar ; ivar++) { min[ivar] = FLT_MAX; max[ivar] = -FLT_MAX; }
453  for (UInt_t ivar=0; ivar<ntgts; ivar++) { tgmin[ivar] = FLT_MAX; tgmax[ivar] = -FLT_MAX; }
454  for (UInt_t ivar=0; ivar<nvis; ivar++) { vmin[ivar] = FLT_MAX; vmax[ivar] = -FLT_MAX; }
455 
456  // perform event loop
457 
458  for (Int_t i=0; i<ds->GetNEvents(); i++) {
459  const Event * ev = ds->GetEvent(i);
460  for (UInt_t ivar=0; ivar<nvar; ivar++) {
461  Double_t v = ev->GetValue(ivar);
462  if (v<min[ivar]) min[ivar] = v;
463  if (v>max[ivar]) max[ivar] = v;
464  }
465  for (UInt_t itgt=0; itgt<ntgts; itgt++) {
466  Double_t v = ev->GetTarget(itgt);
467  if (v<tgmin[itgt]) tgmin[itgt] = v;
468  if (v>tgmax[itgt]) tgmax[itgt] = v;
469  }
470  for (UInt_t ivis=0; ivis<nvis; ivis++) {
471  Double_t v = ev->GetSpectator(ivis);
472  if (v<vmin[ivis]) vmin[ivis] = v;
473  if (v>vmax[ivis]) vmax[ivis] = v;
474  }
475  }
476 
477  for (UInt_t ivar=0; ivar<nvar; ivar++) {
478  dsi.GetVariableInfo(ivar).SetMin(min[ivar]);
479  dsi.GetVariableInfo(ivar).SetMax(max[ivar]);
480  if( TMath::Abs(max[ivar]-min[ivar]) <= FLT_MIN )
481  Log() << kWARNING << Form("Dataset[%s] : ",dsi.GetName()) << "Variable " << dsi.GetVariableInfo(ivar).GetExpression().Data() << " is constant. Please remove the variable." << Endl;
482  }
483  for (UInt_t ivar=0; ivar<ntgts; ivar++) {
484  dsi.GetTargetInfo(ivar).SetMin(tgmin[ivar]);
485  dsi.GetTargetInfo(ivar).SetMax(tgmax[ivar]);
486  if( TMath::Abs(tgmax[ivar]-tgmin[ivar]) <= FLT_MIN )
487  Log() << kFATAL << Form("Dataset[%s] : ",dsi.GetName()) << "Target " << dsi.GetTargetInfo(ivar).GetExpression().Data() << " is constant. Please remove the variable." << Endl;
488  }
489  for (UInt_t ivar=0; ivar<nvis; ivar++) {
490  dsi.GetSpectatorInfo(ivar).SetMin(vmin[ivar]);
491  dsi.GetSpectatorInfo(ivar).SetMax(vmax[ivar]);
492  // if( TMath::Abs(vmax[ivar]-vmin[ivar]) <= FLT_MIN )
493  // Log() << kWARNING << "Spectator variable " << dsi.GetSpectatorInfo(ivar).GetExpression().Data() << " is constant." << Endl;
494  }
495  delete [] min;
496  delete [] max;
497  delete [] tgmin;
498  delete [] tgmax;
499  delete [] vmin;
500  delete [] vmax;
501 }
502 
503 ////////////////////////////////////////////////////////////////////////////////
504 /// computes correlation matrix for variables "theVars" in tree;
505 /// "theType" defines the required event "type"
506 /// ("type" variable must be present in tree)
507 
509 {
510  // first compute variance-covariance
511  TMatrixD* mat = CalcCovarianceMatrix( ds, classNumber );
512 
513  // now the correlation
514  UInt_t nvar = ds->GetNVariables(), ivar, jvar;
515 
516  for (ivar=0; ivar<nvar; ivar++) {
517  for (jvar=0; jvar<nvar; jvar++) {
518  if (ivar != jvar) {
519  Double_t d = (*mat)(ivar, ivar)*(*mat)(jvar, jvar);
520  if (d > 0) (*mat)(ivar, jvar) /= sqrt(d);
521  else {
522  Log() << kWARNING << Form("Dataset[%s] : ",DataSetInfo().GetName())<< "<GetCorrelationMatrix> Zero variances for variables "
523  << "(" << ivar << ", " << jvar << ") = " << d
524  << Endl;
525  (*mat)(ivar, jvar) = 0;
526  }
527  }
528  }
529  }
530 
531  for (ivar=0; ivar<nvar; ivar++) (*mat)(ivar, ivar) = 1.0;
532 
533  return mat;
534 }
535 
536 ////////////////////////////////////////////////////////////////////////////////
537 /// compute covariance matrix
538 
540 {
541  UInt_t nvar = ds->GetNVariables();
542  UInt_t ivar = 0, jvar = 0;
543 
544  TMatrixD* mat = new TMatrixD( nvar, nvar );
545 
546  // init matrices
547  TVectorD vec(nvar);
548  TMatrixD mat2(nvar, nvar);
549  for (ivar=0; ivar<nvar; ivar++) {
550  vec(ivar) = 0;
551  for (jvar=0; jvar<nvar; jvar++) mat2(ivar, jvar) = 0;
552  }
553 
554  // perform event loop
555  Double_t ic = 0;
556  for (Int_t i=0; i<ds->GetNEvents(); i++) {
557 
558  const Event * ev = ds->GetEvent(i);
559  if (ev->GetClass() != classNumber ) continue;
560 
561  Double_t weight = ev->GetWeight();
562  ic += weight; // count used events
563 
564  for (ivar=0; ivar<nvar; ivar++) {
565 
566  Double_t xi = ev->GetValue(ivar);
567  vec(ivar) += xi*weight;
568  mat2(ivar, ivar) += (xi*xi*weight);
569 
570  for (jvar=ivar+1; jvar<nvar; jvar++) {
571  Double_t xj = ev->GetValue(jvar);
572  mat2(ivar, jvar) += (xi*xj*weight);
573  }
574  }
575  }
576 
577  for (ivar=0; ivar<nvar; ivar++)
578  for (jvar=ivar+1; jvar<nvar; jvar++)
579  mat2(jvar, ivar) = mat2(ivar, jvar); // symmetric matrix
580 
581 
582  // variance-covariance
583  for (ivar=0; ivar<nvar; ivar++) {
584  for (jvar=0; jvar<nvar; jvar++) {
585  (*mat)(ivar, jvar) = mat2(ivar, jvar)/ic - vec(ivar)*vec(jvar)/(ic*ic);
586  }
587  }
588 
589  return mat;
590 }
591 
592 // --------------------------------------- new versions
593 
594 ////////////////////////////////////////////////////////////////////////////////
595 /// the dataset splitting
596 
597 void
599  EvtStatsPerClass& nEventRequests,
600  TString& normMode,
601  UInt_t& splitSeed,
602  TString& splitMode,
603  TString& mixMode)
604 {
605  Configurable splitSpecs( dsi.GetSplitOptions() );
606  splitSpecs.SetConfigName("DataSetFactory");
607  splitSpecs.SetConfigDescription( "Configuration options given in the \"PrepareForTrainingAndTesting\" call; these options define the creation of the data sets used for training and expert validation by TMVA" );
608 
609  splitMode = "Random"; // the splitting mode
610  splitSpecs.DeclareOptionRef( splitMode, "SplitMode",
611  "Method of picking training and testing events (default: random)" );
612  splitSpecs.AddPreDefVal(TString("Random"));
613  splitSpecs.AddPreDefVal(TString("Alternate"));
614  splitSpecs.AddPreDefVal(TString("Block"));
615 
616  mixMode = "SameAsSplitMode"; // the splitting mode
617  splitSpecs.DeclareOptionRef( mixMode, "MixMode",
618  "Method of mixing events of different classes into one dataset (default: SameAsSplitMode)" );
619  splitSpecs.AddPreDefVal(TString("SameAsSplitMode"));
620  splitSpecs.AddPreDefVal(TString("Random"));
621  splitSpecs.AddPreDefVal(TString("Alternate"));
622  splitSpecs.AddPreDefVal(TString("Block"));
623 
624  splitSeed = 100;
625  splitSpecs.DeclareOptionRef( splitSeed, "SplitSeed",
626  "Seed for random event shuffling" );
627 
628  normMode = "EqualNumEvents"; // the weight normalisation modes
629  splitSpecs.DeclareOptionRef( normMode, "NormMode",
630  "Overall renormalisation of event-by-event weights used in the training (NumEvents: average weight of 1 per event, independently for signal and background; EqualNumEvents: average weight of 1 per event for signal, and sum of weights for background equal to sum of weights for signal)" );
631  splitSpecs.AddPreDefVal(TString("None"));
632  splitSpecs.AddPreDefVal(TString("NumEvents"));
633  splitSpecs.AddPreDefVal(TString("EqualNumEvents"));
634 
635  splitSpecs.DeclareOptionRef(fScaleWithPreselEff=kFALSE,"ScaleWithPreselEff","Scale the number of requested events by the eff. of the preselection cuts (or not)" );
636 
637  // the number of events
638 
639  // fill in the numbers
640  for (UInt_t cl = 0; cl < dsi.GetNClasses(); cl++) {
641  TString clName = dsi.GetClassInfo(cl)->GetName();
642  TString titleTrain = TString().Format("Number of training events of class %s (default: 0 = all)",clName.Data()).Data();
643  TString titleTest = TString().Format("Number of test events of class %s (default: 0 = all)",clName.Data()).Data();
644  TString titleSplit = TString().Format("Split in training and test events of class %s (default: 0 = deactivated)",clName.Data()).Data();
645 
646  splitSpecs.DeclareOptionRef( nEventRequests.at(cl).nTrainingEventsRequested, TString("nTrain_")+clName, titleTrain );
647  splitSpecs.DeclareOptionRef( nEventRequests.at(cl).nTestingEventsRequested , TString("nTest_")+clName , titleTest );
648  splitSpecs.DeclareOptionRef( nEventRequests.at(cl).TrainTestSplitRequested , TString("TrainTestSplit_")+clName , titleTest );
649  }
650 
651  splitSpecs.DeclareOptionRef( fVerbose, "V", "Verbosity (default: true)" );
652 
653  splitSpecs.DeclareOptionRef( fVerboseLevel=TString("Info"), "VerboseLevel", "VerboseLevel (Debug/Verbose/Info)" );
654  splitSpecs.AddPreDefVal(TString("Debug"));
655  splitSpecs.AddPreDefVal(TString("Verbose"));
656  splitSpecs.AddPreDefVal(TString("Info"));
657 
659  splitSpecs.DeclareOptionRef(fCorrelations, "Correlations", "Boolean to show correlation output (Default: true)");
660 
661  splitSpecs.ParseOptions();
662  splitSpecs.CheckForUnusedOptions();
663 
664  // output logging verbosity
665  if (Verbose()) fLogger->SetMinType( kVERBOSE );
666  if (fVerboseLevel.CompareTo("Debug") ==0) fLogger->SetMinType( kDEBUG );
667  if (fVerboseLevel.CompareTo("Verbose") ==0) fLogger->SetMinType( kVERBOSE );
668  if (fVerboseLevel.CompareTo("Info") ==0) fLogger->SetMinType( kINFO );
669 
670  // put all to upper case
671  splitMode.ToUpper(); mixMode.ToUpper(); normMode.ToUpper();
672  // adjust mixmode if same as splitmode option has been set
673  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
674  << "\tSplitmode is: \"" << splitMode << "\" the mixmode is: \"" << mixMode << "\"" << Endl;
675  if (mixMode=="SAMEASSPLITMODE") mixMode = splitMode;
676  else if (mixMode!=splitMode)
677  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "DataSet splitmode="<<splitMode
678  <<" differs from mixmode="<<mixMode<<Endl;
679 }
680 
681 ////////////////////////////////////////////////////////////////////////////////
682 /// build empty event vectors
683 /// distributes events between kTraining/kTesting/kMaxTreeType
684 
685 void
687  TMVA::DataInputHandler& dataInput,
689  EvtStatsPerClass& eventCounts)
690 {
691  const UInt_t nclasses = dsi.GetNClasses();
692 
693  eventsmap[ Types::kTraining ] = EventVectorOfClasses(nclasses);
694  eventsmap[ Types::kTesting ] = EventVectorOfClasses(nclasses);
695  eventsmap[ Types::kMaxTreeType ] = EventVectorOfClasses(nclasses);
696 
697  // create the type, weight and boostweight branches
698  const UInt_t nvars = dsi.GetNVariables();
699  const UInt_t ntgts = dsi.GetNTargets();
700  const UInt_t nvis = dsi.GetNSpectators();
701 
702  for (size_t i=0; i<nclasses; i++) {
703  eventCounts[i].varAvLength = new Float_t[nvars];
704  for (UInt_t ivar=0; ivar<nvars; ivar++)
705  eventCounts[i].varAvLength[ivar] = 0;
706  }
707 
708  // Bool_t haveArrayVariable = kFALSE;
709  Bool_t *varIsArray = new Bool_t[nvars];
710 
711  // If there are NaNs in the tree:
712  // => warn if used variables/cuts/weights contain nan (no problem if event is cut out)
713  // => fatal if cut value is nan or (event not cut out and nans somewhere)
714  // Count & collect all these warnings/errors and output them at the end.
715  std::map<TString, int> nanInfWarnings;
716  std::map<TString, int> nanInfErrors;
717 
718  // if we work with chains we need to remember the current tree if
719  // the chain jumps to a new tree we have to reset the formulas
720  for (UInt_t cl=0; cl<nclasses; cl++) {
721 
722  //Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "Create training and testing trees -- looping over class \"" << dsi.GetClassInfo(cl)->GetName() << "\" ..." << Endl;
723 
724  EventStats& classEventCounts = eventCounts[cl];
725 
726  // info output for weights
727  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
728  << "\tWeight expression for class \'" << dsi.GetClassInfo(cl)->GetName() << "\': \""
729  << dsi.GetClassInfo(cl)->GetWeight() << "\"" << Endl;
730 
731  // used for chains only
732  TString currentFileName("");
733 
734  std::vector<TreeInfo>::const_iterator treeIt(dataInput.begin(dsi.GetClassInfo(cl)->GetName()));
735  for (;treeIt!=dataInput.end(dsi.GetClassInfo(cl)->GetName()); ++treeIt) {
736 
737  // read first the variables
738  std::vector<Float_t> vars(nvars);
739  std::vector<Float_t> tgts(ntgts);
740  std::vector<Float_t> vis(nvis);
741  TreeInfo currentInfo = *treeIt;
742 
743  Log() << kDEBUG << "Building event vectors " << currentInfo.GetTreeType() << Endl;
744 
745  EventVector& event_v = eventsmap[currentInfo.GetTreeType()].at(cl);
746 
747  Bool_t isChain = (TString("TChain") == currentInfo.GetTree()->ClassName());
748  currentInfo.GetTree()->LoadTree(0);
749  ChangeToNewTree( currentInfo, dsi );
750 
751  // count number of events in tree before cut
752  classEventCounts.nInitialEvents += currentInfo.GetTree()->GetEntries();
753 
754  // loop over events in ntuple
755  const UInt_t nEvts = currentInfo.GetTree()->GetEntries();
756  for (Long64_t evtIdx = 0; evtIdx < nEvts; evtIdx++) {
757  currentInfo.GetTree()->LoadTree(evtIdx);
758 
759  // may need to reload tree in case of chains
760  if (isChain) {
761  if (currentInfo.GetTree()->GetTree()->GetDirectory()->GetFile()->GetName() != currentFileName) {
762  currentFileName = currentInfo.GetTree()->GetTree()->GetDirectory()->GetFile()->GetName();
763  ChangeToNewTree( currentInfo, dsi );
764  }
765  }
766  currentInfo.GetTree()->GetEntry(evtIdx);
767  Int_t sizeOfArrays = 1;
768  Int_t prevArrExpr = 0;
769 
770  // ======= evaluate all formulas =================
771 
772  // first we check if some of the formulas are arrays
773  for (UInt_t ivar=0; ivar<nvars; ivar++) {
774  Int_t ndata = fInputFormulas[ivar]->GetNdata();
775  classEventCounts.varAvLength[ivar] += ndata;
776  if (ndata == 1) continue;
777  // haveArrayVariable = kTRUE;
778  varIsArray[ivar] = kTRUE;
779  if (sizeOfArrays == 1) {
780  sizeOfArrays = ndata;
781  prevArrExpr = ivar;
782  }
783  else if (sizeOfArrays!=ndata) {
784  Log() << kERROR << Form("Dataset[%s] : ",dsi.GetName())<< "ERROR while preparing training and testing trees:" << Endl;
785  Log() << Form("Dataset[%s] : ",dsi.GetName())<< " multiple array-type expressions of different length were encountered" << Endl;
786  Log() << Form("Dataset[%s] : ",dsi.GetName())<< " location of error: event " << evtIdx
787  << " in tree " << currentInfo.GetTree()->GetName()
788  << " of file " << currentInfo.GetTree()->GetCurrentFile()->GetName() << Endl;
789  Log() << Form("Dataset[%s] : ",dsi.GetName())<< " expression " << fInputFormulas[ivar]->GetTitle() << " has "
790  << Form("Dataset[%s] : ",dsi.GetName()) << ndata << " entries, while" << Endl;
791  Log() << Form("Dataset[%s] : ",dsi.GetName())<< " expression " << fInputFormulas[prevArrExpr]->GetTitle() << " has "
792  << Form("Dataset[%s] : ",dsi.GetName())<< fInputFormulas[prevArrExpr]->GetNdata() << " entries" << Endl;
793  Log() << kFATAL << Form("Dataset[%s] : ",dsi.GetName())<< "Need to abort" << Endl;
794  }
795  }
796 
797  // now we read the information
798  for (Int_t idata = 0; idata<sizeOfArrays; idata++) {
799  Bool_t contains_NaN_or_inf = kFALSE;
800 
801  auto checkNanInf = [&](std::map<TString, int> &msgMap, Float_t value, const char *what, const char *formulaTitle) {
802  if (TMath::IsNaN(value)) {
803  contains_NaN_or_inf = kTRUE;
804  ++msgMap[TString::Format("Dataset[%s] : %s expression resolves to indeterminate value (NaN): %s", dsi.GetName(), what, formulaTitle)];
805  } else if (!TMath::Finite(value)) {
806  contains_NaN_or_inf = kTRUE;
807  ++msgMap[TString::Format("Dataset[%s] : %s expression resolves to infinite value (+inf or -inf): %s", dsi.GetName(), what, formulaTitle)];
808  }
809  };
810 
811  TTreeFormula* formula = 0;
812 
813  // the cut expression
814  Double_t cutVal = 1.;
815  formula = fCutFormulas[cl];
816  if (formula) {
817  Int_t ndata = formula->GetNdata();
818  cutVal = (ndata==1 ?
819  formula->EvalInstance(0) :
820  formula->EvalInstance(idata));
821  checkNanInf(nanInfErrors, cutVal, "Cut", formula->GetTitle());
822  }
823 
824  // if event is cut out, add to warnings, else add to errors.
825  auto &nanMessages = cutVal < 0.5 ? nanInfWarnings : nanInfErrors;
826 
827  // the input variable
828  for (UInt_t ivar=0; ivar<nvars; ivar++) {
829  formula = fInputFormulas[ivar];
830  formula->SetQuickLoad(true);
831  Int_t ndata = formula->GetNdata();
832  vars[ivar] = (ndata == 1 ?
833  formula->EvalInstance(0) :
834  formula->EvalInstance(idata));
835  checkNanInf(nanMessages, vars[ivar], "Input", formula->GetTitle());
836  }
837 
838  // the targets
839  for (UInt_t itrgt=0; itrgt<ntgts; itrgt++) {
840  formula = fTargetFormulas[itrgt];
841  Int_t ndata = formula->GetNdata();
842  tgts[itrgt] = (ndata == 1 ?
843  formula->EvalInstance(0) :
844  formula->EvalInstance(idata));
845  checkNanInf(nanMessages, tgts[itrgt], "Target", formula->GetTitle());
846  }
847 
848  // the spectators
849  for (UInt_t itVis=0; itVis<nvis; itVis++) {
850  formula = fSpectatorFormulas[itVis];
851  Int_t ndata = formula->GetNdata();
852  vis[itVis] = (ndata == 1 ?
853  formula->EvalInstance(0) :
854  formula->EvalInstance(idata));
855  checkNanInf(nanMessages, vis[itVis], "Spectator", formula->GetTitle());
856  }
857 
858 
859  // the weight
860  Float_t weight = currentInfo.GetWeight(); // multiply by tree weight
861  formula = fWeightFormula[cl];
862  if (formula!=0) {
863  Int_t ndata = formula->GetNdata();
864  weight *= (ndata == 1 ?
865  formula->EvalInstance() :
866  formula->EvalInstance(idata));
867  checkNanInf(nanMessages, weight, "Weight", formula->GetTitle());
868  }
869 
870  // Count the events before rejection due to cut or NaN
871  // value (weighted and unweighted)
872  classEventCounts.nEvBeforeCut++;
873  if (!TMath::IsNaN(weight))
874  classEventCounts.nWeEvBeforeCut += weight;
875 
876  // apply the cut, skip rest if cut is not fulfilled
877  if (cutVal<0.5) continue;
878 
879  // global flag if negative weights exist -> can be used
880  // by classifiers who may require special data
881  // treatment (also print warning)
882  if (weight < 0) classEventCounts.nNegWeights++;
883 
884  // now read the event-values (variables and regression targets)
885 
886  if (contains_NaN_or_inf) {
887  Log() << kWARNING << Form("Dataset[%s] : ",dsi.GetName())<< "NaN or +-inf in Event " << evtIdx << Endl;
888  if (sizeOfArrays>1) Log() << kWARNING << Form("Dataset[%s] : ",dsi.GetName())<< " rejected" << Endl;
889  continue;
890  }
891 
892  // Count the events after rejection due to cut or NaN value
893  // (weighted and unweighted)
894  classEventCounts.nEvAfterCut++;
895  classEventCounts.nWeEvAfterCut += weight;
896 
897  // event accepted, fill temporary ntuple
898  event_v.push_back(new Event(vars, tgts , vis, cl , weight));
899  }
900  }
901  currentInfo.GetTree()->ResetBranchAddresses();
902  }
903  }
904 
905  if (!nanInfWarnings.empty()) {
906  Log() << kWARNING << "Found events with NaN and/or +-inf values" << Endl;
907  for (const auto &warning : nanInfWarnings) {
908  auto &log = Log() << kWARNING << warning.first;
909  if (warning.second > 1) log << " (" << warning.second << " times)";
910  log << Endl;
911  }
912  Log() << kWARNING << "These NaN and/or +-infs were all removed by the specified cut, continuing." << Endl;
913  Log() << Endl;
914  }
915 
916  if (!nanInfErrors.empty()) {
917  Log() << kWARNING << "Found events with NaN and/or +-inf values (not removed by cut)" << Endl;
918  for (const auto &error : nanInfErrors) {
919  auto &log = Log() << kWARNING << error.first;
920  if (error.second > 1) log << " (" << error.second << " times)";
921  log << Endl;
922  }
923  Log() << kFATAL << "How am I supposed to train a NaN or +-inf?!" << Endl;
924  }
925 
926  // for output format, get the maximum class name length
927  Int_t maxL = dsi.GetClassNameMaxLength();
928 
929  Log() << kHEADER << Form("[%s] : ",dsi.GetName()) << "Number of events in input trees" << Endl;
930  Log() << kDEBUG << "(after possible flattening of arrays):" << Endl;
931 
932 
933  for (UInt_t cl = 0; cl < dsi.GetNClasses(); cl++) {
934  Log() << kDEBUG //<< Form("[%s] : ",dsi.GetName())
935  << " "
936  << setiosflags(ios::left) << std::setw(maxL) << dsi.GetClassInfo(cl)->GetName()
937  << " -- number of events : "
938  << std::setw(5) << eventCounts[cl].nEvBeforeCut
939  << " / sum of weights: " << std::setw(5) << eventCounts[cl].nWeEvBeforeCut << Endl;
940  }
941 
942  for (UInt_t cl = 0; cl < dsi.GetNClasses(); cl++) {
943  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
944  << " " << std::setw(maxL) << dsi.GetClassInfo(cl)->GetName()
945  <<" tree -- total number of entries: "
946  << std::setw(5) << dataInput.GetEntries(dsi.GetClassInfo(cl)->GetName()) << Endl;
947  }
948 
950  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
951  << "\tPreselection: (will affect number of requested training and testing events)" << Endl;
952  else
953  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
954  << "\tPreselection: (will NOT affect number of requested training and testing events)" << Endl;
955 
956  if (dsi.HasCuts()) {
957  for (UInt_t cl = 0; cl< dsi.GetNClasses(); cl++) {
958  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " " << setiosflags(ios::left) << std::setw(maxL) << dsi.GetClassInfo(cl)->GetName()
959  << " requirement: \"" << dsi.GetClassInfo(cl)->GetCut() << "\"" << Endl;
960  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " "
961  << setiosflags(ios::left) << std::setw(maxL) << dsi.GetClassInfo(cl)->GetName()
962  << " -- number of events passed: "
963  << std::setw(5) << eventCounts[cl].nEvAfterCut
964  << " / sum of weights: " << std::setw(5) << eventCounts[cl].nWeEvAfterCut << Endl;
965  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " "
966  << setiosflags(ios::left) << std::setw(maxL) << dsi.GetClassInfo(cl)->GetName()
967  << " -- efficiency : "
968  << std::setw(6) << eventCounts[cl].nWeEvAfterCut/eventCounts[cl].nWeEvBeforeCut << Endl;
969  }
970  }
971  else Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
972  << " No preselection cuts applied on event classes" << Endl;
973 
974  delete[] varIsArray;
975 
976 }
977 
978 ////////////////////////////////////////////////////////////////////////////////
979 /// Select and distribute unassigned events to kTraining and kTesting
980 
983  EventVectorOfClassesOfTreeType& tmpEventVector,
984  EvtStatsPerClass& eventCounts,
985  const TString& splitMode,
986  const TString& mixMode,
987  const TString& normMode,
988  UInt_t splitSeed)
989 {
990  TMVA::RandomGenerator<TRandom3> rndm(splitSeed);
991 
992  // ==== splitting of undefined events to kTraining and kTesting
993 
994  // if splitMode contains "RANDOM", then shuffle the undefined events
995  if (splitMode.Contains( "RANDOM" ) /*&& !emptyUndefined*/ ) {
996  // random shuffle the undefined events of each class
997  for( UInt_t cls = 0; cls < dsi.GetNClasses(); ++cls ){
998  EventVector& unspecifiedEvents = tmpEventVector[Types::kMaxTreeType].at(cls);
999  if( ! unspecifiedEvents.empty() ) {
1000  Log() << kDEBUG << "randomly shuffling "
1001  << unspecifiedEvents.size()
1002  << " events of class " << cls
1003  << " which are not yet associated to testing or training" << Endl;
1004  std::shuffle(unspecifiedEvents.begin(), unspecifiedEvents.end(), rndm);
1005  }
1006  }
1007  }
1008 
1009  // check for each class the number of training and testing events, the requested number and the available number
1010  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "SPLITTING ========" << Endl;
1011  for( UInt_t cls = 0; cls < dsi.GetNClasses(); ++cls ){
1012  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "---- class " << cls << Endl;
1013  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "check number of training/testing events, requested and available number of events and for class " << cls << Endl;
1014 
1015  // check if enough or too many events are already in the training/testing eventvectors of the class cls
1016  EventVector& eventVectorTraining = tmpEventVector[ Types::kTraining ].at(cls);
1017  EventVector& eventVectorTesting = tmpEventVector[ Types::kTesting ].at(cls);
1018  EventVector& eventVectorUndefined = tmpEventVector[ Types::kMaxTreeType ].at(cls);
1019 
1020  Int_t availableTraining = eventVectorTraining.size();
1021  Int_t availableTesting = eventVectorTesting.size();
1022  Int_t availableUndefined = eventVectorUndefined.size();
1023 
1024  Float_t presel_scale;
1025  if (fScaleWithPreselEff) {
1026  presel_scale = eventCounts[cls].cutScaling();
1027  if (presel_scale < 1)
1028  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " you have opted for scaling the number of requested training/testing events\n to be scaled by the preselection efficiency"<< Endl;
1029  }else{
1030  presel_scale = 1.; // this scaling was too confusing to most people, including me! Sorry... (Helge)
1031  if (eventCounts[cls].cutScaling() < 1)
1032  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " you have opted for interpreting the requested number of training/testing events\n to be the number of events AFTER your preselection cuts" << Endl;
1033 
1034  }
1035 
1036  // If TrainTestSplit_<class> is set, set number of requested training events to split*num_all_events
1037  // Requested number of testing events is set to zero and therefore takes all other events
1038  // The option TrainTestSplit_<class> overrides nTrain_<class> or nTest_<class>
1039  if(eventCounts[cls].TrainTestSplitRequested < 1.0 && eventCounts[cls].TrainTestSplitRequested > 0.0){
1040  eventCounts[cls].nTrainingEventsRequested = Int_t(eventCounts[cls].TrainTestSplitRequested*(availableTraining+availableTesting+availableUndefined));
1041  eventCounts[cls].nTestingEventsRequested = Int_t(0);
1042  }
1043  else if(eventCounts[cls].TrainTestSplitRequested != 0.0) Log() << kFATAL << Form("The option TrainTestSplit_<class> has to be in range (0, 1] but is set to %f.",eventCounts[cls].TrainTestSplitRequested) << Endl;
1044  Int_t requestedTraining = Int_t(eventCounts[cls].nTrainingEventsRequested * presel_scale);
1045  Int_t requestedTesting = Int_t(eventCounts[cls].nTestingEventsRequested * presel_scale);
1046 
1047  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "events in training trees : " << availableTraining << Endl;
1048  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "events in testing trees : " << availableTesting << Endl;
1049  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "events in unspecified trees : " << availableUndefined << Endl;
1050  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "requested for training : " << requestedTraining << Endl;;
1051 
1052  if(presel_scale<1)
1053  Log() << " ( " << eventCounts[cls].nTrainingEventsRequested
1054  << " * " << presel_scale << " preselection efficiency)" << Endl;
1055  else
1056  Log() << Endl;
1057  Log() << kDEBUG << "requested for testing : " << requestedTesting;
1058  if(presel_scale<1)
1059  Log() << " ( " << eventCounts[cls].nTestingEventsRequested
1060  << " * " << presel_scale << " preselection efficiency)" << Endl;
1061  else
1062  Log() << Endl;
1063 
1064  // nomenclature r = available training
1065  // s = available testing
1066  // u = available undefined
1067  // R = requested training
1068  // S = requested testing
1069  // nR = to be used to select training events
1070  // nS = to be used to select test events
1071  // we have the constraint: nR + nS < r+s+u,
1072  // since we can not use more events than we have
1073  // free events: Nfree = u-Thet(R-r)-Thet(S-s)
1074  // nomenclature: Thet(x) = x, if x>0 else 0
1075  // nR = max(R,r) + 0.5 * Nfree
1076  // nS = max(S,s) + 0.5 * Nfree
1077  // nR +nS = R+S + u-R+r-S+s = u+r+s= ok! for R>r
1078  // nR +nS = r+S + u-S+s = u+r+s= ok! for r>R
1079 
1080  // three different cases might occur here
1081  //
1082  // Case a
1083  // requestedTraining and requestedTesting >0
1084  // free events: Nfree = u-Thet(R-r)-Thet(S-s)
1085  // nR = Max(R,r) + 0.5 * Nfree
1086  // nS = Max(S,s) + 0.5 * Nfree
1087  //
1088  // Case b
1089  // exactly one of requestedTraining or requestedTesting >0
1090  // assume training R >0
1091  // nR = max(R,r)
1092  // nS = s+u+r-nR
1093  // and s=nS
1094  //
1095  // Case c
1096  // requestedTraining=0, requestedTesting=0
1097  // Nfree = u-|r-s|
1098  // if NFree >=0
1099  // R = Max(r,s) + 0.5 * Nfree = S
1100  // else if r>s
1101  // R = r; S=s+u
1102  // else
1103  // R = r+u; S=s
1104  //
1105  // Next steps:
1106  // Determination of Event numbers R,S, nR, nS
1107  // distribute undefined events according to nR, nS
1108  // finally determine actual sub samples from nR and nS to be used in training / testing
1109  //
1110 
1111  Int_t useForTesting(0),useForTraining(0);
1112  Int_t allAvailable(availableUndefined + availableTraining + availableTesting);
1113 
1114  if( (requestedTraining == 0) && (requestedTesting == 0)){
1115 
1116  // Case C: balance the number of training and testing events
1117 
1118  if ( availableUndefined >= TMath::Abs(availableTraining - availableTesting) ) {
1119  // enough unspecified are available to equal training and testing
1120  useForTraining = useForTesting = allAvailable/2;
1121  } else {
1122  // all unspecified are assigned to the smaller of training / testing
1123  useForTraining = availableTraining;
1124  useForTesting = availableTesting;
1125  if (availableTraining < availableTesting)
1126  useForTraining += availableUndefined;
1127  else
1128  useForTesting += availableUndefined;
1129  }
1130  requestedTraining = useForTraining;
1131  requestedTesting = useForTesting;
1132  }
1133 
1134  else if (requestedTesting == 0){
1135  // case B
1136  useForTraining = TMath::Max(requestedTraining,availableTraining);
1137  if (allAvailable < useForTraining) {
1138  Log() << kFATAL << Form("Dataset[%s] : ",dsi.GetName())<< "More events requested for training ("
1139  << requestedTraining << ") than available ("
1140  << allAvailable << ")!" << Endl;
1141  }
1142  useForTesting = allAvailable - useForTraining; // the rest
1143  requestedTesting = useForTesting;
1144  }
1145 
1146  else if (requestedTraining == 0){ // case B)
1147  useForTesting = TMath::Max(requestedTesting,availableTesting);
1148  if (allAvailable < useForTesting) {
1149  Log() << kFATAL << Form("Dataset[%s] : ",dsi.GetName())<< "More events requested for testing ("
1150  << requestedTesting << ") than available ("
1151  << allAvailable << ")!" << Endl;
1152  }
1153  useForTraining= allAvailable - useForTesting; // the rest
1154  requestedTraining = useForTraining;
1155  }
1156 
1157  else {
1158  // Case A
1159  // requestedTraining R and requestedTesting S >0
1160  // free events: Nfree = u-Thet(R-r)-Thet(S-s)
1161  // nR = Max(R,r) + 0.5 * Nfree
1162  // nS = Max(S,s) + 0.5 * Nfree
1163  Int_t stillNeedForTraining = TMath::Max(requestedTraining-availableTraining,0);
1164  Int_t stillNeedForTesting = TMath::Max(requestedTesting-availableTesting,0);
1165 
1166  int NFree = availableUndefined - stillNeedForTraining - stillNeedForTesting;
1167  if (NFree <0) NFree = 0;
1168  useForTraining = TMath::Max(requestedTraining,availableTraining) + NFree/2;
1169  useForTesting= allAvailable - useForTraining; // the rest
1170  }
1171 
1172  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "determined event sample size to select training sample from="<<useForTraining<<Endl;
1173  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "determined event sample size to select test sample from="<<useForTesting<<Endl;
1174 
1175 
1176 
1177  // associate undefined events
1178  if( splitMode == "ALTERNATE" ){
1179  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "split 'ALTERNATE'" << Endl;
1180  Int_t nTraining = availableTraining;
1181  Int_t nTesting = availableTesting;
1182  for( EventVector::iterator it = eventVectorUndefined.begin(), itEnd = eventVectorUndefined.end(); it != itEnd; ){
1183  ++nTraining;
1184  if( nTraining <= requestedTraining ){
1185  eventVectorTraining.insert( eventVectorTraining.end(), (*it) );
1186  ++it;
1187  }
1188  if( it != itEnd ){
1189  ++nTesting;
1190  eventVectorTesting.insert( eventVectorTesting.end(), (*it) );
1191  ++it;
1192  }
1193  }
1194  } else {
1195  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "split '" << splitMode << "'" << Endl;
1196 
1197  // test if enough events are available
1198  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "availableundefined : " << availableUndefined << Endl;
1199  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "useForTraining : " << useForTraining << Endl;
1200  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "useForTesting : " << useForTesting << Endl;
1201  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "availableTraining : " << availableTraining << Endl;
1202  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "availableTesting : " << availableTesting << Endl;
1203 
1204  if( availableUndefined<(useForTraining-availableTraining) ||
1205  availableUndefined<(useForTesting -availableTesting ) ||
1206  availableUndefined<(useForTraining+useForTesting-availableTraining-availableTesting ) ){
1207  Log() << kFATAL << Form("Dataset[%s] : ",dsi.GetName())<< "More events requested than available!" << Endl;
1208  }
1209 
1210  // select the events
1211  if (useForTraining>availableTraining){
1212  eventVectorTraining.insert( eventVectorTraining.end() , eventVectorUndefined.begin(), eventVectorUndefined.begin()+ useForTraining- availableTraining );
1213  eventVectorUndefined.erase( eventVectorUndefined.begin(), eventVectorUndefined.begin() + useForTraining- availableTraining);
1214  }
1215  if (useForTesting>availableTesting){
1216  eventVectorTesting.insert( eventVectorTesting.end() , eventVectorUndefined.begin(), eventVectorUndefined.begin()+ useForTesting- availableTesting );
1217  }
1218  }
1219  eventVectorUndefined.clear();
1220 
1221  // finally shorten the event vectors to the requested size by removing random events
1222  if (splitMode.Contains( "RANDOM" )){
1223  UInt_t sizeTraining = eventVectorTraining.size();
1224  if( sizeTraining > UInt_t(requestedTraining) ){
1225  std::vector<UInt_t> indicesTraining( sizeTraining );
1226  // make indices
1227  std::generate( indicesTraining.begin(), indicesTraining.end(), TMVA::Increment<UInt_t>(0) );
1228  // shuffle indices
1229  std::shuffle(indicesTraining.begin(), indicesTraining.end(), rndm);
1230  // erase indices of not needed events
1231  indicesTraining.erase( indicesTraining.begin()+sizeTraining-UInt_t(requestedTraining), indicesTraining.end() );
1232  // delete all events with the given indices
1233  for( std::vector<UInt_t>::iterator it = indicesTraining.begin(), itEnd = indicesTraining.end(); it != itEnd; ++it ){
1234  delete eventVectorTraining.at( (*it) ); // delete event
1235  eventVectorTraining.at( (*it) ) = NULL; // set pointer to NULL
1236  }
1237  // now remove and erase all events with pointer==NULL
1238  eventVectorTraining.erase( std::remove( eventVectorTraining.begin(), eventVectorTraining.end(), (void*)NULL ), eventVectorTraining.end() );
1239  }
1240 
1241  UInt_t sizeTesting = eventVectorTesting.size();
1242  if( sizeTesting > UInt_t(requestedTesting) ){
1243  std::vector<UInt_t> indicesTesting( sizeTesting );
1244  // make indices
1245  std::generate( indicesTesting.begin(), indicesTesting.end(), TMVA::Increment<UInt_t>(0) );
1246  // shuffle indices
1247  std::shuffle(indicesTesting.begin(), indicesTesting.end(), rndm);
1248  // erase indices of not needed events
1249  indicesTesting.erase( indicesTesting.begin()+sizeTesting-UInt_t(requestedTesting), indicesTesting.end() );
1250  // delete all events with the given indices
1251  for( std::vector<UInt_t>::iterator it = indicesTesting.begin(), itEnd = indicesTesting.end(); it != itEnd; ++it ){
1252  delete eventVectorTesting.at( (*it) ); // delete event
1253  eventVectorTesting.at( (*it) ) = NULL; // set pointer to NULL
1254  }
1255  // now remove and erase all events with pointer==NULL
1256  eventVectorTesting.erase( std::remove( eventVectorTesting.begin(), eventVectorTesting.end(), (void*)NULL ), eventVectorTesting.end() );
1257  }
1258  }
1259  else { // erase at end
1260  if( eventVectorTraining.size() < UInt_t(requestedTraining) )
1261  Log() << kWARNING << Form("Dataset[%s] : ",dsi.GetName())<< "DataSetFactory/requested number of training samples larger than size of eventVectorTraining.\n"
1262  << "There is probably an issue. Please contact the TMVA developers." << Endl;
1263  std::for_each( eventVectorTraining.begin()+requestedTraining, eventVectorTraining.end(), DeleteFunctor<Event>() );
1264  eventVectorTraining.erase(eventVectorTraining.begin()+requestedTraining,eventVectorTraining.end());
1265 
1266  if( eventVectorTesting.size() < UInt_t(requestedTesting) )
1267  Log() << kWARNING << Form("Dataset[%s] : ",dsi.GetName())<< "DataSetFactory/requested number of testing samples larger than size of eventVectorTesting.\n"
1268  << "There is probably an issue. Please contact the TMVA developers." << Endl;
1269  std::for_each( eventVectorTesting.begin()+requestedTesting, eventVectorTesting.end(), DeleteFunctor<Event>() );
1270  eventVectorTesting.erase(eventVectorTesting.begin()+requestedTesting,eventVectorTesting.end());
1271  }
1272  }
1273 
1274  TMVA::DataSetFactory::RenormEvents( dsi, tmpEventVector, eventCounts, normMode );
1275 
1276  Int_t trainingSize = 0;
1277  Int_t testingSize = 0;
1278 
1279  // sum up number of training and testing events
1280  for( UInt_t cls = 0; cls < dsi.GetNClasses(); ++cls ){
1281  trainingSize += tmpEventVector[Types::kTraining].at(cls).size();
1282  testingSize += tmpEventVector[Types::kTesting].at(cls).size();
1283  }
1284 
1285  // --- collect all training (testing) events into the training (testing) eventvector
1286 
1287  // create event vectors reserve enough space
1288  EventVector* trainingEventVector = new EventVector();
1289  EventVector* testingEventVector = new EventVector();
1290 
1291  trainingEventVector->reserve( trainingSize );
1292  testingEventVector->reserve( testingSize );
1293 
1294 
1295  // collect the events
1296 
1297  // mixing of kTraining and kTesting data sets
1298  Log() << kDEBUG << " MIXING ============= " << Endl;
1299 
1300  if( mixMode == "ALTERNATE" ){
1301  // Inform user if he tries to use alternate mixmode for
1302  // event classes with different number of events, this works but the alternation stops at the last event of the smaller class
1303  for( UInt_t cls = 1; cls < dsi.GetNClasses(); ++cls ){
1304  if (tmpEventVector[Types::kTraining].at(cls).size() != tmpEventVector[Types::kTraining].at(0).size()){
1305  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "Training sample: You are trying to mix events in alternate mode although the classes have different event numbers. This works but the alternation stops at the last event of the smaller class."<<Endl;
1306  }
1307  if (tmpEventVector[Types::kTesting].at(cls).size() != tmpEventVector[Types::kTesting].at(0).size()){
1308  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "Testing sample: You are trying to mix events in alternate mode although the classes have different event numbers. This works but the alternation stops at the last event of the smaller class."<<Endl;
1309  }
1310  }
1311  typedef EventVector::iterator EvtVecIt;
1312  EvtVecIt itEvent, itEventEnd;
1313 
1314  // insert first class
1315  Log() << kDEBUG << "insert class 0 into training and test vector" << Endl;
1316  trainingEventVector->insert( trainingEventVector->end(), tmpEventVector[Types::kTraining].at(0).begin(), tmpEventVector[Types::kTraining].at(0).end() );
1317  testingEventVector->insert( testingEventVector->end(), tmpEventVector[Types::kTesting].at(0).begin(), tmpEventVector[Types::kTesting].at(0).end() );
1318 
1319  // insert other classes
1320  EvtVecIt itTarget;
1321  for( UInt_t cls = 1; cls < dsi.GetNClasses(); ++cls ){
1322  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "insert class " << cls << Endl;
1323  // training vector
1324  itTarget = trainingEventVector->begin() - 1; // start one before begin
1325  // loop over source
1326  for( itEvent = tmpEventVector[Types::kTraining].at(cls).begin(), itEventEnd = tmpEventVector[Types::kTraining].at(cls).end(); itEvent != itEventEnd; ++itEvent ){
1327  // if( std::distance( itTarget, trainingEventVector->end()) < Int_t(cls+1) ) {
1328  if( (trainingEventVector->end() - itTarget) < Int_t(cls+1) ) {
1329  itTarget = trainingEventVector->end();
1330  trainingEventVector->insert( itTarget, itEvent, itEventEnd ); // fill in the rest without mixing
1331  break;
1332  }else{
1333  itTarget += cls+1;
1334  trainingEventVector->insert( itTarget, (*itEvent) ); // fill event
1335  }
1336  }
1337  // testing vector
1338  itTarget = testingEventVector->begin() - 1;
1339  // loop over source
1340  for( itEvent = tmpEventVector[Types::kTesting].at(cls).begin(), itEventEnd = tmpEventVector[Types::kTesting].at(cls).end(); itEvent != itEventEnd; ++itEvent ){
1341  // if( std::distance( itTarget, testingEventVector->end()) < Int_t(cls+1) ) {
1342  if( ( testingEventVector->end() - itTarget ) < Int_t(cls+1) ) {
1343  itTarget = testingEventVector->end();
1344  testingEventVector->insert( itTarget, itEvent, itEventEnd ); // fill in the rest without mixing
1345  break;
1346  }else{
1347  itTarget += cls+1;
1348  testingEventVector->insert( itTarget, (*itEvent) ); // fill event
1349  }
1350  }
1351  }
1352  }else{
1353  for( UInt_t cls = 0; cls < dsi.GetNClasses(); ++cls ){
1354  trainingEventVector->insert( trainingEventVector->end(), tmpEventVector[Types::kTraining].at(cls).begin(), tmpEventVector[Types::kTraining].at(cls).end() );
1355  testingEventVector->insert ( testingEventVector->end(), tmpEventVector[Types::kTesting].at(cls).begin(), tmpEventVector[Types::kTesting].at(cls).end() );
1356  }
1357  }
1358  // delete the tmpEventVector (but not the events therein)
1359  tmpEventVector[Types::kTraining].clear();
1360  tmpEventVector[Types::kTesting].clear();
1361 
1362  tmpEventVector[Types::kMaxTreeType].clear();
1363 
1364  if (mixMode == "RANDOM") {
1365  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "shuffling events"<<Endl;
1366 
1367  std::shuffle(trainingEventVector->begin(), trainingEventVector->end(), rndm);
1368  std::shuffle(testingEventVector->begin(), testingEventVector->end(), rndm);
1369  }
1370 
1371  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "trainingEventVector " << trainingEventVector->size() << Endl;
1372  Log() << kDEBUG << Form("Dataset[%s] : ",dsi.GetName())<< "testingEventVector " << testingEventVector->size() << Endl;
1373 
1374  // create dataset
1375  DataSet* ds = new DataSet(dsi);
1376 
1377  // Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "Create internal training tree" << Endl;
1378  ds->SetEventCollection(trainingEventVector, Types::kTraining );
1379  // Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "Create internal testing tree" << Endl;
1380  ds->SetEventCollection(testingEventVector, Types::kTesting );
1381 
1382 
1383  if (ds->GetNTrainingEvents() < 1){
1384  Log() << kFATAL << "Dataset " << std::string(dsi.GetName()) << " does not have any training events, I better stop here and let you fix that one first " << Endl;
1385  }
1386 
1387  if (ds->GetNTestEvents() < 1) {
1388  Log() << kERROR << "Dataset " << std::string(dsi.GetName()) << " does not have any testing events, guess that will cause problems later..but for now, I continue " << Endl;
1389  }
1390 
1391  delete trainingEventVector;
1392  delete testingEventVector;
1393  return ds;
1394 
1395 }
1396 
1397 ////////////////////////////////////////////////////////////////////////////////
1398 /// renormalisation of the TRAINING event weights
1399 /// - none (kind of obvious) .. use the weights as supplied by the
1400 /// user.. (we store however the relative weight for later use)
1401 /// - numEvents
1402 /// - equalNumEvents reweight the training events such that the sum of all
1403 /// backgr. (class > 0) weights equal that of the signal (class 0)
1404 
1405 void
1407  EventVectorOfClassesOfTreeType& tmpEventVector,
1408  const EvtStatsPerClass& eventCounts,
1409  const TString& normMode )
1410 {
1411 
1412 
1413  // print rescaling info
1414  // ---------------------------------
1415  // compute sizes and sums of weights
1416  Int_t trainingSize = 0;
1417  Int_t testingSize = 0;
1418 
1419  ValuePerClass trainingSumWeightsPerClass( dsi.GetNClasses() );
1420  ValuePerClass testingSumWeightsPerClass( dsi.GetNClasses() );
1421 
1422  NumberPerClass trainingSizePerClass( dsi.GetNClasses() );
1423  NumberPerClass testingSizePerClass( dsi.GetNClasses() );
1424 
1425  Double_t trainingSumSignalWeights = 0;
1426  Double_t trainingSumBackgrWeights = 0; // Backgr. includes all classes that are not signal
1427  Double_t testingSumSignalWeights = 0;
1428  Double_t testingSumBackgrWeights = 0; // Backgr. includes all classes that are not signal
1429 
1430 
1431 
1432  for( UInt_t cls = 0, clsEnd = dsi.GetNClasses(); cls < clsEnd; ++cls ){
1433  trainingSizePerClass.at(cls) = tmpEventVector[Types::kTraining].at(cls).size();
1434  testingSizePerClass.at(cls) = tmpEventVector[Types::kTesting].at(cls).size();
1435 
1436  trainingSize += trainingSizePerClass.back();
1437  testingSize += testingSizePerClass.back();
1438 
1439  // the functional solution
1440  // sum up the weights in Double_t although the individual weights are Float_t to prevent rounding issues in addition of floating points
1441  //
1442  // accumulate --> does what the name says
1443  // begin() and end() denote the range of the vector to be accumulated
1444  // Double_t(0) tells accumulate the type and the starting value
1445  // compose_binary creates a BinaryFunction of ...
1446  // std::plus<Double_t>() knows how to sum up two doubles
1447  // null<Double_t>() leaves the first argument (the running sum) unchanged and returns it
1448  //
1449  // all together sums up all the event-weights of the events in the vector and returns it
1450  trainingSumWeightsPerClass.at(cls) =
1451  std::accumulate(tmpEventVector[Types::kTraining].at(cls).begin(),
1452  tmpEventVector[Types::kTraining].at(cls).end(),
1453  Double_t(0), [](Double_t w, const TMVA::Event *E) { return w + E->GetOriginalWeight(); });
1454 
1455  testingSumWeightsPerClass.at(cls) =
1456  std::accumulate(tmpEventVector[Types::kTesting].at(cls).begin(),
1457  tmpEventVector[Types::kTesting].at(cls).end(),
1458  Double_t(0), [](Double_t w, const TMVA::Event *E) { return w + E->GetOriginalWeight(); });
1459 
1460  if ( cls == dsi.GetSignalClassIndex()){
1461  trainingSumSignalWeights += trainingSumWeightsPerClass.at(cls);
1462  testingSumSignalWeights += testingSumWeightsPerClass.at(cls);
1463  }else{
1464  trainingSumBackgrWeights += trainingSumWeightsPerClass.at(cls);
1465  testingSumBackgrWeights += testingSumWeightsPerClass.at(cls);
1466  }
1467  }
1468 
1469  // ---------------------------------
1470  // compute renormalization factors
1471 
1472  ValuePerClass renormFactor( dsi.GetNClasses() );
1473 
1474 
1475  // for information purposes
1476  dsi.SetNormalization( normMode );
1477  // !! these will be overwritten later by the 'rescaled' ones if
1478  // NormMode != None !!!
1479  dsi.SetTrainingSumSignalWeights(trainingSumSignalWeights);
1480  dsi.SetTrainingSumBackgrWeights(trainingSumBackgrWeights);
1481  dsi.SetTestingSumSignalWeights(testingSumSignalWeights);
1482  dsi.SetTestingSumBackgrWeights(testingSumBackgrWeights);
1483 
1484 
1485  if (normMode == "NONE") {
1486  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "No weight renormalisation applied: use original global and event weights" << Endl;
1487  return;
1488  }
1489  //changed by Helge 27.5.2013 What on earth was done here before? I still remember the idea behind this which apparently was
1490  //NOT understood by the 'programmer' :) .. the idea was to have SAME amount of effective TRAINING data for signal and background.
1491  // Testing events are totally irrelevant for this and might actually skew the whole normalisation!!
1492  else if (normMode == "NUMEVENTS") {
1493  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
1494  << "\tWeight renormalisation mode: \"NumEvents\": renormalises all event classes " << Endl;
1495  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
1496  << " such that the effective (weighted) number of events in each class equals the respective " << Endl;
1497  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
1498  << " number of events (entries) that you demanded in PrepareTrainingAndTestTree(\"\",\"nTrain_Signal=.. )" << Endl;
1499  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
1500  << " ... i.e. such that Sum[i=1..N_j]{w_i} = N_j, j=0,1,2..." << Endl;
1501  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
1502  << " ... (note that N_j is the sum of TRAINING events (nTrain_j...with j=Signal,Background.." << Endl;
1503  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
1504  << " ..... Testing events are not renormalised nor included in the renormalisation factor! )"<< Endl;
1505 
1506  for( UInt_t cls = 0, clsEnd = dsi.GetNClasses(); cls < clsEnd; ++cls ){
1507  // renormFactor.at(cls) = ( (trainingSizePerClass.at(cls) + testingSizePerClass.at(cls))/
1508  // (trainingSumWeightsPerClass.at(cls) + testingSumWeightsPerClass.at(cls)) );
1509  //changed by Helge 27.5.2013
1510  renormFactor.at(cls) = ((Float_t)trainingSizePerClass.at(cls) )/
1511  (trainingSumWeightsPerClass.at(cls)) ;
1512  }
1513  }
1514  else if (normMode == "EQUALNUMEVENTS") {
1515  //changed by Helge 27.5.2013 What on earth was done here before? I still remember the idea behind this which apparently was
1516  //NOT understood by the 'programmer' :) .. the idea was to have SAME amount of effective TRAINING data for signal and background.
1517  //done here was something like having each data source normalized to its number of entries and this even for training+testing together.
1518  // what should this have been good for ???
1519 
1520  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "Weight renormalisation mode: \"EqualNumEvents\": renormalises all event classes ..." << Endl;
1521  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " such that the effective (weighted) number of events in each class is the same " << Endl;
1522  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " (and equals the number of events (entries) given for class=0 )" << Endl;
1523  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "... i.e. such that Sum[i=1..N_j]{w_i} = N_classA, j=classA, classB, ..." << Endl;
1524  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << "... (note that N_j is the sum of TRAINING events" << Endl;
1525  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << " ..... Testing events are not renormalised nor included in the renormalisation factor!)" << Endl;
1526 
1527  // normalize to size of first class
1528  UInt_t referenceClass = 0;
1529  for (UInt_t cls = 0, clsEnd = dsi.GetNClasses(); cls < clsEnd; ++cls ) {
1530  renormFactor.at(cls) = Float_t(trainingSizePerClass.at(referenceClass))/
1531  (trainingSumWeightsPerClass.at(cls));
1532  }
1533  }
1534  else {
1535  Log() << kFATAL << Form("Dataset[%s] : ",dsi.GetName())<< "<PrepareForTrainingAndTesting> Unknown NormMode: " << normMode << Endl;
1536  }
1537 
1538  // ---------------------------------
1539  // now apply the normalization factors
1540  Int_t maxL = dsi.GetClassNameMaxLength();
1541  for (UInt_t cls = 0, clsEnd = dsi.GetNClasses(); cls<clsEnd; ++cls) {
1542  Log() << kDEBUG //<< Form("Dataset[%s] : ",dsi.GetName())
1543  << "--> Rescale " << setiosflags(ios::left) << std::setw(maxL)
1544  << dsi.GetClassInfo(cls)->GetName() << " event weights by factor: " << renormFactor.at(cls) << Endl;
1545  for (EventVector::iterator it = tmpEventVector[Types::kTraining].at(cls).begin(),
1546  itEnd = tmpEventVector[Types::kTraining].at(cls).end(); it != itEnd; ++it){
1547  (*it)->SetWeight ((*it)->GetWeight() * renormFactor.at(cls));
1548  }
1549 
1550  }
1551 
1552 
1553  // print out the result
1554  // (same code as before --> this can be done nicer )
1555  //
1556 
1557  Log() << kINFO //<< Form("Dataset[%s] : ",dsi.GetName())
1558  << "Number of training and testing events" << Endl;
1559  Log() << kDEBUG << "\tafter rescaling:" << Endl;
1560  Log() << kINFO //<< Form("Dataset[%s] : ",dsi.GetName())
1561  << "---------------------------------------------------------------------------" << Endl;
1562 
1563  trainingSumSignalWeights = 0;
1564  trainingSumBackgrWeights = 0; // Backgr. includes all classes that are not signal
1565  testingSumSignalWeights = 0;
1566  testingSumBackgrWeights = 0; // Backgr. includes all classes that are not signal
1567 
1568  for( UInt_t cls = 0, clsEnd = dsi.GetNClasses(); cls < clsEnd; ++cls ){
1569  trainingSumWeightsPerClass.at(cls) =
1570  std::accumulate(tmpEventVector[Types::kTraining].at(cls).begin(),
1571  tmpEventVector[Types::kTraining].at(cls).end(),
1572  Double_t(0), [](Double_t w, const TMVA::Event *E) { return w + E->GetOriginalWeight(); });
1573 
1574  testingSumWeightsPerClass.at(cls) =
1575  std::accumulate(tmpEventVector[Types::kTesting].at(cls).begin(),
1576  tmpEventVector[Types::kTesting].at(cls).end(),
1577  Double_t(0), [](Double_t w, const TMVA::Event *E) { return w + E->GetOriginalWeight(); });
1578 
1579  if ( cls == dsi.GetSignalClassIndex()){
1580  trainingSumSignalWeights += trainingSumWeightsPerClass.at(cls);
1581  testingSumSignalWeights += testingSumWeightsPerClass.at(cls);
1582  }else{
1583  trainingSumBackgrWeights += trainingSumWeightsPerClass.at(cls);
1584  testingSumBackgrWeights += testingSumWeightsPerClass.at(cls);
1585  }
1586 
1587  // output statistics
1588 
1589  Log() << kINFO //<< Form("Dataset[%s] : ",dsi.GetName())
1590  << setiosflags(ios::left) << std::setw(maxL)
1591  << dsi.GetClassInfo(cls)->GetName() << " -- "
1592  << "training events : " << trainingSizePerClass.at(cls) << Endl;
1593  Log() << kDEBUG << "\t(sum of weights: " << trainingSumWeightsPerClass.at(cls) << ")"
1594  << " - requested were " << eventCounts[cls].nTrainingEventsRequested << " events" << Endl;
1595  Log() << kINFO //<< Form("Dataset[%s] : ",dsi.GetName())
1596  << setiosflags(ios::left) << std::setw(maxL)
1597  << dsi.GetClassInfo(cls)->GetName() << " -- "
1598  << "testing events : " << testingSizePerClass.at(cls) << Endl;
1599  Log() << kDEBUG << "\t(sum of weights: " << testingSumWeightsPerClass.at(cls) << ")"
1600  << " - requested were " << eventCounts[cls].nTestingEventsRequested << " events" << Endl;
1601  Log() << kINFO //<< Form("Dataset[%s] : ",dsi.GetName())
1602  << setiosflags(ios::left) << std::setw(maxL)
1603  << dsi.GetClassInfo(cls)->GetName() << " -- "
1604  << "training and testing events: "
1605  << (trainingSizePerClass.at(cls)+testingSizePerClass.at(cls)) << Endl;
1606  Log() << kDEBUG << "\t(sum of weights: "
1607  << (trainingSumWeightsPerClass.at(cls)+testingSumWeightsPerClass.at(cls)) << ")" << Endl;
1608  if(eventCounts[cls].nEvAfterCut<eventCounts[cls].nEvBeforeCut) {
1609  Log() << kINFO << Form("Dataset[%s] : ",dsi.GetName()) << setiosflags(ios::left) << std::setw(maxL)
1610  << dsi.GetClassInfo(cls)->GetName() << " -- "
1611  << "due to the preselection a scaling factor has been applied to the numbers of requested events: "
1612  << eventCounts[cls].cutScaling() << Endl;
1613  }
1614  }
1615  Log() << kINFO << Endl;
1616 
1617  // for information purposes
1618  dsi.SetTrainingSumSignalWeights(trainingSumSignalWeights);
1619  dsi.SetTrainingSumBackgrWeights(trainingSumBackgrWeights);
1620  dsi.SetTestingSumSignalWeights(testingSumSignalWeights);
1621  dsi.SetTestingSumBackgrWeights(testingSumBackgrWeights);
1622 
1623 
1624 }
1625 
TTree * GetTree() const
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:32
UInt_t GetNVariables() const
Definition: DataSetInfo.h:110
std::vector< EventVector > EventVectorOfClasses
void SetTrainingSumBackgrWeights(Double_t trainingSumBackgrWeights)
Definition: DataSetInfo.h:118
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:158
long long Long64_t
Definition: RtypesCore.h:69
const TString & GetInternalName() const
Definition: VariableInfo.h:58
std::vector< VariableInfo > & GetSpectatorInfos()
Definition: DataSetInfo.h:104
std::vector< TTreeFormula * > fInputFormulas
float Float_t
Definition: RtypesCore.h:53
std::vector< TTreeFormula * > fCutFormulas
std::vector< Double_t > ValuePerClass
UInt_t GetNVariables() const
access the number of variables through the datasetinfo
Definition: DataSet.cxx:216
void SetTrainingSumSignalWeights(Double_t trainingSumSignalWeights)
Definition: DataSetInfo.h:117
void SetTestingSumBackgrWeights(Double_t testingSumBackgrWeights)
Definition: DataSetInfo.h:120
void BuildEventVector(DataSetInfo &dsi, DataInputHandler &dataInput, EventVectorOfClassesOfTreeType &eventsmap, EvtStatsPerClass &eventCounts)
build empty event vectors distributes events between kTraining/kTesting/kMaxTreeType ...
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
UInt_t GetNClasses() const
Definition: DataSetInfo.h:136
Bool_t IsNaN(Double_t x)
Definition: TMath.h:891
STL namespace.
std::vector< TreeInfo >::const_iterator end(const TString &className) const
void CalcMinMax(DataSet *, DataSetInfo &dsi)
compute covariance matrix
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
virtual Int_t GetNdim() const
Definition: TFormula.h:237
const TString & GetExpression() const
Definition: VariableInfo.h:57
std::vector< int > NumberPerClass
std::map< Types::ETreeType, EventVectorOfClasses > EventVectorOfClassesOfTreeType
double sqrt(double)
void InitOptions(DataSetInfo &dsi, EvtStatsPerClass &eventsmap, TString &normMode, UInt_t &splitSeed, TString &splitMode, TString &mixMode)
the dataset splitting
DataSet * BuildDynamicDataSet(DataSetInfo &)
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:2286
UInt_t GetNSpectators() const
access the number of targets through the datasetinfo
Definition: DataSet.cxx:232
MsgLogger & Log() const
message logger
std::vector< TTreeFormula * > fWeightFormula
UInt_t GetClass() const
Definition: Event.h:81
void SetTestingSumSignalWeights(Double_t testingSumSignalWeights)
Definition: DataSetInfo.h:119
void PrintCorrelationMatrix(const TString &className)
calculates the correlation matrices for signal and background, prints them to standard output...
void SetMin(Double_t v)
Definition: VariableInfo.h:69
void SetMinType(EMsgType minType)
Definition: MsgLogger.h:72
Int_t Finite(Double_t x)
Check if it is finite with a mask in order to be consistent in presence of fast math.
Definition: TMath.h:770
void RenormEvents(DataSetInfo &dsi, EventVectorOfClassesOfTreeType &eventsmap, const EvtStatsPerClass &eventCounts, const TString &normMode)
renormalisation of the TRAINING event weights
Class that contains all the data information.
Definition: DataSetInfo.h:60
virtual Int_t GetNcodes() const
Definition: TTreeFormula.h:193
Double_t GetWeight() const
return the event weight - depending on whether the flag IgnoreNegWeightsInTraining is or not...
Definition: Event.cxx:382
Long64_t GetNTrainingEvents() const
Definition: DataSet.h:79
Used to pass a selection expression to the Tree drawing routine.
Definition: TTreeFormula.h:58
void ChangeToNewTree(TreeInfo &, const DataSetInfo &)
While the data gets copied into the local training and testing trees, the input tree can change (for ...
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:22
void SetCorrelationMatrix(const TString &className, TMatrixD *matrix)
Class that contains all the data information.
Definition: DataSet.h:69
DataSetFactory()
constructor
TMatrixD * CalcCorrelationMatrix(DataSet *, const UInt_t classNumber)
computes correlation matrix for variables "theVars" in tree; "theType" defines the required event "ty...
Float_t GetTarget(UInt_t itgt) const
Definition: Event.h:97
Int_t LargestCommonDivider(Int_t a, Int_t b)
UInt_t GetNTargets() const
Definition: DataSetInfo.h:111
Types::ETreeType GetTreeType() const
Class that contains all the data information.
SVector< double, 2 > v
Definition: Dict.h:5
ClassInfo * GetClassInfo(Int_t clNum) const
std::vector< TTreeFormula * > fSpectatorFormulas
auto * a
Definition: textangle.C:12
DataSet * CreateDataSet(DataSetInfo &, DataInputHandler &)
steering the creation of a new dataset
VariableInfo & GetTargetInfo(Int_t i)
Definition: DataSetInfo.h:101
std::vector< TreeInfo >::const_iterator begin(const TString &className) const
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
UInt_t GetNSpectators(bool all=kTRUE) const
UInt_t GetSignalClassIndex()
Definition: DataSetInfo.h:139
std::vector< TTreeFormula * > fTargetFormulas
std::vector< Event *> EventVector
constexpr Double_t E()
Base of natural log: .
Definition: TMath.h:97
Long64_t GetNTestEvents() const
Definition: DataSet.h:80
void SetMax(Double_t v)
Definition: VariableInfo.h:70
const Bool_t kFALSE
Definition: RtypesCore.h:88
Float_t GetValue(UInt_t ivar) const
return value of i&#39;th variable
Definition: Event.cxx:237
DataSet * BuildInitialDataSet(DataSetInfo &, TMVA::DataInputHandler &)
if no entries, than create a DataSet with one Event which uses dynamic variables (pointers to variabl...
~DataSetFactory()
destructor
virtual Int_t GetNdata()
Return number of available instances in the formula.
#define d(i)
Definition: RSha256.hxx:102
virtual TLeaf * GetLeaf(Int_t n) const
Return leaf corresponding to serial number n.
double Double_t
Definition: RtypesCore.h:55
const TString & GetClassName() const
VariableInfo & GetSpectatorInfo(Int_t i)
Definition: DataSetInfo.h:106
UInt_t GetEntries(const TString &name) const
void SetEventCollection(std::vector< Event *> *, Types::ETreeType, Bool_t deleteEvents=true)
Sets the event collection (by DataSetFactory)
Definition: DataSet.cxx:250
VariableInfo & GetVariableInfo(Int_t i)
Definition: DataSetInfo.h:96
T EvalInstance(Int_t i=0, const char *stringStack[]=0)
Evaluate this treeformula.
ClassInfo * AddClass(const TString &className)
Int_t GetClassNameMaxLength() const
virtual const char * GetName() const
Returns name of object.
Definition: DataSetInfo.h:67
Long64_t GetNClassEvents(Int_t type, UInt_t classNumber)
Definition: DataSet.cxx:168
ostringstream derivative to redirect and format output
Definition: MsgLogger.h:59
void SetConfigName(const char *n)
Definition: Configurable.h:63
virtual const char * GetTitle() const
Returns title of object.
Definition: TObject.cxx:401
Abstract ClassifierFactory template that handles arbitrary types.
const TCut & GetCut() const
Definition: ClassInfo.h:64
std::vector< EventStats > EvtStatsPerClass
const TString & GetSplitOptions() const
Definition: DataSetInfo.h:167
Bool_t HasCuts() const
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:200
Double_t GetOriginalWeight() const
Definition: Event.h:79
Bool_t CheckTTreeFormula(TTreeFormula *ttf, const TString &expression, Bool_t &hasDollar)
checks a TTreeFormula for problems
void SetNumber(const UInt_t index)
Definition: ClassInfo.h:59
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
Long64_t GetNEvents(Types::ETreeType type=Types::kMaxTreeType) const
Definition: DataSet.h:217
void SetQuickLoad(Bool_t quick)
Definition: TTreeFormula.h:207
const TString & GetWeight() const
Definition: ClassInfo.h:63
DataSet * MixEvents(DataSetInfo &dsi, EventVectorOfClassesOfTreeType &eventsmap, EvtStatsPerClass &eventCounts, const TString &splitMode, const TString &mixMode, const TString &normMode, UInt_t splitSeed)
Select and distribute unassigned events to kTraining and kTesting.
TBranch * GetBranch() const
Definition: TLeaf.h:71
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
virtual Bool_t IsOnTerminalBranch() const
Definition: TLeaf.h:98
Double_t GetWeight() const
UInt_t GetNTargets() const
access the number of targets through the datasetinfo
Definition: DataSet.cxx:224
Float_t GetSpectator(UInt_t ivar) const
return spectator content
Definition: Event.cxx:262
const Bool_t kTRUE
Definition: RtypesCore.h:87
void SetNormalization(const TString &norm)
Definition: DataSetInfo.h:115
TMatrixD * CalcCovarianceMatrix(DataSet *, const UInt_t classNumber)
compute covariance matrix
const Event * GetEvent() const
Definition: DataSet.cxx:202
std::vector< VariableInfo > & GetVariableInfos()
Definition: DataSetInfo.h:94
double log(double)
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
std::vector< TString > * GetClassList() const