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