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