ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MethodCategory.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Joerg Stelzer, Helge Voss, Eckhard von Toerne
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : MethodCompositeBase *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Virtual base class for all MVA method *
12  * *
13  * Authors (alphabetical): *
14  * Andreas Hoecker <Andreas.Hocker@cern.ch> - CERN, Switzerland *
15  * Nadim Sah <Nadim.Sah@cern.ch> - Berlin, Germany *
16  * Peter Speckmayer <Peter.Speckmazer@cern.ch> - CERN, Switzerland *
17  * Joerg Stelzer <Joerg.Stelzer@cern.ch> - MSU East Lansing, USA *
18  * Helge Voss <Helge.Voss@cern.ch> - MPI-K Heidelberg, Germany *
19  * Jan Therhaag <Jan.Therhaag@cern.ch> - U of Bonn, Germany *
20  * Eckhard v. Toerne <evt@uni-bonn.de> - U of Bonn, Germany *
21  * *
22  * Copyright (c) 2005-2011: *
23  * CERN, Switzerland *
24  * MSU East Lansing, USA *
25  * MPI-K Heidelberg, Germany *
26  * U. of Bonn, Germany *
27  * *
28  * Redistribution and use in source and binary forms, with or without *
29  * modification, are permitted according to the terms listed in LICENSE *
30  * (http://tmva.sourceforge.net/LICENSE) *
31  **********************************************************************************/
32 
33 //__________________________________________________________________________
34 //
35 // This class is meant to allow categorisation of the data. For different //
36 // categories, different classifiers may be booked and different variab- //
37 // les may be considered. The aim is to account for the difference that //
38 // is due to different locations/angles. //
39 ////////////////////////////////////////////////////////////////////////////////
40 
41 #include "TMVA/MethodCategory.h"
42 
43 #include <algorithm>
44 #include <iomanip>
45 #include <vector>
46 #include <iostream>
47 
48 #include "Riostream.h"
49 #include "TRandom3.h"
50 #include "TMath.h"
51 #include "TObjString.h"
52 #include "TH1F.h"
53 #include "TGraph.h"
54 #include "TSpline.h"
55 #include "TDirectory.h"
56 #include "TTreeFormula.h"
57 
58 #include "TMVA/ClassifierFactory.h"
59 #include "TMVA/Config.h"
60 #include "TMVA/DataSet.h"
61 #include "TMVA/DataSetInfo.h"
62 #include "TMVA/DataSetManager.h"
63 #include "TMVA/IMethod.h"
64 #include "TMVA/MethodBase.h"
66 #include "TMVA/MsgLogger.h"
67 #include "TMVA/PDF.h"
68 #include "TMVA/Ranking.h"
69 #include "TMVA/Timer.h"
70 #include "TMVA/Tools.h"
71 #include "TMVA/Types.h"
72 #include "TMVA/VariableInfo.h"
74 
75 REGISTER_METHOD(Category)
76 
77 ClassImp(TMVA::MethodCategory)
78 
79 ////////////////////////////////////////////////////////////////////////////////
80 /// standard constructor
81 
82  TMVA::MethodCategory::MethodCategory( const TString& jobName,
83  const TString& methodTitle,
84  DataSetInfo& theData,
85  const TString& theOption,
86  TDirectory* theTargetDir )
87  : TMVA::MethodCompositeBase( jobName, Types::kCategory, methodTitle, theData, theOption, theTargetDir ),
88  fCatTree(0),
89  fDataSetManager(NULL)
90 {
91 }
92 
93 ////////////////////////////////////////////////////////////////////////////////
94 /// constructor from weight file
95 
97  const TString& theWeightFile,
98  TDirectory* theTargetDir )
99  : TMVA::MethodCompositeBase( Types::kCategory, dsi, theWeightFile, theTargetDir ),
100  fCatTree(0),
101  fDataSetManager(NULL)
102 {
103 }
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 /// destructor
107 
109 {
110  std::vector<TTreeFormula*>::iterator formIt = fCatFormulas.begin();
111  std::vector<TTreeFormula*>::iterator lastF = fCatFormulas.end();
112  for(;formIt!=lastF; ++formIt) delete *formIt;
113  delete fCatTree;
114 }
115 
116 ////////////////////////////////////////////////////////////////////////////////
117 /// check whether method category has analysis type
118 /// the method type has to be the same for all sub-methods
119 
121 {
122  std::vector<IMethod*>::iterator itrMethod = fMethods.begin();
123 
124  // iterate over methods and check whether they have the analysis type
125  for(; itrMethod != fMethods.end(); ++itrMethod ) {
126  if ( !(*itrMethod)->HasAnalysisType(type, numberClasses, numberTargets) )
127  return kFALSE;
128  }
129  return kTRUE;
130 }
131 
132 ////////////////////////////////////////////////////////////////////////////////
133 /// options for this method
134 
136 {
137 }
138 
139 ////////////////////////////////////////////////////////////////////////////////
140 /// adds sub-classifier for a category
141 
143  const TString& theVariables,
144  Types::EMVA theMethod ,
145  const TString& theTitle,
146  const TString& theOptions )
147 {
148  std::string addedMethodName = std::string(Types::Instance().GetMethodName(theMethod));
149 
150  Log() << kINFO << "Adding sub-classifier: " << addedMethodName << "::" << theTitle << Endl;
151 
152  DataSetInfo& dsi = CreateCategoryDSI(theCut, theVariables, theTitle);
153 
154  IMethod* addedMethod = ClassifierFactory::Instance().Create(addedMethodName,GetJobName(),theTitle,dsi,theOptions);
155 
156  MethodBase *method = (dynamic_cast<MethodBase*>(addedMethod));
157  if(method==0) return 0;
158 
159  method->SetAnalysisType( fAnalysisType );
160  method->SetupMethod();
161  method->ParseOptions();
162  method->ProcessSetup();
163 
164  // set or create correct method base dir for added method
165  const TString dirName(Form("Method_%s",method->GetMethodTypeName().Data()));
166  TDirectory * dir = BaseDir()->GetDirectory(dirName);
167  if (dir != 0) method->SetMethodBaseDir( dir );
168  else method->SetMethodBaseDir( BaseDir()->mkdir(dirName,Form("Directory for all %s methods", method->GetMethodTypeName().Data())) );
169 
170  // method->SetBaseDir(eigenes base dir, gucken ob Fisher dir existiert, sonst erzeugen )
171 
172  // check-for-unused-options is performed; may be overridden by derived
173  // classes
174  method->CheckSetup();
175 
176  // disable writing of XML files and standalone classes for sub methods
177  method->DisableWriting( kTRUE );
178 
179  // store method, cut and variable names and create cut formula
180  fMethods.push_back(method);
181  fCategoryCuts.push_back(theCut);
182  fVars.push_back(theVariables);
183 
184  DataSetInfo& primaryDSI = DataInfo();
185 
186  UInt_t newSpectatorIndex = primaryDSI.GetSpectatorInfos().size();
187  fCategorySpecIdx.push_back(newSpectatorIndex);
188 
189  primaryDSI.AddSpectator( Form("%s_cat%i:=%s", GetName(),(int)fMethods.size(),theCut.GetTitle()),
190  Form("%s:%s",GetName(),method->GetName()),
191  "pass", 0, 0, 'C' );
192 
193  return method;
194 }
195 
196 ////////////////////////////////////////////////////////////////////////////////
197 /// create a DataSetInfo object for a sub-classifier
198 
200  const TString& theVariables,
201  const TString& theTitle)
202 {
203  // create a new dsi with name: theTitle+"_dsi"
204  TString dsiName=theTitle+"_dsi";
205  DataSetInfo& oldDSI = DataInfo();
206  DataSetInfo* dsi = new DataSetInfo(dsiName);
207 
208  // register the new dsi
209  // DataSetManager::Instance().AddDataSetInfo(*dsi); // DSMTEST replaced by following line
210  fDataSetManager->AddDataSetInfo(*dsi);
211 
212  // copy the targets and spectators from the old dsi to the new dsi
213  std::vector<VariableInfo>::iterator itrVarInfo;
214 
215  for (itrVarInfo = oldDSI.GetTargetInfos().begin(); itrVarInfo != oldDSI.GetTargetInfos().end(); itrVarInfo++)
216  dsi->AddTarget(*itrVarInfo);
217 
218  for (itrVarInfo = oldDSI.GetSpectatorInfos().begin(); itrVarInfo != oldDSI.GetSpectatorInfos().end(); itrVarInfo++)
219  dsi->AddSpectator(*itrVarInfo);
220 
221  // split string that contains the variables into tiny little pieces
222  std::vector<TString> variables = gTools().SplitString(theVariables,':' );
223 
224  // prepare to create varMap
225  std::vector<UInt_t> varMap;
226  Int_t counter=0;
227 
228  // add the variables that were specified in theVariables
229  std::vector<TString>::iterator itrVariables;
230  Bool_t found = kFALSE;
231 
232  // iterate over all variables in 'variables' and add them
233  for (itrVariables = variables.begin(); itrVariables != variables.end(); itrVariables++) {
234  counter=0;
235 
236  // check the variables of the old dsi for the variable that we want to add
237  for (itrVarInfo = oldDSI.GetVariableInfos().begin(); itrVarInfo != oldDSI.GetVariableInfos().end(); itrVarInfo++) {
238  if((*itrVariables==itrVarInfo->GetLabel()) ) { // || (*itrVariables==itrVarInfo->GetExpression())) {
239  // don't compare the expression, since the user might take two times the same expression, but with different labels
240  // and apply different transformations to the variables.
241  dsi->AddVariable(*itrVarInfo);
242  varMap.push_back(counter);
243  found = kTRUE;
244  }
245  counter++;
246  }
247 
248  // check the spectators of the old dsi for the variable that we want to add
249  for (itrVarInfo = oldDSI.GetSpectatorInfos().begin(); itrVarInfo != oldDSI.GetSpectatorInfos().end(); itrVarInfo++) {
250  if((*itrVariables==itrVarInfo->GetLabel()) ) { // || (*itrVariables==itrVarInfo->GetExpression())) {
251  // don't compare the expression, since the user might take two times the same expression, but with different labels
252  // and apply different transformations to the variables.
253  dsi->AddVariable(*itrVarInfo);
254  varMap.push_back(counter);
255  found = kTRUE;
256  }
257  counter++;
258  }
259 
260  // if the variable is neither in the variables nor in the spectators, we abort
261  if (!found) {
262  Log() << kFATAL <<"The variable " << itrVariables->Data() << " was not found and could not be added " << Endl;
263  }
264  found = kFALSE;
265  }
266 
267  // in the case that no variables are specified, add the default-variables from the original dsi
268  if (theVariables=="") {
269  for (UInt_t i=0; i<oldDSI.GetVariableInfos().size(); i++) {
270  dsi->AddVariable(oldDSI.GetVariableInfos()[i]);
271  varMap.push_back(i);
272  }
273  }
274 
275  // add the variable map 'varMap' to the vector of varMaps
276  fVarMaps.push_back(varMap);
277 
278  // set classes and cuts
279  UInt_t nClasses=oldDSI.GetNClasses();
280  TString className;
281 
282  for (UInt_t i=0; i<nClasses; i++) {
283  className = oldDSI.GetClassInfo(i)->GetName();
284  dsi->AddClass(className);
285  dsi->SetCut(oldDSI.GetCut(i),className);
286  dsi->AddCut(theCut,className);
287  dsi->SetWeightExpression(oldDSI.GetWeightExpression(i),className);
288  }
289 
290  // set split options, root dir and normalization for the new dsi
291  dsi->SetSplitOptions(oldDSI.GetSplitOptions());
292  dsi->SetRootDir(oldDSI.GetRootDir());
293  TString norm(oldDSI.GetNormalization().Data());
294  dsi->SetNormalization(norm);
295 
296  DataSetInfo& dsiReference= (*dsi);
297 
298  return dsiReference;
299 }
300 
301 ////////////////////////////////////////////////////////////////////////////////
302 /// initialize the method
303 
305 {
306 }
307 
308 ////////////////////////////////////////////////////////////////////////////////
309 /// initialize the circular tree
310 
312 {
313  delete fCatTree;
314 
315  std::vector<VariableInfo>::const_iterator viIt;
316  const std::vector<VariableInfo>& vars = dsi.GetVariableInfos();
317  const std::vector<VariableInfo>& specs = dsi.GetSpectatorInfos();
318 
319  Bool_t hasAllExternalLinks = kTRUE;
320  for (viIt = vars.begin(); viIt != vars.end(); ++viIt)
321  if( viIt->GetExternalLink() == 0 ) {
322  hasAllExternalLinks = kFALSE;
323  break;
324  }
325  for (viIt = specs.begin(); viIt != specs.end(); ++viIt)
326  if( viIt->GetExternalLink() == 0 ) {
327  hasAllExternalLinks = kFALSE;
328  break;
329  }
330 
331  if(!hasAllExternalLinks) return;
332 
333  {
334  // Rather than having TTree::TTree add to the current directory and then remove it, let
335  // make sure to not add it in the first place.
336  // The add-then-remove can lead to a problem if gDirectory points to the same directory (for example
337  // gROOT) in the current thread and another one (and both try to add to the directory at the same time).
338  TDirectory::TContext ctxt(nullptr);
339  fCatTree = new TTree(Form("Circ%s",GetMethodName().Data()),"Circlar Tree for categorization");
340  fCatTree->SetCircular(1);
341  }
342 
343  for (viIt = vars.begin(); viIt != vars.end(); ++viIt) {
344  const VariableInfo& vi = *viIt;
345  fCatTree->Branch(vi.GetExpression(),(Float_t*)vi.GetExternalLink(), TString(vi.GetExpression())+TString("/F"));
346  }
347  for (viIt = specs.begin(); viIt != specs.end(); ++viIt) {
348  const VariableInfo& vi = *viIt;
349  if(vi.GetVarType()=='C') continue;
350  fCatTree->Branch(vi.GetExpression(),(Float_t*)vi.GetExternalLink(), TString(vi.GetExpression())+TString("/F"));
351  }
352 
353  for(UInt_t cat=0; cat!=fCategoryCuts.size(); ++cat) {
354  fCatFormulas.push_back(new TTreeFormula(Form("Category_%i",cat), fCategoryCuts[cat].GetTitle(), fCatTree));
355  }
356 }
357 
358 ////////////////////////////////////////////////////////////////////////////////
359 /// train all sub-classifiers
360 
362 {
363  // specify the minimum # of training events and set 'classification'
364  const Int_t MinNoTrainingEvents = 10;
365 
366  Types::EAnalysisType analysisType = GetAnalysisType();
367 
368  // start the training
369  Log() << kINFO << "Train all sub-classifiers for "
370  << (analysisType == Types::kRegression ? "Regression" : "Classification") << " ..." << Endl;
371 
372  // don't do anything if no sub-classifier booked
373  if (fMethods.empty()) {
374  Log() << kINFO << "...nothing found to train" << Endl;
375  return;
376  }
377 
378  std::vector<IMethod*>::iterator itrMethod;
379 
380  // iterate over all booked sub-classifiers and train them
381  for (itrMethod = fMethods.begin(); itrMethod != fMethods.end(); ++itrMethod ) {
382 
383  MethodBase* mva = dynamic_cast<MethodBase*>(*itrMethod);
384  if(!mva) continue;
385  mva->SetAnalysisType( analysisType );
386  if (!mva->HasAnalysisType( analysisType,
387  mva->DataInfo().GetNClasses(),
388  mva->DataInfo().GetNTargets() ) ) {
389  Log() << kWARNING << "Method " << mva->GetMethodTypeName() << " is not capable of handling " ;
390  if (analysisType == Types::kRegression)
391  Log() << "regression with " << mva->DataInfo().GetNTargets() << " targets." << Endl;
392  else
393  Log() << "classification with " << mva->DataInfo().GetNClasses() << " classes." << Endl;
394  itrMethod = fMethods.erase( itrMethod );
395  continue;
396  }
397  if (mva->Data()->GetNTrainingEvents() >= MinNoTrainingEvents) {
398 
399  Log() << kINFO << "Train method: " << mva->GetMethodName() << " for "
400  << (analysisType == Types::kRegression ? "Regression" : "Classification") << Endl;
401  mva->TrainMethod();
402  Log() << kINFO << "Training finished" << Endl;
403 
404  } else {
405 
406  Log() << kWARNING << "Method " << mva->GetMethodName()
407  << " not trained (training tree has less entries ["
408  << mva->Data()->GetNTrainingEvents()
409  << "] than required [" << MinNoTrainingEvents << "]" << Endl;
410 
411  Log() << kERROR << " w/o training/test events for that category, I better stop here and let you fix " << Endl;
412  Log() << kFATAL << "that one first, otherwise things get too messy later ... " << Endl;
413 
414  }
415  }
416 
417  if (analysisType != Types::kRegression) {
418 
419  // variable ranking
420  Log() << kINFO << "Begin ranking of input variables..." << Endl;
421  for (itrMethod = fMethods.begin(); itrMethod != fMethods.end(); itrMethod++) {
422  MethodBase* mva = dynamic_cast<MethodBase*>(*itrMethod);
423  if (mva && mva->Data()->GetNTrainingEvents() >= MinNoTrainingEvents) {
424  const Ranking* ranking = (*itrMethod)->CreateRanking();
425  if (ranking != 0)
426  ranking->Print();
427  else
428  Log() << kINFO << "No variable ranking supplied by classifier: "
429  << dynamic_cast<MethodBase*>(*itrMethod)->GetMethodName() << Endl;
430  }
431  }
432  }
433 }
434 
435 ////////////////////////////////////////////////////////////////////////////////
436 /// create XML description of Category classifier
437 
438 void TMVA::MethodCategory::AddWeightsXMLTo( void* parent ) const
439 {
440  void* wght = gTools().AddChild(parent, "Weights");
441  gTools().AddAttr( wght, "NSubMethods", fMethods.size() );
442  void* submethod(0);
443 
444  std::vector<IMethod*>::iterator itrMethod;
445 
446  // iterate over methods and write them to XML file
447  for (UInt_t i=0; i<fMethods.size(); i++) {
448  MethodBase* method = dynamic_cast<MethodBase*>(fMethods[i]);
449  submethod = gTools().AddChild(wght, "SubMethod");
450  gTools().AddAttr(submethod, "Index", i);
451  gTools().AddAttr(submethod, "Method", method->GetMethodTypeName() + "::" + method->GetMethodName());
452  gTools().AddAttr(submethod, "Cut", fCategoryCuts[i]);
453  gTools().AddAttr(submethod, "Variables", fVars[i]);
454  method->WriteStateToXML( submethod );
455  }
456 }
457 
458 ////////////////////////////////////////////////////////////////////////////////
459 /// read weights of sub-classifiers of MethodCategory from xml weight file
460 
462 {
463  UInt_t nSubMethods;
464  TString fullMethodName;
465  TString methodType;
466  TString methodTitle;
467  TString theCutString;
468  TString theVariables;
469  Int_t titleLength;
470  gTools().ReadAttr( wghtnode, "NSubMethods", nSubMethods );
471  void* subMethodNode = gTools().GetChild(wghtnode);
472 
473  Log() << kINFO << "Recreating sub-classifiers from XML-file " << Endl;
474 
475  // recreate all sub-methods from weight file
476  for (UInt_t i=0; i<nSubMethods; i++) {
477  gTools().ReadAttr( subMethodNode, "Method", fullMethodName );
478  gTools().ReadAttr( subMethodNode, "Cut", theCutString );
479  gTools().ReadAttr( subMethodNode, "Variables", theVariables );
480 
481  // determine sub-method type
482  methodType = fullMethodName(0,fullMethodName.Index("::"));
483  if (methodType.Contains(" ")) methodType = methodType(methodType.Last(' ')+1,methodType.Length());
484 
485  // determine sub-method title
486  titleLength = fullMethodName.Length()-fullMethodName.Index("::")-2;
487  methodTitle = fullMethodName(fullMethodName.Index("::")+2,titleLength);
488 
489  // reconstruct dsi for sub-method
490  DataSetInfo& dsi = CreateCategoryDSI(TCut(theCutString), theVariables, methodTitle);
491 
492  // recreate sub-method from weights and add to fMethods
493  MethodBase* method = dynamic_cast<MethodBase*>( ClassifierFactory::Instance().Create( methodType.Data(),
494  dsi, "none" ) );
495  if(method==0)
496  Log() << kFATAL << "Could not create sub-method " << method << " from XML." << Endl;
497 
498  method->SetupMethod();
499  method->ReadStateFromXML(subMethodNode);
500 
501  fMethods.push_back(method);
502  fCategoryCuts.push_back(TCut(theCutString));
503  fVars.push_back(theVariables);
504 
505  DataSetInfo& primaryDSI = DataInfo();
506 
507  UInt_t spectatorIdx = 10000;
508  UInt_t counter=0;
509 
510  // find the spectator index
511  std::vector<VariableInfo>& spectators=primaryDSI.GetSpectatorInfos();
512  std::vector<VariableInfo>::iterator itrVarInfo;
513  TString specName= Form("%s_cat%i", GetName(),(int)fCategorySpecIdx.size()+1);
514 
515  for (itrVarInfo = spectators.begin(); itrVarInfo != spectators.end(); ++itrVarInfo, ++counter) {
516  if((specName==itrVarInfo->GetLabel()) || (specName==itrVarInfo->GetExpression())) {
517  spectatorIdx=counter;
518  fCategorySpecIdx.push_back(spectatorIdx);
519  break;
520  }
521  }
522 
523  subMethodNode = gTools().GetNextChild(subMethodNode);
524  }
525 
526  InitCircularTree(DataInfo());
527 
528 }
529 
530 ////////////////////////////////////////////////////////////////////////////////
531 /// process user options
532 
534 {
535 }
536 
537 ////////////////////////////////////////////////////////////////////////////////
538 /// Get help message text
539 ///
540 /// typical length of text line:
541 /// "|--------------------------------------------------------------|"
542 
544 {
545  Log() << Endl;
546  Log() << gTools().Color("bold") << "--- Short description:" << gTools().Color("reset") << Endl;
547  Log() << Endl;
548  Log() << "This method allows to define different categories of events. The" <<Endl;
549  Log() << "categories are defined via cuts on the variables. For each" << Endl;
550  Log() << "category, a different classifier and set of variables can be" <<Endl;
551  Log() << "specified. The categories which are defined for this method must" << Endl;
552  Log() << "be disjoint." << Endl;
553 }
554 
555 ////////////////////////////////////////////////////////////////////////////////
556 /// no ranking
557 
559 {
560  return 0;
561 }
562 
563 ////////////////////////////////////////////////////////////////////////////////
564 
566 {
567  // if it's not a simple 'spectator' variable (0 or 1) that the categories are defined by
568  // (but rather some 'formula' (i.e. eta>0), then this formulas are stored in fCatTree and that
569  // one will be evaluated.. (the formulae return 'true' or 'false'
570  if (fCatTree) {
571  if (methodIdx>=fCatFormulas.size()) {
572  Log() << kFATAL << "Large method index " << methodIdx << ", number of category formulas = "
573  << fCatFormulas.size() << Endl;
574  }
575  TTreeFormula* f = fCatFormulas[methodIdx];
576  return f->EvalInstance(0) > 0.5;
577  }
578  // otherwise, it simply looks if "variable == true" ("greater 0.5 to be "sure" )
579  else {
580 
581  // checks whether an event lies within a cut
582  if (methodIdx>=fCategorySpecIdx.size()) {
583  Log() << kFATAL << "Unknown method index " << methodIdx << " maximum allowed index="
584  << fCategorySpecIdx.size() << Endl;
585  }
586  UInt_t spectatorIdx = fCategorySpecIdx[methodIdx];
587  Float_t specVal = ev->GetSpectator(spectatorIdx);
588  Bool_t pass = (specVal>0.5);
589  return pass;
590  }
591 }
592 
593 ////////////////////////////////////////////////////////////////////////////////
594 /// returns the mva value of the right sub-classifier
595 
597 {
598  if (fMethods.empty()) return 0;
599 
600  UInt_t methodToUse = 0;
601  const Event* ev = GetEvent();
602 
603  // determine which sub-classifier to use for this event
604  Int_t suitableCutsN = 0;
605 
606  for (UInt_t i=0; i<fMethods.size(); ++i) {
607  if (PassesCut(ev, i)) {
608  ++suitableCutsN;
609  methodToUse=i;
610  }
611  }
612 
613  if (suitableCutsN == 0) {
614  Log() << kWARNING << "Event does not lie within the cut of any sub-classifier." << Endl;
615  return 0;
616  }
617 
618  if (suitableCutsN > 1) {
619  Log() << kFATAL << "The defined categories are not disjoint." << Endl;
620  return 0;
621  }
622 
623  // get mva value from the suitable sub-classifier
624  ev->SetVariableArrangement(&fVarMaps[methodToUse]);
625  Double_t mvaValue = dynamic_cast<MethodBase*>(fMethods[methodToUse])->GetMvaValue(ev,err,errUpper);
626  ev->SetVariableArrangement(0);
627 
628  return mvaValue;
629 }
630 
631 
632 
633 ////////////////////////////////////////////////////////////////////////////////
634 /// returns the mva value of the right sub-classifier
635 
636 const std::vector<Float_t> &TMVA::MethodCategory::GetRegressionValues()
637 {
638  if (fMethods.empty()) return MethodBase::GetRegressionValues();
639 
640  UInt_t methodToUse = 0;
641  const Event* ev = GetEvent();
642 
643  // determine which sub-classifier to use for this event
644  Int_t suitableCutsN = 0;
645 
646  for (UInt_t i=0; i<fMethods.size(); ++i) {
647  if (PassesCut(ev, i)) {
648  ++suitableCutsN;
649  methodToUse=i;
650  }
651  }
652 
653  if (suitableCutsN == 0) {
654  Log() << kWARNING << "Event does not lie within the cut of any sub-classifier." << Endl;
656  }
657 
658  if (suitableCutsN > 1) {
659  Log() << kFATAL << "The defined categories are not disjoint." << Endl;
661  }
662  MethodBase* meth = dynamic_cast<MethodBase*>(fMethods[methodToUse]);
663  if (!meth){
664  Log() << kFATAL << "method not found in Category Regression method" << Endl;
666  }
667  // get mva value from the suitable sub-classifier
668  return meth->GetRegressionValues(ev);
669 }
670 
IMethod * Create(const std::string &name, const TString &job, const TString &title, DataSetInfo &dsi, const TString &option)
creates the method if needed based on the method name using the creator function the factory has stor...
static ClassifierFactory & Instance()
access to the ClassifierFactory singleton creates the instance if needed
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52
void Init()
initialize the method
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition: TString.cxx:851
TMVA::IMethod * AddMethod(const TCut &, const TString &theVariables, Types::EMVA theMethod, const TString &theTitle, const TString &theOptions)
adds sub-classifier for a category
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
VariableInfo & AddTarget(const TString &expression, const TString &title, const TString &unit, Double_t min, Double_t max, Bool_t normalized=kTRUE, void *external=0)
add a variable (can be a complex expression) to the set of variables used in the MV analysis ...
void AddWeightsXMLTo(void *parent) const
create XML description of Category classifier
const TString GetWeightExpression(Int_t i) const
Definition: DataSetInfo.h:161
void ReadStateFromXML(void *parent)
std::vector< VariableInfo > & GetSpectatorInfos()
Definition: DataSetInfo.h:122
Ssiz_t Length() const
Definition: TString.h:390
void variables(TString fin="TMVA.root", TString dirName="InputVariables_Id", TString title="TMVA Input Variables", Bool_t isRegression=kFALSE, Bool_t useTMVAStyle=kTRUE)
Definition: variables.cxx:10
void SetVariableArrangement(std::vector< UInt_t > *const m) const
set the variable arrangement
Definition: Event.cxx:187
float Float_t
Definition: RtypesCore.h:53
void SetCut(const TCut &cut, const TString &className)
set the cut for the classes
void InitCircularTree(const DataSetInfo &dsi)
initialize the circular tree
const TString & GetExpression() const
Definition: VariableInfo.h:62
const char * GetName() const
Returns name of object.
Definition: MethodBase.h:299
UInt_t GetNClasses() const
Definition: DataSetInfo.h:152
void SetMethodBaseDir(TDirectory *methodDir)
Definition: MethodBase.h:337
DataSet * Data() const
Definition: MethodBase.h:363
EAnalysisType
Definition: Types.h:124
UInt_t GetNTargets() const
Definition: DataSetInfo.h:129
const std::vector< Float_t > & GetRegressionValues(const TMVA::Event *const ev)
Definition: MethodBase.h:186
virtual const std::vector< Float_t > & GetRegressionValues()
Definition: MethodBase.h:193
Basic string class.
Definition: TString.h:137
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual ~MethodCategory(void)
destructor
ClassImp(TMVA::MethodCategory) TMVA
standard constructor
Float_t GetSpectator(UInt_t ivar) const
return spectator content
Definition: Event.cxx:256
void AddCut(const TCut &cut, const TString &className)
set the cut for the classes
virtual Bool_t HasAnalysisType(Types::EAnalysisType type, UInt_t numberClasses, UInt_t numberTargets)=0
void AddAttr(void *node, const char *, const T &value, Int_t precision=16)
Definition: Tools.h:308
void * AddChild(void *parent, const char *childname, const char *content=0, bool isRootNode=false)
add child node
Definition: Tools.cxx:1134
void GetHelpMessage() const
Get help message text.
TFile * f
const TString & GetMethodName() const
Definition: MethodBase.h:296
const char * Data() const
Definition: TString.h:349
static Types & Instance()
the the single instance of "Types" if existin already, or create it (Signleton)
Definition: Types.cxx:64
Tools & gTools()
Definition: Tools.cxx:79
void DeclareOptions()
options for this method
void * GetChild(void *parent, const char *childname=0)
get child node
Definition: Tools.cxx:1158
TDirectory * GetRootDir() const
Definition: DataSetInfo.h:187
const TString & GetNormalization() const
Definition: DataSetInfo.h:132
std::vector< std::vector< double > > Data
virtual void ParseOptions()
options parser
void SetupMethod()
setup of methods
Definition: MethodBase.cxx:302
std::vector< VariableInfo > & GetTargetInfos()
Definition: DataSetInfo.h:117
MethodCategory(const TString &jobName, const TString &methodTitle, DataSetInfo &theData, const TString &theOption="", TDirectory *theTargetDir=NULL)
Used to pass a selection expression to the Tree drawing routine.
Definition: TTreeFormula.h:64
A specialized string object used for TTree selections.
Definition: TCut.h:27
const Int_t MinNoTrainingEvents
Definition: Factory.cxx:86
void SetSplitOptions(const TString &so)
Definition: DataSetInfo.h:182
const Ranking * CreateRanking()
no ranking
char GetVarType() const
Definition: VariableInfo.h:67
std::string GetMethodName(TCppMethod_t)
Definition: Cppyy.cxx:706
TString dirName
Definition: demos.C:9
TMVA::DataSetInfo & CreateCategoryDSI(const TCut &, const TString &, const TString &)
create a DataSetInfo object for a sub-classifier
ClassInfo * GetClassInfo(Int_t clNum) const
void SetWeightExpression(const TString &exp, const TString &className="")
set the weight expressions for the classes if class name is specified, set only for this class if cla...
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
void Train(void)
train all sub-classifiers
const TString & GetName() const
Definition: ClassInfo.h:72
void ReadAttr(void *node, const char *, T &value)
Definition: Tools.h:295
const TCut & GetCut(Int_t i) const
Definition: DataSetInfo.h:165
const TString & GetSplitOptions() const
Definition: DataSetInfo.h:183
virtual void CheckSetup()
check may be overridden by derived class (sometimes, eg, fitters are used which can only be implement...
Definition: MethodBase.cxx:327
double Double_t
Definition: RtypesCore.h:55
Describe directory structure in memory.
Definition: TDirectory.h:44
int type
Definition: TGX11.cxx:120
void dir(char *path=0)
Definition: rootalias.C:30
void * GetNextChild(void *prevchild, const char *childname=0)
XML helpers.
Definition: Tools.cxx:1170
DataSetInfo & DataInfo() const
Definition: MethodBase.h:364
virtual Bool_t HasAnalysisType(Types::EAnalysisType type, UInt_t numberClasses, UInt_t)
check whether method category has analysis type the method type has to be the same for all sub-method...
T EvalInstance(Int_t i=0, const char *stringStack[]=0)
Evaluate this treeformula.
ClassInfo * AddClass(const TString &className)
void ProcessSetup()
process all options the "CheckForUnusedOptions" is done in an independent call, since it may be overr...
Definition: MethodBase.cxx:317
const TString & Color(const TString &)
human readable color strings
Definition: Tools.cxx:837
Double_t GetMvaValue(Double_t *err=0, Double_t *errUpper=0)
returns the mva value of the right sub-classifier
void * GetExternalLink() const
Definition: VariableInfo.h:85
#define REGISTER_METHOD(CLASS)
for example
VariableInfo & AddSpectator(const TString &expression, const TString &title, const TString &unit, Double_t min, Double_t max, char type= 'F', Bool_t normalized=kTRUE, void *external=0)
add a spectator (can be a complex expression) to the set of spectator variables used in the MV analys...
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
#define NULL
Definition: Rtypes.h:82
VariableInfo & AddVariable(const TString &expression, const TString &title="", const TString &unit="", Double_t min=0, Double_t max=0, char varType='F', Bool_t normalized=kTRUE, void *external=0)
add a variable (can be a complex expression) to the set of variables used in the MV analysis ...
void DisableWriting(Bool_t setter)
Definition: MethodBase.h:396
Long64_t GetNTrainingEvents() const
Definition: DataSet.h:90
A TTree object has a header with a name and a title.
Definition: TTree.h:98
void ReadWeightsFromXML(void *wghtnode)
read weights of sub-classifiers of MethodCategory from xml weight file
Bool_t PassesCut(const Event *ev, UInt_t methodIdx)
virtual void Print() const
get maximum length of variable names
Definition: Ranking.cxx:111
TString GetMethodTypeName() const
Definition: MethodBase.h:297
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:582
const Bool_t kTRUE
Definition: Rtypes.h:91
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
void WriteStateToXML(void *parent) const
general method used in writing the header of the weight files where the used variables, variable transformation type etc.
void SetNormalization(const TString &norm)
Definition: DataSetInfo.h:133
void SetRootDir(TDirectory *d)
Definition: DataSetInfo.h:186
virtual const std::vector< Float_t > & GetRegressionValues()
returns the mva value of the right sub-classifier
std::vector< TString > SplitString(const TString &theOpt, const char separator) const
splits the option string at 'separator' and fills the list 'splitV' with the primitive strings ...
Definition: Tools.cxx:1207
std::vector< VariableInfo > & GetVariableInfos()
Definition: DataSetInfo.h:112
virtual void SetAnalysisType(Types::EAnalysisType type)
Definition: MethodBase.h:390
Definition: math.cpp:60
void ProcessOptions()
process user options