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