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