#ifndef ROOT_TMVA_MethodBase
#define ROOT_TMVA_MethodBase
#include <iosfwd>
#include <vector>
#ifndef ROOT_TMVA_IMethod
#include "TMVA/IMethod.h"
#endif
#ifndef ROOT_TMVA_Tools
#include "TMVA/Tools.h"
#endif
#ifndef ROOT_TMVA_Configurable
#include "TMVA/Configurable.h"
#endif
#ifndef ROOT_TMVA_Types
#include "TMVA/Types.h"
#endif
#ifndef ROOT_TMVA_VariableTransformBase
#include "TMVA/VariableTransformBase.h"
#endif
#ifndef ROOT_TMVA_DataSet
#include "TMVA/DataSet.h"
#endif
class TGraph;
class TTree;
class TDirectory;
class TSpline;
namespace TMVA {
class Ranking;
class MsgLogger;
class PDF;
class TSpline1;
class MethodCuts;
class MethodBase : public IMethod, public Configurable {
public:
enum EWeightFileType { kROOT=0, kTEXT };
MethodBase( const TString& jobName,
const TString& methodTitle,
DataSet& theData,
const TString& theOption = "",
TDirectory* theBaseDir = 0 );
MethodBase( DataSet& theData,
const TString& weightFile,
TDirectory* theBaseDir = 0 );
virtual ~MethodBase();
virtual void AddClassifierToTestTree( TTree* theTestTree );
void TrainMethod();
virtual void Train() = 0;
virtual void Test( TTree* theTestTree = 0 );
virtual Double_t GetMvaValue() = 0;
virtual Double_t GetProba( Double_t mvaVal, Double_t ap_sig );
virtual Double_t GetRarity( Double_t mvaVal, Types::ESBType reftype = Types::kBackground ) const;
virtual const Ranking* CreateRanking() = 0;
virtual void MakeClass( const TString& classFileName = "" ) const;
void PrintHelpMessage() const;
void WriteStateToFile () const;
void WriteStateToStream ( std::ostream& tf, Bool_t isClass = kFALSE ) const;
void WriteStateToStream ( TFile& rf ) const;
void ReadStateFromFile ();
void ReadStateFromStream( std::istream& tf );
void ReadStateFromStream( TFile& rf );
virtual void WriteWeightsToStream ( std::ostream& tf ) const = 0;
virtual void WriteWeightsToStream ( TFile& ) const {}
virtual void ReadWeightsFromStream( std::istream& tf ) = 0;
virtual void ReadWeightsFromStream( TFile& ) {}
virtual void WriteEvaluationHistosToFile();
virtual void WriteMonitoringHistosToFile() const;
virtual Double_t GetEfficiency( TString, TTree*, Double_t& err );
virtual Double_t GetTrainingEfficiency( TString );
virtual Double_t GetSignificance() const;
virtual Double_t GetMaximumSignificance( Double_t SignalEvents, Double_t BackgroundEvents,
Double_t& optimal_significance_value ) const;
virtual Double_t GetSeparation( TH1*, TH1* ) const;
virtual Double_t GetSeparation( PDF* pdfS = 0, PDF* pdfB = 0 ) const;
const TString& GetJobName () const { return fJobName; }
const TString& GetMethodName () const { return fMethodName; }
const TString& GetMethodTitle() const { return fMethodTitle; }
Types::EMVA GetMethodType () const { return fMethodType; }
const char* GetName () const { return GetMethodName().Data(); }
const TString& GetTestvarName() const { return fTestvar; }
const TString GetProbaName () const { return fTestvar + "_Proba"; }
void SetMethodName ( TString methodName ) { fMethodName = methodName; }
void SetMethodTitle( TString methodTitle ) { fMethodTitle = methodTitle; }
void SetMethodType ( Types::EMVA methodType ) { fMethodType = methodType; }
void SetTestvarPrefix( TString prefix ) { fTestvarPrefix = prefix; }
void SetTestvarName( const TString & v="" ) { fTestvar = (v=="")?(fTestvarPrefix + GetMethodTitle()):v; }
Int_t GetNvar() const { return fNvar; }
const TString& GetInputVar( int i ) const { return Data().GetInternalVarName(i); }
const TString& GetInputExp( int i ) const { return Data().GetExpression(i); }
Double_t GetRMS( Int_t ivar ) const { return GetVarTransform().Variable(ivar).GetRMS(); }
Double_t GetXmin( Int_t ivar ) const { return GetVarTransform().Variable(ivar).GetMin(); }
Double_t GetXmax( Int_t ivar ) const { return GetVarTransform().Variable(ivar).GetMax(); }
Double_t GetSignalReferenceCut() const { return fSignalReferenceCut; }
VariableTransformBase& GetVarTransform() const { return *fVarTransform; }
TDirectory* BaseDir() const;
TDirectory* MethodBaseDir() const;
UInt_t GetTrainingTMVAVersionCode() const { return fTMVATrainingVersion; }
UInt_t GetTrainingROOTVersionCode() const { return fROOTTrainingVersion; }
TString GetTrainingTMVAVersionString() const;
TString GetTrainingROOTVersionString() const;
DataSet& Data() const { return fData; }
Event& GetEvent() const { return GetVarTransform().GetEvent(); }
Bool_t ReadEvent( TTree* tr, UInt_t ievt, Types::ESBType type = Types::kMaxSBType ) const;
Bool_t ReadTrainingEvent( UInt_t ievt, Types::ESBType type = Types::kMaxSBType ) const;
Bool_t ReadTestEvent ( UInt_t ievt, Types::ESBType type = Types::kMaxSBType ) const;
Bool_t IsSignalEvent() const { return GetEvent().IsSignal(); }
Double_t GetEventVal ( Int_t ivar ) const;
Double_t GetEventValNormalised(Int_t ivar) const;
Double_t GetEventWeight() const { return GetEvent().GetWeight(); }
virtual Bool_t IsSignalLike() { return GetMvaValue() > GetSignalReferenceCut() ? kTRUE : kFALSE; }
protected:
TDirectory* LocalTDir() const { return Data().LocalRootDir(); }
void SetWeightFileName( TString );
TString GetWeightFileName() const;
TString GetWeightFileDir() const { return fFileDir; }
void SetWeightFileDir( TString fileDir );
Bool_t IsNormalised() const { return fNormalise; }
void SetNormalised( Bool_t norm ) { fNormalise = norm; }
void SetNvar( Int_t n ) { fNvar = n; }
Types::EVariableTransform GetVariableTransform() const { return fVariableTransform; }
void SetSignalReferenceCut( Double_t cut ) { fSignalReferenceCut = cut; }
Bool_t Verbose() const { return fVerbose; }
Bool_t Help () const { return fHelp; }
const TString& GetInternalVarName( Int_t ivar ) const { return (*fInputVars)[ivar]; }
const TString& GetOriginalVarName( Int_t ivar ) const { return Data().GetExpression(ivar); }
Bool_t HasTrainingTree() const { return Data().GetTrainingTree() != 0; }
TTree* GetTrainingTree() const {
if (GetVariableTransform() != Types::kNone) {
fLogger << kFATAL << "Trying to access correlated Training tree in method "
<< GetMethodName() << Endl;
}
return Data().GetTrainingTree();
}
TTree* GetTestTree() const {
if (GetVariableTransform() != Types::kNone) {
fLogger << kFATAL << "Trying to access correlated Training tree in method "
<< GetMethodName() << Endl;
}
return Data().GetTestTree();
}
virtual void DeclareOptions();
virtual void ProcessOptions();
virtual void MakeClassSpecific( std::ostream&, const TString& = "" ) const {}
virtual void MakeClassSpecificHeader( std::ostream&, const TString& = "" ) const {}
static MethodBase* GetThisBase() { return fgThisBase; }
void Statistics( Types::ETreeType treeType, const TString& theVarName,
Double_t&, Double_t&, Double_t&,
Double_t&, Double_t&, Double_t&, Bool_t norm = kFALSE );
Bool_t CheckSanity( TTree* theTree = 0 );
Bool_t TxtWeightsOnly() const { return fTxtWeightsOnly; }
private:
enum ECutOrientation { kNegative = -1, kPositive = +1 };
ECutOrientation GetCutOrientation() const { return fCutOrientation; }
void ResetThisBase() { fgThisBase = this; }
void Init();
void CreateMVAPdfs();
Bool_t HasMVAPdfs() const { return fHasMVAPdfs; }
static Double_t IGetEffForRoot( Double_t );
Double_t GetEffForRoot ( Double_t );
Bool_t GetLine(std::istream& fin, char * buf );
protected:
Ranking* fRanking;
std::vector<TString>* fInputVars;
Int_t fNbins;
Int_t fNbinsH;
private:
friend class MethodCuts;
DataSet& fData;
Double_t fSignalReferenceCut;
Types::ESBType fVariableTransformType;
TString fJobName;
TString fMethodName;
Types::EMVA fMethodType;
TString fMethodTitle;
TString fTestvar;
TString fTestvarPrefix;
UInt_t fTMVATrainingVersion;
UInt_t fROOTTrainingVersion;
Bool_t fNormalise;
Int_t fNvar;
TDirectory* fBaseDir;
TDirectory* fMethodBaseDir;
TString fFileDir;
TString fWeightFile;
private:
TH1* fHistS_plotbin;
TH1* fHistB_plotbin;
TH1* fHistTrS_plotbin;
TH1* fHistTrB_plotbin;
TH1* fProbaS_plotbin;
TH1* fProbaB_plotbin;
TH1* fRarityS_plotbin;
TH1* fRarityB_plotbin;
TH1* fHistS_highbin;
TH1* fHistB_highbin;
TH1* fEffS;
TH1* fEffB;
TH1* fEffBvsS;
TH1* fRejBvsS;
TH1* finvBeffvsSeff;
TH1* fTrainEffS;
TH1* fTrainEffB;
TH1* fTrainEffBvsS;
TH1* fTrainRejBvsS;
Int_t fNbinsMVAPdf;
Int_t fNsmoothMVAPdf;
PDF* fMVAPdfS;
PDF* fMVAPdfB;
TGraph* fGraphS;
TGraph* fGraphB;
TGraph* fGrapheffBvsS;
PDF* fSplS;
PDF* fSplB;
TSpline* fSpleffBvsS;
TGraph* fGraphTrainS;
TGraph* fGraphTrainB;
TGraph* fGraphTrainEffBvsS;
PDF* fSplTrainS;
PDF* fSplTrainB;
TSpline* fSplTrainEffBvsS;
private:
Double_t fMeanS;
Double_t fMeanB;
Double_t fRmsS;
Double_t fRmsB;
Double_t fXmin;
Double_t fXmax;
Bool_t fUseDecorr;
Types::EVariableTransform fVariableTransform;
VariableTransformBase* fVarTransform;
TString fVarTransformString;
TString fVariableTransformTypeString;
Bool_t fVerbose;
TString fVerbosityLevelString;
EMsgType fVerbosityLevel;
Bool_t fHelp;
Bool_t fHasMVAPdfs;
Bool_t fTxtWeightsOnly;
private:
ECutOrientation fCutOrientation;
TSpline1* fSplRefS;
TSpline1* fSplRefB;
TSpline1* fSplTrainRefS;
TSpline1* fSplTrainRefB;
private:
static MethodBase* fgThisBase;
protected:
mutable MsgLogger fLogger;
ClassDef(MethodBase,0)
};
}
inline Bool_t TMVA::MethodBase::ReadEvent( TTree* tr, UInt_t ievt, Types::ESBType type ) const
{
if (type == Types::kMaxSBType) type = fVariableTransformType;
fVarTransform->ReadEvent(tr, ievt, type);
return kTRUE;
}
inline Bool_t TMVA::MethodBase::ReadTrainingEvent( UInt_t ievt, Types::ESBType type ) const
{
return ReadEvent( Data().GetTrainingTree(), ievt, type );
}
inline Bool_t TMVA::MethodBase::ReadTestEvent( UInt_t ievt, Types::ESBType type ) const
{
return ReadEvent( Data().GetTestTree(), ievt, type );
}
inline Double_t TMVA::MethodBase::GetEventVal( Int_t ivar ) const
{
if (IsNormalised()) return GetEventValNormalised(ivar);
else return GetEvent().GetVal(ivar);
}
inline Double_t TMVA::MethodBase::GetEventValNormalised(Int_t ivar) const
{
return gTools().NormVariable( GetEvent().GetVal(ivar), GetXmin(ivar), GetXmax(ivar) );
}
inline TString TMVA::MethodBase::GetTrainingTMVAVersionString() const
{
UInt_t a = GetTrainingTMVAVersionCode() & 0xff0000; a>>=16;
UInt_t b = GetTrainingTMVAVersionCode() & 0x00ff00; b>>=8;
UInt_t c = GetTrainingTMVAVersionCode() & 0x0000ff;
return TString(Form("%i.%i.%i",a,b,c));
}
inline TString TMVA::MethodBase::GetTrainingROOTVersionString() const
{
UInt_t a = GetTrainingROOTVersionCode() & 0xff0000; a>>=16;
UInt_t b = GetTrainingROOTVersionCode() & 0x00ff00; b>>=8;
UInt_t c = GetTrainingROOTVersionCode() & 0x0000ff;
return TString(Form("%i.%02i/%02i",a,b,c));
}
#endif
Last change: Sat Nov 1 10:21:40 2008
Last generated: 2008-11-01 10:21
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.