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