ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MethodPDEFoam.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Tancredi Carli, Dominik Dannheim, Alexander Voigt
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate Data analysis *
6  * Package: TMVA *
7  * Class : MethodPDEFoam *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementation (see header for description) *
12  * *
13  * Authors (alphabetical): *
14  * Tancredi Carli - CERN, Switzerland *
15  * Dominik Dannheim - CERN, Switzerland *
16  * Alexander Voigt - TU Dresden, Germany *
17  * Peter Speckmayer - CERN, Switzerland *
18  * *
19  * Original author of the TFoam implementation: *
20  * S. Jadach - Institute of Nuclear Physics, Cracow, Poland *
21  * *
22  * Copyright (c) 2008, 2010: *
23  * CERN, Switzerland *
24  * MPI-K Heidelberg, Germany *
25  * *
26  * Redistribution and use in source and binary forms, with or without *
27  * modification, are permitted according to the terms listed in LICENSE *
28  * (http://tmva.sourceforge.net/LICENSE) *
29  **********************************************************************************/
30 
31 //_______________________________________________________________________
32 //
33 // MethodPDEFoam
34 //
35 // The PDEFoam method is an extension of the PDERS method, which
36 // divides the multi-dimensional phase space in a finite number of
37 // hyper-rectangles (cells) of constant event density. This "foam" of
38 // cells is filled with averaged probability-density information
39 // sampled from a training event sample.
40 //
41 // For a given number of cells, the binning algorithm adjusts the size
42 // and position of the cells inside the multidimensional phase space
43 // based on a binary-split algorithm, minimizing the variance of the
44 // event density in the cell.
45 // The binned event density information of the final foam is stored in
46 // binary trees, allowing for a fast and memory-efficient
47 // classification of events.
48 //
49 // The implementation of PDEFoam is based on the Monte-Carlo
50 // integration package TFoam included in the analysis package ROOT.
51 // _______________________________________________________________________
52 
53 #include "TMVA/MethodPDEFoam.h"
54 
55 #include "TMath.h"
56 #include "TH1F.h"
57 #include "TFile.h"
58 
59 #include "TMVA/ClassifierFactory.h"
60 #include "TMVA/Config.h"
61 #include "TMVA/CrossEntropy.h"
62 #include "TMVA/DataSet.h"
63 #include "TMVA/DataSetInfo.h"
64 #include "TMVA/Event.h"
65 #include "TMVA/GiniIndex.h"
68 #include "TMVA/MethodBase.h"
69 #include "TMVA/MsgLogger.h"
70 #include "TMVA/Ranking.h"
71 #include "TMVA/SdivSqrtSplusB.h"
72 #include "TMVA/SeparationBase.h"
73 #include "TMVA/Tools.h"
74 #include "TMVA/Types.h"
75 #include "TMVA/VariableInfo.h"
76 
77 REGISTER_METHOD(PDEFoam)
78 
79 ClassImp(TMVA::MethodPDEFoam)
80 
81 ////////////////////////////////////////////////////////////////////////////////
82 /// init PDEFoam objects
83 
84 TMVA::MethodPDEFoam::MethodPDEFoam( const TString& jobName,
85  const TString& methodTitle,
86  DataSetInfo& dsi,
87  const TString& theOption,
88  TDirectory* theTargetDir ) :
89  MethodBase( jobName, Types::kPDEFoam, methodTitle, dsi, theOption, theTargetDir )
90  , fSigBgSeparated(kFALSE)
91  , fFrac(0.001)
92  , fDiscrErrCut(-1.0)
93  , fVolFrac(1.0/15.0)
94  , fnCells(999)
95  , fnActiveCells(500)
96  , fnSampl(2000)
97  , fnBin(5)
98  , fEvPerBin(10000)
99  , fCompress(kTRUE)
100  , fMultiTargetRegression(kFALSE)
101  , fNmin(100)
102  , fCutNmin(kTRUE)
103  , fMaxDepth(0)
104  , fKernelStr("None")
105  , fKernel(kNone)
106  , fKernelEstimator(NULL)
107  , fTargetSelectionStr("Mean")
108  , fTargetSelection(kMean)
109  , fFillFoamWithOrigWeights(kFALSE)
110  , fUseYesNoCell(kFALSE)
111  , fDTLogic("None")
112  , fDTSeparation(kFoam)
113  , fPeekMax(kTRUE)
114  , fXmin()
115  , fXmax()
116  , fFoam()
117 {
118 }
119 
120 ////////////////////////////////////////////////////////////////////////////////
121 /// constructor from weight file
122 
124  const TString& theWeightFile,
125  TDirectory* theTargetDir ) :
126  MethodBase( Types::kPDEFoam, dsi, theWeightFile, theTargetDir )
127  , fSigBgSeparated(kFALSE)
128  , fFrac(0.001)
129  , fDiscrErrCut(-1.0)
130  , fVolFrac(1.0/15.0)
131  , fnCells(999)
132  , fnActiveCells(500)
133  , fnSampl(2000)
134  , fnBin(5)
135  , fEvPerBin(10000)
136  , fCompress(kTRUE)
137  , fMultiTargetRegression(kFALSE)
138  , fNmin(100)
139  , fCutNmin(kTRUE)
140  , fMaxDepth(0)
141  , fKernelStr("None")
142  , fKernel(kNone)
143  , fKernelEstimator(NULL)
144  , fTargetSelectionStr("Mean")
145  , fTargetSelection(kMean)
146  , fFillFoamWithOrigWeights(kFALSE)
147  , fUseYesNoCell(kFALSE)
148  , fDTLogic("None")
149  , fDTSeparation(kFoam)
150  , fPeekMax(kTRUE)
151  , fXmin()
152  , fXmax()
153  , fFoam()
154 {
155 }
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// PDEFoam can handle classification with multiple classes and regression
159 /// with one or more regression-targets
160 
162 {
163  if (type == Types::kClassification && numberClasses == 2) return kTRUE;
164  if (type == Types::kMulticlass ) return kTRUE;
165  if (type == Types::kRegression) return kTRUE;
166  return kFALSE;
167 }
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 /// default initialization called by all constructors
171 
173 {
174  // init PDEFoam options
175  fSigBgSeparated = kFALSE; // default: unified foam
176  fFrac = 0.001; // fraction of outlier events
177  fDiscrErrCut = -1.; // cut on discriminator error
178  fVolFrac = 1./15.; // range searching box size
179  fnActiveCells = 500; // number of active cells to create
180  fnCells = fnActiveCells*2-1; // total number of cells
181  fnSampl = 2000; // number of sampling points in cell
182  fnBin = 5; // number of bins in edge histogram
183  fEvPerBin = 10000; // number of events per bin
184  fNmin = 100; // minimum number of events in cell
185  fMaxDepth = 0; // cell tree depth (default: unlimited)
186  fFillFoamWithOrigWeights = kFALSE; // fill orig. weights into foam
187  fUseYesNoCell = kFALSE; // return -1 or 1 for bg or signal events
188  fDTLogic = "None"; // decision tree algorithmus
189  fDTSeparation = kFoam; // separation type
190 
191  fKernel = kNone; // default: use no kernel
192  fKernelEstimator= NULL; // kernel estimator used during evaluation
193  fTargetSelection= kMean; // default: use mean for target selection (only multi target regression!)
194 
195  fCompress = kTRUE; // compress ROOT output file
196  fMultiTargetRegression = kFALSE; // multi-target regression
197 
198  DeleteFoams();
199 
200  if (fUseYesNoCell)
201  SetSignalReferenceCut( 0.0 ); // MVA output in [-1, 1]
202  else
203  SetSignalReferenceCut( 0.5 ); // MVA output in [0, 1]
204 }
205 
206 ////////////////////////////////////////////////////////////////////////////////
207 ///
208 /// Declare MethodPDEFoam options
209 ///
210 
212 {
213  DeclareOptionRef( fSigBgSeparated = kFALSE, "SigBgSeparate", "Separate foams for signal and background" );
214  DeclareOptionRef( fFrac = 0.001, "TailCut", "Fraction of outlier events that are excluded from the foam in each dimension" );
215  DeclareOptionRef( fVolFrac = 1./15., "VolFrac", "Size of sampling box, used for density calculation during foam build-up (maximum value: 1.0 is equivalent to volume of entire foam)");
216  DeclareOptionRef( fnActiveCells = 500, "nActiveCells", "Maximum number of active cells to be created by the foam");
217  DeclareOptionRef( fnSampl = 2000, "nSampl", "Number of generated MC events per cell");
218  DeclareOptionRef( fnBin = 5, "nBin", "Number of bins in edge histograms");
219  DeclareOptionRef( fCompress = kTRUE, "Compress", "Compress foam output file");
220  DeclareOptionRef( fMultiTargetRegression = kFALSE, "MultiTargetRegression", "Do regression with multiple targets");
221  DeclareOptionRef( fNmin = 100, "Nmin", "Number of events in cell required to split cell");
222  DeclareOptionRef( fMaxDepth = 0, "MaxDepth", "Maximum depth of cell tree (0=unlimited)");
223  DeclareOptionRef( fFillFoamWithOrigWeights = kFALSE, "FillFoamWithOrigWeights", "Fill foam with original or boost weights");
224  DeclareOptionRef( fUseYesNoCell = kFALSE, "UseYesNoCell", "Return -1 or 1 for bkg or signal like events");
225  DeclareOptionRef( fDTLogic = "None", "DTLogic", "Use decision tree algorithm to split cells");
226  AddPreDefVal(TString("None"));
227  AddPreDefVal(TString("GiniIndex"));
228  AddPreDefVal(TString("MisClassificationError"));
229  AddPreDefVal(TString("CrossEntropy"));
230  AddPreDefVal(TString("GiniIndexWithLaplace"));
231  AddPreDefVal(TString("SdivSqrtSplusB"));
232 
233  DeclareOptionRef( fKernelStr = "None", "Kernel", "Kernel type used");
234  AddPreDefVal(TString("None"));
235  AddPreDefVal(TString("Gauss"));
236  AddPreDefVal(TString("LinNeighbors"));
237  DeclareOptionRef( fTargetSelectionStr = "Mean", "TargetSelection", "Target selection method");
238  AddPreDefVal(TString("Mean"));
239  AddPreDefVal(TString("Mpv"));
240 }
241 
242 
243 ////////////////////////////////////////////////////////////////////////////////
244 /// options that are used ONLY for the READER to ensure backward compatibility
245 
248  DeclareOptionRef(fCutNmin = kTRUE, "CutNmin", "Requirement for minimal number of events in cell");
249  DeclareOptionRef(fPeekMax = kTRUE, "PeekMax", "Peek cell with max. loss for the next split");
250 }
251 
252 ////////////////////////////////////////////////////////////////////////////////
253 /// process user options
254 
256 {
257  if (!(fFrac>=0. && fFrac<=1.)) {
258  Log() << kWARNING << "TailCut not in [0.,1] ==> using 0.001 instead" << Endl;
259  fFrac = 0.001;
260  }
261 
262  if (fnActiveCells < 1) {
263  Log() << kWARNING << "invalid number of active cells specified: "
264  << fnActiveCells << "; setting nActiveCells=2" << Endl;
265  fnActiveCells = 2;
266  }
267  fnCells = fnActiveCells*2-1;
268 
269  // DT logic is only applicable if a single foam is trained
270  if (fSigBgSeparated && fDTLogic != "None") {
271  Log() << kFATAL << "Decision tree logic works only for a single foam (SigBgSeparate=F)" << Endl;
272  }
273 
274  // set separation to use
275  if (fDTLogic == "None")
276  fDTSeparation = kFoam;
277  else if (fDTLogic == "GiniIndex")
278  fDTSeparation = kGiniIndex;
279  else if (fDTLogic == "MisClassificationError")
280  fDTSeparation = kMisClassificationError;
281  else if (fDTLogic == "CrossEntropy")
282  fDTSeparation = kCrossEntropy;
283  else if (fDTLogic == "GiniIndexWithLaplace")
284  fDTSeparation = kGiniIndexWithLaplace;
285  else if (fDTLogic == "SdivSqrtSplusB")
286  fDTSeparation = kSdivSqrtSplusB;
287  else {
288  Log() << kWARNING << "Unknown separation type: " << fDTLogic
289  << ", setting to None" << Endl;
290  fDTLogic = "None";
291  fDTSeparation = kFoam;
292  }
293 
294  if (fKernelStr == "None" ) fKernel = kNone;
295  else if (fKernelStr == "Gauss" ) fKernel = kGaus;
296  else if (fKernelStr == "LinNeighbors") fKernel = kLinN;
297 
298  if (fTargetSelectionStr == "Mean" ) fTargetSelection = kMean;
299  else fTargetSelection = kMpv;
300  // sanity check: number of targets > 1 and MultiTargetRegression=F
301  // makes no sense --> set MultiTargetRegression=T
302  if (DoRegression() && Data()->GetNTargets() > 1 && !fMultiTargetRegression) {
303  Log() << kWARNING << "Warning: number of targets > 1"
304  << " and MultiTargetRegression=F was set, this makes no sense!"
305  << " --> I'm setting MultiTargetRegression=T" << Endl;
306  fMultiTargetRegression = kTRUE;
307  }
308 }
309 
310 ////////////////////////////////////////////////////////////////////////////////
311 /// destructor
312 
314 {
315  DeleteFoams();
316 
317  if (fKernelEstimator != NULL)
318  delete fKernelEstimator;
319 }
320 
321 ////////////////////////////////////////////////////////////////////////////////
322 /// Determine foam range [fXmin, fXmax] for all dimensions, such
323 /// that a fraction of 'fFrac' events lie outside the foam.
324 
326 {
327  fXmin.clear();
328  fXmax.clear();
329  UInt_t kDim = GetNvar(); // == Data()->GetNVariables();
330  UInt_t tDim = Data()->GetNTargets();
331  UInt_t vDim = Data()->GetNVariables();
332  if (fMultiTargetRegression)
333  kDim += tDim;
334 
335  Float_t *xmin = new Float_t[kDim];
336  Float_t *xmax = new Float_t[kDim];
337 
338  // set default values
339  for (UInt_t dim=0; dim<kDim; dim++) {
340  xmin[dim] = FLT_MAX;
341  xmax[dim] = FLT_MIN;
342  }
343 
344  Log() << kDEBUG << "Number of training events: " << Data()->GetNTrainingEvents() << Endl;
345  Int_t nevoutside = (Int_t)((Data()->GetNTrainingEvents())*(fFrac)); // number of events that are outside the range
346  Int_t rangehistbins = 10000; // number of bins in histos
347 
348  // loop over all testing singnal and BG events and clac minimal and
349  // maximal value of every variable
350  for (Long64_t i=0; i<(GetNEvents()); i++) { // events loop
351  const Event* ev = GetEvent(i);
352  for (UInt_t dim=0; dim<kDim; dim++) { // variables loop
353  Float_t val;
354  if (fMultiTargetRegression) {
355  if (dim < vDim)
356  val = ev->GetValue(dim);
357  else
358  val = ev->GetTarget(dim-vDim);
359  }
360  else
361  val = ev->GetValue(dim);
362 
363  if (val<xmin[dim])
364  xmin[dim] = val;
365  if (val>xmax[dim])
366  xmax[dim] = val;
367  }
368  }
369 
370  // Create and fill histograms for each dimension (with same events
371  // as before), to determine range based on number of events outside
372  // the range
373  TH1F **range_h = new TH1F*[kDim];
374  for (UInt_t dim=0; dim<kDim; dim++) {
375  range_h[dim] = new TH1F(Form("range%i", dim), "range", rangehistbins, xmin[dim], xmax[dim]);
376  }
377 
378  // fill all testing events into histos
379  for (Long64_t i=0; i<GetNEvents(); i++) {
380  const Event* ev = GetEvent(i);
381  for (UInt_t dim=0; dim<kDim; dim++) {
382  if (fMultiTargetRegression) {
383  if (dim < vDim)
384  range_h[dim]->Fill(ev->GetValue(dim));
385  else
386  range_h[dim]->Fill(ev->GetTarget(dim-vDim));
387  }
388  else
389  range_h[dim]->Fill(ev->GetValue(dim));
390  }
391  }
392 
393  // calc Xmin, Xmax from Histos
394  for (UInt_t dim=0; dim<kDim; dim++) {
395  for (Int_t i=1; i<(rangehistbins+1); i++) { // loop over bins
396  if (range_h[dim]->Integral(0, i) > nevoutside) { // calc left limit (integral over bins 0..i = nevoutside)
397  xmin[dim]=range_h[dim]->GetBinLowEdge(i);
398  break;
399  }
400  }
401  for (Int_t i=rangehistbins; i>0; i--) { // calc right limit (integral over bins i..max = nevoutside)
402  if (range_h[dim]->Integral(i, (rangehistbins+1)) > nevoutside) {
403  xmax[dim]=range_h[dim]->GetBinLowEdge(i+1);
404  break;
405  }
406  }
407  }
408  // now xmin[] and xmax[] contain upper/lower limits for every dimension
409 
410  // copy xmin[], xmax[] values to the class variable
411  fXmin.clear();
412  fXmax.clear();
413  for (UInt_t dim=0; dim<kDim; dim++) {
414  fXmin.push_back(xmin[dim]);
415  fXmax.push_back(xmax[dim]);
416  }
417 
418 
419  delete[] xmin;
420  delete[] xmax;
421 
422  // delete histos
423  for (UInt_t dim=0; dim<kDim; dim++)
424  delete range_h[dim];
425  delete[] range_h;
426 
427  return;
428 }
429 
430 ////////////////////////////////////////////////////////////////////////////////
431 /// Train PDE-Foam depending on the set options
432 
434 {
435  Log() << kVERBOSE << "Calculate Xmin and Xmax for every dimension" << Endl;
436  CalcXminXmax();
437 
438  // delete foams
439  DeleteFoams();
440 
441  // start training
442  if (DoRegression()) {
443  if (fMultiTargetRegression)
444  TrainMultiTargetRegression();
445  else
446  TrainMonoTargetRegression();
447  }
448  else {
449  if (DoMulticlass())
450  TrainMultiClassification();
451  else {
452  if (DataInfo().GetNormalization() != "EQUALNUMEVENTS" ) {
453  Log() << kINFO << "NormMode=" << DataInfo().GetNormalization()
454  << " chosen. Note that only NormMode=EqualNumEvents"
455  << " ensures that Discriminant values correspond to"
456  << " signal probabilities." << Endl;
457  }
458 
459  Log() << kDEBUG << "N_sig for training events: " << Data()->GetNEvtSigTrain() << Endl;
460  Log() << kDEBUG << "N_bg for training events: " << Data()->GetNEvtBkgdTrain() << Endl;
461  Log() << kDEBUG << "User normalization: " << DataInfo().GetNormalization().Data() << Endl;
462 
463  if (fSigBgSeparated)
464  TrainSeparatedClassification();
465  else
466  TrainUnifiedClassification();
467  }
468  }
469 
470  // delete the binary search tree in order to save memory
471  for(UInt_t i=0; i<fFoam.size(); i++) {
472  if(fFoam.at(i))
473  fFoam.at(i)->DeleteBinarySearchTree();
474  }
475 }
476 
477 ////////////////////////////////////////////////////////////////////////////////
478 /// Creation of 2 separated foams: one for signal events, one for
479 /// backgound events. At the end the foam cells of fFoam[0] will
480 /// contain the average number of signal events and fFoam[1] will
481 /// contain the average number of background events.
482 
484 {
485  TString foamcaption[2];
486  foamcaption[0] = "SignalFoam";
487  foamcaption[1] = "BgFoam";
488 
489  for(int i=0; i<2; i++) {
490  // create 2 PDEFoams
491  fFoam.push_back( InitFoam(foamcaption[i], kSeparate) );
492 
493  Log() << kVERBOSE << "Filling binary search tree of " << foamcaption[i]
494  << " with events" << Endl;
495  // insert event to BinarySearchTree
496  for (Long64_t k=0; k<GetNEvents(); ++k) {
497  const Event* ev = GetEvent(k);
498  if ((i==0 && DataInfo().IsSignal(ev)) || (i==1 && !DataInfo().IsSignal(ev)))
499  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
500  fFoam.back()->FillBinarySearchTree(ev);
501  }
502 
503  Log() << kINFO << "Build up " << foamcaption[i] << Endl;
504  fFoam.back()->Create(); // build foam
505 
506  Log() << kVERBOSE << "Filling foam cells with events" << Endl;
507  // loop over all events -> fill foam cells
508  for (Long64_t k=0; k<GetNEvents(); ++k) {
509  const Event* ev = GetEvent(k);
510  Float_t weight = fFillFoamWithOrigWeights ? ev->GetOriginalWeight() : ev->GetWeight();
511  if ((i==0 && DataInfo().IsSignal(ev)) || (i==1 && !DataInfo().IsSignal(ev)))
512  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
513  fFoam.back()->FillFoamCells(ev, weight);
514  }
515  }
516 }
517 
518 ////////////////////////////////////////////////////////////////////////////////
519 /// Create only one unified foam (fFoam[0]) whose cells contain the
520 /// average discriminator (N_sig)/(N_sig + N_bg)
521 
523 {
524  fFoam.push_back( InitFoam("DiscrFoam", kDiscr, fSignalClass) );
525 
526  Log() << kVERBOSE << "Filling binary search tree of discriminator foam with events" << Endl;
527  // insert event to BinarySearchTree
528  for (Long64_t k=0; k<GetNEvents(); ++k) {
529  const Event* ev = GetEvent(k);
530  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
531  fFoam.back()->FillBinarySearchTree(ev);
532  }
533 
534  Log() << kINFO << "Build up discriminator foam" << Endl;
535  fFoam.back()->Create(); // build foam
536 
537  Log() << kVERBOSE << "Filling foam cells with events" << Endl;
538  // loop over all training events -> fill foam cells with N_sig and N_Bg
539  for (Long64_t k=0; k<GetNEvents(); ++k) {
540  const Event* ev = GetEvent(k);
541  Float_t weight = fFillFoamWithOrigWeights ? ev->GetOriginalWeight() : ev->GetWeight();
542  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
543  fFoam.back()->FillFoamCells(ev, weight);
544  }
545 
546  Log() << kVERBOSE << "Calculate cell discriminator"<< Endl;
547  // calc discriminator (and it's error) for each cell
548  fFoam.back()->Finalize();
549 }
550 
551 ////////////////////////////////////////////////////////////////////////////////
552 /// Create one unified foam (see TrainUnifiedClassification()) for
553 /// each class, where the cells of foam i (fFoam[i]) contain the
554 /// average fraction of events of class i, i.e.
555 ///
556 /// D = number events of class i / total number of events
557 
559 {
560  for (UInt_t iClass=0; iClass<DataInfo().GetNClasses(); ++iClass) {
561 
562  fFoam.push_back( InitFoam(Form("MultiClassFoam%u",iClass), kMultiClass, iClass) );
563 
564  Log() << kVERBOSE << "Filling binary search tree of multiclass foam "
565  << iClass << " with events" << Endl;
566  // insert event to BinarySearchTree
567  for (Long64_t k=0; k<GetNEvents(); ++k) {
568  const Event* ev = GetEvent(k);
569  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
570  fFoam.back()->FillBinarySearchTree(ev);
571  }
572 
573  Log() << kINFO << "Build up multiclass foam " << iClass << Endl;
574  fFoam.back()->Create(); // build foam
575 
576  Log() << kVERBOSE << "Filling foam cells with events" << Endl;
577  // loop over all training events and fill foam cells with signal
578  // and background events
579  for (Long64_t k=0; k<GetNEvents(); ++k) {
580  const Event* ev = GetEvent(k);
581  Float_t weight = fFillFoamWithOrigWeights ? ev->GetOriginalWeight() : ev->GetWeight();
582  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
583  fFoam.back()->FillFoamCells(ev, weight);
584  }
585 
586  Log() << kVERBOSE << "Calculate cell discriminator"<< Endl;
587  // calc discriminator (and it's error) for each cell
588  fFoam.back()->Finalize();
589  }
590 }
591 
592 ////////////////////////////////////////////////////////////////////////////////
593 /// Training one (mono target regression) foam, whose cells contain
594 /// the average 0th target. The dimension of the foam = number of
595 /// non-targets (= number of variables).
596 
598 {
599  if (Data()->GetNTargets() != 1) {
600  Log() << kFATAL << "Can't do mono-target regression with "
601  << Data()->GetNTargets() << " targets!" << Endl;
602  }
603 
604  Log() << kDEBUG << "MethodPDEFoam: number of Targets: " << Data()->GetNTargets() << Endl;
605 
606  fFoam.push_back( InitFoam("MonoTargetRegressionFoam", kMonoTarget) );
607 
608  Log() << kVERBOSE << "Filling binary search tree with events" << Endl;
609  // insert event to BinarySearchTree
610  for (Long64_t k=0; k<GetNEvents(); ++k) {
611  const Event* ev = GetEvent(k);
612  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
613  fFoam.back()->FillBinarySearchTree(ev);
614  }
615 
616  Log() << kINFO << "Build mono target regression foam" << Endl;
617  fFoam.back()->Create(); // build foam
618 
619  Log() << kVERBOSE << "Filling foam cells with events" << Endl;
620  // loop over all events -> fill foam cells with target
621  for (Long64_t k=0; k<GetNEvents(); ++k) {
622  const Event* ev = GetEvent(k);
623  Float_t weight = fFillFoamWithOrigWeights ? ev->GetOriginalWeight() : ev->GetWeight();
624  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
625  fFoam.back()->FillFoamCells(ev, weight);
626  }
627 
628  Log() << kVERBOSE << "Calculate average cell targets"<< Endl;
629  // calc weight (and it's error) for each cell
630  fFoam.back()->Finalize();
631 }
632 
633 ////////////////////////////////////////////////////////////////////////////////
634 /// Training one (multi target regression) foam, whose cells contain
635 /// the average event density. The dimension of the foam = number
636 /// of non-targets + number of targets.
637 
639 {
640  Log() << kDEBUG << "Number of variables: " << Data()->GetNVariables() << Endl;
641  Log() << kDEBUG << "Number of Targets: " << Data()->GetNTargets() << Endl;
642  Log() << kDEBUG << "Dimension of foam: " << Data()->GetNVariables()+Data()->GetNTargets() << Endl;
643  if (fKernel==kLinN)
644  Log() << kFATAL << "LinNeighbors kernel currently not supported"
645  << " for multi target regression" << Endl;
646 
647  fFoam.push_back( InitFoam("MultiTargetRegressionFoam", kMultiTarget) );
648 
649  Log() << kVERBOSE << "Filling binary search tree of multi target regression foam with events"
650  << Endl;
651  // insert event to BinarySearchTree
652  for (Long64_t k=0; k<GetNEvents(); ++k) {
653  Event *ev = new Event(*GetEvent(k));
654  // since in multi-target regression targets are handled like
655  // variables --> remove targets and add them to the event variabels
656  std::vector<Float_t> targets(ev->GetTargets());
657  const UInt_t nVariables = ev->GetValues().size();
658  for (UInt_t i = 0; i < targets.size(); ++i)
659  ev->SetVal(i+nVariables, targets.at(i));
660  ev->GetTargets().clear();
661  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
662  fFoam.back()->FillBinarySearchTree(ev);
663  // since the binary search tree copies the event, one can delete
664  // it
665  delete ev;
666  }
667 
668  Log() << kINFO << "Build multi target regression foam" << Endl;
669  fFoam.back()->Create(); // build foam
670 
671  Log() << kVERBOSE << "Filling foam cells with events" << Endl;
672  // loop over all events -> fill foam cells with number of events
673  for (Long64_t k=0; k<GetNEvents(); ++k) {
674  Event *ev = new Event(*GetEvent(k));
675  // since in multi-target regression targets are handled like
676  // variables --> remove targets and add them to the event variabels
677  std::vector<Float_t> targets = ev->GetTargets();
678  const UInt_t nVariables = ev->GetValues().size();
679  Float_t weight = fFillFoamWithOrigWeights ? ev->GetOriginalWeight() : ev->GetWeight();
680  for (UInt_t i = 0; i < targets.size(); ++i)
681  ev->SetVal(i+nVariables, targets.at(i));
682  ev->GetTargets().clear();
683  if (!(IgnoreEventsWithNegWeightsInTraining() && ev->GetWeight()<=0))
684  fFoam.back()->FillFoamCells(ev, weight);
685  // since the PDEFoam copies the event, one can delete it
686  delete ev;
687  }
688 }
689 
690 ////////////////////////////////////////////////////////////////////////////////
691 /// Return Mva-Value.
692 ///
693 /// In case of 'fSigBgSeparated==false' (one unifiend PDEFoam was
694 /// trained) the function returns the content of the cell, which
695 /// corresponds to the current TMVA::Event, i.e. D =
696 /// N_sig/(N_bg+N_sig).
697 ///
698 /// In case of 'fSigBgSeparated==true' (two separate PDEFoams were
699 /// trained) the function returns
700 ///
701 /// D = Density_sig/(Density_sig+Density_bg)
702 ///
703 /// where 'Density_sig' is the content of the cell in the signal
704 /// PDEFoam (fFoam[0]) and 'Density_bg' is the content of the cell
705 /// in the background PDEFoam (fFoam[1]).
706 ///
707 /// In both cases the error on the discriminant is stored in 'err'
708 /// and 'errUpper'. (Of course err and errUpper must be non-zero
709 /// and point to valid address to make this work.)
710 
712 {
713  const Event* ev = GetEvent();
714  Double_t discr = 0.;
715 
716  if (fSigBgSeparated) {
717  std::vector<Float_t> xvec = ev->GetValues();
718 
719  Double_t density_sig = 0.; // calc signal event density
720  Double_t density_bg = 0.; // calc background event density
721  density_sig = fFoam.at(0)->GetCellValue(xvec, kValueDensity, fKernelEstimator);
722  density_bg = fFoam.at(1)->GetCellValue(xvec, kValueDensity, fKernelEstimator);
723 
724  // calc disciminator (normed!)
725  if ( (density_sig+density_bg) > 0 )
726  discr = density_sig/(density_sig+density_bg);
727  else
728  discr = 0.5; // assume 50% signal probability, if no events found (bad assumption, but can be overruled by cut on error)
729  }
730  else { // Signal and Bg not separated
731  // get discriminator direct from the foam
732  discr = fFoam.at(0)->GetCellValue(ev->GetValues(), kValue, fKernelEstimator);
733  }
734 
735  // calculate the error
736  if (err || errUpper) {
737  const Double_t discr_error = CalculateMVAError();
738  if (err != 0) *err = discr_error;
739  if (errUpper != 0) *errUpper = discr_error;
740  }
741 
742  if (fUseYesNoCell)
743  return (discr < 0.5 ? -1 : 1);
744  else
745  return discr;
746 }
747 
748 ////////////////////////////////////////////////////////////////////////////////
749 /// Calculate the error on the Mva value
750 ///
751 /// If fSigBgSeparated == true the error is calculated from the
752 /// number of events in the signal and background PDEFoam cells.
753 ///
754 /// If fSigBgSeparated == false, the error is taken directly from
755 /// the PDEFoam cell.
756 
758 {
759  const Event* ev = GetEvent(); // current event
760  Double_t mvaError = 0.0; // the error on the Mva value
761 
762  if (fSigBgSeparated) {
763  const std::vector<Float_t>& xvec = ev->GetValues();
764 
765  const Double_t neventsB = fFoam.at(1)->GetCellValue(xvec, kValue, fKernelEstimator);
766  const Double_t neventsS = fFoam.at(0)->GetCellValue(xvec, kValue, fKernelEstimator);
767  const Double_t scaleB = 1.;
768  // estimation of statistical error on counted signal/background events
769  const Double_t errorS = neventsS == 0 ? 1.0 : TMath::Sqrt(neventsS);
770  const Double_t errorB = neventsB == 0 ? 1.0 : TMath::Sqrt(neventsB);
771 
772  if ((neventsS > 1e-10) || (neventsB > 1e-10)) {
773  // eq. (5) in paper T.Carli, B.Koblitz 2002
774  mvaError = TMath::Sqrt(Sqr(scaleB * neventsB / Sqr(neventsS + scaleB * neventsB) * errorS) +
775  Sqr(scaleB * neventsS / Sqr(neventsS + scaleB * neventsB) * errorB));
776  } else {
777  mvaError = 1.0;
778  }
779  } else { // Signal and Bg not separated
780  // get discriminator error direct from the foam
781  mvaError = fFoam.at(0)->GetCellValue(ev->GetValues(), kValueError, fKernelEstimator);
782  }
783 
784  return mvaError;
785 }
786 
787 ////////////////////////////////////////////////////////////////////////////////
788 /// Get the multiclass MVA response for the PDEFoam classifier. The
789 /// returned MVA values are normalized, i.e. their sum equals 1.
790 
791 const std::vector<Float_t>& TMVA::MethodPDEFoam::GetMulticlassValues()
792 {
793  const TMVA::Event *ev = GetEvent();
794  std::vector<Float_t> xvec = ev->GetValues();
795 
796  if (fMulticlassReturnVal == NULL)
797  fMulticlassReturnVal = new std::vector<Float_t>();
798  fMulticlassReturnVal->clear();
799  fMulticlassReturnVal->reserve(DataInfo().GetNClasses());
800 
801  std::vector<Float_t> temp; // temp class. values
802  UInt_t nClasses = DataInfo().GetNClasses();
803  temp.reserve(nClasses);
804  for (UInt_t iClass = 0; iClass < nClasses; ++iClass) {
805  temp.push_back(fFoam.at(iClass)->GetCellValue(xvec, kValue, fKernelEstimator));
806  }
807 
808  for (UInt_t iClass = 0; iClass < nClasses; ++iClass) {
809  Float_t norm = 0.0; // normalization
810  for (UInt_t j = 0; j < nClasses; ++j) {
811  if (iClass != j)
812  norm += exp(temp[j] - temp[iClass]);
813  }
814  fMulticlassReturnVal->push_back(1.0 / (1.0 + norm));
815  }
816 
817  return *fMulticlassReturnVal;
818 }
819 
820 ////////////////////////////////////////////////////////////////////////////////
821 /// Compute ranking of input variables from the number of cuts made
822 /// in each PDEFoam dimension. The PDEFoam dimension (the variable)
823 /// for which the most cuts were done is ranked highest.
824 
826 {
827  // create the ranking object
828  fRanking = new Ranking(GetName(), "Variable Importance");
829  std::vector<Float_t> importance(GetNvar(), 0);
830 
831  // determine variable importances
832  for (UInt_t ifoam = 0; ifoam < fFoam.size(); ++ifoam) {
833  // get the number of cuts made in every dimension of foam
834  PDEFoamCell *root_cell = fFoam.at(ifoam)->GetRootCell();
835  std::vector<UInt_t> nCuts(fFoam.at(ifoam)->GetTotDim(), 0);
836  GetNCuts(root_cell, nCuts);
837 
838  // fill the importance vector (ignoring the target dimensions in
839  // case of a multi-target regression foam)
840  UInt_t sumOfCuts = 0;
841  std::vector<Float_t> tmp_importance;
842  for (UInt_t ivar = 0; ivar < GetNvar(); ++ivar) {
843  sumOfCuts += nCuts.at(ivar);
844  tmp_importance.push_back( nCuts.at(ivar) );
845  }
846  // normalization of the variable importances of this foam: the
847  // sum of all variable importances equals 1 for this foam
848  for (UInt_t ivar = 0; ivar < GetNvar(); ++ivar) {
849  if (sumOfCuts > 0)
850  tmp_importance.at(ivar) /= sumOfCuts;
851  else
852  tmp_importance.at(ivar) = 0;
853  }
854  // the overall variable importance is the average over all foams
855  for (UInt_t ivar = 0; ivar < GetNvar(); ++ivar) {
856  importance.at(ivar) += tmp_importance.at(ivar) / fFoam.size();
857  }
858  }
859 
860  // fill ranking vector
861  for (UInt_t ivar = 0; ivar < GetNvar(); ++ivar) {
862  fRanking->AddRank(Rank(GetInputLabel(ivar), importance.at(ivar)));
863  }
864 
865  return fRanking;
866 }
867 
868 ////////////////////////////////////////////////////////////////////////////////
869 /// Fill in 'nCuts' the number of cuts made in every foam dimension,
870 /// starting at the root cell 'cell'.
871 ///
872 /// Parameters:
873 ///
874 /// - cell - root cell to start the counting from
875 ///
876 /// - nCuts - the number of cuts are saved in this vector
877 
878 void TMVA::MethodPDEFoam::GetNCuts(PDEFoamCell *cell, std::vector<UInt_t> &nCuts)
879 {
880  if (cell == NULL || cell->GetStat() == 1) // cell is active
881  return;
882 
883  nCuts.at(cell->GetBest())++;
884 
885  if (cell->GetDau0() != NULL)
886  GetNCuts(cell->GetDau0(), nCuts);
887  if (cell->GetDau1() != NULL)
888  GetNCuts(cell->GetDau1(), nCuts);
889 }
890 
891 ////////////////////////////////////////////////////////////////////////////////
892 /// Set Xmin, Xmax for every dimension in the given pdefoam object
893 
895 {
896  if (!pdefoam){
897  Log() << kFATAL << "Null pointer given!" << Endl;
898  return;
899  }
900 
901  UInt_t num_vars = GetNvar();
902  if (fMultiTargetRegression)
903  num_vars += Data()->GetNTargets();
904 
905  for (UInt_t idim=0; idim<num_vars; idim++) { // set upper/ lower limit in foam
906  Log()<< kDEBUG << "foam: SetXmin[dim="<<idim<<"]: " << fXmin.at(idim) << Endl;
907  Log()<< kDEBUG << "foam: SetXmax[dim="<<idim<<"]: " << fXmax.at(idim) << Endl;
908  pdefoam->SetXmin(idim, fXmin.at(idim));
909  pdefoam->SetXmax(idim, fXmax.at(idim));
910  }
911 }
912 
913 ////////////////////////////////////////////////////////////////////////////////
914 /// Create a new PDEFoam, set the PDEFoam options (nCells, nBin,
915 /// Xmin, Xmax, etc.) and initialize the PDEFoam by calling
916 /// pdefoam->Initialize().
917 ///
918 /// Parameters:
919 ///
920 /// - foamcaption - name of PDEFoam object
921 ///
922 /// - ft - type of PDEFoam
923 /// Candidates are:
924 /// - kSeparate - creates TMVA::PDEFoamEvent
925 /// - kDiscr - creates TMVA::PDEFoamDiscriminant
926 /// - kMonoTarget - creates TMVA::PDEFoamTarget
927 /// - kMultiTarget - creates TMVA::MultiTarget
928 /// - kMultiClass - creates TMVA::PDEFoamDiscriminant
929 ///
930 /// If 'fDTSeparation != kFoam' then a TMVA::PDEFoamDecisionTree
931 /// is created (the separation type depends on fDTSeparation).
932 ///
933 /// - cls - marked event class (optional, default value = 0)
934 
936 {
937  // number of foam dimensions
938  Int_t dim = 1;
939  if (ft == kMultiTarget)
940  // dimension of foam = number of targets + non-targets
941  dim = Data()->GetNTargets() + Data()->GetNVariables();
942  else
943  dim = GetNvar();
944 
945  // calculate range-searching box
946  std::vector<Double_t> box;
947  for (Int_t idim = 0; idim < dim; ++idim) {
948  box.push_back((fXmax.at(idim) - fXmin.at(idim))* fVolFrac);
949  }
950 
951  // create PDEFoam and PDEFoamDensityBase
952  PDEFoam *pdefoam = NULL;
953  PDEFoamDensityBase *density = NULL;
954  if (fDTSeparation == kFoam) {
955  // use PDEFoam algorithm
956  switch (ft) {
957  case kSeparate:
958  pdefoam = new PDEFoamEvent(foamcaption);
959  density = new PDEFoamEventDensity(box);
960  break;
961  case kMultiTarget:
962  pdefoam = new PDEFoamMultiTarget(foamcaption, fTargetSelection);
963  density = new PDEFoamEventDensity(box);
964  break;
965  case kDiscr:
966  case kMultiClass:
967  pdefoam = new PDEFoamDiscriminant(foamcaption, cls);
968  density = new PDEFoamDiscriminantDensity(box, cls);
969  break;
970  case kMonoTarget:
971  pdefoam = new PDEFoamTarget(foamcaption, 0);
972  density = new PDEFoamTargetDensity(box, 0);
973  break;
974  default:
975  Log() << kFATAL << "Unknown PDEFoam type!" << Endl;
976  break;
977  }
978  } else {
979  // create a decision tree like PDEFoam
980 
981  // create separation type class, which is owned by
982  // PDEFoamDecisionTree (i.e. PDEFoamDecisionTree will delete it)
983  SeparationBase *sepType = NULL;
984  switch (fDTSeparation) {
985  case kGiniIndex:
986  sepType = new GiniIndex();
987  break;
989  sepType = new MisClassificationError();
990  break;
991  case kCrossEntropy:
992  sepType = new CrossEntropy();
993  break;
995  sepType = new GiniIndexWithLaplace();
996  break;
997  case kSdivSqrtSplusB:
998  sepType = new SdivSqrtSplusB();
999  break;
1000  default:
1001  Log() << kFATAL << "Separation type " << fDTSeparation
1002  << " currently not supported" << Endl;
1003  break;
1004  }
1005  switch (ft) {
1006  case kDiscr:
1007  case kMultiClass:
1008  pdefoam = new PDEFoamDecisionTree(foamcaption, sepType, cls);
1009  density = new PDEFoamDecisionTreeDensity(box, cls);
1010  break;
1011  default:
1012  Log() << kFATAL << "Decision tree cell split algorithm is only"
1013  << " available for (multi) classification with a single"
1014  << " PDE-Foam (SigBgSeparate=F)" << Endl;
1015  break;
1016  }
1017  }
1018 
1019  if (pdefoam) pdefoam->SetDensity(density);
1020  else Log() << kFATAL << "PDEFoam pointer not set, exiting.." << Endl;
1021 
1022  // create pdefoam kernel
1023  fKernelEstimator = CreatePDEFoamKernel();
1024 
1025  // set fLogger attributes
1026  pdefoam->Log().SetMinType(this->Log().GetMinType());
1027 
1028  // set PDEFoam parameters
1029  pdefoam->SetDim( dim);
1030  pdefoam->SetnCells( fnCells); // optional
1031  pdefoam->SetnSampl( fnSampl); // optional
1032  pdefoam->SetnBin( fnBin); // optional
1033  pdefoam->SetEvPerBin( fEvPerBin); // optional
1034 
1035  // cuts
1036  pdefoam->SetNmin(fNmin);
1037  pdefoam->SetMaxDepth(fMaxDepth); // maximum cell tree depth
1038 
1039  // Init PDEFoam
1040  pdefoam->Initialize();
1041 
1042  // Set Xmin, Xmax
1043  SetXminXmax(pdefoam);
1044 
1045  return pdefoam;
1046 }
1047 
1048 ////////////////////////////////////////////////////////////////////////////////
1049 /// Return regression values for both multi- and mono-target regression
1050 
1051 const std::vector<Float_t>& TMVA::MethodPDEFoam::GetRegressionValues()
1052 {
1053  if (fRegressionReturnVal == 0) fRegressionReturnVal = new std::vector<Float_t>();
1054  fRegressionReturnVal->clear();
1055  fRegressionReturnVal->reserve(Data()->GetNTargets());
1056 
1057  const Event* ev = GetEvent();
1058  std::vector<Float_t> vals = ev->GetValues(); // get array of event variables (non-targets)
1059 
1060  if (vals.empty()) {
1061  Log() << kWARNING << "<GetRegressionValues> value vector is empty. " << Endl;
1062  }
1063 
1064  if (fMultiTargetRegression) {
1065  // create std::map from event variables
1066  std::map<Int_t, Float_t> xvec;
1067  for (UInt_t i=0; i<vals.size(); ++i)
1068  xvec.insert(std::pair<Int_t, Float_t>(i, vals.at(i)));
1069  // get the targets
1070  std::vector<Float_t> targets = fFoam.at(0)->GetCellValue( xvec, kValue );
1071 
1072  // sanity check
1073  if (targets.size() != Data()->GetNTargets())
1074  Log() << kFATAL << "Something wrong with multi-target regression foam: "
1075  << "number of targets does not match the DataSet()" << Endl;
1076  for(UInt_t i=0; i<targets.size(); i++)
1077  fRegressionReturnVal->push_back(targets.at(i));
1078  }
1079  else {
1080  fRegressionReturnVal->push_back(fFoam.at(0)->GetCellValue(vals, kValue, fKernelEstimator));
1081  }
1082 
1083  // apply inverse transformation to regression values
1084  Event * evT = new Event(*ev);
1085  for (UInt_t itgt = 0; itgt < Data()->GetNTargets(); itgt++) {
1086  evT->SetTarget(itgt, fRegressionReturnVal->at(itgt) );
1087  }
1088  const Event* evT2 = GetTransformationHandler().InverseTransform( evT );
1089  fRegressionReturnVal->clear();
1090  for (UInt_t itgt = 0; itgt < Data()->GetNTargets(); itgt++) {
1091  fRegressionReturnVal->push_back( evT2->GetTarget(itgt) );
1092  }
1093 
1094  delete evT;
1095 
1096  return (*fRegressionReturnVal);
1097 }
1098 
1099 ////////////////////////////////////////////////////////////////////////////////
1100 /// create a pdefoam kernel estimator, depending on the current
1101 /// value of fKernel
1102 
1104 {
1105  switch (fKernel) {
1106  case kNone:
1107  return new PDEFoamKernelTrivial();
1108  case kLinN:
1109  return new PDEFoamKernelLinN();
1110  case kGaus:
1111  return new PDEFoamKernelGauss(fVolFrac/2.0);
1112  default:
1113  Log() << kFATAL << "Kernel: " << fKernel << " not supported!" << Endl;
1114  return NULL;
1115  }
1116  return NULL;
1117 }
1118 
1119 ////////////////////////////////////////////////////////////////////////////////
1120 /// Deletes all trained foams
1121 
1123 {
1124  for (UInt_t i=0; i<fFoam.size(); i++)
1125  if (fFoam.at(i)) delete fFoam.at(i);
1126  fFoam.clear();
1127 }
1128 
1129 ////////////////////////////////////////////////////////////////////////////////
1130 /// reset MethodPDEFoam:
1131 /// - delete all PDEFoams
1132 /// - delete the kernel estimator
1133 
1135 {
1136  DeleteFoams();
1137 
1138  if (fKernelEstimator != NULL) {
1139  delete fKernelEstimator;
1140  fKernelEstimator = NULL;
1141  }
1142 }
1143 
1144 ////////////////////////////////////////////////////////////////////////////////
1145 
1147 {}
1148 
1149 ////////////////////////////////////////////////////////////////////////////////
1150 /// create XML output of PDEFoam method variables
1151 
1152 void TMVA::MethodPDEFoam::AddWeightsXMLTo( void* parent ) const
1153 {
1154  void* wght = gTools().AddChild(parent, "Weights");
1155  gTools().AddAttr( wght, "SigBgSeparated", fSigBgSeparated );
1156  gTools().AddAttr( wght, "Frac", fFrac );
1157  gTools().AddAttr( wght, "DiscrErrCut", fDiscrErrCut );
1158  gTools().AddAttr( wght, "VolFrac", fVolFrac );
1159  gTools().AddAttr( wght, "nCells", fnCells );
1160  gTools().AddAttr( wght, "nSampl", fnSampl );
1161  gTools().AddAttr( wght, "nBin", fnBin );
1162  gTools().AddAttr( wght, "EvPerBin", fEvPerBin );
1163  gTools().AddAttr( wght, "Compress", fCompress );
1164  gTools().AddAttr( wght, "DoRegression", DoRegression() );
1165  gTools().AddAttr( wght, "CutNmin", fNmin>0 );
1166  gTools().AddAttr( wght, "Nmin", fNmin );
1167  gTools().AddAttr( wght, "CutRMSmin", false );
1168  gTools().AddAttr( wght, "RMSmin", 0.0 );
1169  gTools().AddAttr( wght, "Kernel", KernelToUInt(fKernel) );
1170  gTools().AddAttr( wght, "TargetSelection", TargetSelectionToUInt(fTargetSelection) );
1171  gTools().AddAttr( wght, "FillFoamWithOrigWeights", fFillFoamWithOrigWeights );
1172  gTools().AddAttr( wght, "UseYesNoCell", fUseYesNoCell );
1173 
1174  // save foam borders Xmin[i], Xmax[i]
1175  void *xmin_wrap;
1176  for (UInt_t i=0; i<fXmin.size(); i++){
1177  xmin_wrap = gTools().AddChild( wght, "Xmin" );
1178  gTools().AddAttr( xmin_wrap, "Index", i );
1179  gTools().AddAttr( xmin_wrap, "Value", fXmin.at(i) );
1180  }
1181  void *xmax_wrap;
1182  for (UInt_t i=0; i<fXmax.size(); i++){
1183  xmax_wrap = gTools().AddChild( wght, "Xmax" );
1184  gTools().AddAttr( xmax_wrap, "Index", i );
1185  gTools().AddAttr( xmax_wrap, "Value", fXmax.at(i) );
1186  }
1187 
1188  // write foams to xml file
1189  WriteFoamsToFile();
1190 }
1191 
1192 ////////////////////////////////////////////////////////////////////////////////
1193 /// Write PDEFoams to file
1194 
1196 {
1197  // fill variable names into foam
1198  FillVariableNamesToFoam();
1199 
1200  TString rfname( GetWeightFileName() );
1201 
1202  // replace in case of txt weight file
1203  rfname.ReplaceAll( TString(".") + gConfig().GetIONames().fWeightFileExtension + ".txt", ".xml" );
1204 
1205  // add foam indicator to distinguish from main weight file
1206  rfname.ReplaceAll( ".xml", "_foams.root" );
1207 
1208  TFile *rootFile = 0;
1209  if (fCompress) rootFile = new TFile(rfname, "RECREATE", "foamfile", 9);
1210  else rootFile = new TFile(rfname, "RECREATE");
1211 
1212  // write the foams
1213  for (UInt_t i=0; i<fFoam.size(); ++i) {
1214  Log() << "writing foam " << fFoam.at(i)->GetFoamName().Data()
1215  << " to file" << Endl;
1216  fFoam.at(i)->Write(fFoam.at(i)->GetFoamName().Data());
1217  }
1218 
1219  rootFile->Close();
1220  Log() << kINFO << "Foams written to file: "
1221  << gTools().Color("lightblue") << rfname << gTools().Color("reset") << Endl;
1222 }
1223 
1224 ////////////////////////////////////////////////////////////////////////////////
1225 /// read options and internal parameters
1226 
1228 {
1229  istr >> fSigBgSeparated; // Seperate Sig and Bg, or not
1230  istr >> fFrac; // Fraction used for calc of Xmin, Xmax
1231  istr >> fDiscrErrCut; // cut on discrimant error
1232  istr >> fVolFrac; // volume fraction (used for density calculation during buildup)
1233  istr >> fnCells; // Number of Cells (500)
1234  istr >> fnSampl; // Number of MC events per cell in build-up (1000)
1235  istr >> fnBin; // Number of bins in build-up (100)
1236  istr >> fEvPerBin; // Maximum events (equiv.) per bin in buid-up (1000)
1237  istr >> fCompress; // compress output file
1238 
1239  Bool_t regr;
1240  istr >> regr; // regression foam
1241  SetAnalysisType( (regr ? Types::kRegression : Types::kClassification ) );
1242 
1243  Bool_t CutNmin, CutRMSmin; // dummy for backwards compatib.
1244  Float_t RMSmin; // dummy for backwards compatib.
1245  istr >> CutNmin; // cut on minimal number of events in cell
1246  istr >> fNmin;
1247  istr >> CutRMSmin; // cut on minimal RMS in cell
1248  istr >> RMSmin;
1249 
1250  UInt_t ker = 0;
1251  istr >> ker; // used kernel for GetMvaValue()
1252  fKernel = UIntToKernel(ker);
1253 
1254  UInt_t ts = 0;
1255  istr >> ts; // used method for target selection
1256  fTargetSelection = UIntToTargetSelection(ts);
1257 
1258  istr >> fFillFoamWithOrigWeights; // fill foam with original event weights
1259  istr >> fUseYesNoCell; // return -1 or 1 for bg or signal event
1260 
1261  // clear old range and prepare new range
1262  fXmin.clear();
1263  fXmax.clear();
1264  UInt_t kDim = GetNvar();
1265  if (fMultiTargetRegression)
1266  kDim += Data()->GetNTargets();
1267  fXmin.assign(kDim, 0);
1268  fXmax.assign(kDim, 0);
1269 
1270  // read range
1271  for (UInt_t i=0; i<kDim; i++)
1272  istr >> fXmin.at(i);
1273  for (UInt_t i=0; i<kDim; i++)
1274  istr >> fXmax.at(i);
1275 
1276  // read pure foams from file
1277  ReadFoamsFromFile();
1278 }
1279 
1280 ////////////////////////////////////////////////////////////////////////////////
1281 /// read PDEFoam variables from xml weight file
1282 
1284 {
1285  gTools().ReadAttr( wghtnode, "SigBgSeparated", fSigBgSeparated );
1286  gTools().ReadAttr( wghtnode, "Frac", fFrac );
1287  gTools().ReadAttr( wghtnode, "DiscrErrCut", fDiscrErrCut );
1288  gTools().ReadAttr( wghtnode, "VolFrac", fVolFrac );
1289  gTools().ReadAttr( wghtnode, "nCells", fnCells );
1290  gTools().ReadAttr( wghtnode, "nSampl", fnSampl );
1291  gTools().ReadAttr( wghtnode, "nBin", fnBin );
1292  gTools().ReadAttr( wghtnode, "EvPerBin", fEvPerBin );
1293  gTools().ReadAttr( wghtnode, "Compress", fCompress );
1294  Bool_t regr; // dummy for backwards compatib.
1295  gTools().ReadAttr( wghtnode, "DoRegression", regr );
1296  Bool_t CutNmin; // dummy for backwards compatib.
1297  gTools().ReadAttr( wghtnode, "CutNmin", CutNmin );
1298  gTools().ReadAttr( wghtnode, "Nmin", fNmin );
1299  Bool_t CutRMSmin; // dummy for backwards compatib.
1300  Float_t RMSmin; // dummy for backwards compatib.
1301  gTools().ReadAttr( wghtnode, "CutRMSmin", CutRMSmin );
1302  gTools().ReadAttr( wghtnode, "RMSmin", RMSmin );
1303  UInt_t ker = 0;
1304  gTools().ReadAttr( wghtnode, "Kernel", ker );
1305  fKernel = UIntToKernel(ker);
1306  UInt_t ts = 0;
1307  gTools().ReadAttr( wghtnode, "TargetSelection", ts );
1308  fTargetSelection = UIntToTargetSelection(ts);
1309  if (gTools().HasAttr(wghtnode, "FillFoamWithOrigWeights"))
1310  gTools().ReadAttr( wghtnode, "FillFoamWithOrigWeights", fFillFoamWithOrigWeights );
1311  if (gTools().HasAttr(wghtnode, "UseYesNoCell"))
1312  gTools().ReadAttr( wghtnode, "UseYesNoCell", fUseYesNoCell );
1313 
1314  // clear old range [Xmin, Xmax] and prepare new range for reading
1315  fXmin.clear();
1316  fXmax.clear();
1317  UInt_t kDim = GetNvar();
1318  if (fMultiTargetRegression)
1319  kDim += Data()->GetNTargets();
1320  fXmin.assign(kDim, 0);
1321  fXmax.assign(kDim, 0);
1322 
1323  // read foam range
1324  void *xmin_wrap = gTools().GetChild( wghtnode );
1325  for (UInt_t counter=0; counter<kDim; counter++) {
1326  UInt_t i=0;
1327  gTools().ReadAttr( xmin_wrap , "Index", i );
1328  if (i>=kDim)
1329  Log() << kFATAL << "dimension index out of range:" << i << Endl;
1330  gTools().ReadAttr( xmin_wrap , "Value", fXmin.at(i) );
1331  xmin_wrap = gTools().GetNextChild( xmin_wrap );
1332  }
1333 
1334  void *xmax_wrap = xmin_wrap;
1335  for (UInt_t counter=0; counter<kDim; counter++) {
1336  UInt_t i=0;
1337  gTools().ReadAttr( xmax_wrap , "Index", i );
1338  if (i>=kDim)
1339  Log() << kFATAL << "dimension index out of range:" << i << Endl;
1340  gTools().ReadAttr( xmax_wrap , "Value", fXmax.at(i) );
1341  xmax_wrap = gTools().GetNextChild( xmax_wrap );
1342  }
1343 
1344  // if foams exist, delete them
1345  DeleteFoams();
1346 
1347  // read pure foams from file
1348  ReadFoamsFromFile();
1349 
1350  // recreate the pdefoam kernel estimator
1351  if (fKernelEstimator != NULL)
1352  delete fKernelEstimator;
1353  fKernelEstimator = CreatePDEFoamKernel();
1354 }
1355 
1356 ////////////////////////////////////////////////////////////////////////////////
1357 /// Reads a foam with name 'foamname' from file, and returns a clone
1358 /// of the foam. The given ROOT file must be open. (The ROOT file
1359 /// will not be closed in this function.)
1360 ///
1361 /// Parameters:
1362 ///
1363 /// - file - an open ROOT file
1364 ///
1365 /// - foamname - name of foam to load from the file
1366 ///
1367 /// Returns:
1368 ///
1369 /// If a foam with name 'foamname' exists in the file, then it is
1370 /// read from the file, cloned and returned. If a foam with name
1371 /// 'foamname' does not exist in the file or the clone operation
1372 /// does not succeed, then NULL is returned.
1373 
1375 {
1376  if (file == NULL) {
1377  Log() << kWARNING << "<ReadClonedFoamFromFile>: NULL pointer given" << Endl;
1378  return NULL;
1379  }
1380 
1381  // try to load the foam from the file
1382  PDEFoam *foam = (PDEFoam*) file->Get(foamname);
1383  if (foam == NULL) {
1384  return NULL;
1385  }
1386  // try to clone the foam
1387  foam = (PDEFoam*) foam->Clone();
1388  if (foam == NULL) {
1389  Log() << kWARNING << "<ReadClonedFoamFromFile>: " << foamname
1390  << " could not be cloned!" << Endl;
1391  return NULL;
1392  }
1393 
1394  return foam;
1395 }
1396 
1397 ////////////////////////////////////////////////////////////////////////////////
1398 /// read foams from file
1399 
1401 {
1402  TString rfname( GetWeightFileName() );
1403 
1404  // replace in case of txt weight file
1405  rfname.ReplaceAll( TString(".") + gConfig().GetIONames().fWeightFileExtension + ".txt", ".xml" );
1406 
1407  // add foam indicator to distinguish from main weight file
1408  rfname.ReplaceAll( ".xml", "_foams.root" );
1409 
1410  Log() << kINFO << "Read foams from file: " << gTools().Color("lightblue")
1411  << rfname << gTools().Color("reset") << Endl;
1412  TFile *rootFile = new TFile( rfname, "READ" );
1413  if (rootFile->IsZombie()) Log() << kFATAL << "Cannot open file \"" << rfname << "\"" << Endl;
1414 
1415  // read foams from file
1416  if (DoRegression()) {
1417  if (fMultiTargetRegression)
1418  fFoam.push_back(ReadClonedFoamFromFile(rootFile, "MultiTargetRegressionFoam"));
1419  else
1420  fFoam.push_back(ReadClonedFoamFromFile(rootFile, "MonoTargetRegressionFoam"));
1421  } else {
1422  if (fSigBgSeparated) {
1423  fFoam.push_back(ReadClonedFoamFromFile(rootFile, "SignalFoam"));
1424  fFoam.push_back(ReadClonedFoamFromFile(rootFile, "BgFoam"));
1425  } else {
1426  // try to load discriminator foam
1427  PDEFoam *foam = ReadClonedFoamFromFile(rootFile, "DiscrFoam");
1428  if (foam != NULL)
1429  fFoam.push_back(foam);
1430  else {
1431  // load multiclass foams
1432  for (UInt_t iClass=0; iClass<DataInfo().GetNClasses(); ++iClass) {
1433  fFoam.push_back(ReadClonedFoamFromFile(rootFile, Form("MultiClassFoam%u",iClass)));
1434  }
1435  }
1436  }
1437  }
1438 
1439  // Close the root file. Note, that the foams are still present in
1440  // memory!
1441  rootFile->Close();
1442  delete rootFile;
1443 
1444  for (UInt_t i=0; i<fFoam.size(); ++i) {
1445  if (!fFoam.at(0))
1446  Log() << kFATAL << "Could not load foam!" << Endl;
1447  }
1448 }
1449 
1450 ////////////////////////////////////////////////////////////////////////////////
1451 /// convert UInt_t to EKernel (used for reading weight files)
1452 
1454 {
1455  switch(iker) {
1456  case 0: return kNone;
1457  case 1: return kGaus;
1458  case 2: return kLinN;
1459  default:
1460  Log() << kWARNING << "<UIntToKernel>: unknown kernel number: " << iker << Endl;
1461  return kNone;
1462  }
1463  return kNone;
1464 }
1465 
1466 ////////////////////////////////////////////////////////////////////////////////
1467 /// convert UInt_t to ETargetSelection (used for reading weight files)
1468 
1470 {
1471  switch(its) {
1472  case 0: return kMean;
1473  case 1: return kMpv;
1474  default:
1475  Log() << kWARNING << "<UIntToTargetSelection>: unknown method TargetSelection: " << its << Endl;
1476  return kMean;
1477  }
1478  return kMean;
1479 }
1480 
1481 ////////////////////////////////////////////////////////////////////////////////
1482 /// store the variable names in all foams
1483 
1485 {
1486  for (UInt_t ifoam=0; ifoam<fFoam.size(); ifoam++) {
1487  for (Int_t idim=0; idim<fFoam.at(ifoam)->GetTotDim(); idim++) {
1488  if(fMultiTargetRegression && (UInt_t)idim>=DataInfo().GetNVariables())
1489  fFoam.at(ifoam)->AddVariableName(DataInfo().GetTargetInfo(idim-DataInfo().GetNVariables()).GetExpression().Data());
1490  else
1491  fFoam.at(ifoam)->AddVariableName(DataInfo().GetVariableInfo(idim).GetExpression().Data());
1492  }
1493  }
1494 }
1495 
1496 ////////////////////////////////////////////////////////////////////////////////
1497 /// write PDEFoam-specific classifier response
1498 /// NOT IMPLEMENTED YET!
1499 
1500 void TMVA::MethodPDEFoam::MakeClassSpecific( std::ostream& /*fout*/, const TString& /*className*/ ) const
1501 {
1502 }
1503 
1504 ////////////////////////////////////////////////////////////////////////////////
1505 /// provide help message
1506 
1508 {
1509  Log() << Endl;
1510  Log() << gTools().Color("bold") << "--- Short description:" << gTools().Color("reset") << Endl;
1511  Log() << Endl;
1512  Log() << "PDE-Foam is a variation of the PDE-RS method using a self-adapting" << Endl;
1513  Log() << "binning method to divide the multi-dimensional variable space into a" << Endl;
1514  Log() << "finite number of hyper-rectangles (cells). The binning algorithm " << Endl;
1515  Log() << "adjusts the size and position of a predefined number of cells such" << Endl;
1516  Log() << "that the variance of the signal and background densities inside the " << Endl;
1517  Log() << "cells reaches a minimum" << Endl;
1518  Log() << Endl;
1519  Log() << gTools().Color("bold") << "--- Use of booking options:" << gTools().Color("reset") << Endl;
1520  Log() << Endl;
1521  Log() << "The PDEFoam classifier supports two different algorithms: " << Endl;
1522  Log() << Endl;
1523  Log() << " (1) Create one foam, which stores the signal over background" << Endl;
1524  Log() << " probability density. During foam buildup the variance of the" << Endl;
1525  Log() << " discriminant inside the cells is minimised." << Endl;
1526  Log() << Endl;
1527  Log() << " Booking option: SigBgSeparated=F" << Endl;
1528  Log() << Endl;
1529  Log() << " (2) Create two separate foams, one for the signal events and one for" << Endl;
1530  Log() << " background events. During foam buildup the variance of the" << Endl;
1531  Log() << " event density inside the cells is minimised separately for" << Endl;
1532  Log() << " signal and background." << Endl;
1533  Log() << Endl;
1534  Log() << " Booking option: SigBgSeparated=T" << Endl;
1535  Log() << Endl;
1536  Log() << "The following options can be set (the listed values are found to be a" << Endl;
1537  Log() << "good starting point for most applications):" << Endl;
1538  Log() << Endl;
1539  Log() << " SigBgSeparate False Separate Signal and Background" << Endl;
1540  Log() << " TailCut 0.001 Fraction of outlier events that excluded" << Endl;
1541  Log() << " from the foam in each dimension " << Endl;
1542  Log() << " VolFrac 0.0666 Volume fraction (used for density calculation" << Endl;
1543  Log() << " during foam build-up) " << Endl;
1544  Log() << " nActiveCells 500 Maximal number of active cells in final foam " << Endl;
1545  Log() << " nSampl 2000 Number of MC events per cell in foam build-up " << Endl;
1546  Log() << " nBin 5 Number of bins used in foam build-up " << Endl;
1547  Log() << " Nmin 100 Number of events in cell required to split cell" << Endl;
1548  Log() << " Kernel None Kernel type used (possible valuses are: None," << Endl;
1549  Log() << " Gauss)" << Endl;
1550  Log() << " Compress True Compress foam output file " << Endl;
1551  Log() << Endl;
1552  Log() << " Additional regression options:" << Endl;
1553  Log() << Endl;
1554  Log() << "MultiTargetRegression False Do regression with multiple targets " << Endl;
1555  Log() << " TargetSelection Mean Target selection method (possible valuses are: " << Endl;
1556  Log() << " Mean, Mpv)" << Endl;
1557  Log() << Endl;
1558  Log() << gTools().Color("bold") << "--- Performance optimisation:" << gTools().Color("reset") << Endl;
1559  Log() << Endl;
1560  Log() << "The performance of the two implementations was found to be similar for" << Endl;
1561  Log() << "most examples studied. For the same number of cells per foam, the two-" << Endl;
1562  Log() << "foam option approximately doubles the amount of computer memory needed" << Endl;
1563  Log() << "during classification. For special cases where the event-density" << Endl;
1564  Log() << "distribution of signal and background events is very different, the" << Endl;
1565  Log() << "two-foam option was found to perform significantly better than the" << Endl;
1566  Log() << "option with only one foam." << Endl;
1567  Log() << Endl;
1568  Log() << "In order to gain better classification performance we recommend to set" << Endl;
1569  Log() << "the parameter \"nActiveCells\" to a high value." << Endl;
1570  Log() << Endl;
1571  Log() << "The parameter \"VolFrac\" specifies the size of the sampling volume" << Endl;
1572  Log() << "during foam buildup and should be tuned in order to achieve optimal" << Endl;
1573  Log() << "performance. A larger box leads to a reduced statistical uncertainty" << Endl;
1574  Log() << "for small training samples and to smoother sampling. A smaller box on" << Endl;
1575  Log() << "the other hand increases the sensitivity to statistical fluctuations" << Endl;
1576  Log() << "in the training samples, but for sufficiently large training samples" << Endl;
1577  Log() << "it will result in a more precise local estimate of the sampled" << Endl;
1578  Log() << "density. In general, higher dimensional problems require larger box" << Endl;
1579  Log() << "sizes, due to the reduced average number of events per box volume. The" << Endl;
1580  Log() << "default value of 0.0666 was optimised for an example with 5" << Endl;
1581  Log() << "observables and training samples of the order of 50000 signal and" << Endl;
1582  Log() << "background events each." << Endl;
1583  Log() << Endl;
1584  Log() << "Furthermore kernel weighting can be activated, which will lead to an" << Endl;
1585  Log() << "additional performance improvement. Note that Gauss weighting will" << Endl;
1586  Log() << "significantly increase the response time of the method. LinNeighbors" << Endl;
1587  Log() << "weighting performs a linear interpolation with direct neighbor cells" << Endl;
1588  Log() << "for each dimension and is much faster than Gauss weighting." << Endl;
1589  Log() << Endl;
1590  Log() << "The classification results were found to be rather insensitive to the" << Endl;
1591  Log() << "values of the parameters \"nSamples\" and \"nBin\"." << Endl;
1592 }
void Train(void)
Train PDE-Foam depending on the set options.
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:823
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
float xmin
Definition: THbookFile.cxx:93
virtual void Reset()
reset MethodPDEFoam:
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
long long Long64_t
Definition: RtypesCore.h:69
void GetNCuts(PDEFoamCell *cell, std::vector< UInt_t > &nCuts)
Fill in 'nCuts' the number of cuts made in every foam dimension, starting at the root cell 'cell'...
float Float_t
Definition: RtypesCore.h:53
PDEFoam * InitFoam(TString, EFoamType, UInt_t cls=0)
Create a new PDEFoam, set the PDEFoam options (nCells, nBin, Xmin, Xmax, etc.) and initialize the PDE...
MsgLogger & Log() const
Definition: PDEFoam.h:262
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
void PrintCoefficients(void)
Config & gConfig()
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
EAnalysisType
Definition: Types.h:124
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
Bool_t IsZombie() const
Definition: TObject.h:141
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
void SetXmin(Int_t idim, Double_t wmin)
set lower foam bound in dimension idim
Definition: PDEFoam.cxx:276
ClassImp(TMVA::MethodPDEFoam) TMVA
init PDEFoam objects
int Int_t
Definition: RtypesCore.h:41
void TrainUnifiedClassification(void)
Create only one unified foam (fFoam[0]) whose cells contain the average discriminator (N_sig)/(N_sig ...
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
void SetTarget(UInt_t itgt, Float_t value)
set the target value (dimension itgt) to value
Definition: Event.cxx:354
Double_t GetWeight() const
return the event weight - depending on whether the flag IgnoreNegWeightsInTraining is or not...
Definition: Event.cxx:376
Int_t GetStat() const
Definition: PDEFoamCell.h:97
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition: fillpatterns.C:1
void AddAttr(void *node, const char *, const T &value, Int_t precision=16)
Definition: Tools.h:308
void SetVal(UInt_t ivar, Float_t val)
set variable ivar to val
Definition: Event.cxx:335
void * AddChild(void *parent, const char *childname, const char *content=0, bool isRootNode=false)
add child node
Definition: Tools.cxx:1134
Tools & gTools()
Definition: Tools.cxx:79
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TObject.cxx:203
void GetHelpMessage() const
provide help message
virtual Double_t GetBinLowEdge(Int_t bin) const
return bin lower edge for 1D historam Better to use h1.GetXaxis().GetBinLowEdge(bin) ...
Definition: TH1.cxx:8481
void ReadWeightsFromStream(std::istream &i)
read options and internal parameters
virtual ~MethodPDEFoam(void)
destructor
void * GetChild(void *parent, const char *childname=0)
get child node
Definition: Tools.cxx:1158
std::vector< std::vector< double > > Data
EFoamType
Definition: PDEFoam.h:74
void SetMinType(EMsgType minType)
Definition: MsgLogger.h:76
void DeclareOptions()
Declare MethodPDEFoam options.
Double_t GetOriginalWeight() const
Definition: Event.h:84
void SetXminXmax(TMVA::PDEFoam *)
Set Xmin, Xmax for every dimension in the given pdefoam object.
PDEFoam * ReadClonedFoamFromFile(TFile *, const TString &)
Reads a foam with name 'foamname' from file, and returns a clone of the foam.
void SetMaxDepth(UInt_t maxdepth)
Definition: PDEFoam.h:229
void WriteFoamsToFile() const
Write PDEFoams to file.
std::vector< Float_t > & GetTargets()
Definition: Event.h:102
void SetnSampl(Long_t nSampl)
Definition: PDEFoam.h:212
#define None
Definition: TGWin32.h:68
void SetNmin(UInt_t val)
Definition: PDEFoam.h:227
void SetXmax(Int_t idim, Double_t wmax)
set upper foam bound in dimension idim
Definition: PDEFoam.cxx:287
Double_t CalculateMVAError()
Calculate the error on the Mva value.
void Init(void)
default initialization called by all constructors
void CalcXminXmax()
Determine foam range [fXmin, fXmax] for all dimensions, such that a fraction of 'fFrac' events lie ou...
virtual Bool_t HasAnalysisType(Types::EAnalysisType type, UInt_t numberClasses, UInt_t numberTargets)
PDEFoam can handle classification with multiple classes and regression with one or more regression-ta...
Double_t Mean(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:824
EKernel UIntToKernel(UInt_t iker)
convert UInt_t to EKernel (used for reading weight files)
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
const Handle_t kNone
Definition: GuiTypes.h:89
ETargetSelection UIntToTargetSelection(UInt_t its)
convert UInt_t to ETargetSelection (used for reading weight files)
Double_t GetMvaValue(Double_t *err=0, Double_t *errUpper=0)
Return Mva-Value.
void TrainMultiTargetRegression(void)
Training one (multi target regression) foam, whose cells contain the average event density...
void SetDensity(PDEFoamDensityBase *dens)
Definition: PDEFoam.h:216
MethodPDEFoam(const TString &jobName, const TString &methodTitle, DataSetInfo &dsi, const TString &theOption="PDEFoam", TDirectory *theTargetDir=0)
void ReadAttr(void *node, const char *, T &value)
Definition: Tools.h:295
float xmax
Definition: THbookFile.cxx:93
void SetDim(Int_t kDim)
Sets dimension of cubical space.
Definition: PDEFoam.cxx:261
void FillVariableNamesToFoam() const
store the variable names in all foams
void TrainSeparatedClassification(void)
Creation of 2 separated foams: one for signal events, one for backgound events.
Int_t GetBest() const
Definition: PDEFoamCell.h:84
tuple file
Definition: fildir.py:20
void ReadWeightsFromXML(void *wghtnode)
read PDEFoam variables from xml weight file
double Double_t
Definition: RtypesCore.h:55
void TrainMultiClassification()
Create one unified foam (see TrainUnifiedClassification()) for each class, where the cells of foam i ...
Describe directory structure in memory.
Definition: TDirectory.h:44
void SetnCells(Long_t nCells)
Definition: PDEFoam.h:211
int type
Definition: TGX11.cxx:120
void * GetNextChild(void *prevchild, const char *childname=0)
XML helpers.
Definition: Tools.cxx:1170
Float_t GetValue(UInt_t ivar) const
return value of i'th variable
Definition: Event.cxx:231
PDEFoamKernelBase * CreatePDEFoamKernel()
create a pdefoam kernel estimator, depending on the current value of fKernel
void TrainMonoTargetRegression(void)
Training one (mono target regression) foam, whose cells contain the average 0th target.
void Initialize()
Definition: PDEFoam.h:195
const TString & Color(const TString &)
human readable color strings
Definition: Tools.cxx:837
Float_t GetTarget(UInt_t itgt) const
Definition: Event.h:101
const std::vector< Float_t > & GetMulticlassValues()
Get the multiclass MVA response for the PDEFoam classifier.
#define REGISTER_METHOD(CLASS)
for example
const Ranking * CreateRanking()
Compute ranking of input variables from the number of cuts made in each PDEFoam dimension.
void DeleteFoams()
Deletes all trained foams.
std::vector< Float_t > & GetValues()
Definition: Event.h:93
PDEFoamCell * GetDau1() const
Definition: PDEFoamCell.h:101
virtual void DeclareCompatibilityOptions()
options that are used ONLY for the READER to ensure backward compatibility they are hence without any...
Definition: MethodBase.cxx:606
void SetnBin(Int_t nBin)
Definition: PDEFoam.h:213
void SetEvPerBin(Int_t EvPerBin)
Definition: PDEFoam.h:214
#define NULL
Definition: Rtypes.h:82
virtual const std::vector< Float_t > & GetRegressionValues()
Return regression values for both multi- and mono-target regression.
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
void ProcessOptions()
process user options
void MakeClassSpecific(std::ostream &, const TString &) const
write PDEFoam-specific classifier response NOT IMPLEMENTED YET!
double exp(double)
const Bool_t kTRUE
Definition: Rtypes.h:91
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
void AddWeightsXMLTo(void *parent) const
create XML output of PDEFoam method variables
void DeclareCompatibilityOptions()
options that are used ONLY for the READER to ensure backward compatibility
Definition: math.cpp:60
void ReadFoamsFromFile()
read foams from file
virtual void Close(Option_t *option="")
Close a file.
Definition: TFile.cxx:898
PDEFoamCell * GetDau0() const
Definition: PDEFoamCell.h:100