library: libTMVA
#include "MethodCuts.h"

TMVA::MethodCuts


class description - header file - source file
viewCVS header - viewCVS source

class TMVA::MethodCuts: public TMVA::MethodBase

Inheritance Inherited Members Includes Libraries
Class Charts

Function Members (Methods)

Display options:
Show inherited
Show non-public
public:
virtual~MethodCuts()
voidTObject::AbstractMethod(const char* method) const
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTObject::Clear(Option_t* = "")
virtual TObject*TObject::Clone(const char* newname = "") const
virtual Int_tTObject::Compare(const TObject* obj) const
Double_tComputeEstimator(const vector<Double_t>&)
virtual voidTObject::Copy(TObject& object) const
virtual const TMVA::Ranking*CreateRanking()
TMVA::DataSet&TMVA::MethodBase::Data() const
virtual voidDeclareOptions()
virtual voidTObject::Delete(Option_t* option = "")
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() const
virtual TObject*TObject::DrawClone(Option_t* option = "") const
virtual voidTObject::Dump() const
virtual voidTObject::Error(const char* method, const char* msgfmt) const
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
Double_tTMVA::MethodBase::GetEffForRoot(Double_t)
virtual Double_tGetEfficiency(TString, TTree*)
Double_tTMVA::MethodBase::GetEventVal(Int_t ivar) const
Double_tTMVA::MethodBase::GetEventValNormalized(Int_t ivar) const
Double_tTMVA::MethodBase::GetEventWeight() const
virtual const char*TObject::GetIconName() const
const TString&TMVA::MethodBase::GetInputExp(int i) const
const TString&TMVA::MethodBase::GetInputVar(int i) const
virtual const TString&TMVA::MethodBase::GetJobName() const
virtual const TString&TMVA::MethodBase::GetMethodName() const
virtual const TString&TMVA::MethodBase::GetMethodTitle() const
virtual const TMVA::Types::EMVATMVA::MethodBase::GetMethodType() const
virtual Double_tGetmuTransform(TTree*)
virtual Double_tGetMvaValue()
virtual const char*TMVA::MethodBase::GetName() const
Int_tTMVA::MethodBase::GetNvar() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Double_tTMVA::MethodBase::GetOptimalSignificance(Double_t SignalEvents, Double_t BackgroundEvents, Double_t& optimal_significance_value) const
virtual Option_t*TObject::GetOption() const
TStringTMVA::MethodBase::GetOptions() const
virtual TMVA::Types::EPreprocessingMethodTMVA::MethodBase::GetPreprocessingMethod() const
virtual Double_tGetSeparation()
virtual Double_tGetSignificance()
TTree*TMVA::MethodBase::GetTestTree() const
static TMVA::MethodBase*TMVA::MethodBase::GetThisBase()
virtual const char*TObject::GetTitle() const
virtual Double_tGetTrainingEfficiency(TString)
TTree*TMVA::MethodBase::GetTrainingTree() const
virtual UInt_tTObject::GetUniqueID() const
virtual TStringTMVA::MethodBase::GetWeightFileDir() const
virtual TStringTMVA::MethodBase::GetWeightFileExtension() const
TStringTMVA::MethodBase::GetWeightFileName() const
TMVA::MethodBase::EWeightFileTypeTMVA::MethodBase::GetWeightFileType() const
Double_tTMVA::MethodBase::GetXmax(Int_t ivar, TMVA::Types::EPreprocessingMethod corr = Types::kNone) const
Double_tTMVA::MethodBase::GetXmax(const TString& var, TMVA::Types::EPreprocessingMethod corr = Types::kNone) const
Double_tTMVA::MethodBase::GetXmin(Int_t ivar, TMVA::Types::EPreprocessingMethod corr = Types::kNone) const
Double_tTMVA::MethodBase::GetXmin(const TString& var, TMVA::Types::EPreprocessingMethod corr = Types::kNone) const
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTObject::Hash() const
Bool_tTMVA::MethodBase::HasTrainingTree() const
static Double_tTMVA::MethodBase::IGetEffForRoot(Double_t)
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() const
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
virtual Bool_tTMVA::MethodBase::IsOK() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTMVA::MethodBase::IsSignalLike()
virtual Bool_tTObject::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTObject::ls(Option_t* option = "") const
voidTObject::MayNotUse(const char* method) const
TMVA::MethodCutsMethodCuts(TMVA::DataSet& theData, TString theWeightFile, TDirectory* theTargetDir = NULL)
TMVA::MethodCutsMethodCuts(TString jobName, TString methodTitle, TMVA::DataSet& theData, TString theOption = MC:150:10000:, TDirectory* theTargetFile = 0)
Double_tTMVA::MethodBase::Norm(Int_t ivar, Double_t x) const
Double_tTMVA::MethodBase::Norm(TString var, Double_t x) const
virtual Bool_tTObject::Notify()
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
TMVA::IMethod&TMVA::IMethod::operator=(const TMVA::IMethod&)
virtual voidTObject::Paint(Option_t* option = "")
virtual voidTObject::Pop()
virtual voidTMVA::MethodBase::PrepareEvaluationTree(TTree* theTestTree)
virtual voidTObject::Print(Option_t* option = "") const
virtual voidProcessOptions()
virtual Int_tTObject::Read(const char* name)
virtual voidTMVA::MethodBase::ReadStateFromFile()
virtual voidTMVA::MethodBase::ReadStateFromStream(istream& i)
virtual Bool_tTMVA::MethodBase::ReadTestEvent(UInt_t ievt, TMVA::Types::ESBType type = Types::kMaxSBType)
Bool_tTMVA::MethodBase::ReadTrainingEvent(UInt_t ievt, TMVA::Types::ESBType type = Types::kMaxSBType)
virtual voidReadWeightsFromStream(istream& istr)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") const
virtual voidTObject::SavePrimitive(ostream& out, Option_t* option = "")
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
virtual voidTObject::SetDrawOption(Option_t* option = "")
static voidTObject::SetDtorOnly(void* obj)
virtual voidTMVA::MethodBase::SetJobName(TString jobName)
voidTMVA::MethodBase::SetMethodName(TString methodName)
voidTMVA::MethodBase::SetMethodTitle(TString methodTitle)
voidTMVA::MethodBase::SetMethodType(TMVA::Types::EMVA methodType)
voidTMVA::MethodBase::SetNvar(Int_t n)
static voidTObject::SetObjectStat(Bool_t stat)
voidTMVA::MethodBase::SetPreprocessingMethod(TMVA::Types::EPreprocessingMethod m)
voidSetTestSignalEfficiency(Double_t eff)
virtual voidTObject::SetUniqueID(UInt_t uid)
voidTMVA::MethodBase::SetVerbose(Bool_t v = kTRUE)
virtual voidTMVA::MethodBase::SetWeightFileDir(TString fileDir)
virtual voidTMVA::MethodBase::SetWeightFileExtension(TString fileExtension)
voidTMVA::MethodBase::SetWeightFileName(TString)
voidTMVA::MethodBase::SetWeightFileType(TMVA::MethodBase::EWeightFileType w)
voidTMVA::MethodBase::SetXmax(Int_t ivar, Double_t x, TMVA::Types::EPreprocessingMethod corr = Types::kNone)
voidTMVA::MethodBase::SetXmax(const TString& var, Double_t x, TMVA::Types::EPreprocessingMethod corr = Types::kNone)
voidTMVA::MethodBase::SetXmin(Int_t ivar, Double_t x, TMVA::Types::EPreprocessingMethod corr = Types::kNone)
voidTMVA::MethodBase::SetXmin(const TString& var, Double_t x, TMVA::Types::EPreprocessingMethod corr = Types::kNone)
virtual voidShowMembers(TMemberInspector& insp, char* parent)
virtual voidStreamer(TBuffer& b)
voidStreamerNVirtual(TBuffer& b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
virtual voidTest(TTree* theTestTree)
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual voidTMVA::MethodBase::TestInit(TTree* theTestTree = 0)
static TMVA::MethodCuts*ThisCuts()
virtual voidTrain()
voidTMVA::MethodBase::TrainMethod()
virtual voidTObject::UseCurrentStyle()
Bool_tTMVA::MethodBase::Verbose() const
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = "0", Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = "0", Int_t option = 0, Int_t bufsize = 0) const
virtual voidTMVA::MethodBase::WriteEvaluationHistosToFile(TDirectory* targetDir)
virtual voidWriteMonitoringHistosToFile() const
voidTMVA::MethodBase::WriteStateToFile() const
virtual voidTMVA::MethodBase::WriteStateToStream(ostream& o) const
virtual voidWriteWeightsToStream(ostream& o) const
protected:
TDirectory*TMVA::MethodBase::BaseDir() const
Bool_tTMVA::MethodBase::CheckSanity(TTree* theTree = 0)
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
voidTMVA::MethodBase::EnableLooseOptions(Bool_t b = kTRUE)
TMVA::MethodBase::ECutOrientationTMVA::MethodBase::GetCutOrientation() const
TMVA::Types::ESBTypeTMVA::MethodBase::GetPreprocessingType() const
Double_tTMVA::MethodBase::GetSignalReferenceCut() const
const TString&TMVA::MethodBase::GetTestvarName() const
const TString&TMVA::MethodBase::GetTestvarPrefix() const
const TList&TMVA::MethodBase::ListOfOptions() const
TDirectory*TMVA::MethodBase::LocalTDir() const
voidTObject::MakeZombie()
voidTMVA::MethodBase::ParseOptions(Bool_t verbose = kTRUE)
voidTMVA::MethodBase::PrintOptions() const
voidTMVA::MethodBase::ReadOptionsFromStream(istream& istr)
voidTMVA::MethodBase::ResetThisBase()
voidTMVA::MethodBase::SetPreprocessingType(TMVA::Types::ESBType t)
voidTMVA::MethodBase::SetSignalReferenceCut(Double_t cut)
voidTMVA::MethodBase::SetTestvarName()
voidTMVA::MethodBase::SetTestvarName(TString v)
voidTMVA::MethodBase::SetTestvarPrefix(TString prefix)
voidTMVA::MethodBase::Statistics(TMVA::Types::ETreeType treeType, const TString& theVarName, Double_t&, Double_t&, Double_t&, Double_t&, Double_t&, Double_t&, Bool_t norm = kFALSE)
voidTMVA::MethodBase::WriteOptionsToStream(ostream& o) const
private:
voidCreateVariablePDFs()
voidGetEffsfromPDFs(Double_t* cutMin, Double_t* cutMax, Double_t& effS, Double_t& effB)
voidGetEffsfromSelection(Double_t* cutMin, Double_t* cutMax, Double_t& effS, Double_t& effB)
voidInitCuts()
voidMatchCutsToPars(Double_t*, Double_t*, Double_t*)
voidMatchParsToCuts(const vector<Double_t>&, Double_t*, Double_t*)
voidMatchParsToCuts(Double_t*, Double_t*, Double_t*)
Bool_tSanityChecks()

Data Members

public:
enum EConstrainType { kConstrainEffS
kConstrainEffB
};
enum EFitMethodType { kUseMonteCarlo
kUseGeneticAlgorithm
kUseSimulatedAnnealing
};
enum EEffMethod { kUseEventSelection
kUsePDFs
};
enum EFitParameters { kNotEnforced
kForceMin
kForceMax
kForceSmart
kForceVerySmart
};
enum TMVA::MethodBase::EWeightFileType { kROOT
kTEXT
};
enum TMVA::MethodBase::ECutOrientation { kNegative
kPositive
};
enum TObject::EStatusBits { kCanDelete
kMustCleanup
kObjInCanvas
kIsReferenced
kHasUUID
kCannotPick
kNoContextMenu
kInvalidObject
};
enum TObject::[unnamed] { kIsOnHeap
kNotDeleted
kZombie
kBitMask
kSingleKey
kOverwrite
kWriteDelete
};
protected:
TMVA::Ranking*TMVA::MethodBase::fRankingranking
vector<TString>*TMVA::MethodBase::fInputVarsvector of input variables used in MVA
Bool_tTMVA::MethodBase::fIsOKstatus of sanity checks
TH1*TMVA::MethodBase::fHistS_plotbinMVA plots used for graphics representation (signal)
TH1*TMVA::MethodBase::fHistB_plotbinMVA plots used for graphics representation (background)
TH1*TMVA::MethodBase::fHistS_highbinMVA plots used for efficiency calculations (signal)
TH1*TMVA::MethodBase::fHistB_highbinMVA plots used for efficiency calculations (background)
TH1*TMVA::MethodBase::fEffSefficiency plot (signal)
TH1*TMVA::MethodBase::fEffBefficiency plot (background)
TH1*TMVA::MethodBase::fEffBvsSbackground efficiency versus signal efficiency
TH1*TMVA::MethodBase::fRejBvsSbackground rejection (=1-eff.) versus signal efficiency
TH1*TMVA::MethodBase::fHistBhatSworking histograms needed for mu-transform (signal)
TH1*TMVA::MethodBase::fHistBhatBworking histograms needed for mu-transform (background)
TH1*TMVA::MethodBase::fHistMuSmu-transform (signal)
TH1*TMVA::MethodBase::fHistMuBmu-transform (background)
TH1*TMVA::MethodBase::fTrainEffSTraining efficiency plot (signal)
TH1*TMVA::MethodBase::fTrainEffBTraining efficiency plot (background)
TH1*TMVA::MethodBase::fTrainEffBvsSTraining background efficiency versus signal efficiency
TH1*TMVA::MethodBase::fTrainRejBvsSTraining background rejection (=1-eff.) versus signal efficiency
Double_tTMVA::MethodBase::fX
Double_tTMVA::MethodBase::fMode
TGraph*TMVA::MethodBase::fGraphSgraphs used for splines for efficiency (signal)
TGraph*TMVA::MethodBase::fGraphBgraphs used for splines for efficiency (background)
TGraph*TMVA::MethodBase::fGrapheffBvsSgraphs used for splines for signal eff. versus background eff.
TMVA::PDF*TMVA::MethodBase::fSplSPDFs of MVA distribution (signal)
TMVA::PDF*TMVA::MethodBase::fSplBPDFs of MVA distribution (background)
TSpline*TMVA::MethodBase::fSpleffBvsSsplines for signal eff. versus background eff.
TGraph*TMVA::MethodBase::fGraphTrainSgraphs used for splines for training efficiency (signal)
TGraph*TMVA::MethodBase::fGraphTrainBgraphs used for splines for training efficiency (background)
TGraph*TMVA::MethodBase::fGraphTrainEffBvsSgraphs used for splines for training signal eff. versus background eff.
TMVA::PDF*TMVA::MethodBase::fSplTrainSPDFs of training MVA distribution (signal)
TMVA::PDF*TMVA::MethodBase::fSplTrainBPDFs of training MVA distribution (background)
TSpline*TMVA::MethodBase::fSplTrainEffBvsSsplines for training signal eff. versus background eff.
Int_tTMVA::MethodBase::fNbinsnumber of bins in representative histograms
Int_tTMVA::MethodBase::fNbinsHnumber of bins in evaluation histograms
TMVA::MethodBase::ECutOrientationTMVA::MethodBase::fCutOrientation+1 if Sig>Bkg, -1 otherwise
TMVA::TSpline1*TMVA::MethodBase::fSplRefShelper splines for RootFinder (signal)
TMVA::TSpline1*TMVA::MethodBase::fSplRefBhelper splines for RootFinder (background)
TMVA::TSpline1*TMVA::MethodBase::fSplTrainRefShelper splines for RootFinder (signal)
TMVA::TSpline1*TMVA::MethodBase::fSplTrainRefBhelper splines for RootFinder (background)
TMVA::OptionBase*TMVA::MethodBase::fLastDeclaredOptionlast declared option
TListTMVA::MethodBase::fListOfOptionsoption list
TMVA::MsgLoggerTMVA::MethodBase::fLoggermessage logger
private:
TMVA::MethodCuts::EConstrainTypefConstrainType
TStringfFitMethodSchosen fit method (string)
TMVA::MethodCuts::EFitMethodTypefFitMethodchosen fit method
TStringfEffMethodSchosen efficiency calculation method (string)
TMVA::MethodCuts::EEffMethodfEffMethodchosen efficiency calculation method
vector<EFitParameters>*fFitParamsvector for series of fit methods
Double_tfTestSignalEffused to test optimized signal efficiency
Double_tfEffSMinused to test optimized signal efficiency
Double_tfEffSMaxused to test optimized signal efficiency
TMVA::BinarySearchTree*fBinaryTreeS
TMVA::BinarySearchTree*fBinaryTreeB
Int_tfGA_nstepsGA settings: number of steps
Int_tfGA_cyclesGA settings: number of pre-calc steps
Int_tfGA_popSizeGA settings: population size
Int_tfGA_SC_stepsGA settings: SC_steps
Int_tfGA_SC_offstepsGA settings: SC_offsteps
Double_tfGA_SC_factorGA settings: SC_factor
Int_tfSA_MaxCallsmax number of FCN calls
Double_tfSA_TemperatureGradientstarting value for temperature gradient
Bool_tfSA_UseAdaptiveTemperaturecompute temperature steps on the fly
Double_tfSA_InitialTemperatureinitial temperature (depends on FCN)
Double_tfSA_MinTemperatureminimum temperature before SA quit
Double_tfSA_Epsrelative required FCN accuracy at minimum
Int_tfSA_NFunLoopsnumber of FCN loops
Int_tfSA_NEpstest parameter
Int_tfNRandCutsnumber of random cut samplings
Double_t**fCutMinminimum requirement
Double_t**fCutMaxmaximum requirement
Double_t*fTmpCutMintemporary minimum requirement
Double_t*fTmpCutMaxtemporary maximum requirement
TStringfAllVars
TStringfAllVarsI[10]
Int_tfNparnumber of parameters in fit (default: 2*Nvar)
Double_tfEffRefreference efficiency
vector<Int_t>*fRangeSignused to match cuts to fit parameters (and vice versa)
TRandom*fRandomrandom generator for MC optimisation method
vector<Double_t>*fMeanSmeans of variables (signal)
vector<Double_t>*fMeanBmeans of variables (background)
vector<Double_t>*fRmsSRMSs of variables (signal)
vector<Double_t>*fRmsBRMSs of variables (background)
vector<Double_t>*fXminminimum values of variables
vector<Double_t>*fXmaxmaximum values of variables
TH1*fEffBvsSLocalintermediate eff. background versus eff signal histo
vector<TH1*>*fVarHistSreference histograms (signal)
vector<TH1*>*fVarHistBreference histograms (background)
vector<TH1*>*fVarHistS_smoothsmoothed reference histograms (signal)
vector<TH1*>*fVarHistB_smoothsmoothed reference histograms (background)
vector<PDF*>*fVarPdfSreference PDFs (signal)
vector<PDF*>*fVarPdfBreference PDFs (background)
static TMVA::MethodCuts*fgThisCutsused for function reference (GA)

Class Description


/*
  Multivariate optimisation of signal efficiency for given background  
  efficiency, applying rectangular minimum and maximum requirements.

  

Also implemented is a "decorrelate/diagonlized cuts approach", which improves over the uncorrelated cuts ansatz by transforming linearly the input variables into a diagonal space, using the square-root of the covariance matrix.

Other optimisation criteria, such as maximising the signal significance- squared, S^2/(S+B), with S and B being the signal and background yields, correspond to a particular point in the optimised background rejection versus signal efficiency curve. This working point requires the knowledge of the expected yields, which is not the case in general. Note also that for rare signals, Poissonian statistics should be used, which modifies the significance criterion.

The rectangular cut of a volume in the variable space is performed using a binary tree to sort the training events. This provides a significant reduction in computing time (up to several orders of magnitudes, depending on the complexity of the problem at hand).

Technically, optimisation is achieved in TMVA by two methods:

  1. Monte Carlo generation using uniform priors for the lower cut value, and the cut width, thrown within the variable ranges.
  2. A Genetic Algorithm (GA) searches for the optimal ("fittest") cut sample. The GA is configurable by many external settings through the option string. For difficult cases (such as many variables), some tuning may be necessary to achieve satisfying results

Attempts to use Minuit fits (Simplex ot Migrad) instead have not shown superior results, and often failed due to convergence at local minima.

The tests we have performed so far showed that in generic applications, the GA is superior to MC sampling, and hence GA is the default method. It is worthwhile trying both anyway. Decorrelated (or "diagonalized") Cuts

See class description for Method Likelihood for a detailed explanation. */

MethodCuts( TString jobName, TString methodTitle, DataSet& theData, TString theOption, TDirectory* theTargetDir )
 standard constructor
 ---------------------------------------------------------------------------------- 
 format of option string: "OptMethod:EffMethod:Option_var1:...:Option_varn:Decorr"
 "OptMethod" can be:
     - "GA"    : Genetic Algorithm (recommended)
     - "SA"    : Simulated Annealing
     - "MC"    : Monte-Carlo optimization 
 "EffMethod" can be:
     - "EffSel": compute efficiency by event counting
     - "EffPDF": compute efficiency from PDFs
 === For "GA" method ======
 "Option_var1++" are (see GA for explanation of parameters):
     - fGA_nsteps        
     - fGA_cycles        
     - fGA_popSize
     - fGA_SC_steps        
     - fGA_SC_offsteps 
     - fGA_SC_factor   
 === For "SA" method ======
 "Option_var1++" are (see SA for explanation of parameters):
     - fSA_MaxCalls                
     - fSA_TemperatureGradient      
     - fSA_UseAdaptiveTemperature    
     - fSA_InitialTemperature        
     - fSA_MinTemperature        
     - fSA_Eps                       
     - fSA_NFunLoops                 
     - fSA_NEps                      
 === For "MC" method ======
 "Option_var1" is number of random samples
 "Option_var2++" can be 
     - "FMax"  : ForceMax   (the max cut is fixed to maximum of variable i)
     - "FMin"  : ForceMin   (the min cut is fixed to minimum of variable i)
     - "FSmart": ForceSmart (the min or max cut is fixed to min/max, based on mean value)
     - Adding "All" to "option_vari", eg, "AllFSmart" will use this option for all variables
     - if "option_vari" is empty (== ""), no assumptions on cut min/max are made
 "Decorr" can be:
     - omitted : Decorrelation not used
     - "D"     : Decorrelates variables, evaluation events decorrelated with signal decorrelation matrix
     - "DS"    : Decorrelates variables, evaluation events decorrelated with signal decorrelation matrix
     - "DB"    : Decorrelates variables, evaluation events decorrelated with background decorrelation matrix
 ---------------------------------------------------------------------------------- 
MethodCuts( DataSet& theData, TString theWeightFile, TDirectory* theTargetDir )
 construction from weight file
void InitCuts( void )
 default initialisation called by all constructors
~MethodCuts( void )
 destructor
void DeclareOptions()
 define the options (their key words) that can be set in the option string 
 know options:
 Method             <string> Minimization method
    available values are:        MC Monte Carlo <default>
                                 GA Genetic Algorithm
                                 SA Simulated annealing

 EffMethod          <string> Efficiency selection method
    available values are:        EffSel <default>
                                 EffPDF

 MC_NRandCuts       <int>    Number of random cuts to estimate the efficiency for the MC method
 MC_AllVarProp      <string> Property of all variables for the MC method
    available values are:        AllNotEnforced <default>
                                 AllFMax
                                 AllFMin
                                 AllFSmart
                                 AllFVerySmart
 MC_Var1Prop        <string> Property of variable 1 for the MC method (taking precedence over the
    globale setting. The same values as for the global option are available. Variables 1..10 can be
    set this way


 GA_nsteps           <int>   Number of steps for the genetic algorithm
 GA_cycles           <int>   Number of generations for the genetic algorithm
 GA_popSize          <int>   Size of the population for the genetic algorithm
 GA_SC_steps         <int>   Number of steps for the genetic algorithm
 GA_SC_offsteps      <int>    for the genetic algorithm
 GA_SC_factor        <float>  for the genetic algorithm


 SA_MaxCalls                 <int>      maximum number of calls for simulated annealing
 SA_TemperatureGradient      <float>    temperature gradient for simulated annealing
 SA_UseAdaptiveTemperature   <bool>     use of adaptive temperature for simulated annealing
 SA_InitialTemperature       <float>    initial temperature for simulated annealing
 SA_MinTemperature           <float>    minimum temperature for simulated annealing 
 SA_Eps                      <int>      number of epochs for simulated annealing
 SA_NFunLoops                <int>      number of loops for simulated annealing      
 SA_NEps                     <int>      number of epochs for simulated annealing
void ProcessOptions()
 process user options
Double_t GetMvaValue()
 cut evaluation: returns 1.0 if event passed, 0.0 otherwise
void Train( void )
 training method: here the cuts are optimised for the training sample
void Test( TTree* )
 not used 
Double_t ComputeEstimator( const std::vector<Double_t>& par )
 returns estimator for "cut fitness" used by GA
 there are two requirements:
 1) the signal efficiency must be equal to the required one in the 
    efficiency scan
 2) the background efficiency must be as small as possible
 the requirement 1) has priority over 2)
void MatchParsToCuts( const std::vector<Double_t> & par, Double_t* cutMin, Double_t* cutMax )
 translates parameters into cuts
void MatchCutsToPars( Double_t* par, Double_t* cutMin, Double_t* cutMax )
 translates cuts into parameters
void GetEffsfromPDFs( Double_t* cutMin, Double_t* cutMax, Double_t& effS, Double_t& effB )
 compute signal and background efficiencies from PDFs 
 for given cut sample
void GetEffsfromSelection( Double_t* cutMin, Double_t* cutMax, Double_t& effS, Double_t& effB)
 compute signal and background efficiencies from event counting 
 for given cut sample
void CreateVariablePDFs( void )
 for PDF method: create efficiency reference histograms and PDFs
Bool_t SanityChecks( void )
 basic checks to ensure that assumptions on variable order are satisfied
void WriteWeightsToStream( ostream & o )
 first the dimensions
void ReadWeightsFromStream( istream& istr )
 read the cuts from stream
void WriteMonitoringHistosToFile( void )
 write histograms and PDFs to file for monitoring purposes
Double_t GetTrainingEfficiency( TString theString)
 - overloaded function to create background efficiency (rejection) versus
   signal efficiency plot (first call of this function)
 - the function returns the signal efficiency at background efficiency
   indicated in theString

 "theString" must have two entries:
 [0]: "Efficiency"
 [1]: the value of background efficiency at which the signal efficiency 
      is to be returned
Double_t GetEfficiency( TString theString, TTree* /*theTree*/ )
 - overloaded function to create background efficiency (rejection) versus
   signal efficiency plot (first call of this function)
 - the function returns the signal efficiency at background efficiency
   indicated in theString

 "theString" must have two entries:
 [0]: "Efficiency"
 [1]: the value of background efficiency at which the signal efficiency 
      is to be returned
Double_t GetSignificance( void )
 also overwrite:
{ return 0; }
void SetTestSignalEfficiency( Double_t eff )
{ fTestSignalEff = eff; }
MethodCuts* ThisCuts( void )
 static pointer to this object
{ return fgThisCuts; }
const Ranking* CreateRanking()
 ranking of input variables
{ return 0; }
void MatchParsToCuts( const std::vector<Double_t> &, Double_t*, Double_t* )
 the definition of fit parameters can be different from the actual 
 cut requirements; these functions provide the matching

Author: Andreas Hoecker, Matt Jachowski, Peter Speckmayer, Helge Voss, Kai Voss
Last update: root/tmva $Id: MethodCuts.cxx,v 1.10 2006/11/20 15:35:28 brun Exp $
Copyright (c) 2005: *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.