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