48 Bool_t UseOffsetMethod = kTRUE;
50 void TMVAClassificationCategory()
55 std::cout << std::endl <<
"==> Start TMVAClassificationCategory" << std::endl;
58 TMVA::Tools::Instance();
60 bool batchMode =
false;
63 TString outfileName(
"TMVA.root" );
64 TFile* outputFile = TFile::Open( outfileName,
"RECREATE" );
68 std::string factoryOptions(
"!V:!Silent:Transformations=I;D;P;G,D" );
69 if (batchMode) factoryOptions +=
":!Color:!DrawProgressBar";
71 TMVA::Factory *factory =
new TMVA::Factory(
"TMVAClassificationCategory", outputFile, factoryOptions );
74 TMVA::DataLoader *dataloader=
new TMVA::DataLoader(
"dataset");
77 dataloader->AddVariable(
"var1",
'F' );
78 dataloader->AddVariable(
"var2",
'F' );
79 dataloader->AddVariable(
"var3",
'F' );
80 dataloader->AddVariable(
"var4",
'F' );
85 dataloader->AddSpectator(
"eta" );
89 TString fname = TString(gSystem->DirName(__FILE__) ) +
"/data/";
90 if (gSystem->AccessPathName( fname +
"toy_sigbkg_categ_offset.root")) {
92 fname = gROOT->GetTutorialDir() +
"/tmva/data/";
94 if (UseOffsetMethod) fname +=
"toy_sigbkg_categ_offset.root";
95 else fname +=
"toy_sigbkg_categ_varoff.root";
96 if (!gSystem->AccessPathName( fname )) {
98 std::cout <<
"--- TMVAClassificationCategory: Accessing " << fname << std::endl;
99 input = TFile::Open( fname );
103 std::cout <<
"ERROR: could not open data file: " << fname << std::endl;
107 TTree *signalTree = (TTree*)input->Get(
"TreeS");
108 TTree *background = (TTree*)input->Get(
"TreeB");
111 Double_t signalWeight = 1.0;
112 Double_t backgroundWeight = 1.0;
115 dataloader->AddSignalTree ( signalTree, signalWeight );
116 dataloader->AddBackgroundTree( background, backgroundWeight );
123 dataloader->PrepareTrainingAndTestTree( mycuts, mycutb,
124 "nTrain_Signal=0:nTrain_Background=0:SplitMode=Random:NormMode=NumEvents:!V" );
129 factory->BookMethod( dataloader, TMVA::Types::kFisher,
"Fisher",
"!H:!V:Fisher" );
132 factory->BookMethod( dataloader, TMVA::Types::kLikelihood,
"Likelihood",
133 "!H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmoothBkg[1]=10:NSmooth=1:NAvEvtPerBin=50" );
136 TMVA::MethodCategory* mcat = 0;
139 TString theCat1Vars =
"var1:var2:var3:var4";
140 TString theCat2Vars = (UseOffsetMethod ?
"var1:var2:var3:var4" :
"var1:var2:var3");
143 TMVA::MethodBase* fiCat = factory->BookMethod( dataloader, TMVA::Types::kCategory,
"FisherCat",
"" );
144 mcat =
dynamic_cast<TMVA::MethodCategory*
>(fiCat);
145 mcat->AddMethod(
"abs(eta)<=1.3", theCat1Vars, TMVA::Types::kFisher,
"Category_Fisher_1",
"!H:!V:Fisher" );
146 mcat->AddMethod(
"abs(eta)>1.3", theCat2Vars, TMVA::Types::kFisher,
"Category_Fisher_2",
"!H:!V:Fisher" );
149 TMVA::MethodBase* liCat = factory->BookMethod( dataloader, TMVA::Types::kCategory,
"LikelihoodCat",
"" );
150 mcat =
dynamic_cast<TMVA::MethodCategory*
>(liCat);
151 mcat->AddMethod(
"abs(eta)<=1.3",theCat1Vars, TMVA::Types::kLikelihood,
152 "Category_Likelihood_1",
"!H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmoothBkg[1]=10:NSmooth=1:NAvEvtPerBin=50" );
153 mcat->AddMethod(
"abs(eta)>1.3", theCat2Vars, TMVA::Types::kLikelihood,
154 "Category_Likelihood_2",
"!H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmoothBkg[1]=10:NSmooth=1:NAvEvtPerBin=50" );
159 factory->TrainAllMethods();
162 factory->TestAllMethods();
165 factory->EvaluateAllMethods();
172 std::cout <<
"==> Wrote root file: " << outputFile->GetName() << std::endl;
173 std::cout <<
"==> TMVAClassificationCategory is done!" << std::endl;
180 if (!gROOT->IsBatch()) TMVA::TMVAGui( outfileName );
182 int main(
int argc,
char** argv )
184 TMVAClassificationCategory();