Logo ROOT   6.07/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 
46 
47 
48 #include "TMVA/Event.h"
49 #include "TMVA/MsgLogger.h"
50 #include "TMVA/PDEFoam.h"
51 #include "TMVA/PDEFoamCell.h"
52 #include "TMVA/PDEFoamKernelBase.h"
53 #include "TMVA/Types.h"
54 
55 #include "TDirectory.h"
56 #include "TMath.h"
57 #include "TH2D.h"
58 
59 #include <climits>
60 
62 
63 ////////////////////////////////////////////////////////////////////////////////
64 /// Default constructor for streamer, user should not use it.
65 
67 : PDEFoam()
68  , fClass(0)
69 {
70 }
71 
72 ////////////////////////////////////////////////////////////////////////////////
73 
75  : PDEFoam(name)
76  , fClass(cls)
77 {}
78 
79 ////////////////////////////////////////////////////////////////////////////////
80 /// Copy Constructor NOT IMPLEMENTED (NEVER USED)
81 
83  : PDEFoam(from)
84  , fClass(from.fClass)
85 {
86  Log() << kFATAL << "COPY CONSTRUCTOR NOT IMPLEMENTED" << Endl;
87 }
88 
89 ////////////////////////////////////////////////////////////////////////////////
90 /// This function fills an event into the discriminant PDEFoam. The
91 /// event weight 'wt' is filled into cell element 0 if the event is
92 /// of class fClass, and filled into cell element 1 otherwise.
93 
95 {
96  // find corresponding foam cell
97  std::vector<Float_t> values = ev->GetValues();
98  std::vector<Float_t> tvalues = VarTransform(values);
99  PDEFoamCell *cell = FindCell(tvalues);
100 
101  // 0. Element: Number of signal events (even class == fClass)
102  // 1. Element: Number of background events times normalization
103  if (ev->GetClass() == fClass)
104  SetCellElement(cell, 0, GetCellElement(cell, 0) + wt);
105  else
106  SetCellElement(cell, 1, GetCellElement(cell, 1) + wt);
107 }
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 /// Calc discriminator and its error for every cell and save it to
111 /// the cell.
112 
114 {
115  // loop over cells
116  for (Long_t iCell = 0; iCell <= fLastCe; iCell++) {
117  if (!(fCells[iCell]->GetStat()))
118  continue;
119 
120  Double_t n_sig = GetCellElement(fCells[iCell], 0); // get number of signal events
121  Double_t n_bg = GetCellElement(fCells[iCell], 1); // get number of bg events
122 
123  if (n_sig < 0.) {
124  Log() << kWARNING << "Negative number of signal events in cell " << iCell
125  << ": " << n_sig << ". Set to 0." << Endl;
126  n_sig = 0.;
127  }
128  if (n_bg < 0.) {
129  Log() << kWARNING << "Negative number of background events in cell " << iCell
130  << ": " << n_bg << ". Set to 0." << Endl;
131  n_bg = 0.;
132  }
133 
134  // calculate discriminant
135  if (n_sig + n_bg > 0) {
136  // discriminant
137  SetCellElement(fCells[iCell], 0, n_sig / (n_sig + n_bg));
138  // discriminant error
139  SetCellElement(fCells[iCell], 1, TMath::Sqrt(Sqr(n_sig / Sqr(n_sig + n_bg))*n_sig +
140  Sqr(n_bg / Sqr(n_sig + n_bg))*n_bg));
141 
142  } else {
143  SetCellElement(fCells[iCell], 0, 0.5); // set discriminator
144  SetCellElement(fCells[iCell], 1, 1.); // set discriminator error
145  }
146  }
147 }
148 
149 ////////////////////////////////////////////////////////////////////////////////
150 /// Project foam variable idim1 and variable idim2 to histogram.
151 /// The projection algorithm is modified such that the z axis range
152 /// of the returned histogram is [0, 1], as necessary for the
153 /// interpretation as a discriminator. This is done by weighting
154 /// the cell values (in case of cell_value = kValue) by the cell
155 /// volume in all dimensions, excluding 'idim1' and 'idim2'.
156 ///
157 /// Parameters:
158 ///
159 /// - idim1, idim2 - dimensions to project to
160 ///
161 /// - cell_value - the cell value to draw
162 ///
163 /// - kernel - a PDEFoam kernel (optional). If NULL is given, the
164 /// kernel is ignored and the pure cell values are
165 /// plotted.
166 ///
167 /// - nbin - number of bins in x and y direction of result histogram
168 /// (optional, default is 50).
169 ///
170 /// Returns:
171 /// a 2-dimensional histogram
172 
174 {
175  // avoid plotting of wrong dimensions
176  if ((idim1 >= GetTotDim()) || (idim1 < 0) ||
177  (idim2 >= GetTotDim()) || (idim2 < 0) ||
178  (idim1 == idim2))
179  Log() << kFATAL << "<Project2>: wrong dimensions given: "
180  << idim1 << ", " << idim2 << Endl;
181 
182  // root can not handle too many bins in one histogram --> catch this
183  // Furthermore, to have more than 1000 bins in the histogram doesn't make
184  // sense.
185  if (nbin > 1000) {
186  Log() << kWARNING << "Warning: number of bins too big: " << nbin
187  << " Using 1000 bins for each dimension instead." << Endl;
188  nbin = 1000;
189  } else if (nbin < 1) {
190  Log() << kWARNING << "Wrong bin number: " << nbin
191  << "; set nbin=50" << Endl;
192  nbin = 50;
193  }
194 
195  // create result histogram
196  TString hname(Form("h_%d_vs_%d", idim1, idim2));
197 
198  // if histogram with this name already exists, delete it
199  TH2D* h1 = (TH2D*)gDirectory->Get(hname.Data());
200  if (h1) delete h1;
201  h1 = new TH2D(hname.Data(), Form("var%d vs var%d", idim1, idim2), nbin, fXmin[idim1], fXmax[idim1], nbin, fXmin[idim2], fXmax[idim2]);
202 
203  if (!h1) Log() << kFATAL << "ERROR: Can not create histo" << hname << Endl;
204  if (cell_value == kValue)
207 
208  // ============== start projection algorithm ================
209  // loop over all histogram bins (2-dim)
210  for (Int_t xbin = 1; xbin <= h1->GetNbinsX(); ++xbin) {
211  for (Int_t ybin = 1; ybin <= h1->GetNbinsY(); ++ybin) {
212  // calculate the phase space point, which corresponds to this
213  // bin combination
214  std::map<Int_t, Float_t> txvec;
215  txvec[idim1] = VarTransform(idim1, h1->GetXaxis()->GetBinCenter(xbin));
216  txvec[idim2] = VarTransform(idim2, h1->GetYaxis()->GetBinCenter(ybin));
217 
218  // find the cells, which corresponds to this phase space
219  // point
220  std::vector<TMVA::PDEFoamCell*> cells = FindCells(txvec);
221 
222  // loop over cells and fill the histogram with the cell
223  // values
224  Float_t sum_cv = 0; // sum of the cell values
225  for (std::vector<TMVA::PDEFoamCell*>::const_iterator it = cells.begin();
226  it != cells.end(); ++it) {
227  // get cell position and size
228  PDEFoamVect cellPosi(GetTotDim()), cellSize(GetTotDim());
229  (*it)->GetHcub(cellPosi, cellSize);
230  // Create complete event vector from txvec. The missing
231  // coordinates of txvec are set to the cell center.
232  std::vector<Float_t> tvec;
233  for (Int_t i = 0; i < GetTotDim(); ++i) {
234  if (i != idim1 && i != idim2)
235  tvec.push_back(cellPosi[i] + 0.5 * cellSize[i]);
236  else
237  tvec.push_back(txvec[i]);
238  }
239  // get the cell value using the kernel
240  Float_t cv = 0;
241  if (kernel != NULL) {
242  cv = kernel->Estimate(this, tvec, cell_value);
243  } else {
244  cv = GetCellValue(FindCell(tvec), cell_value);
245  }
246  if (cell_value == kValue) {
247  // calculate cell volume in other dimensions (not
248  // including idim1 and idim2)
249  Float_t area_cell = 1.;
250  for (Int_t d1 = 0; d1 < GetTotDim(); ++d1) {
251  if ((d1 != idim1) && (d1 != idim2))
252  area_cell *= cellSize[d1];
253  }
254  // calc discriminator * (cell area times foam area)
255  // foam is normalized -> length of foam = 1.0
256  cv *= area_cell;
257  }
258  sum_cv += cv;
259  }
260 
261  // fill the bin content
262  h1->SetBinContent(xbin, ybin, sum_cv + h1->GetBinContent(xbin, ybin));
263  }
264  }
265 
266  return h1;
267 }
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
PDEFoamDiscriminant()
Default constructor for streamer, user should not use it.
float Float_t
Definition: RtypesCore.h:53
MsgLogger & Log() const
Definition: PDEFoam.h:262
Double_t * fXmax
Definition: PDEFoam.h:129
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:301
void SetCellElement(PDEFoamCell *cell, UInt_t i, Double_t value)
Set cell element i of cell to value.
Definition: PDEFoam.cxx:1432
std::vector< TMVA::PDEFoamCell * > FindCells(const std::vector< Float_t > &) const
Find all cells, that contain txvec.
Definition: PDEFoam.cxx:1165
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:925
virtual Float_t Estimate(PDEFoam *, std::vector< Float_t > &, ECellValue)=0
Int_t GetTotDim() const
Definition: PDEFoam.h:219
TH1F * h1
Definition: legend1.C:5
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:90
Double_t GetCellElement(const PDEFoamCell *cell, UInt_t i) const
Returns cell element i of cell &#39;cell&#39;.
Definition: PDEFoam.cxx:1416
T Sqr(T x) const
Definition: PDEFoam.h:183
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
TAxis * GetYaxis()
Definition: TH1.h:325
Float_t VarTransform(Int_t idim, Float_t x) const
Definition: PDEFoam.h:303
Double_t * fXmin
Definition: PDEFoam.h:128
REAL epsilon
Definition: triangle.c:617
long Long_t
Definition: RtypesCore.h:50
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
ECellValue
Definition: PDEFoam.h:85
PDEFoamCell ** fCells
Definition: PDEFoam.h:118
THist< 2, double, THistStatContent, THistStatUncertainty > TH2D
Definition: THist.hxx:307
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 &#39;xvec&#39; and re...
Definition: PDEFoam.cxx:1016
UInt_t GetClass() const
Definition: Event.h:89
virtual void Finalize()
Calc discriminator and its error for every cell and save it to the cell.
TAxis * GetZaxis()
Definition: TH1.h:326
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:464
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:302
Abstract ClassifierFactory template that handles arbitrary types.
std::vector< Float_t > & GetValues()
Definition: Event.h:96
#define NULL
Definition: Rtypes.h:82
Int_t fLastCe
Definition: PDEFoam.h:117
PDEFoamCell * FindCell(const std::vector< Float_t > &) const
Find cell that contains &#39;xvec&#39; (in foam coordinates [0,1]).
Definition: PDEFoam.cxx:1080
#define gDirectory
Definition: TDirectory.h:221
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2475
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
char name[80]
Definition: TGX11.cxx:109
TAxis * GetXaxis()
Definition: TH1.h:324
tomato 2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:296