ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MethodRuleFit.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Fredrik Tegenfeldt
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : MethodRuleFit *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * Implementation (see header file for description) *
12  * *
13  * Authors (alphabetical): *
14  * Fredrik Tegenfeldt <Fredrik.Tegenfeldt@cern.ch> - Iowa State U., USA *
15  * *
16  * Copyright (c) 2005: *
17  * CERN, Switzerland *
18  * Iowa State U. *
19  * MPI-K Heidelberg, Germany *
20  * *
21  * Redistribution and use in source and binary forms, with or without *
22  * modification, are permitted according to the terms listed in LICENSE *
23  * (http://tmva.sourceforge.net/LICENSE) *
24  **********************************************************************************/
25 
26 //_______________________________________________________________________
27 //
28 // J Friedman's RuleFit method
29 //_______________________________________________________________________
30 
31 #include "TMVA/MethodRuleFit.h"
32 
33 #include <algorithm>
34 #include <list>
35 
36 #include "Riostream.h"
37 #include "TRandom3.h"
38 #include "TMath.h"
39 #include "TMatrix.h"
40 #include "TDirectory.h"
41 
42 #include "TMVA/ClassifierFactory.h"
43 #include "TMVA/Config.h"
44 #include "TMVA/CrossEntropy.h"
45 #include "TMVA/DataSet.h"
46 #include "TMVA/DecisionTree.h"
47 #include "TMVA/GiniIndex.h"
49 #include "TMVA/MsgLogger.h"
50 #include "TMVA/Ranking.h"
51 #include "TMVA/RuleFitAPI.h"
52 #include "TMVA/SdivSqrtSplusB.h"
53 #include "TMVA/SeparationBase.h"
54 #include "TMVA/Timer.h"
55 #include "TMVA/Tools.h"
56 #include "TMVA/Types.h"
57 
58 using std::min;
59 
60 REGISTER_METHOD(RuleFit)
61 
62 ClassImp(TMVA::MethodRuleFit)
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 /// standard constructor
66 
67 TMVA::MethodRuleFit::MethodRuleFit( const TString& jobName,
68  const TString& methodTitle,
69  DataSetInfo& theData,
70  const TString& theOption,
71  TDirectory* theTargetDir ) :
72  MethodBase( jobName, Types::kRuleFit, methodTitle, theData, theOption, theTargetDir )
73  , fSignalFraction(0)
74  , fNTImportance(0)
75  , fNTCoefficient(0)
76  , fNTSupport(0)
77  , fNTNcuts(0)
78  , fNTNvars(0)
79  , fNTPtag(0)
80  , fNTPss(0)
81  , fNTPsb(0)
82  , fNTPbs(0)
83  , fNTPbb(0)
84  , fNTSSB(0)
85  , fNTType(0)
86  , fUseRuleFitJF(kFALSE)
87  , fRFNrules(0)
88  , fRFNendnodes(0)
89  , fNTrees(0)
90  , fTreeEveFrac(0)
91  , fSepType(0)
92  , fMinFracNEve(0)
93  , fMaxFracNEve(0)
94  , fNCuts(0)
95  , fPruneMethod(TMVA::DecisionTree::kCostComplexityPruning)
96  , fPruneStrength(0)
97  , fUseBoost(kFALSE)
98  , fGDPathEveFrac(0)
99  , fGDValidEveFrac(0)
100  , fGDTau(0)
101  , fGDTauPrec(0)
102  , fGDTauMin(0)
103  , fGDTauMax(0)
104  , fGDTauScan(0)
105  , fGDPathStep(0)
106  , fGDNPathSteps(0)
107  , fGDErrScale(0)
108  , fMinimp(0)
109  , fRuleMinDist(0)
110  , fLinQuantile(0)
111 {
112 }
113 
114 ////////////////////////////////////////////////////////////////////////////////
115 /// constructor from weight file
116 
118  const TString& theWeightFile,
119  TDirectory* theTargetDir ) :
120  MethodBase( Types::kRuleFit, theData, theWeightFile, theTargetDir )
121  , fSignalFraction(0)
122  , fNTImportance(0)
123  , fNTCoefficient(0)
124  , fNTSupport(0)
125  , fNTNcuts(0)
126  , fNTNvars(0)
127  , fNTPtag(0)
128  , fNTPss(0)
129  , fNTPsb(0)
130  , fNTPbs(0)
131  , fNTPbb(0)
132  , fNTSSB(0)
133  , fNTType(0)
134  , fUseRuleFitJF(kFALSE)
135  , fRFNrules(0)
136  , fRFNendnodes(0)
137  , fNTrees(0)
138  , fTreeEveFrac(0)
139  , fSepType(0)
140  , fMinFracNEve(0)
141  , fMaxFracNEve(0)
142  , fNCuts(0)
143  , fPruneMethod(TMVA::DecisionTree::kCostComplexityPruning)
144  , fPruneStrength(0)
145  , fUseBoost(kFALSE)
146  , fGDPathEveFrac(0)
147  , fGDValidEveFrac(0)
148  , fGDTau(0)
149  , fGDTauPrec(0)
150  , fGDTauMin(0)
151  , fGDTauMax(0)
152  , fGDTauScan(0)
153  , fGDPathStep(0)
154  , fGDNPathSteps(0)
155  , fGDErrScale(0)
156  , fMinimp(0)
157  , fRuleMinDist(0)
158  , fLinQuantile(0)
159 {
160 }
161 
162 ////////////////////////////////////////////////////////////////////////////////
163 /// destructor
164 
166 {
167  for (UInt_t i=0; i<fEventSample.size(); i++) delete fEventSample[i];
168  for (UInt_t i=0; i<fForest.size(); i++) delete fForest[i];
169 }
170 
171 ////////////////////////////////////////////////////////////////////////////////
172 /// RuleFit can handle classification with 2 classes
173 
175 {
176  if (type == Types::kClassification && numberClasses == 2) return kTRUE;
177  return kFALSE;
178 }
179 
180 ////////////////////////////////////////////////////////////////////////////////
181 /// define the options (their key words) that can be set in the option string
182 /// know options.
183 ///---------
184 /// general
185 ///---------
186 /// RuleFitModule <string>
187 /// available values are: RFTMVA - use TMVA implementation
188 /// RFFriedman - use Friedmans original implementation
189 ///----------------------
190 /// Path search (fitting)
191 ///----------------------
192 /// GDTau <float> gradient-directed path: fit threshhold, default
193 /// GDTauPrec <float> gradient-directed path: precision of estimated tau
194 /// GDStep <float> gradient-directed path: step size
195 /// GDNSteps <float> gradient-directed path: number of steps
196 /// GDErrScale <float> stop scan when error>scale*errmin
197 ///-----------------
198 /// Tree generation
199 ///-----------------
200 /// fEventsMin <float> minimum fraction of events in a splittable node
201 /// fEventsMax <float> maximum fraction of events in a splittable node
202 /// nTrees <float> number of trees in forest.
203 /// ForestType <string>
204 /// available values are: Random - create forest using random subsample and only random variables subset at each node
205 /// AdaBoost - create forest with boosted events
206 ///
207 ///-----------------
208 /// Model creation
209 ///-----------------
210 /// RuleMinDist <float> min distance allowed between rules
211 /// MinImp <float> minimum rule importance accepted
212 /// Model <string> model to be used
213 /// available values are: ModRuleLinear <default>
214 /// ModRule
215 /// ModLinear
216 ///
217 ///-----------------
218 /// Friedmans module
219 ///-----------------
220 /// RFWorkDir <string> directory where Friedmans module (rf_go.exe) is installed
221 /// RFNrules <int> maximum number of rules allowed
222 /// RFNendnodes <int> average number of end nodes in the forest of trees
223 ///
224 
226 {
227  DeclareOptionRef(fGDTau=-1, "GDTau", "Gradient-directed (GD) path: default fit cut-off");
228  DeclareOptionRef(fGDTauPrec=0.01, "GDTauPrec", "GD path: precision of tau");
229  DeclareOptionRef(fGDPathStep=0.01, "GDStep", "GD path: step size");
230  DeclareOptionRef(fGDNPathSteps=10000, "GDNSteps", "GD path: number of steps");
231  DeclareOptionRef(fGDErrScale=1.1, "GDErrScale", "Stop scan when error > scale*errmin");
232  DeclareOptionRef(fLinQuantile, "LinQuantile", "Quantile of linear terms (removes outliers)");
233  DeclareOptionRef(fGDPathEveFrac=0.5, "GDPathEveFrac", "Fraction of events used for the path search");
234  DeclareOptionRef(fGDValidEveFrac=0.5, "GDValidEveFrac", "Fraction of events used for the validation");
235  // tree options
236  DeclareOptionRef(fMinFracNEve=0.1, "fEventsMin", "Minimum fraction of events in a splittable node");
237  DeclareOptionRef(fMaxFracNEve=0.9, "fEventsMax", "Maximum fraction of events in a splittable node");
238  DeclareOptionRef(fNTrees=20, "nTrees", "Number of trees in forest.");
239 
240  DeclareOptionRef(fForestTypeS="AdaBoost", "ForestType", "Method to use for forest generation (AdaBoost or RandomForest)");
241  AddPreDefVal(TString("AdaBoost"));
242  AddPreDefVal(TString("Random"));
243  // rule cleanup options
244  DeclareOptionRef(fRuleMinDist=0.001, "RuleMinDist", "Minimum distance between rules");
245  DeclareOptionRef(fMinimp=0.01, "MinImp", "Minimum rule importance accepted");
246  // rule model option
247  DeclareOptionRef(fModelTypeS="ModRuleLinear", "Model", "Model to be used");
248  AddPreDefVal(TString("ModRule"));
249  AddPreDefVal(TString("ModRuleLinear"));
250  AddPreDefVal(TString("ModLinear"));
251  DeclareOptionRef(fRuleFitModuleS="RFTMVA", "RuleFitModule","Which RuleFit module to use");
252  AddPreDefVal(TString("RFTMVA"));
253  AddPreDefVal(TString("RFFriedman"));
254 
255  DeclareOptionRef(fRFWorkDir="./rulefit", "RFWorkDir", "Friedman\'s RuleFit module (RFF): working dir");
256  DeclareOptionRef(fRFNrules=2000, "RFNrules", "RFF: Mximum number of rules");
257  DeclareOptionRef(fRFNendnodes=4, "RFNendnodes", "RFF: Average number of end nodes");
258 }
259 
260 ////////////////////////////////////////////////////////////////////////////////
261 /// process the options specified by the user
262 
264 {
265  if (IgnoreEventsWithNegWeightsInTraining()) {
266  Log() << kFATAL << "Mechanism to ignore events with negative weights in training not yet available for method: "
267  << GetMethodTypeName()
268  << " --> please remove \"IgnoreNegWeightsInTraining\" option from booking string."
269  << Endl;
270  }
271 
272  fRuleFitModuleS.ToLower();
273  if (fRuleFitModuleS == "rftmva") fUseRuleFitJF = kFALSE;
274  else if (fRuleFitModuleS == "rffriedman") fUseRuleFitJF = kTRUE;
275  else fUseRuleFitJF = kTRUE;
276 
277  fSepTypeS.ToLower();
278  if (fSepTypeS == "misclassificationerror") fSepType = new MisClassificationError();
279  else if (fSepTypeS == "giniindex") fSepType = new GiniIndex();
280  else if (fSepTypeS == "crossentropy") fSepType = new CrossEntropy();
281  else fSepType = new SdivSqrtSplusB();
282 
283  fModelTypeS.ToLower();
284  if (fModelTypeS == "modlinear" ) fRuleFit.SetModelLinear();
285  else if (fModelTypeS == "modrule" ) fRuleFit.SetModelRules();
286  else fRuleFit.SetModelFull();
287 
288  fPruneMethodS.ToLower();
289  if (fPruneMethodS == "expectederror" ) fPruneMethod = DecisionTree::kExpectedErrorPruning;
290  else if (fPruneMethodS == "costcomplexity" ) fPruneMethod = DecisionTree::kCostComplexityPruning;
291  else fPruneMethod = DecisionTree::kNoPruning;
292 
293  fForestTypeS.ToLower();
294  if (fForestTypeS == "random" ) fUseBoost = kFALSE;
295  else if (fForestTypeS == "adaboost" ) fUseBoost = kTRUE;
296  else fUseBoost = kTRUE;
297  //
298  // if creating the forest by boosting the events
299  // the full training sample is used per tree
300  // -> only true for the TMVA version of RuleFit.
301  if (fUseBoost && (!fUseRuleFitJF)) fTreeEveFrac = 1.0;
302 
303  // check event fraction for tree generation
304  // if <0 set to automatic number
305  if (fTreeEveFrac<=0) {
306  Int_t nevents = Data()->GetNTrainingEvents();
307  Double_t n = static_cast<Double_t>(nevents);
308  fTreeEveFrac = min( 0.5, (100.0 +6.0*sqrt(n))/n);
309  }
310  // verify ranges of options
311  VerifyRange(Log(), "nTrees", fNTrees,0,100000,20);
312  VerifyRange(Log(), "MinImp", fMinimp,0.0,1.0,0.0);
313  VerifyRange(Log(), "GDTauPrec", fGDTauPrec,1e-5,5e-1);
314  VerifyRange(Log(), "GDTauMin", fGDTauMin,0.0,1.0);
315  VerifyRange(Log(), "GDTauMax", fGDTauMax,fGDTauMin,1.0);
316  VerifyRange(Log(), "GDPathStep", fGDPathStep,0.0,100.0,0.01);
317  VerifyRange(Log(), "GDErrScale", fGDErrScale,1.0,100.0,1.1);
318  VerifyRange(Log(), "GDPathEveFrac", fGDPathEveFrac,0.01,0.9,0.5);
319  VerifyRange(Log(), "GDValidEveFrac",fGDValidEveFrac,0.01,1.0-fGDPathEveFrac,1.0-fGDPathEveFrac);
320  VerifyRange(Log(), "fEventsMin", fMinFracNEve,0.0,1.0);
321  VerifyRange(Log(), "fEventsMax", fMaxFracNEve,fMinFracNEve,1.0);
322 
323  fRuleFit.GetRuleEnsemblePtr()->SetLinQuantile(fLinQuantile);
324  fRuleFit.GetRuleFitParamsPtr()->SetGDTauRange(fGDTauMin,fGDTauMax);
325  fRuleFit.GetRuleFitParamsPtr()->SetGDTau(fGDTau);
326  fRuleFit.GetRuleFitParamsPtr()->SetGDTauPrec(fGDTauPrec);
327  fRuleFit.GetRuleFitParamsPtr()->SetGDTauScan(fGDTauScan);
328  fRuleFit.GetRuleFitParamsPtr()->SetGDPathStep(fGDPathStep);
329  fRuleFit.GetRuleFitParamsPtr()->SetGDNPathSteps(fGDNPathSteps);
330  fRuleFit.GetRuleFitParamsPtr()->SetGDErrScale(fGDErrScale);
331  fRuleFit.SetImportanceCut(fMinimp);
332  fRuleFit.SetRuleMinDist(fRuleMinDist);
333 
334 
335  // check if Friedmans module is used.
336  // print a message concerning the options.
337  if (fUseRuleFitJF) {
338  Log() << kINFO << "" << Endl;
339  Log() << kINFO << "--------------------------------------" <<Endl;
340  Log() << kINFO << "Friedmans RuleFit module is selected." << Endl;
341  Log() << kINFO << "Only the following options are used:" << Endl;
342  Log() << kINFO << Endl;
343  Log() << kINFO << gTools().Color("bold") << " Model" << gTools().Color("reset") << Endl;
344  Log() << kINFO << gTools().Color("bold") << " RFWorkDir" << gTools().Color("reset") << Endl;
345  Log() << kINFO << gTools().Color("bold") << " RFNrules" << gTools().Color("reset") << Endl;
346  Log() << kINFO << gTools().Color("bold") << " RFNendnodes" << gTools().Color("reset") << Endl;
347  Log() << kINFO << gTools().Color("bold") << " GDNPathSteps" << gTools().Color("reset") << Endl;
348  Log() << kINFO << gTools().Color("bold") << " GDPathStep" << gTools().Color("reset") << Endl;
349  Log() << kINFO << gTools().Color("bold") << " GDErrScale" << gTools().Color("reset") << Endl;
350  Log() << kINFO << "--------------------------------------" <<Endl;
351  Log() << kINFO << Endl;
352  }
353 
354  // Select what weight to use in the 'importance' rule visualisation plots.
355  // Note that if UseCoefficientsVisHists() is selected, the following weight is used:
356  // w = rule coefficient * rule support
357  // The support is a positive number which is 0 if no events are accepted by the rule.
358  // Normally the importance gives more useful information.
359  //
360  //fRuleFit.UseCoefficientsVisHists();
361  fRuleFit.UseImportanceVisHists();
362 
363  fRuleFit.SetMsgType( Log().GetMinType() );
364 
365  if (HasTrainingTree()) InitEventSample();
366 
367 }
368 
369 ////////////////////////////////////////////////////////////////////////////////
370 /// initialize the monitoring ntuple
371 
373 {
374  BaseDir()->cd();
375  fMonitorNtuple= new TTree("MonitorNtuple_RuleFit","RuleFit variables");
376  fMonitorNtuple->Branch("importance",&fNTImportance,"importance/D");
377  fMonitorNtuple->Branch("support",&fNTSupport,"support/D");
378  fMonitorNtuple->Branch("coefficient",&fNTCoefficient,"coefficient/D");
379  fMonitorNtuple->Branch("ncuts",&fNTNcuts,"ncuts/I");
380  fMonitorNtuple->Branch("nvars",&fNTNvars,"nvars/I");
381  fMonitorNtuple->Branch("type",&fNTType,"type/I");
382  fMonitorNtuple->Branch("ptag",&fNTPtag,"ptag/D");
383  fMonitorNtuple->Branch("pss",&fNTPss,"pss/D");
384  fMonitorNtuple->Branch("psb",&fNTPsb,"psb/D");
385  fMonitorNtuple->Branch("pbs",&fNTPbs,"pbs/D");
386  fMonitorNtuple->Branch("pbb",&fNTPbb,"pbb/D");
387  fMonitorNtuple->Branch("soversb",&fNTSSB,"soversb/D");
388 }
389 
390 ////////////////////////////////////////////////////////////////////////////////
391 /// default initialization
392 
394 {
395  // the minimum requirement to declare an event signal-like
396  SetSignalReferenceCut( 0.0 );
397 
398  // set variables that used to be options
399  // any modifications are then made in ProcessOptions()
400  fLinQuantile = 0.025; // Quantile of linear terms (remove outliers)
401  fTreeEveFrac = -1.0; // Fraction of events used to train each tree
402  fNCuts = 20; // Number of steps during node cut optimisation
403  fSepTypeS = "GiniIndex"; // Separation criterion for node splitting; see BDT
404  fPruneMethodS = "NONE"; // Pruning method; see BDT
405  fPruneStrength = 3.5; // Pruning strength; see BDT
406  fGDTauMin = 0.0; // Gradient-directed path: min fit threshold (tau)
407  fGDTauMax = 1.0; // Gradient-directed path: max fit threshold (tau)
408  fGDTauScan = 1000; // Gradient-directed path: number of points scanning for best tau
409 
410 }
411 
412 ////////////////////////////////////////////////////////////////////////////////
413 /// write all Events from the Tree into a vector of Events, that are
414 /// more easily manipulated.
415 /// This method should never be called without existing trainingTree, as it
416 /// the vector of events from the ROOT training tree
417 
419 {
420  if (Data()->GetNEvents()==0) Log() << kFATAL << "<Init> Data().TrainingTree() is zero pointer" << Endl;
421 
422  Int_t nevents = Data()->GetNEvents();
423  for (Int_t ievt=0; ievt<nevents; ievt++){
424  const Event * ev = GetEvent(ievt);
425  fEventSample.push_back( new Event(*ev));
426  }
427  if (fTreeEveFrac<=0) {
428  Double_t n = static_cast<Double_t>(nevents);
429  fTreeEveFrac = min( 0.5, (100.0 +6.0*sqrt(n))/n);
430  }
431  if (fTreeEveFrac>1.0) fTreeEveFrac=1.0;
432  //
433  std::random_shuffle(fEventSample.begin(), fEventSample.end());
434  //
435  Log() << kDEBUG << "Set sub-sample fraction to " << fTreeEveFrac << Endl;
436 }
437 
438 ////////////////////////////////////////////////////////////////////////////////
439 
441 {
443  // training of rules
444 
445  InitMonitorNtuple();
446 
447  // fill the STL Vector with the event sample
448  this->InitEventSample();
449 
450  if (fUseRuleFitJF) {
451  TrainJFRuleFit();
452  }
453  else {
454  TrainTMVARuleFit();
455  }
456  fRuleFit.GetRuleEnsemblePtr()->ClearRuleMap();
458 }
459 
460 ////////////////////////////////////////////////////////////////////////////////
461 /// training of rules using TMVA implementation
462 
464 {
465  if (IsNormalised()) Log() << kFATAL << "\"Normalise\" option cannot be used with RuleFit; "
466  << "please remove the optoin from the configuration string, or "
467  << "use \"!Normalise\""
468  << Endl;
469 
470  // timer
471  Timer timer( 1, GetName() );
472 
473  // test tree nmin cut -> for debug purposes
474  // the routine will generate trees with stopping cut on N(eve) given by
475  // a fraction between [20,N(eve)-1].
476  //
477  // MakeForestRnd();
478  // exit(1);
479  //
480 
481  // Init RuleFit object and create rule ensemble
482  // + make forest & rules
483  fRuleFit.Initialize( this );
484 
485  // Make forest of decision trees
486  // if (fRuleFit.GetRuleEnsemble().DoRules()) fRuleFit.MakeForest();
487 
488  // Fit the rules
489  Log() << kDEBUG << "Fitting rule coefficients ..." << Endl;
490  fRuleFit.FitCoefficients();
491 
492  // Calculate importance
493  Log() << kDEBUG << "Computing rule and variable importance" << Endl;
494  fRuleFit.CalcImportance();
495 
496  // Output results and fill monitor ntuple
497  fRuleFit.GetRuleEnsemblePtr()->Print();
498  //
499  Log() << kDEBUG << "Filling rule ntuple" << Endl;
500  UInt_t nrules = fRuleFit.GetRuleEnsemble().GetRulesConst().size();
501  const Rule *rule;
502  for (UInt_t i=0; i<nrules; i++ ) {
503  rule = fRuleFit.GetRuleEnsemble().GetRulesConst(i);
504  fNTImportance = rule->GetRelImportance();
505  fNTSupport = rule->GetSupport();
506  fNTCoefficient = rule->GetCoefficient();
507  fNTType = (rule->IsSignalRule() ? 1:-1 );
508  fNTNvars = rule->GetRuleCut()->GetNvars();
509  fNTNcuts = rule->GetRuleCut()->GetNcuts();
510  fNTPtag = fRuleFit.GetRuleEnsemble().GetRulePTag(i); // should be identical with support
511  fNTPss = fRuleFit.GetRuleEnsemble().GetRulePSS(i);
512  fNTPsb = fRuleFit.GetRuleEnsemble().GetRulePSB(i);
513  fNTPbs = fRuleFit.GetRuleEnsemble().GetRulePBS(i);
514  fNTPbb = fRuleFit.GetRuleEnsemble().GetRulePBB(i);
515  fNTSSB = rule->GetSSB();
516  fMonitorNtuple->Fill();
517  }
518  Log() << kDEBUG << "Training done" << Endl;
519 
520  fRuleFit.MakeVisHists();
521 
522  fRuleFit.MakeDebugHists();
523 }
524 
525 ////////////////////////////////////////////////////////////////////////////////
526 /// training of rules using Jerome Friedmans implementation
527 
529 {
530  fRuleFit.InitPtrs( this );
531  Data()->SetCurrentType(Types::kTraining);
532  UInt_t nevents = Data()->GetNTrainingEvents();
533  std::vector<const TMVA::Event*> tmp;
534  for (Long64_t ievt=0; ievt<nevents; ievt++) {
535  const Event *event = GetEvent(ievt);
536  tmp.push_back(event);
537  }
538  fRuleFit.SetTrainingEvents( tmp );
539 
540  RuleFitAPI *rfAPI = new RuleFitAPI( this, &fRuleFit, Log().GetMinType() );
541 
542  rfAPI->WelcomeMessage();
543 
544  // timer
545  Timer timer( 1, GetName() );
546 
547  Log() << kINFO << "Training ..." << Endl;
548  rfAPI->TrainRuleFit();
549 
550  Log() << kDEBUG << "reading model summary from rf_go.exe output" << Endl;
551  rfAPI->ReadModelSum();
552 
553  // fRuleFit.GetRuleEnsemblePtr()->MakeRuleMap();
554 
555  Log() << kDEBUG << "calculating rule and variable importance" << Endl;
556  fRuleFit.CalcImportance();
557 
558  // Output results and fill monitor ntuple
559  fRuleFit.GetRuleEnsemblePtr()->Print();
560  //
561  fRuleFit.MakeVisHists();
562 
563  delete rfAPI;
564 
565  Log() << kDEBUG << "done training" << Endl;
566 }
567 
568 ////////////////////////////////////////////////////////////////////////////////
569 /// computes ranking of input variables
570 
572 {
573  // create the ranking object
574  fRanking = new Ranking( GetName(), "Importance" );
575 
576  for (UInt_t ivar=0; ivar<GetNvar(); ivar++) {
577  fRanking->AddRank( Rank( GetInputLabel(ivar), fRuleFit.GetRuleEnsemble().GetVarImportance(ivar) ) );
578  }
579 
580  return fRanking;
581 }
582 
583 ////////////////////////////////////////////////////////////////////////////////
584 /// add the rules to XML node
585 
586 void TMVA::MethodRuleFit::AddWeightsXMLTo( void* parent ) const
587 {
588  fRuleFit.GetRuleEnsemble().AddXMLTo( parent );
589 }
590 
591 ////////////////////////////////////////////////////////////////////////////////
592 /// read rules from an std::istream
593 
595 {
596  fRuleFit.GetRuleEnsemblePtr()->ReadRaw( istr );
597 }
598 
599 ////////////////////////////////////////////////////////////////////////////////
600 /// read rules from XML node
601 
603 {
604  fRuleFit.GetRuleEnsemblePtr()->ReadFromXML( wghtnode );
605 }
606 
607 ////////////////////////////////////////////////////////////////////////////////
608 /// returns MVA value for given event
609 
611 {
612  // cannot determine error
613  NoErrorCalc(err, errUpper);
614 
615  return fRuleFit.EvalEvent( *GetEvent() );
616 }
617 
618 ////////////////////////////////////////////////////////////////////////////////
619 /// write special monitoring histograms to file (here ntuple)
620 
622 {
623  BaseDir()->cd();
624  Log() << kINFO << "Write monitoring ntuple to file: " << BaseDir()->GetPath() << Endl;
625  fMonitorNtuple->Write();
626 }
627 
628 ////////////////////////////////////////////////////////////////////////////////
629 /// write specific classifier response
630 
631 void TMVA::MethodRuleFit::MakeClassSpecific( std::ostream& fout, const TString& className ) const
632 {
633  Int_t dp = fout.precision();
634  fout << " // not implemented for class: \"" << className << "\"" << std::endl;
635  fout << "};" << std::endl;
636  fout << "void " << className << "::Initialize(){}" << std::endl;
637  fout << "void " << className << "::Clear(){}" << std::endl;
638  fout << "double " << className << "::GetMvaValue__( const std::vector<double>& inputValues ) const {" << std::endl;
639  fout << " double rval=" << std::setprecision(10) << fRuleFit.GetRuleEnsemble().GetOffset() << ";" << std::endl;
640  MakeClassRuleCuts(fout);
641  MakeClassLinear(fout);
642  fout << " return rval;" << std::endl;
643  fout << "}" << std::endl;
644  fout << std::setprecision(dp);
645 }
646 
647 ////////////////////////////////////////////////////////////////////////////////
648 /// print out the rule cuts
649 
650 void TMVA::MethodRuleFit::MakeClassRuleCuts( std::ostream& fout ) const
651 {
652  Int_t dp = fout.precision();
653  if (!fRuleFit.GetRuleEnsemble().DoRules()) {
654  fout << " //" << std::endl;
655  fout << " // ==> MODEL CONTAINS NO RULES <==" << std::endl;
656  fout << " //" << std::endl;
657  return;
658  }
659  const RuleEnsemble *rens = &(fRuleFit.GetRuleEnsemble());
660  const std::vector< Rule* > *rules = &(rens->GetRulesConst());
661  const RuleCut *ruleCut;
662  //
663  std::list< std::pair<Double_t,Int_t> > sortedRules;
664  for (UInt_t ir=0; ir<rules->size(); ir++) {
665  sortedRules.push_back( std::pair<Double_t,Int_t>( (*rules)[ir]->GetImportance()/rens->GetImportanceRef(),ir ) );
666  }
667  sortedRules.sort();
668  //
669  fout << " //" << std::endl;
670  fout << " // here follows all rules ordered in importance (most important first)" << std::endl;
671  fout << " // at the end of each line, the relative importance of the rule is given" << std::endl;
672  fout << " //" << std::endl;
673  //
674  for ( std::list< std::pair<double,int> >::reverse_iterator itpair = sortedRules.rbegin();
675  itpair != sortedRules.rend(); itpair++ ) {
676  UInt_t ir = itpair->second;
677  Double_t impr = itpair->first;
678  ruleCut = (*rules)[ir]->GetRuleCut();
679  if (impr<rens->GetImportanceCut()) fout << " //" << std::endl;
680  fout << " if (" << std::flush;
681  for (UInt_t ic=0; ic<ruleCut->GetNvars(); ic++) {
682  Double_t sel = ruleCut->GetSelector(ic);
683  Double_t valmin = ruleCut->GetCutMin(ic);
684  Double_t valmax = ruleCut->GetCutMax(ic);
685  Bool_t domin = ruleCut->GetCutDoMin(ic);
686  Bool_t domax = ruleCut->GetCutDoMax(ic);
687  //
688  if (ic>0) fout << "&&" << std::flush;
689  if (domin) {
690  fout << "(" << std::setprecision(10) << valmin << std::flush;
691  fout << "<inputValues[" << sel << "])" << std::flush;
692  }
693  if (domax) {
694  if (domin) fout << "&&" << std::flush;
695  fout << "(inputValues[" << sel << "]" << std::flush;
696  fout << "<" << std::setprecision(10) << valmax << ")" <<std::flush;
697  }
698  }
699  fout << ") rval+=" << std::setprecision(10) << (*rules)[ir]->GetCoefficient() << ";" << std::flush;
700  fout << " // importance = " << Form("%3.3f",impr) << std::endl;
701  }
702  fout << std::setprecision(dp);
703 }
704 
705 ////////////////////////////////////////////////////////////////////////////////
706 /// print out the linear terms
707 
708 void TMVA::MethodRuleFit::MakeClassLinear( std::ostream& fout ) const
709 {
710  if (!fRuleFit.GetRuleEnsemble().DoLinear()) {
711  fout << " //" << std::endl;
712  fout << " // ==> MODEL CONTAINS NO LINEAR TERMS <==" << std::endl;
713  fout << " //" << std::endl;
714  return;
715  }
716  fout << " //" << std::endl;
717  fout << " // here follows all linear terms" << std::endl;
718  fout << " // at the end of each line, the relative importance of the term is given" << std::endl;
719  fout << " //" << std::endl;
720  const RuleEnsemble *rens = &(fRuleFit.GetRuleEnsemble());
721  UInt_t nlin = rens->GetNLinear();
722  for (UInt_t il=0; il<nlin; il++) {
723  if (rens->IsLinTermOK(il)) {
724  Double_t norm = rens->GetLinNorm(il);
725  Double_t imp = rens->GetLinImportance(il)/rens->GetImportanceRef();
726  fout << " rval+="
727  // << std::setprecision(10) << rens->GetLinCoefficients(il)*norm << "*std::min(" << setprecision(10) << rens->GetLinDP(il)
728  // << ", std::max( inputValues[" << il << "]," << std::setprecision(10) << rens->GetLinDM(il) << "));"
729  << std::setprecision(10) << rens->GetLinCoefficients(il)*norm
730  << "*std::min( double(" << std::setprecision(10) << rens->GetLinDP(il)
731  << "), std::max( double(inputValues[" << il << "]), double(" << std::setprecision(10) << rens->GetLinDM(il) << ")));"
732  << std::flush;
733  fout << " // importance = " << Form("%3.3f",imp) << std::endl;
734  }
735  }
736 }
737 
738 ////////////////////////////////////////////////////////////////////////////////
739 /// get help message text
740 ///
741 /// typical length of text line:
742 /// "|--------------------------------------------------------------|"
743 
745 {
746  TString col = gConfig().WriteOptionsReference() ? TString() : gTools().Color("bold");
747  TString colres = gConfig().WriteOptionsReference() ? TString() : gTools().Color("reset");
748  TString brk = gConfig().WriteOptionsReference() ? "<br>" : "";
749 
750  Log() << Endl;
751  Log() << col << "--- Short description:" << colres << Endl;
752  Log() << Endl;
753  Log() << "This method uses a collection of so called rules to create a" << Endl;
754  Log() << "discriminating scoring function. Each rule consists of a series" << Endl;
755  Log() << "of cuts in parameter space. The ensemble of rules are created" << Endl;
756  Log() << "from a forest of decision trees, trained using the training data." << Endl;
757  Log() << "Each node (apart from the root) corresponds to one rule." << Endl;
758  Log() << "The scoring function is then obtained by linearly combining" << Endl;
759  Log() << "the rules. A fitting procedure is applied to find the optimum" << Endl;
760  Log() << "set of coefficients. The goal is to find a model with few rules" << Endl;
761  Log() << "but with a strong discriminating power." << Endl;
762  Log() << Endl;
763  Log() << col << "--- Performance optimisation:" << colres << Endl;
764  Log() << Endl;
765  Log() << "There are two important considerations to make when optimising:" << Endl;
766  Log() << Endl;
767  Log() << " 1. Topology of the decision tree forest" << brk << Endl;
768  Log() << " 2. Fitting of the coefficients" << Endl;
769  Log() << Endl;
770  Log() << "The maximum complexity of the rules is defined by the size of" << Endl;
771  Log() << "the trees. Large trees will yield many complex rules and capture" << Endl;
772  Log() << "higher order correlations. On the other hand, small trees will" << Endl;
773  Log() << "lead to a smaller ensemble with simple rules, only capable of" << Endl;
774  Log() << "modeling simple structures." << Endl;
775  Log() << "Several parameters exists for controlling the complexity of the" << Endl;
776  Log() << "rule ensemble." << Endl;
777  Log() << Endl;
778  Log() << "The fitting procedure searches for a minimum using a gradient" << Endl;
779  Log() << "directed path. Apart from step size and number of steps, the" << Endl;
780  Log() << "evolution of the path is defined by a cut-off parameter, tau." << Endl;
781  Log() << "This parameter is unknown and depends on the training data." << Endl;
782  Log() << "A large value will tend to give large weights to a few rules." << Endl;
783  Log() << "Similarily, a small value will lead to a large set of rules" << Endl;
784  Log() << "with similar weights." << Endl;
785  Log() << Endl;
786  Log() << "A final point is the model used; rules and/or linear terms." << Endl;
787  Log() << "For a given training sample, the result may improve by adding" << Endl;
788  Log() << "linear terms. If best performance is optained using only linear" << Endl;
789  Log() << "terms, it is very likely that the Fisher discriminant would be" << Endl;
790  Log() << "a better choice. Ideally the fitting procedure should be able to" << Endl;
791  Log() << "make this choice by giving appropriate weights for either terms." << Endl;
792  Log() << Endl;
793  Log() << col << "--- Performance tuning via configuration options:" << colres << Endl;
794  Log() << Endl;
795  Log() << "I. TUNING OF RULE ENSEMBLE:" << Endl;
796  Log() << Endl;
797  Log() << " " << col << "ForestType " << colres
798  << ": Recomended is to use the default \"AdaBoost\"." << brk << Endl;
799  Log() << " " << col << "nTrees " << colres
800  << ": More trees leads to more rules but also slow" << Endl;
801  Log() << " performance. With too few trees the risk is" << Endl;
802  Log() << " that the rule ensemble becomes too simple." << brk << Endl;
803  Log() << " " << col << "fEventsMin " << colres << brk << Endl;
804  Log() << " " << col << "fEventsMax " << colres
805  << ": With a lower min, more large trees will be generated" << Endl;
806  Log() << " leading to more complex rules." << Endl;
807  Log() << " With a higher max, more small trees will be" << Endl;
808  Log() << " generated leading to more simple rules." << Endl;
809  Log() << " By changing this range, the average complexity" << Endl;
810  Log() << " of the rule ensemble can be controlled." << brk << Endl;
811  Log() << " " << col << "RuleMinDist " << colres
812  << ": By increasing the minimum distance between" << Endl;
813  Log() << " rules, fewer and more diverse rules will remain." << Endl;
814  Log() << " Initially it is a good idea to keep this small" << Endl;
815  Log() << " or zero and let the fitting do the selection of" << Endl;
816  Log() << " rules. In order to reduce the ensemble size," << Endl;
817  Log() << " the value can then be increased." << Endl;
818  Log() << Endl;
819  // "|--------------------------------------------------------------|"
820  Log() << "II. TUNING OF THE FITTING:" << Endl;
821  Log() << Endl;
822  Log() << " " << col << "GDPathEveFrac " << colres
823  << ": fraction of events in path evaluation" << Endl;
824  Log() << " Increasing this fraction will improve the path" << Endl;
825  Log() << " finding. However, a too high value will give few" << Endl;
826  Log() << " unique events available for error estimation." << Endl;
827  Log() << " It is recomended to usethe default = 0.5." << brk << Endl;
828  Log() << " " << col << "GDTau " << colres
829  << ": cutoff parameter tau" << Endl;
830  Log() << " By default this value is set to -1.0." << Endl;
831  // "|----------------|---------------------------------------------|"
832  Log() << " This means that the cut off parameter is" << Endl;
833  Log() << " automatically estimated. In most cases" << Endl;
834  Log() << " this should be fine. However, you may want" << Endl;
835  Log() << " to fix this value if you already know it" << Endl;
836  Log() << " and want to reduce on training time." << brk << Endl;
837  Log() << " " << col << "GDTauPrec " << colres
838  << ": precision of estimated tau" << Endl;
839  Log() << " Increase this precision to find a more" << Endl;
840  Log() << " optimum cut-off parameter." << brk << Endl;
841  Log() << " " << col << "GDNStep " << colres
842  << ": number of steps in path search" << Endl;
843  Log() << " If the number of steps is too small, then" << Endl;
844  Log() << " the program will give a warning message." << Endl;
845  Log() << Endl;
846  Log() << "III. WARNING MESSAGES" << Endl;
847  Log() << Endl;
848  Log() << col << "Risk(i+1)>=Risk(i) in path" << colres << brk << Endl;
849  Log() << col << "Chaotic behaviour of risk evolution." << colres << Endl;
850  // "|----------------|---------------------------------------------|"
851  Log() << " The error rate was still decreasing at the end" << Endl;
852  Log() << " By construction the Risk should always decrease." << Endl;
853  Log() << " However, if the training sample is too small or" << Endl;
854  Log() << " the model is overtrained, such warnings can" << Endl;
855  Log() << " occur." << Endl;
856  Log() << " The warnings can safely be ignored if only a" << Endl;
857  Log() << " few (<3) occur. If more warnings are generated," << Endl;
858  Log() << " the fitting fails." << Endl;
859  Log() << " A remedy may be to increase the value" << brk << Endl;
860  Log() << " "
861  << col << "GDValidEveFrac" << colres
862  << " to 1.0 (or a larger value)." << brk << Endl;
863  Log() << " In addition, if "
864  << col << "GDPathEveFrac" << colres
865  << " is too high" << Endl;
866  Log() << " the same warnings may occur since the events" << Endl;
867  Log() << " used for error estimation are also used for" << Endl;
868  Log() << " path estimation." << Endl;
869  Log() << " Another possibility is to modify the model - " << Endl;
870  Log() << " See above on tuning the rule ensemble." << Endl;
871  Log() << Endl;
872  Log() << col << "The error rate was still decreasing at the end of the path"
873  << colres << Endl;
874  Log() << " Too few steps in path! Increase "
875  << col << "GDNSteps" << colres << "." << Endl;
876  Log() << Endl;
877  Log() << col << "Reached minimum early in the search" << colres << Endl;
878 
879  Log() << " Minimum was found early in the fitting. This" << Endl;
880  Log() << " may indicate that the used step size "
881  << col << "GDStep" << colres << "." << Endl;
882  Log() << " was too large. Reduce it and rerun." << Endl;
883  Log() << " If the results still are not OK, modify the" << Endl;
884  Log() << " model either by modifying the rule ensemble" << Endl;
885  Log() << " or add/remove linear terms" << Endl;
886 }
void DeclareOptions()
define the options (their key words) that can be set in the option string know options.
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:823
void Init(void)
default initialization
void WelcomeMessage()
welcome message
Definition: RuleFitAPI.cxx:78
void ReadWeightsFromXML(void *wghtnode)
read rules from XML node
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
void TrainRuleFit()
Definition: RuleFitAPI.h:201
long long Long64_t
Definition: RtypesCore.h:69
void ReadWeightsFromStream(std::istream &istr)
read rules from an std::istream
const std::vector< Double_t > & GetLinCoefficients() const
Definition: RuleEnsemble.h:279
Double_t GetLinDP(int i) const
Definition: RuleEnsemble.h:294
Config & gConfig()
void InitMonitorNtuple()
initialize the monitoring ntuple
EAnalysisType
Definition: Types.h:124
void WriteMonitoringHistosToFile(void) const
write special monitoring histograms to file (here ntuple)
Basic string class.
Definition: TString.h:137
Bool_t IsSignalRule() const
Definition: Rule.h:125
Double_t GetCutMax(Int_t is) const
Definition: RuleCut.h:75
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
const Ranking * CreateRanking()
computes ranking of input variables
void TrainJFRuleFit()
training of rules using Jerome Friedmans implementation
Double_t GetImportanceRef() const
Definition: RuleEnsemble.h:274
UInt_t GetSelector(Int_t is) const
Definition: RuleCut.h:73
Double_t GetMvaValue(Double_t *err=0, Double_t *errUpper=0)
returns MVA value for given event
const std::vector< Double_t > & GetLinImportance() const
Definition: RuleEnsemble.h:281
double sqrt(double)
Tools & gTools()
Definition: Tools.cxx:79
Char_t GetCutDoMax(Int_t is) const
Definition: RuleCut.h:77
TStopwatch timer
Definition: pirndm.C:37
Double_t GetCutMin(Int_t is) const
Definition: RuleCut.h:74
void ProcessOptions()
process the options specified by the user
void MakeClassLinear(std::ostream &) const
print out the linear terms
Bool_t IsLinTermOK(int i) const
Definition: RuleEnsemble.h:303
std::vector< std::vector< double > > Data
Double_t GetSSB() const
Definition: Rule.h:123
Bool_t ReadModelSum()
read model from rulefit.sum
Definition: RuleFitAPI.cxx:546
void TrainTMVARuleFit()
training of rules using TMVA implementation
const RuleCut * GetRuleCut() const
Definition: Rule.h:145
UInt_t GetNvars() const
Definition: RuleCut.h:72
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:594
Double_t GetLinDM(int i) const
Definition: RuleEnsemble.h:293
void AddWeightsXMLTo(void *parent) const
add the rules to XML node
const RuleEnsemble * GetRuleEnsemble() const
Definition: Rule.h:146
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
UInt_t GetNcuts() const
get number of cuts
Definition: RuleCut.cxx:158
const std::vector< TMVA::Rule * > & GetRulesConst() const
Definition: RuleEnsemble.h:277
void GetHelpMessage() const
get help message text
UInt_t GetNLinear() const
Definition: RuleEnsemble.h:283
double Double_t
Definition: RtypesCore.h:55
Describe directory structure in memory.
Definition: TDirectory.h:44
int type
Definition: TGX11.cxx:120
Double_t GetRelImportance() const
Definition: Rule.h:108
Char_t GetCutDoMin(Int_t is) const
Definition: RuleCut.h:76
Double_t GetCoefficient() const
Definition: Rule.h:147
const TString & Color(const TString &)
human readable color strings
Definition: Tools.cxx:837
#define REGISTER_METHOD(CLASS)
for example
ClassImp(TMVA::MethodRuleFit) TMVA
standard constructor
Double_t GetSupport() const
Definition: Rule.h:148
Bool_t WriteOptionsReference() const
Definition: Config.h:66
void MakeClassSpecific(std::ostream &, const TString &) const
write specific classifier response
A TTree object has a header with a name and a title.
Definition: TTree.h:98
const std::vector< Double_t > & GetLinNorm() const
Definition: RuleEnsemble.h:280
virtual Bool_t HasAnalysisType(Types::EAnalysisType type, UInt_t numberClasses, UInt_t)
RuleFit can handle classification with 2 classes.
virtual ~MethodRuleFit(void)
destructor
const Bool_t kTRUE
Definition: Rtypes.h:91
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
const Int_t n
Definition: legend1.C:16
Definition: math.cpp:60
MethodRuleFit(const TString &jobName, const TString &methodTitle, DataSetInfo &theData, const TString &theOption="", TDirectory *theTargetDir=0)
void InitEventSample(void)
write all Events from the Tree into a vector of Events, that are more easily manipulated.
void MakeClassRuleCuts(std::ostream &) const
print out the rule cuts