ROOT  6.06/09
Reference Guide
PDEFoamDiscriminant.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Tancredi Carli, Dominik Dannheim, Alexander Voigt
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Classes: PDEFoamDiscriminant *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementation. *
12  * *
13  * Authors (alphabetical): *
14  * Tancredi Carli - CERN, Switzerland *
15  * Dominik Dannheim - CERN, Switzerland *
16  * S. Jadach - Institute of Nuclear Physics, Cracow, Poland *
17  * Alexander Voigt - TU Dresden, Germany *
18  * Peter Speckmayer - CERN, Switzerland *
19  * *
20  * Copyright (c) 2008, 2010: *
21  * CERN, Switzerland *
22  * MPI-K Heidelberg, Germany *
23  * *
24  * Redistribution and use in source and binary forms, with or without *
25  * modification, are permitted according to the terms listed in LICENSE *
26  * (http://tmva.sourceforge.net/LICENSE) *
27  **********************************************************************************/
28 
29 //_____________________________________________________________________
30 //
31 // PDEFoamDiscriminant
32 //
33 // This PDEFoam variant stores in every cell the discriminant
34 //
35 // D = #events with given class / total number of events
36 //
37 // as well as the statistical error on the discriminant. It therefore
38 // acts as a discriminant estimator. It should be booked together
39 // with the PDEFoamDiscriminantDensity density estimator, which
40 // returns the discriminant density at a given phase space point
41 // during the foam build-up.
42 //
43 //_____________________________________________________________________
44 
45 #include <climits>
46 
47 #ifndef ROOT_TMath
48 #include "TMath.h"
49 #endif
50 
51 #ifndef ROOT_TMVA_PDEFoamDiscriminant
53 #endif
54 
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 /// Default constructor for streamer, user should not use it.
59 
60 TMVA::PDEFoamDiscriminant::PDEFoamDiscriminant()
61  : PDEFoam()
62  , fClass(0)
63 {
64 }
65 
66 ////////////////////////////////////////////////////////////////////////////////
67 
69  : PDEFoam(name)
70  , fClass(cls)
71 {}
72 
73 ////////////////////////////////////////////////////////////////////////////////
74 /// Copy Constructor NOT IMPLEMENTED (NEVER USED)
75 
77  : PDEFoam(from)
78  , fClass(from.fClass)
79 {
80  Log() << kFATAL << "COPY CONSTRUCTOR NOT IMPLEMENTED" << Endl;
81 }
82 
83 ////////////////////////////////////////////////////////////////////////////////
84 /// This function fills an event into the discriminant PDEFoam. The
85 /// event weight 'wt' is filled into cell element 0 if the event is
86 /// of class fClass, and filled into cell element 1 otherwise.
87 
89 {
90  // find corresponding foam cell
91  std::vector<Float_t> values = ev->GetValues();
92  std::vector<Float_t> tvalues = VarTransform(values);
93  PDEFoamCell *cell = FindCell(tvalues);
94 
95  // 0. Element: Number of signal events (even class == fClass)
96  // 1. Element: Number of background events times normalization
97  if (ev->GetClass() == fClass)
98  SetCellElement(cell, 0, GetCellElement(cell, 0) + wt);
99  else
100  SetCellElement(cell, 1, GetCellElement(cell, 1) + wt);
101 }
102 
103 ////////////////////////////////////////////////////////////////////////////////
104 /// Calc discriminator and its error for every cell and save it to
105 /// the cell.
106 
108 {
109  // loop over cells
110  for (Long_t iCell = 0; iCell <= fLastCe; iCell++) {
111  if (!(fCells[iCell]->GetStat()))
112  continue;
113 
114  Double_t n_sig = GetCellElement(fCells[iCell], 0); // get number of signal events
115  Double_t n_bg = GetCellElement(fCells[iCell], 1); // get number of bg events
116 
117  if (n_sig < 0.) {
118  Log() << kWARNING << "Negative number of signal events in cell " << iCell
119  << ": " << n_sig << ". Set to 0." << Endl;
120  n_sig = 0.;
121  }
122  if (n_bg < 0.) {
123  Log() << kWARNING << "Negative number of background events in cell " << iCell
124  << ": " << n_bg << ". Set to 0." << Endl;
125  n_bg = 0.;
126  }
127 
128  // calculate discriminant
129  if (n_sig + n_bg > 0) {
130  // discriminant
131  SetCellElement(fCells[iCell], 0, n_sig / (n_sig + n_bg));
132  // discriminant error
133  SetCellElement(fCells[iCell], 1, TMath::Sqrt(Sqr(n_sig / Sqr(n_sig + n_bg))*n_sig +
134  Sqr(n_bg / Sqr(n_sig + n_bg))*n_bg));
135 
136  } else {
137  SetCellElement(fCells[iCell], 0, 0.5); // set discriminator
138  SetCellElement(fCells[iCell], 1, 1.); // set discriminator error
139  }
140  }
141 }
142 
143 ////////////////////////////////////////////////////////////////////////////////
144 /// Project foam variable idim1 and variable idim2 to histogram.
145 /// The projection algorithm is modified such that the z axis range
146 /// of the returned histogram is [0, 1], as necessary for the
147 /// interpretation as a discriminator. This is done by weighting
148 /// the cell values (in case of cell_value = kValue) by the cell
149 /// volume in all dimensions, excluding 'idim1' and 'idim2'.
150 ///
151 /// Parameters:
152 ///
153 /// - idim1, idim2 - dimensions to project to
154 ///
155 /// - cell_value - the cell value to draw
156 ///
157 /// - kernel - a PDEFoam kernel (optional). If NULL is given, the
158 /// kernel is ignored and the pure cell values are
159 /// plotted.
160 ///
161 /// - nbin - number of bins in x and y direction of result histogram
162 /// (optional, default is 50).
163 ///
164 /// Returns:
165 /// a 2-dimensional histogram
166 
168 {
169  // avoid plotting of wrong dimensions
170  if ((idim1 >= GetTotDim()) || (idim1 < 0) ||
171  (idim2 >= GetTotDim()) || (idim2 < 0) ||
172  (idim1 == idim2))
173  Log() << kFATAL << "<Project2>: wrong dimensions given: "
174  << idim1 << ", " << idim2 << Endl;
175 
176  // root can not handle too many bins in one histogram --> catch this
177  // Furthermore, to have more than 1000 bins in the histogram doesn't make
178  // sense.
179  if (nbin > 1000) {
180  Log() << kWARNING << "Warning: number of bins too big: " << nbin
181  << " Using 1000 bins for each dimension instead." << Endl;
182  nbin = 1000;
183  } else if (nbin < 1) {
184  Log() << kWARNING << "Wrong bin number: " << nbin
185  << "; set nbin=50" << Endl;
186  nbin = 50;
187  }
188 
189  // create result histogram
190  TString hname(Form("h_%d_vs_%d", idim1, idim2));
191 
192  // if histogram with this name already exists, delete it
193  TH2D* h1 = (TH2D*)gDirectory->Get(hname.Data());
194  if (h1) delete h1;
195  h1 = new TH2D(hname.Data(), Form("var%d vs var%d", idim1, idim2), nbin, fXmin[idim1], fXmax[idim1], nbin, fXmin[idim2], fXmax[idim2]);
196 
197  if (!h1) Log() << kFATAL << "ERROR: Can not create histo" << hname << Endl;
198  if (cell_value == kValue)
201 
202  // ============== start projection algorithm ================
203  // loop over all histogram bins (2-dim)
204  for (Int_t xbin = 1; xbin <= h1->GetNbinsX(); ++xbin) {
205  for (Int_t ybin = 1; ybin <= h1->GetNbinsY(); ++ybin) {
206  // calculate the phase space point, which corresponds to this
207  // bin combination
208  std::map<Int_t, Float_t> txvec;
209  txvec[idim1] = VarTransform(idim1, h1->GetXaxis()->GetBinCenter(xbin));
210  txvec[idim2] = VarTransform(idim2, h1->GetYaxis()->GetBinCenter(ybin));
211 
212  // find the cells, which corresponds to this phase space
213  // point
214  std::vector<TMVA::PDEFoamCell*> cells = FindCells(txvec);
215 
216  // loop over cells and fill the histogram with the cell
217  // values
218  Float_t sum_cv = 0; // sum of the cell values
219  for (std::vector<TMVA::PDEFoamCell*>::const_iterator it = cells.begin();
220  it != cells.end(); ++it) {
221  // get cell position and size
222  PDEFoamVect cellPosi(GetTotDim()), cellSize(GetTotDim());
223  (*it)->GetHcub(cellPosi, cellSize);
224  // Create complete event vector from txvec. The missing
225  // coordinates of txvec are set to the cell center.
226  std::vector<Float_t> tvec;
227  for (Int_t i = 0; i < GetTotDim(); ++i) {
228  if (i != idim1 && i != idim2)
229  tvec.push_back(cellPosi[i] + 0.5 * cellSize[i]);
230  else
231  tvec.push_back(txvec[i]);
232  }
233  // get the cell value using the kernel
234  Float_t cv = 0;
235  if (kernel != NULL) {
236  cv = kernel->Estimate(this, tvec, cell_value);
237  } else {
238  cv = GetCellValue(FindCell(tvec), cell_value);
239  }
240  if (cell_value == kValue) {
241  // calculate cell volume in other dimensions (not
242  // including idim1 and idim2)
243  Float_t area_cell = 1.;
244  for (Int_t d1 = 0; d1 < GetTotDim(); ++d1) {
245  if ((d1 != idim1) && (d1 != idim2))
246  area_cell *= cellSize[d1];
247  }
248  // calc discriminator * (cell area times foam area)
249  // foam is normalized -> length of foam = 1.0
250  cv *= area_cell;
251  }
252  sum_cv += cv;
253  }
254 
255  // fill the bin content
256  h1->SetBinContent(xbin, ybin, sum_cv + h1->GetBinContent(xbin, ybin));
257  }
258  }
259 
260  return h1;
261 }
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
THist< 2, double > TH2D
Definition: THist.h:320
float Float_t
Definition: RtypesCore.h:53
ClassImp(TMVA::PDEFoamDiscriminant) TMVA
Default constructor for streamer, user should not use it.
MsgLogger & Log() const
Definition: PDEFoam.h:265
#define gDirectory
Definition: TDirectory.h:218
Basic string class.
Definition: TString.h:137
int Int_t
Definition: RtypesCore.h:41
virtual TH2D * Project2(Int_t, Int_t, ECellValue, PDEFoamKernelBase *, UInt_t)
Project foam variable idim1 and variable idim2 to histogram.
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
const char * Data() const
Definition: TString.h:349
virtual void SetRangeUser(Double_t ufirst, Double_t ulast)
Set the viewing range for the axis from ufirst to ulast (in user coordinates).
Definition: TAxis.cxx:869
virtual Float_t Estimate(PDEFoam *, std::vector< Float_t > &, ECellValue)=0
TClass * fClass
pointer to the foreign object
TH1F * h1
Definition: legend1.C:5
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:90
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
TAxis * GetYaxis()
Definition: TH1.h:320
REAL epsilon
Definition: triangle.c:617
long Long_t
Definition: RtypesCore.h:50
double Double_t
Definition: RtypesCore.h:55
ECellValue
Definition: PDEFoam.h:85
UInt_t GetClass() const
Definition: Event.h:86
virtual void Finalize()
Calc discriminator and its error for every cell and save it to the cell.
#define name(a, b)
Definition: linkTestLib0.cpp:5
TAxis * GetZaxis()
Definition: TH1.h:321
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:449
virtual void FillFoamCells(const Event *ev, Float_t wt)
This function fills an event into the discriminant PDEFoam.
virtual Int_t GetNbinsY() const
Definition: TH1.h:297
Abstract ClassifierFactory template that handles arbitrary types.
std::vector< Float_t > & GetValues()
Definition: Event.h:93
#define NULL
Definition: Rtypes.h:82
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2689
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
Definition: math.cpp:60
TAxis * GetXaxis()
Definition: TH1.h:319
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:297