ROOT  6.06/09
Reference Guide
PDEFoam.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: S.Jadach, Tancredi Carli, Dominik Dannheim, Alexander Voigt
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Classes: PDEFoam *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementations *
12  * *
13  * Authors (alphabetical): *
14  * Tancredi Carli - CERN, Switzerland *
15  * Dominik Dannheim - CERN, Switzerland *
16  * S. Jadach - Institute of Nuclear Physics, Cracow, Poland *
17  * Alexander Voigt - TU Dresden, Germany *
18  * Peter Speckmayer - CERN, Switzerland *
19  * *
20  * Copyright (c) 2008, 2010: *
21  * CERN, Switzerland *
22  * MPI-K Heidelberg, Germany *
23  * *
24  * Redistribution and use in source and binary forms, with or without *
25  * modification, are permitted according to the terms listed in LICENSE *
26  * (http://tmva.sourceforge.net/LICENSE) *
27  **********************************************************************************/
28 
29 //_____________________________________________________________________
30 //
31 // Implementation of PDEFoam
32 //
33 // The PDEFoam method is an extension of the PDERS method, which uses
34 // self-adapting binning to divide the multi-dimensional phase space
35 // in a finite number of hyper-rectangles (boxes).
36 //
37 // For a given number of boxes, the binning algorithm adjusts the size
38 // and position of the boxes inside the multidimensional phase space,
39 // minimizing the variance of the signal and background densities inside
40 // the boxes. The binned density information is stored in binary trees,
41 // allowing for a very fast and memory-efficient classification of
42 // events.
43 //
44 // The implementation of the PDEFoam is based on the monte-carlo
45 // integration package TFoam included in the analysis package ROOT.
46 //
47 // The class TMVA::PDEFoam defines the default interface for the
48 // PDEFoam variants:
49 //
50 // - PDEFoamEvent
51 // - PDEFoamDiscriminant
52 // - PDEFoamTarget
53 // - PDEFoamMultiTarget
54 // - PDEFoamDecisionTree
55 //
56 // Per default PDEFoam stores in the cells the number of events (event
57 // weights) and therefore acts as an event density estimator.
58 // However, the above listed derived classes override this behaviour
59 // to implement certain PDEFoam variations.
60 //
61 // In order to use PDEFoam the user has to set the density estimator
62 // of the type TMVA::PDEFoamDensityBase, which is used to during the foam
63 // build-up. The default PDEFoam should be used with
64 // PDEFoamEventDensity.
65 // _____________________________________________________________________
66 
67 
68 #include <iostream>
69 #include <iomanip>
70 #include <fstream>
71 #include <sstream>
72 #include <cassert>
73 #include <limits>
74 
75 #include "TMVA/Event.h"
76 #include "TMVA/Tools.h"
77 #include "TMVA/PDEFoam.h"
78 #include "TMVA/MsgLogger.h"
79 #include "TMVA/Types.h"
80 
81 #ifndef ROOT_TStyle
82 #include "TStyle.h"
83 #endif
84 #ifndef ROOT_TObject
85 #include "TObject.h"
86 #endif
87 #ifndef ROOT_TH1D
88 #include "TH1D.h"
89 #endif
90 #ifndef ROOT_TMath
91 #include "TMath.h"
92 #endif
93 #ifndef ROOT_TVectorT
94 #include "TVectorT.h"
95 #endif
96 #ifndef ROOT_TRandom3
97 #include "TRandom3.h"
98 #endif
99 #ifndef ROOT_TColor
100 #include "TColor.h"
101 #endif
102 
104 
105 static const Float_t kHigh= FLT_MAX;
106 static const Float_t kVlow=-FLT_MAX;
107 
108 using namespace std;
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 /// Default constructor for streamer, user should not use it.
112 
114  fName("PDEFoam"),
115  fDim(0),
116  fNCells(0),
117  fNBin(5),
118  fNSampl(2000),
119  fEvPerBin(0),
120  fMaskDiv(0),
121  fInhiDiv(0),
122  fNoAct(1),
123  fLastCe(-1),
124  fCells(0),
125  fHistEdg(0),
126  fRvec(0),
127  fPseRan(new TRandom3(4356)),
128  fAlpha(0),
129  fFoamType(kSeparate),
130  fXmin(0),
131  fXmax(0),
132  fNElements(0),
133  fNmin(100),
134  fMaxDepth(0),
135  fVolFrac(1.0/15.0),
136  fFillFoamWithOrigWeights(kFALSE),
137  fDTSeparation(kFoam),
138  fPeekMax(kTRUE),
139  fDistr(NULL),
140  fTimer(new Timer(0, "PDEFoam", kTRUE)),
141  fVariableNames(new TObjArray()),
142  fLogger(new MsgLogger("PDEFoam"))
143 {
144  // fVariableNames may delete it's heap-based content
145  if (fVariableNames)
146  fVariableNames->SetOwner(kTRUE);
147 }
148 
149 ////////////////////////////////////////////////////////////////////////////////
150 /// User constructor, to be employed by the user
151 
153  fName(name),
154  fDim(0),
155  fNCells(1000),
156  fNBin(5),
157  fNSampl(2000),
158  fEvPerBin(0),
159  fMaskDiv(0),
160  fInhiDiv(0),
161  fNoAct(1),
162  fLastCe(-1),
163  fCells(0),
164  fHistEdg(0),
165  fRvec(0),
166  fPseRan(new TRandom3(4356)),
167  fAlpha(0),
168  fFoamType(kSeparate),
169  fXmin(0),
170  fXmax(0),
171  fNElements(0),
172  fNmin(100),
173  fMaxDepth(0),
174  fVolFrac(1.0/15.0),
175  fFillFoamWithOrigWeights(kFALSE),
176  fDTSeparation(kFoam),
177  fPeekMax(kTRUE),
178  fDistr(NULL),
179  fTimer(new Timer(1, "PDEFoam", kTRUE)),
180  fVariableNames(new TObjArray()),
181  fLogger(new MsgLogger("PDEFoam"))
182 {
183  if(strlen(name) > 128)
184  Log() << kFATAL << "Name too long " << name.Data() << Endl;
185 
186  // fVariableNames may delete it's heap-based content
187  if (fVariableNames)
189 }
190 
191 ////////////////////////////////////////////////////////////////////////////////
192 /// Default destructor
193 
195 {
196  delete fVariableNames;
197  delete fTimer;
198  if (fDistr) delete fDistr;
199  if (fPseRan) delete fPseRan;
200  if (fXmin) { delete [] fXmin; fXmin=0; }
201  if (fXmax) { delete [] fXmax; fXmax=0; }
202 
203  ResetCellElements();
204  if(fCells!= 0) {
205  for(Int_t i=0; i<fNCells; i++) delete fCells[i]; // PDEFoamCell*[]
206  delete [] fCells;
207  }
208  delete [] fRvec; //double[]
209  delete [] fAlpha; //double[]
210  delete [] fMaskDiv; //int[]
211  delete [] fInhiDiv; //int[]
212 
213  delete fLogger;
214 }
215 
216 ////////////////////////////////////////////////////////////////////////////////
217 /// Copy Constructor NOT IMPLEMENTED (NEVER USED)
218 
220  TObject(from)
221  , fDim(0)
222  , fNCells(0)
223  , fNBin(0)
224  , fNSampl(0)
225  , fEvPerBin(0)
226  , fMaskDiv(0)
227  , fInhiDiv(0)
228  , fNoAct(0)
229  , fLastCe(0)
230  , fCells(0)
231  , fHistEdg(0)
232  , fRvec(0)
233  , fPseRan(0)
234  , fAlpha(0)
235  , fFoamType(kSeparate)
236  , fXmin(0)
237  , fXmax(0)
238  , fNElements(0)
239  , fNmin(0)
240  , fMaxDepth(0)
241  , fVolFrac(1.0/15.0)
242  , fFillFoamWithOrigWeights(kFALSE)
243  , fDTSeparation(kFoam)
244  , fPeekMax(kTRUE)
245  , fDistr(0)
246  , fTimer(0)
247  , fVariableNames(0)
248  , fLogger(new MsgLogger(*from.fLogger))
249 {
250  Log() << kFATAL << "COPY CONSTRUCTOR NOT IMPLEMENTED" << Endl;
251 
252  // fVariableNames may delete it's heap-based content
253  if (fVariableNames)
255 }
256 
257 ////////////////////////////////////////////////////////////////////////////////
258 /// Sets dimension of cubical space
259 
261 {
262  if (kDim < 1)
263  Log() << kFATAL << "<SetDim>: Dimension is zero or negative!" << Endl;
264 
265  fDim = kDim;
266  if (fXmin) delete [] fXmin;
267  if (fXmax) delete [] fXmax;
268  fXmin = new Double_t[GetTotDim()];
269  fXmax = new Double_t[GetTotDim()];
270 }
271 
272 ////////////////////////////////////////////////////////////////////////////////
273 /// set lower foam bound in dimension idim
274 
276 {
277  if (idim<0 || idim>=GetTotDim())
278  Log() << kFATAL << "<SetXmin>: Dimension out of bounds!" << Endl;
279 
280  fXmin[idim]=wmin;
281 }
282 
283 ////////////////////////////////////////////////////////////////////////////////
284 /// set upper foam bound in dimension idim
285 
287 {
288  if (idim<0 || idim>=GetTotDim())
289  Log() << kFATAL << "<SetXmax>: Dimension out of bounds!" << Endl;
290 
291  fXmax[idim]=wmax;
292 }
293 
294 ////////////////////////////////////////////////////////////////////////////////
295 /// Basic initialization of FOAM invoked by the user.
296 /// IMPORTANT: Random number generator and the distribution object has to be
297 /// provided using SetPseRan and SetRho prior to invoking this initializator!
298 ///
299 /// After the foam is grown, space for 2 variables is reserved in
300 /// every cell. They are used for filling the foam cells.
301 
303 {
304  Bool_t addStatus = TH1::AddDirectoryStatus();
306 
307  if(fPseRan==0) Log() << kFATAL << "Random number generator not set" << Endl;
308  if(fDistr==0) Log() << kFATAL << "Distribution function not set" << Endl;
309  if(fDim==0) Log() << kFATAL << "Zero dimension not allowed" << Endl;
310 
311  /////////////////////////////////////////////////////////////////////////
312  // ALLOCATE SMALL LISTS //
313  // it is done globally, not for each cell, to save on allocation time //
314  /////////////////////////////////////////////////////////////////////////
315  fRvec = new Double_t[fDim]; // Vector of random numbers
316  if(fRvec==0) Log() << kFATAL << "Cannot initialize buffer fRvec" << Endl;
317 
318  if(fDim>0){
319  fAlpha = new Double_t[fDim]; // sum<1 for internal parametrization of the simplex
320  if(fAlpha==0) Log() << kFATAL << "Cannot initialize buffer fAlpha" << Endl;
321  }
322 
323  //====== List of directions inhibited for division
324  if(fInhiDiv == 0){
325  fInhiDiv = new Int_t[fDim];
326  for(Int_t i=0; i<fDim; i++) fInhiDiv[i]=0;
327  }
328  //====== Dynamic mask used in Explore for edge determination
329  if(fMaskDiv == 0){
330  fMaskDiv = new Int_t[fDim];
331  for(Int_t i=0; i<fDim; i++) fMaskDiv[i]=1;
332  }
333  //====== Initialize list of histograms
334  fHistEdg = new TObjArray(fDim); // Initialize list of histograms
335  for(Int_t i=0; i<fDim; i++){
336  TString hname, htitle;
337  hname = fName+TString("_HistEdge_");
338  hname += i;
339  htitle = TString("Edge Histogram No. ");
340  htitle += i;
341  (*fHistEdg)[i] = new TH1D(hname.Data(),htitle.Data(),fNBin,0.0, 1.0); // Initialize histogram for each edge
342  ((TH1D*)(*fHistEdg)[i])->Sumw2();
343  }
344 
345  // ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| //
346  // BUILD-UP of the FOAM //
347  // ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| //
348 
349  // prepare PDEFoam for growing
350  ResetCellElements(); // reset all cell elements
351 
352  // Define and explore root cell(s)
353  InitCells();
354  Grow();
355 
356  TH1::AddDirectory(addStatus);
357 
358  // prepare PDEFoam for the filling with events
359  ResetCellElements(); // reset all cell elements
360 } // Create
361 
362 ////////////////////////////////////////////////////////////////////////////////
363 /// Internal subprogram used by Create.
364 /// It initializes "root part" of the FOAM of the tree of cells.
365 
367 {
368  fLastCe =-1; // Index of the last cell
369  if(fCells!= 0) {
370  for(Int_t i=0; i<fNCells; i++) delete fCells[i];
371  delete [] fCells;
372  }
373 
374  fCells = new(nothrow) PDEFoamCell*[fNCells];
375  if (!fCells) {
376  Log() << kFATAL << "not enough memory to create " << fNCells
377  << " cells" << Endl;
378  }
379  for(Int_t i=0; i<fNCells; i++){
380  fCells[i]= new PDEFoamCell(fDim); // Allocate BIG list of cells
381  fCells[i]->SetSerial(i);
382  }
383 
384  /////////////////////////////////////////////////////////////////////////////
385  // Single Root Hypercube //
386  /////////////////////////////////////////////////////////////////////////////
387  CellFill(1, 0); // 0-th cell ACTIVE
388 
389  // Exploration of the root cell(s)
390  for(Long_t iCell=0; iCell<=fLastCe; iCell++){
391  Explore( fCells[iCell] ); // Exploration of root cell(s)
392  }
393 }//InitCells
394 
395 ////////////////////////////////////////////////////////////////////////////////
396 /// Internal subprogram used by Create.
397 /// It initializes content of the newly allocated active cell.
398 
400 {
401  PDEFoamCell *cell;
402  if (fLastCe==fNCells){
403  Log() << kFATAL << "Too many cells" << Endl;
404  }
405  fLastCe++; // 0-th cell is the first
406 
407  cell = fCells[fLastCe];
408 
409  cell->Fill(status, parent, 0, 0);
410 
411  cell->SetBest( -1); // pointer for planning division of the cell
412  cell->SetXdiv(0.5); // factor for division
413  Double_t xInt2,xDri2;
414  if(parent!=0){
415  xInt2 = 0.5*parent->GetIntg();
416  xDri2 = 0.5*parent->GetDriv();
417  cell->SetIntg(xInt2);
418  cell->SetDriv(xDri2);
419  }else{
420  cell->SetIntg(0.0);
421  cell->SetDriv(0.0);
422  }
423  return fLastCe;
424 }
425 
426 ////////////////////////////////////////////////////////////////////////////////
427 /// Internal subprogram used by Create.
428 /// It explores newly defined cell with help of special short MC sampling.
429 /// As a result, estimates of kTRUE and drive volume is defined/determined
430 /// Average and dispersion of the weight distribution will is found along
431 /// each edge and the best edge (minimum dispersion, best maximum weight)
432 /// is memorized for future use.
433 /// The optimal division point for eventual future cell division is
434 /// determined/recorded. Recorded are also minimum and maximum weight etc.
435 /// The volume estimate in all (inactive) parent cells is updated.
436 /// Note that links to parents and initial volume = 1/2 parent has to be
437 /// already defined prior to calling this routine.
438 ///
439 /// If fNmin > 0 then the total number of (training) events found in
440 /// the cell during the exploration is stored in the cell. This
441 /// information is used withing PeekMax() to avoid splitting cells
442 /// which contain less than fNmin events.
443 
445 {
446  Double_t wt, dx, xBest=0, yBest;
447  Double_t intOld, driOld;
448 
449  Long_t iev;
450  Double_t nevMC;
451  Int_t i, j, k;
452  Int_t nProj, kBest;
453  Double_t ceSum[5], xproj;
454 
455  Double_t event_density = 0;
456  Double_t totevents = 0;
457  Double_t toteventsOld = 0;
458 
459  PDEFoamVect cellSize(fDim);
460  PDEFoamVect cellPosi(fDim);
461 
462  cell->GetHcub(cellPosi,cellSize);
463 
464  PDEFoamCell *parent;
465 
466  Double_t *xRand = new Double_t[fDim];
467 
468  Double_t *volPart=0;
469 
470  // calculate volume scale
471  Double_t vol_scale = 1.0;
472  for (Int_t idim = 0; idim < fDim; ++idim)
473  vol_scale *= fXmax[idim] - fXmin[idim];
474 
475  cell->CalcVolume();
476  dx = cell->GetVolume() * vol_scale;
477  intOld = cell->GetIntg(); //memorize old values,
478  driOld = cell->GetDriv(); //will be needed for correcting parent cells
479  toteventsOld = GetCellElement(cell, 0);
480 
481  /////////////////////////////////////////////////////
482  // Special Short MC sampling to probe cell //
483  /////////////////////////////////////////////////////
484  ceSum[0]=0;
485  ceSum[1]=0;
486  ceSum[2]=0;
487  ceSum[3]=kHigh; //wtmin
488  ceSum[4]=kVlow; //wtmax
489 
490  for (i=0;i<fDim;i++) ((TH1D *)(*fHistEdg)[i])->Reset(); // Reset histograms
491 
492  Double_t nevEff=0.;
493  // ||||||||||||||||||||||||||BEGIN MC LOOP|||||||||||||||||||||||||||||
494  for (iev=0;iev<fNSampl;iev++){
495  MakeAlpha(); // generate uniformly vector inside hypercube
496 
497  if (fDim>0) for (j=0; j<fDim; j++) xRand[j]= cellPosi[j] +fAlpha[j]*(cellSize[j]);
498 
499  wt = dx*Eval(xRand, event_density);
500  totevents += event_density;
501 
502  nProj = 0;
503  if (fDim>0) {
504  for (k=0; k<fDim; k++) {
505  xproj =fAlpha[k];
506  ((TH1D *)(*fHistEdg)[nProj])->Fill(xproj,wt);
507  nProj++;
508  }
509  }
510 
511  ceSum[0] += wt; // sum of weights
512  ceSum[1] += wt*wt; // sum of weights squared
513  ceSum[2]++; // sum of 1
514  if (ceSum[3]>wt) ceSum[3]=wt; // minimum weight;
515  if (ceSum[4]<wt) ceSum[4]=wt; // maximum weight
516  // test MC loop exit condition
517  if (ceSum[1]>0) nevEff = ceSum[0]*ceSum[0]/ceSum[1];
518  else nevEff = 0;
519  if ( nevEff >= fNBin*fEvPerBin) break;
520  } // ||||||||||||||||||||||||||END MC LOOP|||||||||||||||||||||||||||||
521  totevents *= dx;
522 
523  if (fNSampl>0) totevents /= fNSampl;
524 
525  // make shure that, if root cell is explored, more than zero
526  // events were found.
527  if (cell==fCells[0] && ceSum[0]<=0.0){
528  if (ceSum[0]==0.0)
529  Log() << kFATAL << "No events were found during exploration of "
530  << "root cell. Please check PDEFoam parameters nSampl "
531  << "and VolFrac." << Endl;
532  else
533  Log() << kWARNING << "Negative number of events found during "
534  << "exploration of root cell" << Endl;
535  }
536 
537  //------------------------------------------------------------------
538  //--- predefine logics of searching for the best division edge ---
539  for (k=0; k<fDim;k++){
540  fMaskDiv[k] =1; // default is all
541  if ( fInhiDiv[k]==1) fMaskDiv[k] =0; // inhibit some...
542  }
543  kBest=-1;
544 
545  /////////////////////////////////////////////////////////////////////////////
546 
547  nevMC = ceSum[2];
548  Double_t intTrue = ceSum[0]/(nevMC+0.000001);
549  Double_t intDriv=0.;
550 
551  if (kBest == -1) Varedu(ceSum,kBest,xBest,yBest); // determine the best edge,
552  intDriv =sqrt(ceSum[1]/nevMC) -intTrue; // Foam build-up, sqrt(<w**2>) -<w>
553 
554  //=================================================================================
555  cell->SetBest(kBest);
556  cell->SetXdiv(xBest);
557  cell->SetIntg(intTrue);
558  cell->SetDriv(intDriv);
559  SetCellElement(cell, 0, totevents);
560 
561  // correct/update integrals in all parent cells to the top of the tree
562  Double_t parIntg, parDriv;
563  for (parent = cell->GetPare(); parent!=0; parent = parent->GetPare()){
564  parIntg = parent->GetIntg();
565  parDriv = parent->GetDriv();
566  parent->SetIntg( parIntg +intTrue -intOld );
567  parent->SetDriv( parDriv +intDriv -driOld );
568  SetCellElement( parent, 0, GetCellElement(parent, 0) + totevents - toteventsOld);
569  }
570  delete [] volPart;
571  delete [] xRand;
572 }
573 
574 ////////////////////////////////////////////////////////////////////////////////
575 /// Internal subrogram used by Create.
576 /// In determines the best edge candidate and the position of the cell division plane
577 /// in case of the variance reduction for future cell division,
578 /// using results of the MC exploration run stored in fHistEdg
579 
580 void TMVA::PDEFoam::Varedu(Double_t ceSum[5], Int_t &kBest, Double_t &xBest, Double_t &yBest)
581 {
582  Double_t nent = ceSum[2];
583  // Double_t swAll = ceSum[0];
584  Double_t sswAll = ceSum[1];
585  Double_t ssw = sqrt(sswAll)/sqrt(nent);
586  //
587  Double_t sswIn,sswOut,xLo,xUp;
588  kBest =-1;
589  xBest =0.5;
590  yBest =1.0;
591  Double_t maxGain=0.0;
592  // Now go over all projections kProj
593  for(Int_t kProj=0; kProj<fDim; kProj++) {
594  if( fMaskDiv[kProj]) {
595  // initialize search over bins
596  // Double_t sigmIn =0.0; Double_t sigmOut =0.0;
597  Double_t sswtBest = kHigh;
598  Double_t gain =0.0;
599  Double_t xMin=0.0; Double_t xMax=0.0;
600  // Double loop over all pairs jLo<jUp
601  for(Int_t jLo=1; jLo<=fNBin; jLo++) {
602  Double_t aswIn=0; Double_t asswIn=0;
603  for(Int_t jUp=jLo; jUp<=fNBin;jUp++) {
604  aswIn += ((TH1D *)(*fHistEdg)[kProj])->GetBinContent(jUp);
605  asswIn += Sqr(((TH1D *)(*fHistEdg)[kProj])->GetBinError( jUp));
606  xLo=(jLo-1.0)/fNBin;
607  xUp=(jUp*1.0)/fNBin;
608  // swIn = aswIn/nent;
609  // swOut = (swAll-aswIn)/nent;
610  if ( (xUp-xLo) < std::numeric_limits<double>::epsilon()) sswIn=0.;
611  else sswIn = sqrt(asswIn) /sqrt(nent*(xUp-xLo)) *(xUp-xLo);
612  if ( (1.0-xUp+xLo) < std::numeric_limits<double>::epsilon()) sswOut=0.;
613  else if ( sswAll-asswIn < std::numeric_limits<double>::epsilon()) sswOut=0.;
614  else sswOut= sqrt(sswAll-asswIn)/sqrt(nent*(1.0-xUp+xLo)) *(1.0-xUp+xLo);
615  if( (sswIn+sswOut) < sswtBest) {
616  sswtBest = sswIn+sswOut;
617  gain = ssw-sswtBest;
618  // sigmIn = sswIn -swIn; // Debug
619  // sigmOut = sswOut-swOut; // Debug
620  xMin = xLo;
621  xMax = xUp;
622  }
623  }//jUp
624  }//jLo
625  Int_t iLo = (Int_t) (fNBin*xMin);
626  Int_t iUp = (Int_t) (fNBin*xMax);
627 
628  if(gain>=maxGain) {
629  maxGain=gain;
630  kBest=kProj; // <--- !!!!! The best edge
631  xBest=xMin;
632  yBest=xMax;
633  if(iLo == 0) xBest=yBest; // The best division point
634  if(iUp == fNBin) yBest=xBest; // this is not really used
635  }
636  }
637  } //kProj
638 
639  if( (kBest >= fDim) || (kBest<0) )
640  Log() << kFATAL << "Something wrong with kBest" << Endl;
641 } //PDEFoam::Varedu
642 
643 ////////////////////////////////////////////////////////////////////////////////
644 /// Internal subrogram used by Create.
645 /// Provides random vector Alpha 0< Alpha(i) < 1
646 
648 {
649  // simply generate and load kDim uniform random numbers
650  fPseRan->RndmArray(fDim,fRvec); // kDim random numbers needed
651  for(Int_t k=0; k<fDim; k++) fAlpha[k] = fRvec[k];
652 } //MakeAlpha
653 
654 ////////////////////////////////////////////////////////////////////////////////
655 /// Internal subprogram used by Create. It finds cell with maximal
656 /// driver integral for the purpose of the division. This function
657 /// is overridden by the PDEFoam Class to apply cuts on the number
658 /// of events in the cell (fNmin) and the cell tree depth
659 /// (GetMaxDepth() > 0) during cell buildup.
660 
662 {
663  Long_t iCell = -1;
664 
665  Long_t i;
666  Double_t drivMax, driv, xDiv;
667  Bool_t bCutNmin = kTRUE;
668  Bool_t bCutMaxDepth = kTRUE;
669  // drivMax = kVlow;
670  drivMax = 0; // only split cells if gain>0 (this also avoids splitting at cell boundary)
671  for(i=0; i<=fLastCe; i++) {//without root
672  if( fCells[i]->GetStat() == 1 ) {
673  // if driver integral < numeric limit, skip cell
674  driv = fCells[i]->GetDriv();
676  continue;
677 
678  // do not split cell at the edges
679  xDiv = TMath::Abs(fCells[i]->GetXdiv());
682  continue;
683 
684  // apply cut on depth
685  if (GetMaxDepth() > 0)
686  bCutMaxDepth = fCells[i]->GetDepth() < GetMaxDepth();
687 
688  // apply Nmin-cut
689  if (GetNmin() > 0)
690  bCutNmin = GetCellElement(fCells[i], 0) > GetNmin();
691 
692  // choose cell
693  if(driv > drivMax && bCutNmin && bCutMaxDepth) {
694  drivMax = driv;
695  iCell = i;
696  }
697  }
698  }
699 
700  if (iCell == -1){
701  if (!bCutNmin)
702  Log() << kVERBOSE << "Warning: No cell with more than "
703  << GetNmin() << " events found!" << Endl;
704  else if (!bCutMaxDepth)
705  Log() << kVERBOSE << "Warning: Maximum depth reached: "
706  << GetMaxDepth() << Endl;
707  else
708  Log() << kWARNING << "<PDEFoam::PeekMax>: no more candidate cells (drivMax>0) found for further splitting." << Endl;
709  }
710 
711  return(iCell);
712 }
713 
714 ////////////////////////////////////////////////////////////////////////////////
715 /// Internal subrogram used by Create.
716 /// It divides cell iCell into two daughter cells.
717 /// The iCell is retained and tagged as inactive, daughter cells are appended
718 /// at the end of the buffer.
719 /// New vertex is added to list of vertices.
720 /// List of active cells is updated, iCell removed, two daughters added
721 /// and their properties set with help of MC sampling (PDEFoam_Explore)
722 /// Returns Code RC=-1 of buffer limit is reached, fLastCe=fnBuf.
723 
725 {
726  // Double_t xdiv;
727  Int_t kBest;
728 
729  if(fLastCe+1 >= fNCells) Log() << kFATAL << "Buffer limit is reached, fLastCe=fnBuf" << Endl;
730 
731  cell->SetStat(0); // reset cell as inactive
732  fNoAct++;
733 
734  // xdiv = cell->GetXdiv();
735  kBest = cell->GetBest();
736  if( kBest<0 || kBest>=fDim ) Log() << kFATAL << "Wrong kBest" << Endl;
737 
738  //////////////////////////////////////////////////////////////////
739  // define two daughter cells (active) //
740  //////////////////////////////////////////////////////////////////
741 
742  Int_t d1 = CellFill(1, cell);
743  Int_t d2 = CellFill(1, cell);
744  cell->SetDau0((fCells[d1]));
745  cell->SetDau1((fCells[d2]));
746 
747  Explore( (fCells[d1]) );
748  Explore( (fCells[d2]) );
749 
750  return 1;
751 } // PDEFoam_Divide
752 
753 ////////////////////////////////////////////////////////////////////////////////
754 /// Internal subprogram.
755 /// Evaluates (training) distribution.
756 
758 {
759  // Transform variable xRand, since Foam boundaries are [0,1] and
760  // fDistr is filled with events which range in [fXmin,fXmax]
761  //
762  // Transformation: [0, 1] --> [xmin, xmax]
763  std::vector<Double_t> xvec;
764  xvec.reserve(GetTotDim());
765  for (Int_t idim = 0; idim < GetTotDim(); ++idim)
766  xvec.push_back( VarTransformInvers(idim, xRand[idim]) );
767 
768  return GetDistr()->Density(xvec, event_density);
769 }
770 
771 ////////////////////////////////////////////////////////////////////////////////
772 /// Internal subrogram used by Create.
773 /// It grow new cells by the binary division process.
774 /// This function is overridden by the PDEFoam class to stop the foam buildup process
775 /// if one of the cut conditions stop the cell split.
776 
778 {
779  fTimer->Init(fNCells);
780 
781  Long_t iCell;
782  PDEFoamCell* newCell;
783 
784  while ( (fLastCe+2) < fNCells ) { // this condition also checked inside Divide
785  iCell = PeekMax(); // peek up cell with maximum driver integral
786 
787  if ( (iCell<0) || (iCell>fLastCe) ) {
788  Log() << kVERBOSE << "Break: "<< fLastCe+1 << " cells created" << Endl;
789  // remove remaining empty cells
790  for (Long_t jCell=fLastCe+1; jCell<fNCells; jCell++)
791  delete fCells[jCell];
792  fNCells = fLastCe+1;
793  break;
794  }
795  newCell = fCells[iCell];
796 
797  OutputGrow();
798 
799  if ( Divide( newCell )==0) break; // and divide it into two
800  }
801  OutputGrow( kTRUE );
802  CheckAll(1); // set arg=1 for more info
803 
804  Log() << kVERBOSE << GetNActiveCells() << " active cells created" << Endl;
805 }// Grow
806 
807 ////////////////////////////////////////////////////////////////////////////////
808 /// This can be called before Create, after setting kDim
809 /// It defines which variables are excluded in the process of the cell division.
810 /// For example 'FoamX->SetInhiDiv(1, 1);' inhibits division of y-variable.
811 
813 {
814  if(fDim==0) Log() << kFATAL << "SetInhiDiv: fDim=0" << Endl;
815  if(fInhiDiv == 0) {
816  fInhiDiv = new Int_t[ fDim ];
817  for(Int_t i=0; i<fDim; i++) fInhiDiv[i]=0;
818  }
819  //
820  if( ( 0<=iDim) && (iDim<fDim)) {
821  fInhiDiv[iDim] = inhiDiv;
822  } else
823  Log() << kFATAL << "Wrong iDim" << Endl;
824 }//SetInhiDiv
825 
826 ////////////////////////////////////////////////////////////////////////////////
827 /// User utility, miscellaneous and debug.
828 /// Checks all pointers in the tree of cells. This is useful autodiagnostic.
829 /// level=0, no printout, failures causes STOP
830 /// level=1, printout, failures lead to WARNINGS only
831 
833 {
834  Int_t errors, warnings;
835  PDEFoamCell *cell;
836  Long_t iCell;
837 
838  errors = 0; warnings = 0;
839  if (level==1) Log() << kVERBOSE << "Performing consistency checks for created foam" << Endl;
840  for(iCell=1; iCell<=fLastCe; iCell++) {
841  cell = fCells[iCell];
842  // checking general rules
843  if( ((cell->GetDau0()==0) && (cell->GetDau1()!=0) ) ||
844  ((cell->GetDau1()==0) && (cell->GetDau0()!=0) ) ) {
845  errors++;
846  if (level==1) Log() << kFATAL << "ERROR: Cell's no %d has only one daughter " << iCell << Endl;
847  }
848  if( (cell->GetDau0()==0) && (cell->GetDau1()==0) && (cell->GetStat()==0) ) {
849  errors++;
850  if (level==1) Log() << kFATAL << "ERROR: Cell's no %d has no daughter and is inactive " << iCell << Endl;
851  }
852  if( (cell->GetDau0()!=0) && (cell->GetDau1()!=0) && (cell->GetStat()==1) ) {
853  errors++;
854  if (level==1) Log() << kFATAL << "ERROR: Cell's no %d has two daughters and is active " << iCell << Endl;
855  }
856 
857  // checking parents
858  if( (cell->GetPare())!=fCells[0] ) { // not child of the root
859  if ( (cell != cell->GetPare()->GetDau0()) && (cell != cell->GetPare()->GetDau1()) ) {
860  errors++;
861  if (level==1) Log() << kFATAL << "ERROR: Cell's no %d parent not pointing to this cell " << iCell << Endl;
862  }
863  }
864 
865  // checking daughters
866  if(cell->GetDau0()!=0) {
867  if(cell != (cell->GetDau0())->GetPare()) {
868  errors++;
869  if (level==1) Log() << kFATAL << "ERROR: Cell's no %d daughter 0 not pointing to this cell " << iCell << Endl;
870  }
871  }
872  if(cell->GetDau1()!=0) {
873  if(cell != (cell->GetDau1())->GetPare()) {
874  errors++;
875  if (level==1) Log() << kFATAL << "ERROR: Cell's no %d daughter 1 not pointing to this cell " << iCell << Endl;
876  }
877  }
878  if(cell->GetVolume()<1E-50) {
879  errors++;
880  if(level==1) Log() << kFATAL << "ERROR: Cell no. " << iCell << " has Volume of <1E-50" << Endl;
881  }
882  }// loop after cells;
883 
884  // Check for cells with Volume=0
885  for(iCell=0; iCell<=fLastCe; iCell++) {
886  cell = fCells[iCell];
887  if( (cell->GetStat()==1) && (cell->GetVolume()<1E-11) ) {
888  errors++;
889  if(level==1) Log() << kFATAL << "ERROR: Cell no. " << iCell << " is active but Volume is 0 " << Endl;
890  }
891  }
892  // summary
893  if(level==1){
894  Log() << kVERBOSE << "Check has found " << errors << " errors and " << warnings << " warnings." << Endl;
895  }
896  if(errors>0){
897  Info("CheckAll","Check - found total %d errors \n",errors);
898  }
899 } // Check
900 
901 ////////////////////////////////////////////////////////////////////////////////
902 /// Prints geometry of and elements of 'iCell', as well as relations
903 /// to parent and daughter cells.
904 
906 {
907  if (iCell < 0 || iCell > fLastCe) {
908  Log() << kWARNING << "<PrintCell(iCell=" << iCell
909  << ")>: cell number " << iCell << " out of bounds!"
910  << Endl;
911  return;
912  }
913 
914  PDEFoamVect cellPosi(fDim), cellSize(fDim);
915  fCells[iCell]->GetHcub(cellPosi,cellSize);
916  Int_t kBest = fCells[iCell]->GetBest();
917  Double_t xBest = fCells[iCell]->GetXdiv();
918 
919  Log() << "Cell[" << iCell << "]={ ";
920  Log() << " " << fCells[iCell] << " " << Endl; // extra DEBUG
921  Log() << " Xdiv[abs. coord.]="
922  << VarTransformInvers(kBest,cellPosi[kBest] + xBest*cellSize[kBest])
923  << Endl;
924  Log() << " Abs. coord. = (";
925  for (Int_t idim=0; idim<fDim; idim++) {
926  Log() << "dim[" << idim << "]={"
927  << VarTransformInvers(idim,cellPosi[idim]) << ","
928  << VarTransformInvers(idim,cellPosi[idim] + cellSize[idim])
929  << "}";
930  if (idim < fDim-1)
931  Log() << ", ";
932  }
933  Log() << ")" << Endl;
934  fCells[iCell]->Print("1");
935  // print the cell elements
936  Log() << "Elements: [";
937  TVectorD *vec = (TVectorD*)fCells[iCell]->GetElement();
938  if (vec != NULL){
939  for (Int_t i=0; i<vec->GetNrows(); i++){
940  if (i>0) Log() << ", ";
941  Log() << GetCellElement(fCells[iCell], i);
942  }
943  } else
944  Log() << "not set";
945  Log() << "]" << Endl;
946  Log()<<"}"<<Endl;
947 }
948 
949 ////////////////////////////////////////////////////////////////////////////////
950 /// Prints geometry of ALL cells of the FOAM
951 
953 {
954  for(Long_t iCell=0; iCell<=fLastCe; iCell++)
955  PrintCell(iCell);
956 }
957 
958 ////////////////////////////////////////////////////////////////////////////////
959 /// This function fills a weight 'wt' into the PDEFoam cell, which
960 /// corresponds to the given event 'ev'. Per default cell element 0
961 /// is filled with the weight 'wt', and cell element 1 is filled
962 /// with the squared weight. This function can be overridden by a
963 /// subclass in order to change the values stored in the foam cells.
964 
966 {
967  // find corresponding foam cell
968  std::vector<Float_t> values = ev->GetValues();
969  std::vector<Float_t> tvalues = VarTransform(values);
970  PDEFoamCell *cell = FindCell(tvalues);
971 
972  // 0. Element: Sum of weights 'wt'
973  // 1. Element: Sum of weights 'wt' squared
974  SetCellElement(cell, 0, GetCellElement(cell, 0) + wt);
975  SetCellElement(cell, 1, GetCellElement(cell, 1) + wt*wt);
976 }
977 
978 ////////////////////////////////////////////////////////////////////////////////
979 /// Remove the cell elements from all cells.
980 
982 {
983  if (!fCells) return;
984 
985  Log() << kVERBOSE << "Delete cell elements" << Endl;
986  for (Long_t iCell = 0; iCell < fNCells; ++iCell) {
987  TObject* elements = fCells[iCell]->GetElement();
988  if (elements) {
989  delete elements;
990  fCells[iCell]->SetElement(NULL);
991  }
992  }
993 }
994 
995 ////////////////////////////////////////////////////////////////////////////////
996 /// Returns true, if the value of the given cell is undefined.
997 /// Default value: kFALSE. This function can be overridden by
998 /// sub-classes.
999 
1001 {
1002  return kFALSE;
1003 }
1004 
1005 ////////////////////////////////////////////////////////////////////////////////
1006 /// This function finds the cell, which corresponds to the given
1007 /// untransformed event vector 'xvec' and return its value, which is
1008 /// given by the parameter 'cv'. If kernel != NULL, then
1009 /// PDEFoamKernelBase::Estimate() is called on the transformed event
1010 /// variables.
1011 ///
1012 /// Parameters:
1013 ///
1014 /// - xvec - event vector (untransformed, [fXmin,fXmax])
1015 ///
1016 /// - cv - the cell value to return
1017 ///
1018 /// - kernel - PDEFoam kernel estimator. If NULL is given, than the
1019 /// pure cell value is returned
1020 ///
1021 /// Return:
1022 ///
1023 /// The cell value, corresponding to 'xvec', estimated by the given
1024 /// kernel.
1025 
1026 Float_t TMVA::PDEFoam::GetCellValue(const std::vector<Float_t> &xvec, ECellValue cv, PDEFoamKernelBase *kernel)
1027 {
1028  std::vector<Float_t> txvec(VarTransform(xvec));
1029  if (kernel == NULL)
1030  return GetCellValue(FindCell(txvec), cv);
1031  else
1032  return kernel->Estimate(this, txvec, cv);
1033 }
1034 
1035 ////////////////////////////////////////////////////////////////////////////////
1036 /// This function finds all cells, which corresponds to the given
1037 /// (incomplete) untransformed event vector 'xvec' and returns the
1038 /// cell values, according to the parameter 'cv'.
1039 ///
1040 /// Parameters:
1041 ///
1042 /// - xvec - map for the untransformed vector. The key (Int_t) is
1043 /// the dimension, and the value (Float_t) is the event
1044 /// coordinate. Note that not all coordinates have to be
1045 /// specified.
1046 ///
1047 /// - cv - cell values to return
1048 ///
1049 /// Return:
1050 ///
1051 /// cell values from all cells that were found
1052 
1053 std::vector<Float_t> TMVA::PDEFoam::GetCellValue( const std::map<Int_t,Float_t>& xvec, ECellValue cv )
1054 {
1055  // transformed event
1056  std::map<Int_t,Float_t> txvec;
1057  for (std::map<Int_t,Float_t>::const_iterator it=xvec.begin(); it!=xvec.end(); ++it)
1058  txvec.insert(std::pair<Int_t, Float_t>(it->first, VarTransform(it->first, it->second)));
1059 
1060  // find all cells, which correspond to the transformed event
1061  std::vector<PDEFoamCell*> cells = FindCells(txvec);
1062 
1063  // get the cell values
1064  std::vector<Float_t> cell_values;
1065  cell_values.reserve(cells.size());
1066  for (std::vector<PDEFoamCell*>::const_iterator cell_it=cells.begin();
1067  cell_it != cells.end(); ++cell_it)
1068  cell_values.push_back(GetCellValue(*cell_it, cv));
1069 
1070  return cell_values;
1071 }
1072 
1073 ////////////////////////////////////////////////////////////////////////////////
1074 /// Find cell that contains 'xvec' (in foam coordinates [0,1]).
1075 ///
1076 /// Loop to find cell that contains 'xvec' starting at root cell,
1077 /// and traversing binary tree to find the cell quickly. Note, that
1078 /// if 'xvec' lies outside the foam, the cell which is nearest to
1079 /// 'xvec' is returned. (The returned pointer should never be
1080 /// NULL.)
1081 ///
1082 /// Parameters:
1083 ///
1084 /// - xvec - event vector (in foam coordinates [0,1])
1085 ///
1086 /// Return:
1087 ///
1088 /// PDEFoam cell corresponding to 'xvec'
1089 
1090 TMVA::PDEFoamCell* TMVA::PDEFoam::FindCell( const std::vector<Float_t> &xvec ) const
1091 {
1092  PDEFoamVect cellPosi0(GetTotDim()), cellSize0(GetTotDim());
1093  PDEFoamCell *cell, *cell0;
1094 
1095  cell=fCells[0]; // start with root cell
1096  Int_t idim=0;
1097  while (cell->GetStat()!=1) { //go down binary tree until cell is found
1098  idim=cell->GetBest(); // dimension that changed
1099  cell0=cell->GetDau0();
1100  cell0->GetHcub(cellPosi0,cellSize0);
1101 
1102  if (xvec.at(idim)<=cellPosi0[idim]+cellSize0[idim])
1103  cell=cell0;
1104  else
1105  cell=(cell->GetDau1());
1106  }
1107  return cell;
1108 }
1109 
1110 ////////////////////////////////////////////////////////////////////////////////
1111 /// This is a helper function for std::vector<PDEFoamCell*>
1112 /// FindCells(...) and a generalisation of PDEFoamCell* FindCell().
1113 /// It saves in 'cells' all cells, which contain the coordinates
1114 /// specifies in 'txvec'. Note, that not all coordinates have to be
1115 /// specified in 'txvec'.
1116 ///
1117 /// Parameters:
1118 ///
1119 /// - txvec - event vector in foam coordinates [0,1]. The key is
1120 /// the dimension and the value is the event coordinate. Note,
1121 /// that not all coordinates have to be specified.
1122 ///
1123 /// - cell - cell to start searching with (usually root cell
1124 /// fCells[0])
1125 ///
1126 /// - cells - list of cells that were found
1127 
1128 void TMVA::PDEFoam::FindCells(const std::map<Int_t, Float_t> &txvec, PDEFoamCell* cell, std::vector<PDEFoamCell*> &cells) const
1129 {
1130  PDEFoamVect cellPosi0(GetTotDim()), cellSize0(GetTotDim());
1131  PDEFoamCell *cell0;
1132  Int_t idim=0;
1133 
1134  while (cell->GetStat()!=1) { //go down binary tree until cell is found
1135  idim=cell->GetBest(); // dimension that changed
1136 
1137  // check if dimension 'idim' is specified in 'txvec'
1138  map<Int_t, Float_t>::const_iterator it = txvec.find(idim);
1139 
1140  if (it != txvec.end()){
1141  // case 1: cell is splitten in a dimension which is specified
1142  // in txvec
1143  cell0=cell->GetDau0();
1144  cell0->GetHcub(cellPosi0,cellSize0);
1145  // check, whether left daughter cell contains txvec
1146  if (it->second <= cellPosi0[idim] + cellSize0[idim])
1147  cell=cell0;
1148  else
1149  cell=cell->GetDau1();
1150  } else {
1151  // case 2: cell is splitten in target dimension
1152  FindCells(txvec, cell->GetDau0(), cells);
1153  FindCells(txvec, cell->GetDau1(), cells);
1154  return;
1155  }
1156  }
1157  cells.push_back(cell);
1158 }
1159 
1160 ////////////////////////////////////////////////////////////////////////////////
1161 /// Find all cells, that contain txvec. This function can be used,
1162 /// when the dimension of the foam is greater than the dimension of
1163 /// txvec. E.g. this is the case for multi-target regression.
1164 ///
1165 /// Parameters:
1166 ///
1167 /// - txvec - event vector of variables, transformed into foam
1168 /// coordinates [0,1]. The size of txvec can be smaller than the
1169 /// dimension of the foam.
1170 ///
1171 /// Return value:
1172 ///
1173 /// - vector of cells, that fit txvec
1174 
1175 std::vector<TMVA::PDEFoamCell*> TMVA::PDEFoam::FindCells(const std::vector<Float_t> &txvec) const
1176 {
1177  // copy the coordinates from 'txvec' into a map
1178  std::map<Int_t, Float_t> txvec_map;
1179  for (UInt_t i=0; i<txvec.size(); ++i)
1180  txvec_map.insert(std::pair<Int_t, Float_t>(i, txvec.at(i)));
1181 
1182  // the cells found
1183  std::vector<PDEFoamCell*> cells(0);
1184 
1185  // loop over all target dimensions
1186  FindCells(txvec_map, fCells[0], cells);
1187 
1188  return cells;
1189 }
1190 
1191 ////////////////////////////////////////////////////////////////////////////////
1192 /// Find all cells, that contain the coordinates specified in txvec.
1193 /// The key in 'txvec' is the dimension, and the corresponding value
1194 /// is the coordinate. Note, that not all coordinates have to be
1195 /// specified in txvec.
1196 ///
1197 /// Parameters:
1198 ///
1199 /// - txvec - map of coordinates (transformed into foam coordinates
1200 /// [0,1])
1201 ///
1202 /// Return value:
1203 ///
1204 /// - vector of cells, that fit txvec
1205 
1206 std::vector<TMVA::PDEFoamCell*> TMVA::PDEFoam::FindCells(const std::map<Int_t, Float_t> &txvec) const
1207 {
1208  // the cells found
1209  std::vector<PDEFoamCell*> cells(0);
1210 
1211  // loop over all target dimensions
1212  FindCells(txvec, fCells[0], cells);
1213 
1214  return cells;
1215 }
1216 
1217 ////////////////////////////////////////////////////////////////////////////////
1218 /// Draws 1-dimensional foam (= histogram)
1219 ///
1220 /// Parameters:
1221 ///
1222 /// - cell_value - the cell value to draw
1223 ///
1224 /// - nbin - number of bins of result histogram
1225 ///
1226 /// - kernel - a PDEFoam kernel.
1227 
1229 {
1230  // avoid plotting of wrong dimensions
1231  if ( GetTotDim()!=1 )
1232  Log() << kFATAL << "<Draw1Dim>: function can only be used for 1-dimensional foams!"
1233  << Endl;
1234 
1235  TString hname("h_1dim");
1236  TH1D* h1=(TH1D*)gDirectory->Get(hname);
1237  if (h1) delete h1;
1238  h1= new TH1D(hname, "1-dimensional Foam", nbin, fXmin[0], fXmax[0]);
1239 
1240  if (!h1) Log() << kFATAL << "ERROR: Can not create histo" << hname << Endl;
1241 
1242  // loop over all bins
1243  for (Int_t ibinx=1; ibinx<=h1->GetNbinsX(); ++ibinx) {
1244  // get event vector corresponding to bin
1245  std::vector<Float_t> txvec;
1246  txvec.push_back( VarTransform(0, h1->GetBinCenter(ibinx)) );
1247  Float_t val = 0;
1248  if (kernel != NULL) {
1249  // get cell value using the kernel
1250  val = kernel->Estimate(this, txvec, cell_value);
1251  } else {
1252  val = GetCellValue(FindCell(txvec), cell_value);
1253  }
1254  // fill value to histogram
1255  h1->SetBinContent(ibinx, val + h1->GetBinContent(ibinx));
1256  }
1257 
1258  return h1;
1259 }
1260 
1261 ////////////////////////////////////////////////////////////////////////////////
1262 /// Project foam variable idim1 and variable idim2 to histogram.
1263 ///
1264 /// Parameters:
1265 ///
1266 /// - idim1, idim2 - dimensions to project to
1267 ///
1268 /// - cell_value - the cell value to draw
1269 ///
1270 /// - kernel - a PDEFoam kernel (optional). If NULL is given, the
1271 /// kernel is ignored and the pure cell values are
1272 /// plotted.
1273 ///
1274 /// - nbin - number of bins in x and y direction of result histogram
1275 /// (optional, default is 50).
1276 ///
1277 /// Returns:
1278 /// a 2-dimensional histogram
1279 
1280 TH2D* TMVA::PDEFoam::Project2( Int_t idim1, Int_t idim2, ECellValue cell_value, PDEFoamKernelBase *kernel, UInt_t nbin )
1281 {
1282  // avoid plotting of wrong dimensions
1283  if ((idim1>=GetTotDim()) || (idim1<0) ||
1284  (idim2>=GetTotDim()) || (idim2<0) ||
1285  (idim1==idim2) )
1286  Log() << kFATAL << "<Project2>: wrong dimensions given: "
1287  << idim1 << ", " << idim2 << Endl;
1288 
1289  // root can not handle too many bins in one histogram --> catch this
1290  // Furthermore, to have more than 1000 bins in the histogram doesn't make
1291  // sense.
1292  if (nbin>1000){
1293  Log() << kWARNING << "Warning: number of bins too big: " << nbin
1294  << " Using 1000 bins for each dimension instead." << Endl;
1295  nbin = 1000;
1296  } else if (nbin<1) {
1297  Log() << kWARNING << "Wrong bin number: " << nbin
1298  << "; set nbin=50" << Endl;
1299  nbin = 50;
1300  }
1301 
1302  // create result histogram
1303  TString hname(Form("h_%d_vs_%d",idim1,idim2));
1304 
1305  // if histogram with this name already exists, delete it
1306  TH2D* h1=(TH2D*)gDirectory->Get(hname.Data());
1307  if (h1) delete h1;
1308  h1= new TH2D(hname.Data(), Form("var%d vs var%d",idim1,idim2), nbin, fXmin[idim1], fXmax[idim1], nbin, fXmin[idim2], fXmax[idim2]);
1309 
1310  if (!h1) Log() << kFATAL << "ERROR: Can not create histo" << hname << Endl;
1311 
1312  // ============== start projection algorithm ================
1313  // loop over all histogram bins (2-dim)
1314  for (Int_t xbin = 1; xbin <= h1->GetNbinsX(); ++xbin) {
1315  for (Int_t ybin = 1; ybin <= h1->GetNbinsY(); ++ybin) {
1316  // calculate the phase space point, which corresponds to this
1317  // bin combination
1318  std::map<Int_t, Float_t> txvec;
1319  txvec[idim1] = VarTransform(idim1, h1->GetXaxis()->GetBinCenter(xbin));
1320  txvec[idim2] = VarTransform(idim2, h1->GetYaxis()->GetBinCenter(ybin));
1321 
1322  // find the cells, which corresponds to this phase space
1323  // point
1324  std::vector<TMVA::PDEFoamCell*> cells = FindCells(txvec);
1325 
1326  // loop over cells and fill the histogram with the cell
1327  // values
1328  Float_t sum_cv = 0; // sum of the cell values
1329  for (std::vector<TMVA::PDEFoamCell*>::const_iterator it = cells.begin();
1330  it != cells.end(); ++it) {
1331  // get cell position and size
1332  PDEFoamVect cellPosi(GetTotDim()), cellSize(GetTotDim());
1333  (*it)->GetHcub(cellPosi,cellSize);
1334  // Create complete event vector from txvec. The missing
1335  // coordinates of txvec are set to the cell center.
1336  std::vector<Float_t> tvec;
1337  for (Int_t i=0; i<GetTotDim(); ++i) {
1338  if ( i != idim1 && i != idim2 )
1339  tvec.push_back(cellPosi[i] + 0.5*cellSize[i]);
1340  else
1341  tvec.push_back(txvec[i]);
1342  }
1343  if (kernel != NULL) {
1344  // get the cell value using the kernel
1345  sum_cv += kernel->Estimate(this, tvec, cell_value);
1346  } else {
1347  sum_cv += GetCellValue(FindCell(tvec), cell_value);
1348  }
1349  }
1350 
1351  // fill the bin content
1352  h1->SetBinContent(xbin, ybin, sum_cv + h1->GetBinContent(xbin, ybin));
1353  }
1354  }
1355 
1356  return h1;
1357 }
1358 
1359 ////////////////////////////////////////////////////////////////////////////////
1360 /// Returns the cell value of 'cell' corresponding to the given
1361 /// option 'cv'. This function should be overridden by the subclass
1362 /// in order to specify which cell elements to return for a given
1363 /// cell value 'cv'. By default kValue returns cell element 0, and
1364 /// kValueError returns cell element 1.
1365 
1367 {
1368  // calculate cell value (depending on the given option 'cv')
1369  switch (cv) {
1370 
1371  case kValue:
1372  return GetCellElement(cell, 0);
1373 
1374  case kValueError:
1375  return GetCellElement(cell, 1);
1376 
1377  case kValueDensity: {
1378 
1379  Double_t volume = cell->GetVolume();
1380  if (volume > numeric_limits<double>::epsilon()) {
1381  return GetCellValue(cell, kValue)/volume;
1382  } else {
1383  if (volume<=0){
1384  cell->Print("1"); // debug output
1385  Log() << kWARNING << "<GetCellDensity(cell)>: ERROR: cell volume"
1386  << " negative or zero!"
1387  << " ==> return cell density 0!"
1388  << " cell volume=" << volume
1389  << " cell entries=" << GetCellValue(cell, kValue) << Endl;
1390  } else {
1391  Log() << kWARNING << "<GetCellDensity(cell)>: WARNING: cell volume"
1392  << " close to zero!"
1393  << " cell volume: " << volume << Endl;
1394  }
1395  }
1396  }
1397  return 0;
1398 
1399  case kMeanValue:
1400  return cell->GetIntg();
1401 
1402  case kRms:
1403  return cell->GetDriv();
1404 
1405  case kRmsOvMean:
1406  if (cell->GetIntg() != 0)
1407  return cell->GetDriv()/cell->GetIntg();
1408  else
1409  return 0;
1410 
1411  case kCellVolume:
1412  return cell->GetVolume();
1413 
1414  default:
1415  Log() << kFATAL << "<GetCellValue>: unknown cell value" << Endl;
1416  return 0;
1417  }
1418 
1419  return 0;
1420 }
1421 
1422 ////////////////////////////////////////////////////////////////////////////////
1423 /// Returns cell element i of cell 'cell'. If the cell has no
1424 /// elements or the index 'i' is out of range, than 0 is returned.
1425 
1427 {
1428  // dynamic_cast doesn't seem to work here ?!
1429  TVectorD *vec = (TVectorD*)cell->GetElement();
1430 
1431  // if vec is not set or index out of range, return 0
1432  if (!vec || i >= (UInt_t) vec->GetNrows())
1433  return 0;
1434 
1435  return (*vec)(i);
1436 }
1437 
1438 ////////////////////////////////////////////////////////////////////////////////
1439 /// Set cell element i of cell to value. If the cell element i does
1440 /// not exist, it is created.
1441 
1443 {
1444  TVectorD *vec = NULL;
1445 
1446  // if no cell elements are set, create TVectorD with i+1 entries,
1447  // ranging from [0,i]
1448  if (cell->GetElement() == NULL) {
1449  vec = new TVectorD(i+1);
1450  vec->Zero(); // set all values to zero
1451  (*vec)(i) = value; // set element i to value
1452  cell->SetElement(vec);
1453  } else {
1454  // dynamic_cast doesn't seem to work here ?!
1455  vec = (TVectorD*)cell->GetElement();
1456  if (!vec)
1457  Log() << kFATAL << "<SetCellElement> ERROR: cell element is not a TVectorD*" << Endl;
1458  // check vector size and resize if necessary
1459  if (i >= (UInt_t) vec->GetNrows())
1460  vec->ResizeTo(0,i);
1461  // set element i to value
1462  (*vec)(i) = value;
1463  }
1464 }
1465 
1466 ////////////////////////////////////////////////////////////////////////////////
1467 /// Overridden function of PDEFoam to avoid native foam output.
1468 /// Draw TMVA-process bar instead.
1469 
1471 {
1472  if (finished) {
1473  Log() << kINFO << "Elapsed time: " + fTimer->GetElapsedTime()
1474  << " " << Endl;
1475  return;
1476  }
1477 
1478  Int_t modulo = 1;
1479 
1480  if (fNCells >= 100) modulo = Int_t(fNCells/100);
1481  if (fLastCe%modulo == 0) fTimer->DrawProgressBar( fLastCe );
1482 }
1483 
1484 ////////////////////////////////////////////////////////////////////////////////
1485 /// Debugging tool which plots the cells of a 2-dimensional PDEFoam
1486 /// as rectangles in C++ format readable for ROOT.
1487 ///
1488 /// Parameters:
1489 /// - filename - filename of ouput root macro
1490 ///
1491 /// - opt - cell_value, rms, rms_ov_mean
1492 /// If cell_value is set, the following values will be filled into
1493 /// the result histogram:
1494 /// - number of events - in case of classification with 2 separate
1495 /// foams or multi-target regression
1496 /// - discriminator - in case of classification with one
1497 /// unified foam
1498 /// - target - in case of mono-target regression
1499 /// If none of {cell_value, rms, rms_ov_mean} is given, the cells
1500 /// will not be filled.
1501 /// If 'opt' contains the string 'cellnumber', the index of
1502 /// each cell is draw in addition.
1503 ///
1504 /// - createCanvas - whether to create a new canvas or not
1505 ///
1506 /// - colors - whether to fill cells with colors or shades of grey
1507 ///
1508 /// Example:
1509 ///
1510 /// The following commands load a mono-target regression foam from
1511 /// file 'foam.root' and create a ROOT macro 'output.C', which
1512 /// draws all PDEFoam cells with little boxes. The latter are
1513 /// filled with colors according to the target value stored in the
1514 /// cell. Also the cell number is drawn.
1515 ///
1516 /// TFile file("foam.root");
1517 /// TMVA::PDEFoam *foam = (TMVA::PDEFoam*) gDirectory->Get("MonoTargetRegressionFoam");
1518 /// foam->RootPlot2dim("output.C","cell_value,cellnumber");
1519 /// gROOT->Macro("output.C");
1520 
1522  Bool_t createCanvas, Bool_t colors )
1523 {
1524  if (GetTotDim() != 2)
1525  Log() << kFATAL << "RootPlot2dim() can only be used with "
1526  << "two-dimensional foams!" << Endl;
1527 
1528  // select value to plot
1529  ECellValue cell_value = kValue;
1530  Bool_t plotcellnumber = kFALSE;
1531  Bool_t fillcells = kTRUE;
1532  if (opt.Contains("cell_value")){
1533  cell_value = kValue;
1534  } else if (opt.Contains("rms_ov_mean")){
1535  cell_value = kRmsOvMean;
1536  } else if (opt.Contains("rms")){
1537  cell_value = kRms;
1538  } else {
1539  fillcells = kFALSE;
1540  }
1541  if (opt.Contains("cellnumber"))
1542  plotcellnumber = kTRUE;
1543 
1544  // open file (root macro)
1545  std::ofstream outfile(filename, std::ios::out);
1546 
1547  outfile<<"{" << std::endl;
1548 
1549  // declare boxes and set the fill styles
1550  if (!colors) { // define grayscale colors from light to dark,
1551  // starting from color index 1000
1552  outfile << "TColor *graycolors[100];" << std::endl;
1553  outfile << "for (Int_t i=0.; i<100; i++)" << std::endl;
1554  outfile << " graycolors[i]=new TColor(1000+i, 1-(Float_t)i/100.,1-(Float_t)i/100.,1-(Float_t)i/100.);"<< std::endl;
1555  }
1556  if (createCanvas)
1557  outfile << "cMap = new TCanvas(\"" << fName << "\",\"Cell Map for "
1558  << fName << "\",600,600);" << std::endl;
1559 
1560  outfile<<"TBox*a=new TBox();"<<std::endl;
1561  outfile<<"a->SetFillStyle(0);"<<std::endl; // big frame
1562  outfile<<"a->SetLineWidth(4);"<<std::endl;
1563  outfile<<"TBox *b1=new TBox();"<<std::endl; // single cell
1564  outfile<<"TText*t=new TText();"<<std::endl; // text for numbering
1565  if (fillcells) {
1566  outfile << (colors ? "gStyle->SetPalette(1, 0);" : "gStyle->SetPalette(0);")
1567  << std::endl;
1568  outfile <<"b1->SetFillStyle(1001);"<<std::endl;
1569  outfile<<"TBox *b2=new TBox();"<<std::endl; // single cell
1570  outfile <<"b2->SetFillStyle(0);"<<std::endl;
1571  }
1572  else {
1573  outfile <<"b1->SetFillStyle(0);"<<std::endl;
1574  }
1575 
1576  if (fillcells)
1577  (colors ? gStyle->SetPalette(1, 0) : gStyle->SetPalette(0) );
1578 
1579  Float_t zmin = 1E8; // minimal value (for color calculation)
1580  Float_t zmax = -1E8; // maximal value (for color calculation)
1581 
1582  // if cells shall be filled, calculate minimal and maximal plot
1583  // value --> store in zmin and zmax
1584  if (fillcells) {
1585  for (Long_t iCell=1; iCell<=fLastCe; iCell++) {
1586  if ( fCells[iCell]->GetStat() == 1) {
1587  Float_t value = GetCellValue(fCells[iCell], cell_value);
1588  if (value<zmin)
1589  zmin=value;
1590  if (value>zmax)
1591  zmax=value;
1592  }
1593  }
1594  outfile << "// observed minimum and maximum of distribution: " << std::endl;
1595  outfile << "// Float_t zmin = "<< zmin << ";" << std::endl;
1596  outfile << "// Float_t zmax = "<< zmax << ";" << std::endl;
1597  }
1598 
1599  outfile << "// used minimum and maximum of distribution (taking into account log scale if applicable): " << std::endl;
1600  outfile << "Float_t zmin = "<< zmin << ";" << std::endl;
1601  outfile << "Float_t zmax = "<< zmax << ";" << std::endl;
1602 
1603  Float_t x1,y1,x2,y2,x,y; // box and text coordintates
1604  Float_t offs = 0.01;
1605  Float_t lpag = 1-2*offs;
1606  Int_t ncolors = colors ? gStyle->GetNumberOfColors() : 100;
1607  Float_t scale = (ncolors-1)/(zmax - zmin);
1608  PDEFoamVect cellPosi(GetTotDim()), cellSize(GetTotDim());
1609 
1610  // loop over cells and draw a box for every cell (and maybe the
1611  // cell number as well)
1612  outfile << "// =========== Rectangular cells ==========="<< std::endl;
1613  for (Long_t iCell=1; iCell<=fLastCe; iCell++) {
1614  if ( fCells[iCell]->GetStat() == 1) {
1615  fCells[iCell]->GetHcub(cellPosi,cellSize);
1616  x1 = offs+lpag*(cellPosi[0]);
1617  y1 = offs+lpag*(cellPosi[1]);
1618  x2 = offs+lpag*(cellPosi[0]+cellSize[0]);
1619  y2 = offs+lpag*(cellPosi[1]+cellSize[1]);
1620 
1621  if (fillcells) {
1622  // get cell value
1623  Float_t value = GetCellValue(fCells[iCell], cell_value);
1624 
1625  // calculate fill color
1626  Int_t color;
1627  if (colors)
1628  color = gStyle->GetColorPalette(Int_t((value-zmin)*scale));
1629  else
1630  color = 1000+(Int_t((value-zmin)*scale));
1631 
1632  // set fill color of box b1
1633  outfile << "b1->SetFillColor(" << color << ");" << std::endl;
1634  }
1635 
1636  // cell rectangle
1637  outfile<<"b1->DrawBox("<<x1<<","<<y1<<","<<x2<<","<<y2<<");"<<std::endl;
1638  if (fillcells)
1639  outfile<<"b2->DrawBox("<<x1<<","<<y1<<","<<x2<<","<<y2<<");"<<std::endl;
1640 
1641  // cell number
1642  if (plotcellnumber) {
1643  outfile<<"t->SetTextColor(4);"<<std::endl;
1644  if(fLastCe<51)
1645  outfile<<"t->SetTextSize(0.025);"<<std::endl; // text for numbering
1646  else if(fLastCe<251)
1647  outfile<<"t->SetTextSize(0.015);"<<std::endl;
1648  else
1649  outfile<<"t->SetTextSize(0.008);"<<std::endl;
1650  x = offs+lpag*(cellPosi[0]+0.5*cellSize[0]);
1651  y = offs+lpag*(cellPosi[1]+0.5*cellSize[1]);
1652  outfile<<"t->DrawText("<<x<<","<<y<<","<<"\""<<iCell<<"\""<<");"<<std::endl;
1653  }
1654  }
1655  }
1656  outfile<<"// ============== End Rectangles ==========="<< std::endl;
1657 
1658  outfile << "}" << std::endl;
1659  outfile.flush();
1660  outfile.close();
1661 }
1662 
1663 ////////////////////////////////////////////////////////////////////////////////
1664 /// Insert event to internal foam's density estimator
1665 /// PDEFoamDensityBase.
1666 
1668 {
1669  GetDistr()->FillBinarySearchTree(ev);
1670 }
1671 
1672 ////////////////////////////////////////////////////////////////////////////////
1673 /// Delete the foam's density estimator, which contains the binary
1674 /// search tree.
1675 
1677 {
1678  if(fDistr) delete fDistr;
1679  fDistr = NULL;
1680 }
virtual void FillFoamCells(const Event *ev, Float_t wt)
This function fills a weight 'wt' into the PDEFoam cell, which corresponds to the given event 'ev'...
Definition: PDEFoam.cxx:965
An array of TObjects.
Definition: TObjArray.h:39
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition: TVectorT.cxx:291
Random number generator class based on M.
Definition: TRandom3.h:29
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
void SetElement(TObject *fobj)
Definition: PDEFoamCell.h:112
THist< 2, double > TH2D
Definition: THist.h:320
Long_t PeekMax()
Internal subprogram used by Create.
Definition: PDEFoam.cxx:661
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4629
TObjArray * fVariableNames
timer for graphical output
Definition: PDEFoam.h:142
float Float_t
Definition: RtypesCore.h:53
void OutputGrow(Bool_t finished=false)
message logger
Definition: PDEFoam.cxx:1470
MsgLogger & Log() const
Definition: PDEFoam.h:265
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
#define gDirectory
Definition: TDirectory.h:218
Int_t GetNrows() const
Definition: TVectorT.h:81
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
void Fill(Int_t, PDEFoamCell *, PDEFoamCell *, PDEFoamCell *)
Fills in certain data into newly allocated cell.
void Grow()
Internal subrogram used by Create.
Definition: PDEFoam.cxx:777
static const char * filename()
Basic string class.
Definition: TString.h:137
static Bool_t AddDirectoryStatus()
static function: cannot be inlined on Windows/NT
Definition: TH1.cxx:709
void SetXmin(Int_t idim, Double_t wmin)
set lower foam bound in dimension idim
Definition: PDEFoam.cxx:275
TObject * GetElement() const
Definition: PDEFoamCell.h:113
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Int_t GetStat() const
Definition: PDEFoamCell.h:97
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
void SetCellElement(PDEFoamCell *cell, UInt_t i, Double_t value)
Set cell element i of cell to value.
Definition: PDEFoam.cxx:1442
Int_t Divide(PDEFoamCell *)
Internal subrogram used by Create.
Definition: PDEFoam.cxx:724
STL namespace.
void SetBest(Int_t Best)
Definition: PDEFoamCell.h:85
TH1D * Draw1Dim(ECellValue cell_value, Int_t nbin, PDEFoamKernelBase *kernel=NULL)
Draws 1-dimensional foam (= histogram)
Definition: PDEFoam.cxx:1228
void CheckAll(Int_t)
User utility, miscellaneous and debug.
Definition: PDEFoam.cxx:832
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Double_t GetVolume() const
Definition: PDEFoamCell.h:91
std::vector< TMVA::PDEFoamCell * > FindCells(const std::vector< Float_t > &) const
Find all cells, that contain txvec.
Definition: PDEFoam.cxx:1175
THist< 1, double > TH1D
Definition: THist.h:314
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition: TH1.cxx:1231
Int_t CellFill(Int_t, PDEFoamCell *)
Internal subprogram used by Create.
Definition: PDEFoam.cxx:399
const char * Data() const
Definition: TString.h:349
void MakeAlpha()
Internal subrogram used by Create.
Definition: PDEFoam.cxx:647
double sqrt(double)
Double_t GetIntg() const
Definition: PDEFoamCell.h:92
virtual Bool_t CellValueIsUndefined(PDEFoamCell *)
Returns true, if the value of the given cell is undefined.
Definition: PDEFoam.cxx:1000
static const double x2[5]
virtual Float_t Estimate(PDEFoam *, std::vector< Float_t > &, ECellValue)=0
Double_t x[n]
Definition: legend1.C:17
void RootPlot2dim(const TString &filename, TString opt, Bool_t createCanvas=kTRUE, Bool_t colors=kTRUE)
Debugging tool which plots the cells of a 2-dimensional PDEFoam as rectangles in C++ format readable ...
Definition: PDEFoam.cxx:1521
void SetDau1(PDEFoamCell *Daug)
Definition: PDEFoamCell.h:103
void Info(const char *location, const char *msgfmt,...)
Int_t GetNumberOfColors() const
Return number of colors in the color palette.
Definition: TStyle.cxx:796
TH1F * h1
Definition: legend1.C:5
void ResetCellElements()
Remove the cell elements from all cells.
Definition: PDEFoam.cxx:981
virtual void Explore(PDEFoamCell *Cell)
Internal subprogram used by Create.
Definition: PDEFoam.cxx:444
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:24
void Print(Option_t *option) const
Printout of the cell geometry parameters for the debug purpose.
char * out
Definition: TBase64.cxx:29
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:90
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:594
void InitCells()
Internal subprogram used by Create.
Definition: PDEFoam.cxx:366
void PrintCells()
Prints geometry of ALL cells of the FOAM.
Definition: PDEFoam.cxx:952
virtual Double_t GetBinCenter(Int_t bin) const
return bin center for 1D historam Better to use h1.GetXaxis().GetBinCenter(bin)
Definition: TH1.cxx:8470
void Varedu(Double_t[], Int_t &, Double_t &, Double_t &)
Internal subrogram used by Create.
Definition: PDEFoam.cxx:580
void SetXdiv(Double_t Xdiv)
Definition: PDEFoamCell.h:86
void SetXmax(Int_t idim, Double_t wmax)
set upper foam bound in dimension idim
Definition: PDEFoam.cxx:286
Double_t GetCellElement(const PDEFoamCell *cell, UInt_t i) const
Returns cell element i of cell 'cell'.
Definition: PDEFoam.cxx:1426
PDEFoamCell * GetPare() const
Definition: PDEFoamCell.h:99
void CalcVolume()
Calculates volume of the cell using size params which are calculated.
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8543
void SetStat(Int_t Stat)
Definition: PDEFoamCell.h:98
TVectorT< Element > & Zero()
Set vector elements to zero.
Definition: TVectorT.cxx:450
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
Double_t E()
Definition: TMath.h:54
void SetDriv(Double_t Driv)
Definition: PDEFoamCell.h:95
Double_t GetDriv() const
Definition: PDEFoamCell.h:93
TAxis * GetYaxis()
Definition: TH1.h:320
void SetDim(Int_t kDim)
Sets dimension of cubical space.
Definition: PDEFoam.cxx:260
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
void PrintCell(Long_t iCell=0)
Prints geometry of and elements of 'iCell', as well as relations to parent and daughter cells...
Definition: PDEFoam.cxx:905
REAL epsilon
Definition: triangle.c:617
Int_t GetBest() const
Definition: PDEFoamCell.h:84
Int_t GetColorPalette(Int_t i) const
Return color number i in current palette.
Definition: TStyle.cxx:730
long Long_t
Definition: RtypesCore.h:50
Color * colors
Definition: X3DBuffer.c:19
Double_t Eval(Double_t *xRand, Double_t &event_density)
Internal subprogram.
Definition: PDEFoam.cxx:757
static const double x1[5]
double Double_t
Definition: RtypesCore.h:55
ECellValue
Definition: PDEFoam.h:85
void Create()
Basic initialization of FOAM invoked by the user.
Definition: PDEFoam.cxx:302
Double_t y[n]
Definition: legend1.C:17
virtual TH2D * Project2(Int_t idim1, Int_t idim2, ECellValue cell_value=kValue, PDEFoamKernelBase *kernel=NULL, UInt_t nbin=50)
Project foam variable idim1 and variable idim2 to histogram.
Definition: PDEFoam.cxx:1280
ClassImp(TMVA::PDEFoam) static const Float_t kHigh
virtual Float_t GetCellValue(const std::vector< Float_t > &xvec, ECellValue cv, PDEFoamKernelBase *)
This function finds the cell, which corresponds to the given untransformed event vector 'xvec' and re...
Definition: PDEFoam.cxx:1026
#define name(a, b)
Definition: linkTestLib0.cpp:5
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:449
Mother of all ROOT objects.
Definition: TObject.h:58
virtual Int_t GetNbinsY() const
Definition: TH1.h:297
void SetIntg(Double_t Intg)
Definition: PDEFoamCell.h:94
Abstract ClassifierFactory template that handles arbitrary types.
void DeleteBinarySearchTree()
Delete the foam's density estimator, which contains the binary search tree.
Definition: PDEFoam.cxx:1676
std::vector< Float_t > & GetValues()
Definition: Event.h:93
void FillBinarySearchTree(const Event *ev)
Insert event to internal foam's density estimator PDEFoamDensityBase.
Definition: PDEFoam.cxx:1667
PDEFoamCell * GetDau1() const
Definition: PDEFoamCell.h:101
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
virtual ~PDEFoam()
Default destructor.
Definition: PDEFoam.cxx:194
#define NULL
Definition: Rtypes.h:82
PDEFoamCell * FindCell(const std::vector< Float_t > &) const
Find cell that contains 'xvec' (in foam coordinates [0,1]).
Definition: PDEFoam.cxx:1090
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2689
void SetDau0(PDEFoamCell *Daug)
Definition: PDEFoamCell.h:102
static const Float_t kVlow
Definition: PDEFoam.cxx:106
const Bool_t kTRUE
Definition: Rtypes.h:91
void GetHcub(PDEFoamVect &, PDEFoamVect &) const
Provides size and position of the cell These parameter are calculated by analyzing information in all...
float value
Definition: math.cpp:443
void SetPalette(Int_t ncolors=kBird, Int_t *colors=0, Float_t alpha=1.)
See TColor::SetPalette.
Definition: TStyle.cxx:1445
Definition: math.cpp:60
TAxis * GetXaxis()
Definition: TH1.h:319
void SetInhiDiv(Int_t, Int_t)
This can be called before Create, after setting kDim It defines which variables are excluded in the p...
Definition: PDEFoam.cxx:812
PDEFoam()
Default constructor for streamer, user should not use it.
Definition: PDEFoam.cxx:113
PDEFoamCell * GetDau0() const
Definition: PDEFoamCell.h:100
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:297