PDEFoamDistr.cxx

Go to the documentation of this file.
00001 
00002 /**********************************************************************************
00003  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis       *
00004  * Package: TMVA                                                                  *
00005  * Classes: PDEFoamDistr                                                          *
00006  * Web    : http://tmva.sourceforge.net                                           *
00007  *                                                                                *
00008  * Description:                                                                   *
00009  *      The TFDSITR class provides an interface between the Binary search tree    *
00010  *      and the PDEFoam object.  In order to build-up the foam one needs to       *
00011  *      calculate the density of events at a given point (sampling during         *
00012  *      Foam build-up).  The function PDEFoamDistr::Density() does this job.  It  *
00013  *      uses a binary search tree, filled with training events, in order to       *
00014  *      provide this density.                                                     *
00015  *                                                                                *
00016  * Authors (alphabetical):                                                        *
00017  *      Tancredi Carli   - CERN, Switzerland                                      *
00018  *      Dominik Dannheim - CERN, Switzerland                                      *
00019  *      S. Jadach        - Institute of Nuclear Physics, Cracow, Poland           *
00020  *      Alexander Voigt  - CERN, Switzerland                                      *
00021  *      Peter Speckmayer - CERN, Switzerland                                      *
00022  *                                                                                *
00023  * Copyright (c) 2008:                                                            *
00024  *      CERN, Switzerland                                                         *
00025  *      MPI-K Heidelberg, Germany                                                 *
00026  *                                                                                *
00027  * Redistribution and use in source and binary forms, with or without             *
00028  * modification, are permitted according to the terms listed in LICENSE           *
00029  * (http://tmva.sourceforge.net/LICENSE)                                          *
00030  **********************************************************************************/
00031 
00032 #include <cmath>
00033 #include <limits>
00034 
00035 #ifndef ROOT_TMath
00036 #include "TMath.h"
00037 #endif
00038 
00039 #ifndef ROOT_TMVA_PDEFoamDistr
00040 #include "TMVA/PDEFoamDistr.h"
00041 #endif
00042 
00043 ClassImp(TMVA::PDEFoamDistr)
00044 
00045 //_____________________________________________________________________
00046 TMVA::PDEFoamDistr::PDEFoamDistr() 
00047    : TObject(),
00048      fPDEFoam(NULL),
00049      fBst(NULL),
00050      fDensityCalc(kEVENT_DENSITY), // default: fill event density to BinarySearchTree
00051      fLogger( new MsgLogger("PDEFoamDistr"))
00052 {}
00053 
00054 //_____________________________________________________________________
00055 TMVA::PDEFoamDistr::~PDEFoamDistr() 
00056 {
00057    if (fBst)  delete fBst;
00058    delete fLogger;
00059 }
00060 
00061 //_____________________________________________________________________
00062 TMVA::PDEFoamDistr::PDEFoamDistr(const PDEFoamDistr &distr)
00063    : TObject(),
00064      fPDEFoam         (distr.fPDEFoam),
00065      fBst             (distr.fBst),
00066      fDensityCalc     (kEVENT_DENSITY), // default: fill event density to BinarySearchTree
00067      fLogger( new MsgLogger("PDEFoamDistr"))
00068 {
00069    // Copy constructor
00070    Log() << kFATAL << "COPY CONSTRUCTOR NOT IMPLEMENTED" << Endl;
00071 }
00072 
00073 //_____________________________________________________________________
00074 void TMVA::PDEFoamDistr::Initialize()
00075 {
00076    // Initialisation of binary search tree.  
00077    // Set dimension and create new BinarySearchTree.
00078 
00079    if (!GetPDEFoam())
00080       Log() << kFATAL << "<PDEFoamDistr::Initialize()> Pointer to owner not set!" << Endl;
00081 
00082    if (fBst) delete fBst;
00083    fBst = new TMVA::BinarySearchTree();
00084 
00085    if (!fBst){
00086       Log() << kFATAL << "<PDEFoamDistr::Initialize> "
00087             << "ERROR: an not create binary tree !" << Endl;
00088    }
00089 
00090    // set periode (number of variables)
00091    fBst->SetPeriode(GetPDEFoam()->GetTotDim());
00092 }
00093 
00094 //_____________________________________________________________________
00095 void TMVA::PDEFoamDistr::FillBinarySearchTree( const Event* ev, EFoamType ft, Bool_t NoNegWeights )
00096 {
00097    // This method creates an TMVA::Event and inserts it into the
00098    // binary search tree.
00099    //
00100    // If 'NoNegWeights' is true, an event with negative weight will
00101    // not be filled into the foam.  (Default value: false)
00102 
00103    if((NoNegWeights && ev->GetWeight()<=0) || ev->GetOriginalWeight()==0)
00104       return;
00105 
00106    TMVA::Event *event = new TMVA::Event(*ev);
00107  
00108    // set event variables in case of multi-target regression
00109    if (ft==kMultiTarget){
00110       // since in multi target regression targets are handled like
00111       // variables, remove targets and add them to the event variabels
00112       std::vector<Float_t> targets = ev->GetTargets();
00113       for (UInt_t i = 0; i < targets.size(); i++)
00114          event->SetVal(i+ev->GetValues().size(), targets.at(i));
00115       event->GetTargets().clear();
00116    }
00117    fBst->Insert(event);
00118 
00119    delete event;
00120 }
00121 
00122 //_____________________________________________________________________
00123 Double_t TMVA::PDEFoamDistr::Density( Double_t *Xarg, Double_t &event_density )
00124 {
00125    // This function is needed during the foam buildup.
00126    // It return a certain density depending on the selected classification
00127    // or regression options:
00128    //
00129    // In case of separated foams (classification) or multi target regression:
00130    //  - returns event density within volume (specified by VolFrac)
00131    // In case of unified foams: (classification)
00132    //  - returns discriminator (N_sig)/(N_sig + N_bg) divided by volume
00133    //    (specified by VolFrac)
00134    // In case of mono target regression:
00135    //  - returns average target value within volume divided by volume
00136    //    (specified by VolFrac)
00137 
00138    if (!GetPDEFoam())
00139       Log() << kFATAL << "<PDEFoamDistr::Density()> Pointer to owner not set!" << Endl;
00140 
00141    if (!fBst)
00142       Log() << kFATAL << "<PDEFoamDistr::Density()> Binary tree not found!"<< Endl;
00143 
00144    // get PDEFoam properties
00145    Int_t Dim       = GetPDEFoam()->GetTotDim(); // dimension of foam
00146    Float_t VolFrac = GetPDEFoam()->GetVolumeFraction(); // get fVolFrac
00147 
00148    // make the variable Xarg transform, since Foam only knows about x=[0,1]
00149    // transformation [0, 1] --> [xmin, xmax]
00150    for (Int_t idim=0; idim<Dim; idim++)
00151       Xarg[idim] = GetPDEFoam()->VarTransformInvers(idim, Xarg[idim]);
00152 
00153    //create volume around point to be found
00154    std::vector<Double_t> lb(Dim);
00155    std::vector<Double_t> ub(Dim);
00156 
00157    // probevolume relative to hypercube with edge length 1:
00158    const Double_t probevolume_inv = std::pow((VolFrac/2), Dim);
00159 
00160    // set upper and lower bound for search volume
00161    for (Int_t idim = 0; idim < Dim; idim++) {
00162       Double_t volsize=(GetPDEFoam()->GetXmax(idim) 
00163                         - GetPDEFoam()->GetXmin(idim)) / VolFrac;
00164       lb[idim] = Xarg[idim] - volsize;
00165       ub[idim] = Xarg[idim] + volsize;
00166    }
00167 
00168    TMVA::Volume volume(&lb, &ub);                        // volume to search in
00169    std::vector<const TMVA::BinarySearchTreeNode*> nodes; // BST nodes found
00170 
00171    // do range searching
00172    fBst->SearchVolume(&volume, &nodes);
00173 
00174    // normalized density: (number of counted events) / volume / (total
00175    // number of events) should be ~1 on average
00176    const UInt_t count = nodes.size(); // number of events found
00177 
00178    // store density based on total number of events
00179    event_density = count * probevolume_inv;
00180 
00181    Double_t weighted_count = 0.; // number of events found (sum of weights!)
00182    for (UInt_t j=0; j<nodes.size(); j++)
00183       weighted_count += (nodes.at(j))->GetWeight();
00184 
00185    if (FillDiscriminator()){ // calc number of signal events in nodes
00186       Double_t N_sig = 0;    // number of signal events found
00187       // now sum over all nodes->IsSignal;
00188       for (UInt_t j=0; j<count; j++){
00189          if (nodes.at(j)->IsSignal()) N_sig += nodes.at(j)->GetWeight();
00190       }
00191       return (N_sig/(weighted_count+0.1))*probevolume_inv; // return:  (N_sig/N_total) / (cell_volume)
00192    }
00193    else if (FillTarget0()){ // calc sum of weighted target values
00194       Double_t N_tar = 0;   // number of target events found
00195       // now sum over all nodes->GetTarget(0);
00196       for (UInt_t j=0; j<count; j++) {
00197          N_tar += ((nodes.at(j))->GetTargets()).at(0) * ((nodes.at(j))->GetWeight());
00198       }
00199       return (N_tar/(weighted_count+0.1))*probevolume_inv; // return:  (N_tar/N_total) / (cell_volume)
00200    }
00201 
00202    return ((weighted_count+0.1)*probevolume_inv); // return:  N_total(weighted) / cell_volume
00203 }
00204 
00205 //_____________________________________________________________________
00206 void TMVA::PDEFoamDistr::FillHist(PDEFoamCell* cell, std::vector<TH1F*> &hsig, std::vector<TH1F*> &hbkg, std::vector<TH1F*> &hsig_unw, std::vector<TH1F*> &hbkg_unw)
00207 {
00208    // fill the given histograms with signal and background events,
00209    // which are located in the given cell
00210 
00211    if (!GetPDEFoam())
00212       Log() << kFATAL << "<PDEFoamDistr::FillHist> Pointer to owner not set!" << Endl;
00213 
00214    // get PDEFoam properties
00215    Int_t Dim = GetPDEFoam()->GetTotDim(); // dimension of foam
00216 
00217    // sanity check
00218    if (!cell)
00219       Log() << kFATAL << "<PDEFoamDistr::FillHist> Null pointer for cell given!" << Endl;
00220    if (Int_t(hsig.size()) != Dim || Int_t(hbkg.size()) != Dim || 
00221        Int_t(hsig_unw.size()) != Dim || Int_t(hbkg_unw.size()) != Dim)
00222       Log() << kFATAL << "<PDEFoamDistr::FillHist> Edge histograms have wrong size!" << Endl;
00223 
00224    // check histograms
00225    for (Int_t idim=0; idim<Dim; idim++) {
00226       if (!hsig.at(idim) || !hbkg.at(idim) || 
00227           !hsig_unw.at(idim) || !hbkg_unw.at(idim))
00228          Log() << kFATAL << "<PDEFoamDistr::FillHist> Histogram not initialized!" << Endl;
00229    }
00230 
00231    // get cell position and size
00232    PDEFoamVect  cellSize(Dim);
00233    PDEFoamVect  cellPosi(Dim);
00234    cell->GetHcub(cellPosi, cellSize);
00235 
00236    // determine lower and upper cell bound
00237    std::vector<Double_t> lb(Dim); // lower bound
00238    std::vector<Double_t> ub(Dim); // upper bound
00239    for (Int_t idim = 0; idim < Dim; idim++) {
00240       lb[idim] = GetPDEFoam()->VarTransformInvers(idim, cellPosi[idim] - std::numeric_limits<float>::epsilon());
00241       ub[idim] = GetPDEFoam()->VarTransformInvers(idim, cellPosi[idim] + cellSize[idim] + std::numeric_limits<float>::epsilon());
00242    }
00243 
00244    // create TMVA::Volume object needed for searching within the BST
00245    TMVA::Volume volume(&lb, &ub); // volume to search in
00246    std::vector<const TMVA::BinarySearchTreeNode*> nodes; // BST nodes found
00247 
00248    // do range searching
00249    fBst->SearchVolume(&volume, &nodes);
00250 
00251    // calc xmin and xmax of events found in cell
00252    std::vector<Float_t> xmin(Dim, std::numeric_limits<float>::max());
00253    std::vector<Float_t> xmax(Dim, -std::numeric_limits<float>::max());
00254    for (UInt_t iev=0; iev<nodes.size(); iev++) {
00255       std::vector<Float_t> ev = nodes.at(iev)->GetEventV();
00256       for (Int_t idim=0; idim<Dim; idim++) {
00257          if (ev.at(idim) < xmin.at(idim))  xmin.at(idim) = ev.at(idim);
00258          if (ev.at(idim) > xmax.at(idim))  xmax.at(idim) = ev.at(idim);
00259       }
00260    }
00261 
00262    // reset histogram ranges
00263    for (Int_t idim=0; idim<Dim; idim++) {
00264       hsig.at(idim)->GetXaxis()->SetLimits(GetPDEFoam()->VarTransform(idim,xmin.at(idim)), 
00265                                            GetPDEFoam()->VarTransform(idim,xmax.at(idim)));
00266       hbkg.at(idim)->GetXaxis()->SetLimits(GetPDEFoam()->VarTransform(idim,xmin.at(idim)), 
00267                                            GetPDEFoam()->VarTransform(idim,xmax.at(idim)));
00268       hsig_unw.at(idim)->GetXaxis()->SetLimits(GetPDEFoam()->VarTransform(idim,xmin.at(idim)), 
00269                                                GetPDEFoam()->VarTransform(idim,xmax.at(idim)));
00270       hbkg_unw.at(idim)->GetXaxis()->SetLimits(GetPDEFoam()->VarTransform(idim,xmin.at(idim)), 
00271                                                GetPDEFoam()->VarTransform(idim,xmax.at(idim)));
00272       hsig.at(idim)->Reset();
00273       hbkg.at(idim)->Reset();
00274       hsig_unw.at(idim)->Reset();
00275       hbkg_unw.at(idim)->Reset();
00276    }
00277 
00278    // fill histograms
00279    for (UInt_t iev=0; iev<nodes.size(); iev++) {
00280       std::vector<Float_t> ev = nodes.at(iev)->GetEventV();
00281       Float_t              wt = nodes.at(iev)->GetWeight();
00282       Bool_t           signal = nodes.at(iev)->IsSignal();
00283       for (Int_t idim=0; idim<Dim; idim++) {
00284          if (signal) {
00285             hsig.at(idim)->Fill(GetPDEFoam()->VarTransform(idim,ev.at(idim)), wt);
00286             hsig_unw.at(idim)->Fill(GetPDEFoam()->VarTransform(idim,ev.at(idim)), 1);
00287          } else {
00288             hbkg.at(idim)->Fill(GetPDEFoam()->VarTransform(idim,ev.at(idim)), wt);
00289             hbkg_unw.at(idim)->Fill(GetPDEFoam()->VarTransform(idim,ev.at(idim)), 1);
00290          }
00291       }
00292    }
00293 }

Generated on Tue Jul 5 15:25:13 2011 for ROOT_528-00b_version by  doxygen 1.5.1