Logo ROOT  
Reference Guide
PDEFoamEvent.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: PDEFoamEvent *
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 /*! \class TMVA::PDEFoamEvent
30 \ingroup TMVA
31 This PDEFoam variant stores in every cell the sum of event weights
32 and the sum of the squared event weights. It therefore acts as
33 event density estimator. It should be booked together with the
34 PDEFoamEventDensity density estimator, which returns the event
35 weight density at a given phase space point during the foam
36 build-up.
37 */
38 
39 #include "TMVA/PDEFoamEvent.h"
40 
41 #include "TMVA/Event.h"
42 #include "TMVA/MsgLogger.h"
43 #include "TMVA/PDEFoam.h"
44 #include "TMVA/Types.h"
45 
46 #include "Rtypes.h"
47 
48 namespace TMVA {
49  class PDEFoamCell;
50 }
51 class TString;
52 
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 /// Default constructor for streamer, user should not use it.
57 
59 : PDEFoam()
60 {
61 }
62 
63 ////////////////////////////////////////////////////////////////////////////////
64 
66  : PDEFoam(name)
67 {}
68 
69 ////////////////////////////////////////////////////////////////////////////////
70 /// Copy Constructor NOT IMPLEMENTED (NEVER USED)
71 
73  : PDEFoam(from)
74 {
75  Log() << kFATAL << "COPY CONSTRUCTOR NOT IMPLEMENTED" << Endl;
76 }
77 
78 ////////////////////////////////////////////////////////////////////////////////
79 /// This function fills an event weight 'wt' into the PDEFoam. Cell
80 /// element 0 is filled with the weight 'wt', and element 1 is
81 /// filled with the squared weight.
82 
84 {
85  // find corresponding foam cell
86  std::vector<Float_t> values = ev->GetValues();
87  std::vector<Float_t> tvalues = VarTransform(values);
88  PDEFoamCell *cell = FindCell(tvalues);
89 
90  // 0. Element: Sum of event weights 'wt'
91  // 1. Element: Sum of squared event weights 'wt'
92  SetCellElement(cell, 0, GetCellElement(cell, 0) + wt);
93  SetCellElement(cell, 1, GetCellElement(cell, 1) + wt * wt);
94 }
ClassImp
#define ClassImp(name)
Definition: Rtypes.h:364
TMVA::PDEFoamEvent::FillFoamCells
virtual void FillFoamCells(const Event *ev, Float_t wt)
This function fills an event weight 'wt' into the PDEFoam.
Definition: PDEFoamEvent.cxx:83
Float_t
float Float_t
Definition: RtypesCore.h:57
TString
Basic string class.
Definition: TString.h:136
TMVA::PDEFoam
Implementation of PDEFoam.
Definition: PDEFoam.h:79
MsgLogger.h
TMVA::PDEFoamEvent
This PDEFoam variant stores in every cell the sum of event weights and the sum of the squared event w...
Definition: PDEFoamEvent.h:39
TMVA::Event::GetValues
std::vector< Float_t > & GetValues()
Definition: Event.h:94
PDEFoamEvent.h
Event.h
PDEFoam.h
Types.h
TMVA::Endl
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:158
TMVA::PDEFoam::Log
MsgLogger & Log() const
Definition: PDEFoam.h:240
TMVA::PDEFoamEvent::PDEFoamEvent
PDEFoamEvent()
Default constructor for streamer, user should not use it.
Definition: PDEFoamEvent.cxx:58
TMVA::Event
Definition: Event.h:51
name
char name[80]
Definition: TGX11.cxx:110
Rtypes.h
TMVA::PDEFoamCell
Definition: PDEFoamCell.h:41
TMVA
create variable transformations
Definition: GeneticMinimizer.h:22