ROOT  6.06/09
Reference Guide
RuleFit.cxx
Go to the documentation of this file.
1 // @(#)root/tmva $Id$
2 // Author: Andreas Hoecker, Joerg Stelzer, Fredrik Tegenfeldt, Helge Voss
3 
4 /**********************************************************************************
5  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis *
6  * Package: TMVA *
7  * Class : Rule *
8  * Web : http://tmva.sourceforge.net *
9  * *
10  * Description: *
11  * A class describung a 'rule' *
12  * Each internal node of a tree defines a rule from all the parental nodes. *
13  * A rule with 0 or 1 nodes in the list is a root rule -> corresponds to a0. *
14  * Input: a decision tree (in the constructor) *
15  * its coefficient *
16  * *
17  * *
18  * Authors (alphabetical): *
19  * Fredrik Tegenfeldt <Fredrik.Tegenfeldt@cern.ch> - Iowa State U., USA *
20  * *
21  * Copyright (c) 2005: *
22  * CERN, Switzerland *
23  * Iowa State U. *
24  * MPI-K Heidelberg, Germany *
25  * *
26  * Redistribution and use in source and binary forms, with or without *
27  * modification, are permitted according to the terms listed in LICENSE *
28  * (http://tmva.sourceforge.net/LICENSE) *
29  **********************************************************************************/
30 
31 #include <algorithm>
32 
33 #include "TKey.h"
34 #include "TRandom3.h"
35 
36 #include "TMVA/SeparationBase.h"
37 #include "TMVA/GiniIndex.h"
38 #include "TMVA/RuleFit.h"
39 #include "TMVA/MethodRuleFit.h"
40 #include "TMVA/Timer.h"
41 #include "TMVA/Tools.h"
42 #include "TMVA/Factory.h" // for root base dir
43 
45 
46 ////////////////////////////////////////////////////////////////////////////////
47 /// constructor
48 
49 TMVA::RuleFit::RuleFit( const MethodBase *rfbase )
50 : fVisHistsUseImp( kTRUE ),
51  fLogger( new MsgLogger("RuleFit") )
52 {
53  Initialize( rfbase );
54  std::srand( randSEED ); // initialize random number generator used by std::random_shuffle
55 }
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 /// default constructor
59 
61  : fNTreeSample(0)
62  , fNEveEffTrain(0)
63  , fMethodRuleFit(0)
64  , fMethodBase(0)
65  , fVisHistsUseImp( kTRUE )
66  , fLogger( new MsgLogger("RuleFit") )
67 {
68  std::srand( randSEED ); // initialize random number generator used by std::random_shuffle
69 }
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// destructor
73 
75 {
76  delete fLogger;
77 }
78 
79 ////////////////////////////////////////////////////////////////////////////////
80 /// init effective number of events (using event weights)
81 
83 {
84  UInt_t neve = fTrainingEvents.size();
85  if (neve==0) return;
86  //
87  fNEveEffTrain = CalcWeightSum( &fTrainingEvents );
88  //
89 }
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// initialize pointers
93 
94 void TMVA::RuleFit::InitPtrs( const MethodBase *rfbase )
95 {
96  this->SetMethodBase(rfbase);
97  fRuleEnsemble.Initialize( this );
98  fRuleFitParams.SetRuleFit( this );
99 }
100 
101 ////////////////////////////////////////////////////////////////////////////////
102 /// initialize the parameters of the RuleFit method and make rules
103 
105 {
106  InitPtrs(rfbase);
107 
108  if (fMethodRuleFit){
109  fMethodRuleFit->Data()->SetCurrentType(Types::kTraining);
110  UInt_t nevents = fMethodRuleFit->Data()->GetNTrainingEvents();
111  std::vector<const TMVA::Event*> tmp;
112  for (Long64_t ievt=0; ievt<nevents; ievt++) {
113  const Event *event = fMethodRuleFit->GetEvent(ievt);
114  tmp.push_back(event);
115  }
116  SetTrainingEvents( tmp );
117  }
118  // SetTrainingEvents( fMethodRuleFit->GetTrainingEvents() );
119 
120  InitNEveEff();
121 
122  MakeForest();
123 
124  // Make the model - Rule + Linear (if fDoLinear is true)
125  fRuleEnsemble.MakeModel();
126 
127  // init rulefit params
128  fRuleFitParams.Init();
129 
130 }
131 
132 ////////////////////////////////////////////////////////////////////////////////
133 /// set MethodBase
134 
136 {
137  fMethodBase = rfbase;
138  fMethodRuleFit = dynamic_cast<const MethodRuleFit *>(rfbase);
139 }
140 
141 ////////////////////////////////////////////////////////////////////////////////
142 /// copy method
143 
144 void TMVA::RuleFit::Copy( const RuleFit& other )
145 {
146  if(this != &other) {
147  fMethodRuleFit = other.GetMethodRuleFit();
148  fMethodBase = other.GetMethodBase();
149  fTrainingEvents = other.GetTrainingEvents();
150  // fSubsampleEvents = other.GetSubsampleEvents();
151 
152  fForest = other.GetForest();
153  fRuleEnsemble = other.GetRuleEnsemble();
154  }
155 }
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// calculate the sum of weights
159 
160 Double_t TMVA::RuleFit::CalcWeightSum( const std::vector<const Event *> *events, UInt_t neve )
161 {
162  if (events==0) return 0.0;
163  if (neve==0) neve=events->size();
164  //
165  Double_t sumw=0;
166  for (UInt_t ie=0; ie<neve; ie++) {
167  sumw += ((*events)[ie])->GetWeight();
168  }
169  return sumw;
170 }
171 
172 ////////////////////////////////////////////////////////////////////////////////
173 /// set the current message type to that of mlog for this class and all other subtools
174 
176 {
177  fLogger->SetMinType(t);
178  fRuleEnsemble.SetMsgType(t);
179  fRuleFitParams.SetMsgType(t);
180 }
181 
182 ////////////////////////////////////////////////////////////////////////////////
183 /// build the decision tree using fNTreeSample events from fTrainingEventsRndm
184 
186 {
187  if (dt==0) return;
188  if (fMethodRuleFit==0) {
189  Log() << kFATAL << "RuleFit::BuildTree() - Attempting to build a tree NOT from a MethodRuleFit" << Endl;
190  }
191  std::vector<const Event *> evevec;
192  for (UInt_t ie=0; ie<fNTreeSample; ie++) {
193  evevec.push_back(fTrainingEventsRndm[ie]);
194  }
195  dt->BuildTree(evevec);
196  if (fMethodRuleFit->GetPruneMethod() != DecisionTree::kNoPruning) {
197  dt->SetPruneMethod(fMethodRuleFit->GetPruneMethod());
198  dt->SetPruneStrength(fMethodRuleFit->GetPruneStrength());
199  dt->PruneTree();
200  }
201 }
202 
203 ////////////////////////////////////////////////////////////////////////////////
204 /// make a forest of decisiontrees
205 
207 {
208  if (fMethodRuleFit==0) {
209  Log() << kFATAL << "RuleFit::BuildTree() - Attempting to build a tree NOT from a MethodRuleFit" << Endl;
210  }
211  Log() << kDEBUG << "Creating a forest with " << fMethodRuleFit->GetNTrees() << " decision trees" << Endl;
212  Log() << kDEBUG << "Each tree is built using a random subsample with " << fNTreeSample << " events" << Endl;
213  //
214  Timer timer( fMethodRuleFit->GetNTrees(), "RuleFit" );
215 
216  // Double_t fsig;
217  Int_t nsig,nbkg;
218  //
219  TRandom3 rndGen;
220  //
221  //
222  // First save all event weights.
223  // Weights are modifed by the boosting.
224  // Those weights we do not want for the later fitting.
225  //
226  Bool_t useBoost = fMethodRuleFit->UseBoost(); // (AdaBoost (True) or RandomForest/Tree (False)
227 
228  if (useBoost) SaveEventWeights();
229 
230  for (Int_t i=0; i<fMethodRuleFit->GetNTrees(); i++) {
231  // timer.DrawProgressBar(i);
232  if (!useBoost) ReshuffleEvents();
233  nsig=0;
234  nbkg=0;
235  for (UInt_t ie = 0; ie<fNTreeSample; ie++) {
236  if (fMethodBase->DataInfo().IsSignal(fTrainingEventsRndm[ie])) nsig++; // ignore weights here
237  else nbkg++;
238  }
239  // fsig = Double_t(nsig)/Double_t(nsig+nbkg);
240  // do not implement the above in this release...just set it to default
241 
242  DecisionTree *dt=nullptr;
243  Bool_t tryAgain=kTRUE;
244  Int_t ntries=0;
245  const Int_t ntriesMax=10;
246  Double_t frnd = 0.;
247  while (tryAgain) {
248  frnd = 100*rndGen.Uniform( fMethodRuleFit->GetMinFracNEve(), 0.5*fMethodRuleFit->GetMaxFracNEve() );
249  Int_t iclass = 0; // event class being treated as signal during training
250  Bool_t useRandomisedTree = !useBoost;
251  dt = new DecisionTree( fMethodRuleFit->GetSeparationBase(), frnd, fMethodRuleFit->GetNCuts(), &(fMethodRuleFit->DataInfo()), iclass, useRandomisedTree);
252  dt->SetNVars(fMethodBase->GetNvar());
253 
254  BuildTree(dt); // reads fNTreeSample events from fTrainingEventsRndm
255  if (dt->GetNNodes()<3) {
256  delete dt;
257  dt=0;
258  }
259  ntries++;
260  tryAgain = ((dt==0) && (ntries<ntriesMax));
261  }
262  if (dt) {
263  fForest.push_back(dt);
264  if (useBoost) Boost(dt);
265 
266  } else {
267 
268  Log() << kWARNING << "------------------------------------------------------------------" << Endl;
269  Log() << kWARNING << " Failed growing a tree even after " << ntriesMax << " trials" << Endl;
270  Log() << kWARNING << " Possible solutions: " << Endl;
271  Log() << kWARNING << " 1. increase the number of training events" << Endl;
272  Log() << kWARNING << " 2. set a lower min fraction cut (fEventsMin)" << Endl;
273  Log() << kWARNING << " 3. maybe also decrease the max fraction cut (fEventsMax)" << Endl;
274  Log() << kWARNING << " If the above warning occurs rarely only, it can be ignored" << Endl;
275  Log() << kWARNING << "------------------------------------------------------------------" << Endl;
276  }
277 
278  Log() << kDEBUG << "Built tree with minimum cut at N = " << frnd <<"% events"
279  << " => N(nodes) = " << fForest.back()->GetNNodes()
280  << " ; n(tries) = " << ntries
281  << Endl;
282  }
283 
284  // Now restore event weights
285  if (useBoost) RestoreEventWeights();
286 
287  // print statistics on the forest created
288  ForestStatistics();
289 }
290 
291 ////////////////////////////////////////////////////////////////////////////////
292 /// save event weights - must be done before making the forest
293 
295 {
296  fEventWeights.clear();
297  for (std::vector<const Event*>::iterator e=fTrainingEvents.begin(); e!=fTrainingEvents.end(); e++) {
298  Double_t w = (*e)->GetBoostWeight();
299  fEventWeights.push_back(w);
300  }
301 }
302 
303 ////////////////////////////////////////////////////////////////////////////////
304 /// save event weights - must be done before making the forest
305 
307 {
308  UInt_t ie=0;
309  if (fEventWeights.size() != fTrainingEvents.size()) {
310  Log() << kERROR << "RuleFit::RestoreEventWeights() called without having called SaveEventWeights() before!" << Endl;
311  return;
312  }
313  for (std::vector<const Event*>::iterator e=fTrainingEvents.begin(); e!=fTrainingEvents.end(); e++) {
314  (*e)->SetBoostWeight(fEventWeights[ie]);
315  ie++;
316  }
317 }
318 
319 ////////////////////////////////////////////////////////////////////////////////
320 /// Boost the events. The algorithm below is the called AdaBoost.
321 /// See MethodBDT for details.
322 /// Actually, this is a more or less copy of MethodBDT::AdaBoost().
323 
325 {
326  Double_t sumw=0; // sum of initial weights - all events
327  Double_t sumwfalse=0; // idem, only missclassified events
328  //
329  std::vector<Char_t> correctSelected; // <--- boolean stored
330  //
331  for (std::vector<const Event*>::iterator e=fTrainingEvents.begin(); e!=fTrainingEvents.end(); e++) {
332  Bool_t isSignalType = (dt->CheckEvent(*e,kTRUE) > 0.5 );
333  Double_t w = (*e)->GetWeight();
334  sumw += w;
335  //
336  if (isSignalType == fMethodBase->DataInfo().IsSignal(*e)) { // correctly classified
337  correctSelected.push_back(kTRUE);
338  }
339  else { // missclassified
340  sumwfalse+= w;
341  correctSelected.push_back(kFALSE);
342  }
343  }
344  // missclassification error
345  Double_t err = sumwfalse/sumw;
346  // calculate boost weight for missclassified events
347  // use for now the exponent = 1.0
348  // one could have w = ((1-err)/err)^beta
349  Double_t boostWeight = (err>0 ? (1.0-err)/err : 1000.0);
350  Double_t newSumw=0.0;
351  UInt_t ie=0;
352  // set new weight to missclassified events
353  for (std::vector<const Event*>::iterator e=fTrainingEvents.begin(); e!=fTrainingEvents.end(); e++) {
354  if (!correctSelected[ie])
355  (*e)->SetBoostWeight( (*e)->GetBoostWeight() * boostWeight);
356  newSumw+=(*e)->GetWeight();
357  ie++;
358  }
359  // reweight all events
360  Double_t scale = sumw/newSumw;
361  for (std::vector<const Event*>::iterator e=fTrainingEvents.begin(); e!=fTrainingEvents.end(); e++) {
362  (*e)->SetBoostWeight( (*e)->GetBoostWeight() * scale);
363  }
364  Log() << kDEBUG << "boostWeight = " << boostWeight << " scale = " << scale << Endl;
365 }
366 
367 ////////////////////////////////////////////////////////////////////////////////
368 /// summary of statistics of all trees
369 /// * end-nodes: average and spread
370 
372 {
373  UInt_t ntrees = fForest.size();
374  if (ntrees==0) return;
375  const DecisionTree *tree;
376  Double_t sumn2 = 0;
377  Double_t sumn = 0;
378  Double_t nd;
379  for (UInt_t i=0; i<ntrees; i++) {
380  tree = fForest[i];
381  nd = Double_t(tree->GetNNodes());
382  sumn += nd;
383  sumn2 += nd*nd;
384  }
385  Double_t sig = TMath::Sqrt( gTools().ComputeVariance( sumn2, sumn, ntrees ));
386  Log() << kVERBOSE << "Nodes in trees: average & std dev = " << sumn/ntrees << " , " << sig << Endl;
387 }
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 ///
391 /// Fit the coefficients for the rule ensemble
392 ///
393 
395 {
396  Log() << kVERBOSE << "Fitting rule/linear terms" << Endl;
397  fRuleFitParams.MakeGDPath();
398 }
399 
400 ////////////////////////////////////////////////////////////////////////////////
401 /// calculates the importance of each rule
402 
404 {
405  Log() << kVERBOSE << "Calculating importance" << Endl;
406  fRuleEnsemble.CalcImportance();
407  fRuleEnsemble.CleanupRules();
408  fRuleEnsemble.CleanupLinear();
409  fRuleEnsemble.CalcVarImportance();
410  Log() << kVERBOSE << "Filling rule statistics" << Endl;
411  fRuleEnsemble.RuleResponseStats();
412 }
413 
414 ////////////////////////////////////////////////////////////////////////////////
415 /// evaluate single event
416 
418 {
419  return fRuleEnsemble.EvalEvent( e );
420 }
421 
422 ////////////////////////////////////////////////////////////////////////////////
423 /// set the training events randomly
424 
425 void TMVA::RuleFit::SetTrainingEvents( const std::vector<const Event *>& el )
426 {
427  if (fMethodRuleFit==0) Log() << kFATAL << "RuleFit::SetTrainingEvents - MethodRuleFit not initialized" << Endl;
428  UInt_t neve = el.size();
429  if (neve==0) Log() << kWARNING << "An empty sample of training events was given" << Endl;
430 
431  // copy vector
432  fTrainingEvents.clear();
433  fTrainingEventsRndm.clear();
434  for (UInt_t i=0; i<neve; i++) {
435  fTrainingEvents.push_back(static_cast< const Event *>(el[i]));
436  fTrainingEventsRndm.push_back(static_cast< const Event *>(el[i]));
437  }
438 
439  // Re-shuffle the vector, ie, recreate it in a random order
440  std::random_shuffle( fTrainingEventsRndm.begin(), fTrainingEventsRndm.end() );
441 
442  // fraction events per tree
443  fNTreeSample = static_cast<UInt_t>(neve*fMethodRuleFit->GetTreeEveFrac());
444  Log() << kDEBUG << "Number of events per tree : " << fNTreeSample
445  << " ( N(events) = " << neve << " )"
446  << " randomly drawn without replacement" << Endl;
447 }
448 
449 ////////////////////////////////////////////////////////////////////////////////
450 /// draw a random subsample of the training events without replacement
451 
452 void TMVA::RuleFit::GetRndmSampleEvents(std::vector< const Event * > & evevec, UInt_t nevents)
453 {
454  ReshuffleEvents();
455  if ((nevents<fTrainingEventsRndm.size()) && (nevents>0)) {
456  evevec.resize(nevents);
457  for (UInt_t ie=0; ie<nevents; ie++) {
458  evevec[ie] = fTrainingEventsRndm[ie];
459  }
460  }
461  else {
462  Log() << kWARNING << "GetRndmSampleEvents() : requested sub sample size larger than total size (BUG!).";
463  }
464 }
465 ////////////////////////////////////////////////////////////////////////////////
466 /// normalize rule importance hists
467 ///
468 /// if all weights are positive, the scale will be 1/maxweight
469 /// if minimum weight < 0, then the scale will be 1/max(maxweight,abs(minweight))
470 ///
471 
472 void TMVA::RuleFit::NormVisHists(std::vector<TH2F *> & hlist)
473 {
474  if (hlist.empty()) return;
475  //
476  Double_t wmin=0;
477  Double_t wmax=0;
478  Double_t w,wm;
479  Double_t awmin;
480  Double_t scale;
481  for (UInt_t i=0; i<hlist.size(); i++) {
482  TH2F *hs = hlist[i];
483  w = hs->GetMaximum();
484  wm = hs->GetMinimum();
485  if (i==0) {
486  wmin=wm;
487  wmax=w;
488  }
489  else {
490  if (w>wmax) wmax=w;
491  if (wm<wmin) wmin=wm;
492  }
493  }
494  awmin = TMath::Abs(wmin);
495  Double_t usemin,usemax;
496  if (awmin>wmax) {
497  scale = 1.0/awmin;
498  usemin = -1.0;
499  usemax = scale*wmax;
500  }
501  else {
502  scale = 1.0/wmax;
503  usemin = scale*wmin;
504  usemax = 1.0;
505  }
506 
507  //
508  for (UInt_t i=0; i<hlist.size(); i++) {
509  TH2F *hs = hlist[i];
510  hs->Scale(scale);
511  hs->SetMinimum(usemin);
512  hs->SetMaximum(usemax);
513  }
514 }
515 
516 ////////////////////////////////////////////////////////////////////////////////
517 /// Fill cut
518 
519 void TMVA::RuleFit::FillCut(TH2F* h2, const Rule *rule, Int_t vind)
520 {
521  if (rule==0) return;
522  if (h2==0) return;
523  //
524  Double_t rmin, rmax;
525  Bool_t dormin,dormax;
526  Bool_t ruleHasVar = rule->GetRuleCut()->GetCutRange(vind,rmin,rmax,dormin,dormax);
527  if (!ruleHasVar) return;
528  //
529  Int_t firstbin = h2->GetBin(1,1,1);
530  if(firstbin<0) firstbin=0;
531  Int_t lastbin = h2->GetBin(h2->GetNbinsX(),1,1);
532  Int_t binmin=(dormin ? h2->FindBin(rmin,0.5):firstbin);
533  Int_t binmax=(dormax ? h2->FindBin(rmax,0.5):lastbin);
534  Int_t fbin;
535  Double_t xbinw = h2->GetXaxis()->GetBinWidth(firstbin);
536  Double_t fbmin = h2->GetXaxis()->GetBinLowEdge(binmin-firstbin+1);
537  Double_t lbmax = h2->GetXaxis()->GetBinLowEdge(binmax-firstbin+1)+xbinw;
538  Double_t fbfrac = (dormin ? ((fbmin+xbinw-rmin)/xbinw):1.0);
539  Double_t lbfrac = (dormax ? ((rmax-lbmax+xbinw)/xbinw):1.0);
540  Double_t f;
541  Double_t xc;
542  Double_t val;
543 
544  for (Int_t bin = binmin; bin<binmax+1; bin++) {
545  fbin = bin-firstbin+1;
546  if (bin==binmin) {
547  f = fbfrac;
548  }
549  else if (bin==binmax) {
550  f = lbfrac;
551  }
552  else {
553  f = 1.0;
554  }
555  xc = h2->GetXaxis()->GetBinCenter(fbin);
556  //
557  if (fVisHistsUseImp) {
558  val = rule->GetImportance();
559  }
560  else {
561  val = rule->GetCoefficient()*rule->GetSupport();
562  }
563  h2->Fill(xc,0.5,val*f);
564  }
565 }
566 
567 ////////////////////////////////////////////////////////////////////////////////
568 /// fill lin
569 
571 {
572  if (h2==0) return;
573  if (!fRuleEnsemble.DoLinear()) return;
574  //
575  Int_t firstbin = 1;
576  Int_t lastbin = h2->GetNbinsX();
577  Double_t xc;
578  Double_t val;
579  if (fVisHistsUseImp) {
580  val = fRuleEnsemble.GetLinImportance(vind);
581  }
582  else {
583  val = fRuleEnsemble.GetLinCoefficients(vind);
584  }
585  for (Int_t bin = firstbin; bin<lastbin+1; bin++) {
586  xc = h2->GetXaxis()->GetBinCenter(bin);
587  h2->Fill(xc,0.5,val);
588  }
589 }
590 
591 ////////////////////////////////////////////////////////////////////////////////
592 /// fill rule correlation between vx and vy, weighted with either the importance or the coefficient
593 
594 void TMVA::RuleFit::FillCorr(TH2F* h2,const Rule *rule,Int_t vx, Int_t vy)
595 {
596  if (rule==0) return;
597  if (h2==0) return;
598  Double_t val;
599  if (fVisHistsUseImp) {
600  val = rule->GetImportance();
601  }
602  else {
603  val = rule->GetCoefficient()*rule->GetSupport();
604  }
605  //
606  Double_t rxmin, rxmax, rymin, rymax;
607  Bool_t dorxmin, dorxmax, dorymin, dorymax;
608  //
609  // Get range in rule for X and Y
610  //
611  Bool_t ruleHasVarX = rule->GetRuleCut()->GetCutRange(vx,rxmin,rxmax,dorxmin,dorxmax);
612  Bool_t ruleHasVarY = rule->GetRuleCut()->GetCutRange(vy,rymin,rymax,dorymin,dorymax);
613  if (!(ruleHasVarX || ruleHasVarY)) return;
614  // min max of varX and varY in hist
615  Double_t vxmin = (dorxmin ? rxmin:h2->GetXaxis()->GetXmin());
616  Double_t vxmax = (dorxmax ? rxmax:h2->GetXaxis()->GetXmax());
617  Double_t vymin = (dorymin ? rymin:h2->GetYaxis()->GetXmin());
618  Double_t vymax = (dorymax ? rymax:h2->GetYaxis()->GetXmax());
619  // min max bin in X and Y
620  Int_t binxmin = h2->GetXaxis()->FindBin(vxmin);
621  Int_t binxmax = h2->GetXaxis()->FindBin(vxmax);
622  Int_t binymin = h2->GetYaxis()->FindBin(vymin);
623  Int_t binymax = h2->GetYaxis()->FindBin(vymax);
624  // bin widths
625  Double_t xbinw = h2->GetXaxis()->GetBinWidth(binxmin);
626  Double_t ybinw = h2->GetYaxis()->GetBinWidth(binxmin);
627  Double_t xbinmin = h2->GetXaxis()->GetBinLowEdge(binxmin);
628  Double_t xbinmax = h2->GetXaxis()->GetBinLowEdge(binxmax)+xbinw;
629  Double_t ybinmin = h2->GetYaxis()->GetBinLowEdge(binymin);
630  Double_t ybinmax = h2->GetYaxis()->GetBinLowEdge(binymax)+ybinw;
631  // fraction of edges
632  Double_t fxbinmin = (dorxmin ? ((xbinmin+xbinw-vxmin)/xbinw):1.0);
633  Double_t fxbinmax = (dorxmax ? ((vxmax-xbinmax+xbinw)/xbinw):1.0);
634  Double_t fybinmin = (dorymin ? ((ybinmin+ybinw-vymin)/ybinw):1.0);
635  Double_t fybinmax = (dorymax ? ((vymax-ybinmax+ybinw)/ybinw):1.0);
636  //
637  Double_t fx,fy;
638  Double_t xc,yc;
639  // fill histo
640  for (Int_t binx = binxmin; binx<binxmax+1; binx++) {
641  if (binx==binxmin) {
642  fx = fxbinmin;
643  }
644  else if (binx==binxmax) {
645  fx = fxbinmax;
646  }
647  else {
648  fx = 1.0;
649  }
650  xc = h2->GetXaxis()->GetBinCenter(binx);
651  for (Int_t biny = binymin; biny<binymax+1; biny++) {
652  if (biny==binymin) {
653  fy = fybinmin;
654  }
655  else if (biny==binymax) {
656  fy = fybinmax;
657  }
658  else {
659  fy = 1.0;
660  }
661  yc = h2->GetYaxis()->GetBinCenter(biny);
662  h2->Fill(xc,yc,val*fx*fy);
663  }
664  }
665 }
666 
667 ////////////////////////////////////////////////////////////////////////////////
668 /// help routine to MakeVisHists() - fills for all variables
669 
670 void TMVA::RuleFit::FillVisHistCut(const Rule* rule, std::vector<TH2F *> & hlist)
671 {
672  Int_t nhists = hlist.size();
673  Int_t nvar = fMethodBase->GetNvar();
674  if (nhists!=nvar) Log() << kFATAL << "BUG TRAP: number of hists is not equal the number of variables!" << Endl;
675  //
676  std::vector<Int_t> vindex;
677  TString hstr;
678  // not a nice way to do a check...
679  for (Int_t ih=0; ih<nhists; ih++) {
680  hstr = hlist[ih]->GetTitle();
681  for (Int_t iv=0; iv<nvar; iv++) {
682  if (fMethodBase->GetInputTitle(iv) == hstr)
683  vindex.push_back(iv);
684  }
685  }
686  //
687  for (Int_t iv=0; iv<nvar; iv++) {
688  if (rule) {
689  if (rule->ContainsVariable(vindex[iv])) {
690  FillCut(hlist[iv],rule,vindex[iv]);
691  }
692  }
693  else {
694  FillLin(hlist[iv],vindex[iv]);
695  }
696  }
697 }
698 ////////////////////////////////////////////////////////////////////////////////
699 /// help routine to MakeVisHists() - fills for all correlation plots
700 
701 void TMVA::RuleFit::FillVisHistCorr(const Rule * rule, std::vector<TH2F *> & hlist)
702 {
703  if (rule==0) return;
704  Double_t ruleimp = rule->GetImportance();
705  if (!(ruleimp>0)) return;
706  if (ruleimp<fRuleEnsemble.GetImportanceCut()) return;
707  //
708  Int_t nhists = hlist.size();
709  Int_t nvar = fMethodBase->GetNvar();
710  Int_t ncorr = (nvar*(nvar+1)/2)-nvar;
711  if (nhists!=ncorr) Log() << kERROR << "BUG TRAP: number of corr hists is not correct! ncorr = "
712  << ncorr << " nvar = " << nvar << " nhists = " << nhists << Endl;
713  //
714  std::vector< std::pair<Int_t,Int_t> > vindex;
715  TString hstr, var1, var2;
716  Int_t iv1=0,iv2=0;
717  // not a nice way to do a check...
718  for (Int_t ih=0; ih<nhists; ih++) {
719  hstr = hlist[ih]->GetName();
720  if (GetCorrVars( hstr, var1, var2 )) {
721  iv1 = fMethodBase->DataInfo().FindVarIndex( var1 );
722  iv2 = fMethodBase->DataInfo().FindVarIndex( var2 );
723  vindex.push_back( std::pair<Int_t,Int_t>(iv2,iv1) ); // pair X, Y
724  }
725  else {
726  Log() << kERROR << "BUG TRAP: should not be here - failed getting var1 and var2" << Endl;
727  }
728  }
729  //
730  for (Int_t ih=0; ih<nhists; ih++) {
731  if ( (rule->ContainsVariable(vindex[ih].first)) ||
732  (rule->ContainsVariable(vindex[ih].second)) ) {
733  FillCorr(hlist[ih],rule,vindex[ih].first,vindex[ih].second);
734  }
735  }
736 }
737 ////////////////////////////////////////////////////////////////////////////////
738 /// get first and second variables from title
739 
741 {
742  var1="";
743  var2="";
744  if(!title.BeginsWith("scat_")) return kFALSE;
745 
746  TString titleCopy = title(5,title.Length());
747  if(titleCopy.Index("_RF2D")>=0) titleCopy.Remove(titleCopy.Index("_RF2D"));
748 
749  Int_t splitPos = titleCopy.Index("_vs_");
750  if(splitPos>=0) { // there is a _vs_ in the string
751  var1 = titleCopy(0,splitPos);
752  var2 = titleCopy(splitPos+4, titleCopy.Length());
753  return kTRUE;
754  }
755  else {
756  var1 = titleCopy;
757  return kFALSE;
758  }
759 }
760 ////////////////////////////////////////////////////////////////////////////////
761 /// this will create histograms visualizing the rule ensemble
762 
764 {
765  const TString directories[5] = { "InputVariables_Id",
766  "InputVariables_Deco",
767  "InputVariables_PCA",
768  "InputVariables_Gauss",
769  "InputVariables_Gauss_Deco" };
770 
771  const TString corrDirName = "CorrelationPlots";
772 
773  TDirectory* rootDir = Factory::RootBaseDir();
774  TDirectory* varDir = 0;
775  TDirectory* corrDir = 0;
776 
777  TDirectory* methodDir = fMethodBase->BaseDir();
778  TString varDirName;
779  //
780  Bool_t done=(rootDir==0);
781  Int_t type=0;
782  if (done) {
783  Log() << kWARNING << "No basedir - BUG??" << Endl;
784  return;
785  }
786  while (!done) {
787  varDir = (TDirectory*)rootDir->Get( directories[type] );
788  type++;
789  done = ((varDir!=0) || (type>4));
790  }
791  if (varDir==0) {
792  Log() << kWARNING << "No input variable directory found - BUG?" << Endl;
793  return;
794  }
795  corrDir = (TDirectory*)varDir->Get( corrDirName );
796  if (corrDir==0) {
797  Log() << kWARNING << "No correlation directory found" << Endl;
798  Log() << kWARNING << "Check for other warnings related to correlation histograms" << Endl;
799  return;
800  }
801  if (methodDir==0) {
802  Log() << kWARNING << "No rulefit method directory found - BUG?" << Endl;
803  return;
804  }
805 
806  varDirName = varDir->GetName();
807  varDir->cd();
808  //
809  // get correlation plot directory
810  corrDir = (TDirectory *)varDir->Get(corrDirName);
811  if (corrDir==0) {
812  Log() << kWARNING << "No correlation directory found : " << corrDirName << Endl;
813  return;
814  }
815 
816  // how many plots are in the var directory?
817  Int_t noPlots = ((varDir->GetListOfKeys())->GetEntries()) / 2;
818  Log() << kDEBUG << "Got number of plots = " << noPlots << Endl;
819 
820  // loop over all objects in directory
821  std::vector<TH2F *> h1Vector;
822  std::vector<TH2F *> h2CorrVector;
823  TIter next(varDir->GetListOfKeys());
824  TKey *key;
825  while ((key = (TKey*)next())) {
826  // make sure, that we only look at histograms
827  TClass *cl = gROOT->GetClass(key->GetClassName());
828  if (!cl->InheritsFrom(TH1F::Class())) continue;
829  TH1F *sig = (TH1F*)key->ReadObj();
830  TString hname= sig->GetName();
831  Log() << kDEBUG << "Got histogram : " << hname << Endl;
832 
833  // check for all signal histograms
834  if (hname.Contains("__S")){ // found a new signal plot
835  TString htitle = sig->GetTitle();
836  htitle.ReplaceAll("signal","");
837  TString newname = hname;
838  newname.ReplaceAll("__Signal","__RF");
839  newname.ReplaceAll("__S","__RF");
840 
841  methodDir->cd();
842  TH2F *newhist = new TH2F(newname,htitle,sig->GetNbinsX(),sig->GetXaxis()->GetXmin(),sig->GetXaxis()->GetXmax(),
843  1,sig->GetYaxis()->GetXmin(),sig->GetYaxis()->GetXmax());
844  varDir->cd();
845  h1Vector.push_back( newhist );
846  }
847  }
848  //
849  corrDir->cd();
850  TString var1,var2;
851  TIter nextCorr(corrDir->GetListOfKeys());
852  while ((key = (TKey*)nextCorr())) {
853  // make sure, that we only look at histograms
854  TClass *cl = gROOT->GetClass(key->GetClassName());
855  if (!cl->InheritsFrom(TH2F::Class())) continue;
856  TH2F *sig = (TH2F*)key->ReadObj();
857  TString hname= sig->GetName();
858 
859  // check for all signal histograms
860  if ((hname.Contains("scat_")) && (hname.Contains("_Signal"))) {
861  Log() << kDEBUG << "Got histogram (2D) : " << hname << Endl;
862  TString htitle = sig->GetTitle();
863  htitle.ReplaceAll("(Signal)","");
864  TString newname = hname;
865  newname.ReplaceAll("_Signal","_RF2D");
866 
867  methodDir->cd();
868  const Int_t rebin=2;
869  TH2F *newhist = new TH2F(newname,htitle,
870  sig->GetNbinsX()/rebin,sig->GetXaxis()->GetXmin(),sig->GetXaxis()->GetXmax(),
871  sig->GetNbinsY()/rebin,sig->GetYaxis()->GetXmin(),sig->GetYaxis()->GetXmax());
872  if (GetCorrVars( newname, var1, var2 )) {
873  Int_t iv1 = fMethodBase->DataInfo().FindVarIndex(var1);
874  Int_t iv2 = fMethodBase->DataInfo().FindVarIndex(var2);
875  if (iv1<0) {
876  sig->GetYaxis()->SetTitle(var1);
877  }
878  else {
879  sig->GetYaxis()->SetTitle(fMethodBase->GetInputTitle(iv1));
880  }
881  if (iv2<0) {
882  sig->GetXaxis()->SetTitle(var2);
883  }
884  else {
885  sig->GetXaxis()->SetTitle(fMethodBase->GetInputTitle(iv2));
886  }
887  }
888  corrDir->cd();
889  h2CorrVector.push_back( newhist );
890  }
891  }
892 
893 
894  varDir->cd();
895  // fill rules
896  UInt_t nrules = fRuleEnsemble.GetNRules();
897  const Rule *rule;
898  for (UInt_t i=0; i<nrules; i++) {
899  rule = fRuleEnsemble.GetRulesConst(i);
900  FillVisHistCut(rule, h1Vector);
901  }
902  // fill linear terms and normalise hists
903  FillVisHistCut(0, h1Vector);
904  NormVisHists(h1Vector);
905 
906  //
907  corrDir->cd();
908  // fill rules
909  for (UInt_t i=0; i<nrules; i++) {
910  rule = fRuleEnsemble.GetRulesConst(i);
911  FillVisHistCorr(rule, h2CorrVector);
912  }
913  NormVisHists(h2CorrVector);
914 
915  // write histograms to file
916  methodDir->cd();
917  for (UInt_t i=0; i<h1Vector.size(); i++) h1Vector[i]->Write();
918  for (UInt_t i=0; i<h2CorrVector.size(); i++) h2CorrVector[i]->Write();
919 }
920 
921 ////////////////////////////////////////////////////////////////////////////////
922 /// this will create a histograms intended rather for debugging or for the curious user
923 
925 {
926  TDirectory* methodDir = fMethodBase->BaseDir();
927  if (methodDir==0) {
928  Log() << kWARNING << "<MakeDebugHists> No rulefit method directory found - bug?" << Endl;
929  return;
930  }
931  //
932  methodDir->cd();
933  std::vector<Double_t> distances;
934  std::vector<Double_t> fncuts;
935  std::vector<Double_t> fnvars;
936  const Rule *ruleA;
937  const Rule *ruleB;
938  Double_t dABmin=1000000.0;
939  Double_t dABmax=-1.0;
940  UInt_t nrules = fRuleEnsemble.GetNRules();
941  for (UInt_t i=0; i<nrules; i++) {
942  ruleA = fRuleEnsemble.GetRulesConst(i);
943  for (UInt_t j=i+1; j<nrules; j++) {
944  ruleB = fRuleEnsemble.GetRulesConst(j);
945  Double_t dAB = ruleA->RuleDist( *ruleB, kTRUE );
946  if (dAB>-0.5) {
947  UInt_t nc = ruleA->GetNcuts();
948  UInt_t nv = ruleA->GetNumVarsUsed();
949  distances.push_back(dAB);
950  fncuts.push_back(static_cast<Double_t>(nc));
951  fnvars.push_back(static_cast<Double_t>(nv));
952  if (dAB<dABmin) dABmin=dAB;
953  if (dAB>dABmax) dABmax=dAB;
954  }
955  }
956  }
957  //
958  TH1F *histDist = new TH1F("RuleDist","Rule distances",100,dABmin,dABmax);
959  TTree *distNtuple = new TTree("RuleDistNtuple","RuleDist ntuple");
960  Double_t ntDist;
961  Double_t ntNcuts;
962  Double_t ntNvars;
963  distNtuple->Branch("dist", &ntDist, "dist/D");
964  distNtuple->Branch("ncuts",&ntNcuts, "ncuts/D");
965  distNtuple->Branch("nvars",&ntNvars, "nvars/D");
966  //
967  for (UInt_t i=0; i<distances.size(); i++) {
968  histDist->Fill(distances[i]);
969  ntDist = distances[i];
970  ntNcuts = fncuts[i];
971  ntNvars = fnvars[i];
972  distNtuple->Fill();
973  }
974  distNtuple->Write();
975 }
void ForestStatistics()
summary of statistics of all trees
Definition: RuleFit.cxx:371
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52
void SetPruneMethod(EPruneMethod m=kCostComplexityPruning)
Definition: DecisionTree.h:148
void MakeForest()
make a forest of decisiontrees
Definition: RuleFit.cxx:206
virtual Int_t FindBin(Double_t x, Double_t y=0, Double_t z=0)
Return Global bin number corresponding to x,y,z.
Definition: TH1.cxx:3478
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6174
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
static TDirectory * RootBaseDir()
Definition: Factory.h:228
Random number generator class based on M.
Definition: TRandom3.h:29
MsgLogger & Endl(MsgLogger &ml)
Definition: MsgLogger.h:162
long long Long64_t
Definition: RtypesCore.h:69
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:394
const std::vector< const TMVA::Event * > & GetTrainingEvents() const
Definition: RuleFit.h:141
void CalcImportance()
calculates the importance of each rule
Definition: RuleFit.cxx:403
virtual TList * GetListOfKeys() const
Definition: TDirectory.h:155
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
Definition: TDirectory.cxx:727
ClassImp(TSeqCollection) Int_t TSeqCollection TIter next(this)
Return index of object in collection.
Ssiz_t Length() const
Definition: TString.h:390
THist< 2, float > TH2F
Definition: THist.h:321
void FillVisHistCorr(const Rule *rule, std::vector< TH2F * > &hlist)
help routine to MakeVisHists() - fills for all correlation plots
Definition: RuleFit.cxx:701
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
virtual Int_t Fill()
Fill all branches.
Definition: TTree.cxx:4328
void SetMsgType(EMsgType t)
set the current message type to that of mlog for this class and all other subtools ...
Definition: RuleFit.cxx:175
THist< 1, float > TH1F
Definition: THist.h:315
Bool_t GetCorrVars(TString &title, TString &var1, TString &var2)
get first and second variables from title
Definition: RuleFit.cxx:740
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:489
void InitNEveEff()
init effective number of events (using event weights)
Definition: RuleFit.cxx:82
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:395
#define gROOT
Definition: TROOT.h:340
void FitCoefficients()
Fit the coefficients for the rule ensemble.
Definition: RuleFit.cxx:394
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition: TAxis.cxx:511
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
Definition: TString.h:558
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
ClassImp(TMVA::RuleFit) TMVA
constructor
Definition: RuleFit.cxx:44
Double_t RuleDist(const Rule &other, Bool_t useCutValue) const
Returns: -1.0 : rules are NOT equal, i.e, variables and/or cut directions are wrong >=0: rules are eq...
Definition: Rule.cxx:183
UInt_t GetNcuts() const
Definition: Rule.h:139
Tools & gTools()
Definition: Tools.cxx:79
RuleFit(void)
default constructor
Definition: RuleFit.cxx:60
TStopwatch timer
Definition: pirndm.C:37
void BuildTree(TMVA::DecisionTree *dt)
build the decision tree using fNTreeSample events from fTrainingEventsRndm
Definition: RuleFit.cxx:185
void Class()
Definition: Class.C:29
const MethodBase * GetMethodBase() const
Definition: RuleFit.h:153
void SetTrainingEvents(const std::vector< const TMVA::Event * > &el)
set the training events randomly
Definition: RuleFit.cxx:425
void GetRndmSampleEvents(std::vector< const TMVA::Event * > &evevec, UInt_t nevents)
draw a random subsample of the training events without replacement
Definition: RuleFit.cxx:452
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition: TKey.h:30
virtual ~RuleFit(void)
destructor
Definition: RuleFit.cxx:74
ROOT::Math::KDTree< _DataPoint > * BuildTree(const std::vector< const _DataPoint * > &vDataPoints, const unsigned int iBucketSize)
Double_t GetXmin() const
Definition: TAxis.h:137
const RuleCut * GetRuleCut() const
Definition: Rule.h:145
void SetMethodBase(const MethodBase *rfbase)
set MethodBase
Definition: RuleFit.cxx:135
Double_t CheckEvent(const TMVA::Event *, Bool_t UseYesNoLeaf=kFALSE) const
the event e is put into the decision tree (starting at the root node) and the output is NodeType (sig...
const std::vector< const TMVA::DecisionTree * > & GetForest() const
Definition: RuleFit.h:147
void SetNVars(Int_t n)
Definition: DecisionTree.h:202
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TTree.cxx:8811
void RestoreEventWeights()
save event weights - must be done before making the forest
Definition: RuleFit.cxx:306
void FillVisHistCut(const Rule *rule, std::vector< TH2F * > &hlist)
help routine to MakeVisHists() - fills for all variables
Definition: RuleFit.cxx:670
void Copy(const RuleFit &other)
copy method
Definition: RuleFit.cxx:144
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition: tmvaglob.cxx:176
Bool_t ContainsVariable(UInt_t iv) const
check if variable in node
Definition: Rule.cxx:131
void FillCorr(TH2F *h2, const TMVA::Rule *rule, Int_t v1, Int_t v2)
fill rule correlation between vx and vy, weighted with either the importance or the coefficient ...
Definition: RuleFit.cxx:594
void MakeDebugHists()
this will create a histograms intended rather for debugging or for the curious user ...
Definition: RuleFit.cxx:924
static const Int_t randSEED
Definition: RuleFit.h:179
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:256
Double_t GetImportance() const
Definition: Rule.h:151
const RuleEnsemble & GetRuleEnsemble() const
Definition: RuleFit.h:148
EMsgType
Definition: Types.h:61
void SetPruneStrength(Double_t p)
Definition: DecisionTree.h:154
unsigned int UInt_t
Definition: RtypesCore.h:42
void FillLin(TH2F *h2, Int_t vind)
fill lin
Definition: RuleFit.cxx:570
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
Double_t GetWeight(Double_t x) const
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:81
TAxis * GetYaxis()
Definition: TH1.h:320
Bool_t GetCutRange(Int_t sel, Double_t &rmin, Double_t &rmax, Bool_t &dormin, Bool_t &dormax) const
get cut range for a given selector
Definition: RuleCut.cxx:170
const MethodRuleFit * GetMethodRuleFit() const
Definition: RuleFit.h:152
void Boost(TMVA::DecisionTree *dt)
Boost the events.
Definition: RuleFit.cxx:324
virtual Int_t GetBin(Int_t binx, Int_t biny, Int_t binz=0) const
Return Global bin number corresponding to binx,y,z.
Definition: TH2.cxx:956
virtual Int_t FindBin(Double_t x)
Find bin number corresponding to abscissa x.
Definition: TAxis.cxx:264
TString & Remove(Ssiz_t pos)
Definition: TString.h:616
void SaveEventWeights()
save event weights - must be done before making the forest
Definition: RuleFit.cxx:294
double f(double x)
double Double_t
Definition: RtypesCore.h:55
Describe directory structure in memory.
Definition: TDirectory.h:41
int type
Definition: TGX11.cxx:120
void MakeVisHists()
this will create histograms visualizing the rule ensemble
Definition: RuleFit.cxx:763
UInt_t GetNNodes() const
Definition: BinaryTree.h:92
Double_t GetXmax() const
Definition: TAxis.h:138
UInt_t GetNumVarsUsed() const
Definition: Rule.h:136
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:606
void InitPtrs(const TMVA::MethodBase *rfbase)
initialize pointers
Definition: RuleFit.cxx:94
Double_t GetCoefficient() const
Definition: Rule.h:147
Double_t CalcWeightSum(const std::vector< const TMVA::Event * > *events, UInt_t neve=0)
calculate the sum of weights
Definition: RuleFit.cxx:160
void NormVisHists(std::vector< TH2F * > &hlist)
normalize rule importance hists
Definition: RuleFit.cxx:472
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:449
virtual Int_t Branch(TCollection *list, Int_t bufsize=32000, Int_t splitlevel=99, const char *name="")
Create one branch for each element in the collection.
Definition: TTree.cxx:1624
virtual Int_t GetNbinsY() const
Definition: TH1.h:297
Double_t PruneTree(const EventConstList *validationSample=NULL)
prune (get rid of internal nodes) the Decision tree to avoid overtraining serveral different pruning ...
Abstract ClassifierFactory template that handles arbitrary types.
void FillCut(TH2F *h2, const TMVA::Rule *rule, Int_t vind)
Fill cut.
Definition: RuleFit.cxx:519
virtual Bool_t cd(const char *path=0)
Change current directory to "this" directory.
Definition: TDirectory.cxx:433
Double_t GetSupport() const
Definition: Rule.h:148
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
UInt_t BuildTree(const EventConstList &eventSample, DecisionTreeNode *node=NULL)
building the decision tree by recursively calling the splitting of one (root-) node into two daughter...
Double_t EvalEvent(const Event &e)
evaluate single event
Definition: RuleFit.cxx:417
A TTree object has a header with a name and a title.
Definition: TTree.h:94
Bool_t InheritsFrom(const char *cl) const
Return kTRUE if this class inherits from a class with name "classname".
Definition: TClass.cxx:4579
void Initialize(const TMVA::MethodBase *rfbase)
initialize the parameters of the RuleFit method and make rules
Definition: RuleFit.cxx:104
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:582
virtual Double_t GetMaximum(Double_t maxval=FLT_MAX) const
Return maximum value smaller than maxval of bins in the range, unless the value has been overridden b...
Definition: TH1.cxx:7921
const Bool_t kTRUE
Definition: Rtypes.h:91
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:285
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
virtual Double_t GetMinimum(Double_t minval=-FLT_MAX) const
Return minimum value larger than minval of bins in the range, unless the value has been overridden by...
Definition: TH1.cxx:8006
Definition: math.cpp:60
TAxis * GetXaxis()
Definition: TH1.h:319