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