source: trunk/MagicSoft/Mars/mranforest/MRanForest.cc@ 7535

Last change on this file since 7535 was 7535, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 17.8 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Hengstebeck 3/2003 <mailto:hengsteb@physik.hu-berlin.de>
19!
20! Copyright: MAGIC Software Development, 2000-2005
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MRanForest
28//
29// ParameterContainer for Forest structure
30//
31// A random forest can be grown by calling GrowForest.
32// In advance SetupGrow must be called in order to initialize arrays and
33// do some preprocessing.
34// GrowForest() provides the training data for a single tree (bootstrap
35// aggregate procedure)
36//
37// Essentially two random elements serve to provide a "random" forest,
38// namely bootstrap aggregating (which is done in GrowForest()) and random
39// split selection (which is subject to MRanTree::GrowTree())
40//
41/////////////////////////////////////////////////////////////////////////////
42#include "MRanForest.h"
43
44#include <TVector.h>
45#include <TRandom.h>
46
47#include "MHMatrix.h"
48#include "MRanTree.h"
49#include "MData.h"
50#include "MDataArray.h"
51#include "MParList.h"
52
53#include "MLog.h"
54#include "MLogManip.h"
55
56ClassImp(MRanForest);
57
58using namespace std;
59
60// --------------------------------------------------------------------------
61//
62// Default constructor.
63//
64MRanForest::MRanForest(const char *name, const char *title)
65 : fClassify(kTRUE), fNumTrees(100), fNumTry(0), fNdSize(1),
66 fRanTree(NULL), fRules(NULL), fMatrix(NULL), fUserVal(-1)
67{
68 fName = name ? name : "MRanForest";
69 fTitle = title ? title : "Storage container for Random Forest";
70
71 fForest=new TObjArray();
72 fForest->SetOwner(kTRUE);
73}
74
75MRanForest::MRanForest(const MRanForest &rf)
76{
77 // Copy constructor
78 fName = rf.fName;
79 fTitle = rf.fTitle;
80
81 fClassify = rf.fClassify;
82 fNumTrees = rf.fNumTrees;
83 fNumTry = rf.fNumTry;
84 fNdSize = rf.fNdSize;
85 fTreeNo = rf.fTreeNo;
86 fRanTree = NULL;
87
88 fRules=new MDataArray();
89 fRules->Reset();
90
91 MDataArray *newrules=rf.fRules;
92
93 for(Int_t i=0;i<newrules->GetNumEntries();i++)
94 {
95 MData &data=(*newrules)[i];
96 fRules->AddEntry(data.GetRule());
97 }
98
99 // trees
100 fForest=new TObjArray();
101 fForest->SetOwner(kTRUE);
102
103 TObjArray *newforest=rf.fForest;
104 for(Int_t i=0;i<newforest->GetEntries();i++)
105 {
106 MRanTree *rantree=(MRanTree*)newforest->At(i);
107
108 MRanTree *newtree=new MRanTree(*rantree);
109 fForest->Add(newtree);
110 }
111
112 fHadTrue = rf.fHadTrue;
113 fHadEst = rf.fHadEst;
114 fDataSort = rf.fDataSort;
115 fDataRang = rf.fDataRang;
116 fClassPop = rf.fClassPop;
117 fWeight = rf.fWeight;
118 fTreeHad = rf.fTreeHad;
119
120 fNTimesOutBag = rf.fNTimesOutBag;
121
122}
123
124// --------------------------------------------------------------------------
125// Destructor.
126MRanForest::~MRanForest()
127{
128 delete fForest;
129 if (fMatrix)
130 delete fMatrix;
131 if (fRules)
132 delete fRules;
133}
134
135void MRanForest::Print(Option_t *o) const
136{
137 *fLog << inf << GetDescriptor() << ": " << endl;
138 MRanTree *t = GetTree(0);
139 if (t)
140 {
141 *fLog << "Setting up RF for training on target:" << endl;
142 *fLog << " " << t->GetTitle() << endl;
143 }
144 if (fRules)
145 {
146 *fLog << "Following rules are used as input to RF:" << endl;
147 for (Int_t i=0;i<fRules->GetNumEntries();i++)
148 *fLog << " " << i << ") " << (*fRules)[i].GetRule() << endl;
149 }
150 *fLog << "Random forest parameters:" << endl;
151 if (t)
152 {
153 *fLog << " - " << (t->IsClassify()?"classification":"regression") << " tree" << endl;
154 *fLog << " - Number of trys: " << t->GetNumTry() << endl;
155 *fLog << " - Node size: " << t->GetNdSize() << endl;
156 }
157 *fLog << " - Number of trees: " << fNumTrees << endl;
158 *fLog << " - User value: " << fUserVal << endl;
159 *fLog << endl;
160}
161
162void MRanForest::SetNumTrees(Int_t n)
163{
164 //at least 1 tree
165 fNumTrees=TMath::Max(n,1);
166}
167
168void MRanForest::SetNumTry(Int_t n)
169{
170 fNumTry=TMath::Max(n,0);
171}
172
173void MRanForest::SetNdSize(Int_t n)
174{
175 fNdSize=TMath::Max(n,1);
176}
177
178void MRanForest::SetWeights(const TArrayF &weights)
179{
180 fWeight=weights;
181}
182
183void MRanForest::SetGrid(const TArrayD &grid)
184{
185 const int n=grid.GetSize();
186
187 for(int i=0;i<n-1;i++)
188 if(grid[i]>=grid[i+1])
189 {
190 *fLog<<warn<<"Grid points must be in increasing order! Ignoring grid."<<endl;
191 return;
192 }
193
194 fGrid=grid;
195
196 //*fLog<<inf<<"Following "<<n<<" grid points are used:"<<endl;
197 //for(int i=0;i<n;i++)
198 // *fLog<<inf<<" "<<i<<") "<<fGrid[i]<<endl;
199}
200
201MRanTree *MRanForest::GetTree(Int_t i) const
202{
203 return static_cast<MRanTree*>(fForest->UncheckedAt(i));
204}
205
206Int_t MRanForest::GetNumDim() const
207{
208 return fMatrix ? fMatrix->GetNcols() : 0;
209}
210
211Int_t MRanForest::GetNumData() const
212{
213 return fMatrix ? fMatrix->GetNrows() : 0;
214}
215
216Int_t MRanForest::GetNclass() const
217{
218 int maxidx = TMath::LocMax(fClass.GetSize(),fClass.GetArray());
219
220 return int(fClass[maxidx])+1;
221}
222
223void MRanForest::PrepareClasses()
224{
225 const int numdata=fHadTrue.GetSize();
226
227 if(fGrid.GetSize()>0)
228 {
229 // classes given by grid
230 const int ngrid=fGrid.GetSize();
231
232 for(int j=0;j<numdata;j++)
233 {
234 // Array is supposed to be sorted prior to this call.
235 // If match is found, function returns position of element.
236 // If no match found, function gives nearest element smaller
237 // than value.
238 int k=TMath::BinarySearch(ngrid, fGrid.GetArray(), fHadTrue[j]);
239
240 fClass[j] = k;
241 }
242
243 int minidx = TMath::LocMin(fClass.GetSize(),fClass.GetArray());
244 int min = fClass[minidx];
245 if(min!=0) for(int j=0;j<numdata;j++)fClass[j]-=min;
246
247 }else{
248 // classes directly given
249 for (Int_t j=0;j<numdata;j++)
250 fClass[j] = TMath::Nint(fHadTrue[j]);
251 }
252}
253
254Double_t MRanForest::CalcHadroness()
255{
256 TVector event;
257 *fRules >> event;
258
259 return CalcHadroness(event);
260}
261
262Double_t MRanForest::CalcHadroness(const TVector &event)
263{
264 fTreeHad.Set(fNumTrees);
265
266 Double_t hadroness=0;
267 Int_t ntree =0;
268
269 TIter Next(fForest);
270
271 MRanTree *tree;
272 while ((tree=(MRanTree*)Next()))
273 hadroness += (fTreeHad[ntree++]=tree->TreeHad(event));
274
275 return hadroness/ntree;
276}
277
278Bool_t MRanForest::AddTree(MRanTree *rantree=NULL)
279{
280 fRanTree = rantree ? rantree : fRanTree;
281
282 if (!fRanTree) return kFALSE;
283
284 MRanTree *newtree=new MRanTree(*fRanTree);
285 fForest->Add(newtree);
286
287 return kTRUE;
288}
289
290Bool_t MRanForest::SetupGrow(MHMatrix *mat,MParList *plist)
291{
292 //-------------------------------------------------------------------
293 // access matrix, copy last column (target) preliminarily
294 // into fHadTrue
295 if (fMatrix)
296 delete fMatrix;
297 fMatrix = new TMatrix(mat->GetM());
298
299 int dim = fMatrix->GetNcols()-1;
300 int numdata = fMatrix->GetNrows();
301
302 fHadTrue.Set(numdata);
303 fHadTrue.Reset(0);
304
305 for (Int_t j=0;j<numdata;j++)
306 fHadTrue[j] = (*fMatrix)(j,dim);
307
308 // remove last col
309 fMatrix->ResizeTo(numdata,dim);
310
311 //-------------------------------------------------------------------
312 // setup labels for classification/regression
313 fClass.Set(numdata);
314 fClass.Reset(0);
315
316 if (fClassify)
317 PrepareClasses();
318
319 //-------------------------------------------------------------------
320 // allocating and initializing arrays
321 fHadEst.Set(numdata);
322 fHadEst.Reset(0);
323
324 fNTimesOutBag.Set(numdata);
325 fNTimesOutBag.Reset(0);
326
327 fDataSort.Set(dim*numdata);
328 fDataSort.Reset(0);
329
330 fDataRang.Set(dim*numdata);
331 fDataRang.Reset(0);
332
333 if(fWeight.GetSize()!=numdata)
334 {
335 fWeight.Set(numdata);
336 fWeight.Reset(1.);
337 *fLog << inf <<"Setting weights to 1 (no weighting)"<< endl;
338 }
339
340 //-------------------------------------------------------------------
341 // setup rules to be used for classification/regression
342 const MDataArray *allrules=(MDataArray*)mat->GetColumns();
343 if (allrules==NULL)
344 {
345 *fLog << err <<"Rules of matrix == null, exiting"<< endl;
346 return kFALSE;
347 }
348
349 if (allrules->GetNumEntries()!=dim+1)
350 {
351 *fLog << err <<"Rules of matrix " << allrules->GetNumEntries() << " mismatch dimension+1 " << dim+1 << "...exiting."<< endl;
352 return kFALSE;
353 }
354
355 if (fRules)
356 delete fRules;
357
358 fRules = new MDataArray();
359 fRules->Reset();
360
361 const TString target_rule = (*allrules)[dim].GetRule();
362 for (Int_t i=0;i<dim;i++)
363 fRules->AddEntry((*allrules)[i].GetRule());
364
365 *fLog << inf << endl;
366 *fLog << "Setting up RF for training on target:" << endl;
367 *fLog << " " << target_rule.Data() << endl;
368 *fLog << "Following rules are used as input to RF:" << endl;
369 for (Int_t i=0;i<dim;i++)
370 *fLog << " " << i << ") " << (*fRules)[i].GetRule() << endl;
371 *fLog << endl;
372
373 //-------------------------------------------------------------------
374 // prepare (sort) data for fast optimization algorithm
375 if (!CreateDataSort())
376 return kFALSE;
377
378 //-------------------------------------------------------------------
379 // access and init tree container
380 fRanTree = (MRanTree*)plist->FindCreateObj("MRanTree");
381 if(!fRanTree)
382 {
383 *fLog << err << dbginf << "MRanForest, fRanTree not initialized... aborting." << endl;
384 return kFALSE;
385 }
386 //fRanTree->SetName(target_rule); // Is not stored anyhow
387
388 const Int_t tryest = TMath::Nint(TMath::Sqrt(dim));
389
390 *fLog << inf << endl;
391 *fLog << "Following input for the tree growing are used:"<<endl;
392 *fLog << " Number of Trees : "<<fNumTrees<<endl;
393 *fLog << " Number of Trials: "<<(fNumTry==0?tryest:fNumTry)<<(fNumTry==0?" (auto)":"")<<endl;
394 *fLog << " Final Node size : "<<fNdSize<<endl;
395 *fLog << " Using Grid : "<<(fGrid.GetSize()>0?"Yes":"No")<<endl;
396 *fLog << " Number of Events: "<<numdata<<endl;
397 *fLog << " Number of Params: "<<dim<<endl;
398
399 if(fNumTry==0)
400 {
401 fNumTry=tryest;
402 *fLog << inf << endl;
403 *fLog << "Set no. of trials to the recommended value of round(";
404 *fLog << TMath::Sqrt(dim) << ") = " << fNumTry << endl;
405 }
406
407 fRanTree->SetNumTry(fNumTry);
408 fRanTree->SetClassify(fClassify);
409 fRanTree->SetNdSize(fNdSize);
410
411 fTreeNo=0;
412
413 return kTRUE;
414}
415
416Bool_t MRanForest::GrowForest()
417{
418 if(!gRandom)
419 {
420 *fLog << err << dbginf << "gRandom not initialized... aborting." << endl;
421 return kFALSE;
422 }
423
424 fTreeNo++;
425
426 //-------------------------------------------------------------------
427 // initialize running output
428
429 float minfloat=fHadTrue[TMath::LocMin(fHadTrue.GetSize(),fHadTrue.GetArray())];
430 Bool_t calcResolution=(minfloat>0.001);
431
432 if (fTreeNo==1)
433 {
434 *fLog << inf << endl << underline;
435
436 if(calcResolution)
437 *fLog << "no. of tree no. of nodes resolution in % (from oob-data -> overest. of error)" << endl;
438 else
439 *fLog << "no. of tree no. of nodes rms in % (from oob-data -> overest. of error)" << endl;
440 // 12345678901234567890123456789012345678901234567890
441 }
442
443 const Int_t numdata = GetNumData();
444 const Int_t nclass = GetNclass();
445
446 //-------------------------------------------------------------------
447 // bootstrap aggregating (bagging) -> sampling with replacement:
448
449 TArrayF classpopw(nclass);
450 TArrayI jinbag(numdata); // Initialization includes filling with 0
451 TArrayF winbag(numdata); // Initialization includes filling with 0
452
453 float square=0;
454 float mean=0;
455
456 for (Int_t n=0; n<numdata; n++)
457 {
458 // The integer k is randomly (uniformly) chosen from the set
459 // {0,1,...,numdata-1}, which is the set of the index numbers of
460 // all events in the training sample
461
462 const Int_t k = Int_t(gRandom->Rndm()*numdata);
463
464 if(fClassify)
465 classpopw[fClass[k]]+=fWeight[k];
466 else
467 classpopw[0]+=fWeight[k];
468
469 mean +=fHadTrue[k]*fWeight[k];
470 square+=fHadTrue[k]*fHadTrue[k]*fWeight[k];
471
472 winbag[k]+=fWeight[k];
473 jinbag[k]=1;
474
475 }
476
477 //-------------------------------------------------------------------
478 // modifying sorted-data array for in-bag data:
479
480 // In bagging procedure ca. 2/3 of all elements in the original
481 // training sample are used to build the in-bag data
482 TArrayI datsortinbag=fDataSort;
483 Int_t ninbag=0;
484
485 ModifyDataSort(datsortinbag, ninbag, jinbag);
486
487 fRanTree->GrowTree(fMatrix,fHadTrue,fClass,datsortinbag,fDataRang,classpopw,mean,square,
488 jinbag,winbag,nclass);
489
490 //-------------------------------------------------------------------
491 // error-estimates from out-of-bag data (oob data):
492 //
493 // For a single tree the events not(!) contained in the bootstrap sample of
494 // this tree can be used to obtain estimates for the classification error of
495 // this tree.
496 // If you take a certain event, it is contained in the oob-data of 1/3 of
497 // the trees (see comment to ModifyData). This means that the classification error
498 // determined from oob-data is underestimated, but can still be taken as upper limit.
499
500 for (Int_t ievt=0;ievt<numdata;ievt++)
501 {
502 if (jinbag[ievt]>0)
503 continue;
504
505 fHadEst[ievt] +=fRanTree->TreeHad((*fMatrix), ievt);
506 fNTimesOutBag[ievt]++;
507
508 }
509
510 Int_t n=0;
511 Float_t ferr=0;
512
513 for (Int_t ievt=0;ievt<numdata;ievt++)
514 {
515 if(fNTimesOutBag[ievt]!=0)
516 {
517 float val = fHadEst[ievt]/float(fNTimesOutBag[ievt])-fHadTrue[ievt];
518 if(calcResolution) val/=fHadTrue[ievt];
519
520 ferr += val*val;
521 n++;
522 }
523 }
524 ferr = TMath::Sqrt(ferr/n);
525
526 //-------------------------------------------------------------------
527 // give running output
528 *fLog << setw(5) << fTreeNo;
529 *fLog << setw(18) << fRanTree->GetNumEndNodes();
530 *fLog << Form("%18.2f", ferr*100.);
531 *fLog << endl;
532
533 fRanTree->SetError(ferr);
534
535 // adding tree to forest
536 AddTree();
537
538 return fTreeNo<fNumTrees;
539}
540
541Bool_t MRanForest::CreateDataSort()
542{
543 // fDataSort(m,n) is the event number in which fMatrix(m,n) occurs.
544 // fDataRang(m,n) is the rang of fMatrix(m,n), i.e. if rang = r:
545 // fMatrix(m,n) is the r-th highest value of all fMatrix(m,.).
546 //
547 // There may be more then 1 event with rang r (due to bagging).
548
549 const Int_t numdata = GetNumData();
550 const Int_t dim = GetNumDim();
551
552 TArrayF v(numdata);
553 TArrayI isort(numdata);
554
555
556 for (Int_t mvar=0;mvar<dim;mvar++)
557 {
558
559 for(Int_t n=0;n<numdata;n++)
560 {
561 v[n]=(*fMatrix)(n,mvar);
562 isort[n]=n;
563
564 if(TMath::IsNaN(v[n]))
565 {
566 *fLog << err <<"Event no. "<<n<<", matrix column no. "<<mvar;
567 *fLog << err <<" has the value NaN."<<endl;
568 return kFALSE;
569 }
570 }
571
572 TMath::Sort(numdata,v.GetArray(),isort.GetArray(),kFALSE);
573
574 // this sorts the v[n] in ascending order. isort[n] is the event number
575 // of that v[n], which is the n-th from the lowest (assume the original
576 // event numbers are 0,1,...).
577
578 // control sorting
579 for(int n=1;n<numdata;n++)
580 if(v[isort[n-1]]>v[isort[n]])
581 {
582 *fLog << err <<"Event no. "<<n<<", matrix column no. "<<mvar;
583 *fLog << err <<" not at correct sorting position."<<endl;
584 return kFALSE;
585 }
586
587 for(Int_t n=0;n<numdata-1;n++)
588 {
589 const Int_t n1=isort[n];
590 const Int_t n2=isort[n+1];
591
592 fDataSort[mvar*numdata+n]=n1;
593 if(n==0) fDataRang[mvar*numdata+n1]=0;
594 if(v[n1]<v[n2])
595 {
596 fDataRang[mvar*numdata+n2]=fDataRang[mvar*numdata+n1]+1;
597 }else{
598 fDataRang[mvar*numdata+n2]=fDataRang[mvar*numdata+n1];
599 }
600 }
601 fDataSort[(mvar+1)*numdata-1]=isort[numdata-1];
602 }
603 return kTRUE;
604}
605
606void MRanForest::ModifyDataSort(TArrayI &datsortinbag, Int_t ninbag, const TArrayI &jinbag)
607{
608 const Int_t numdim=GetNumDim();
609 const Int_t numdata=GetNumData();
610
611 ninbag=0;
612 for (Int_t n=0;n<numdata;n++)
613 if(jinbag[n]==1) ninbag++;
614
615 for(Int_t m=0;m<numdim;m++)
616 {
617 Int_t k=0;
618 Int_t nt=0;
619 for(Int_t n=0;n<numdata;n++)
620 {
621 if(jinbag[datsortinbag[m*numdata+k]]==1)
622 {
623 datsortinbag[m*numdata+nt]=datsortinbag[m*numdata+k];
624 k++;
625 }else{
626 for(Int_t j=1;j<numdata-k;j++)
627 {
628 if(jinbag[datsortinbag[m*numdata+k+j]]==1)
629 {
630 datsortinbag[m*numdata+nt]=datsortinbag[m*numdata+k+j];
631 k+=j+1;
632 break;
633 }
634 }
635 }
636 nt++;
637 if(nt>=ninbag) break;
638 }
639 }
640}
641
642Bool_t MRanForest::AsciiWrite(ostream &out) const
643{
644 Int_t n=0;
645 MRanTree *tree;
646 TIter forest(fForest);
647
648 while ((tree=(MRanTree*)forest.Next()))
649 {
650 tree->AsciiWrite(out);
651 n++;
652 }
653
654 return n==fNumTrees;
655}
Note: See TracBrowser for help on using the repository browser.