source: trunk/MagicSoft/Mars/mranforest/MRFEnergyEst.cc@ 6979

Last change on this file since 6979 was 6932, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 11.6 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 2/2005 <mailto:hengsteb@physik.hu-berlin.de>
19!
20! Copyright: MAGIC Software Development, 2000-2005
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MRFEnergyEst
28//
29//
30////////////////////////////////////////////////////////////////////////////
31#include "MRFEnergyEst.h"
32
33#include <TFile.h>
34#include <TList.h>
35
36#include <TH1.h>
37#include <TH2.h>
38#include <TStyle.h>
39#include <TCanvas.h>
40#include <TMath.h>
41#include <TVector.h>
42
43#include "MHMatrix.h"
44
45#include "MLog.h"
46#include "MLogManip.h"
47
48#include "MParList.h"
49#include "MTaskList.h"
50#include "MEvtLoop.h"
51
52#include "MRanTree.h"
53#include "MRanForest.h"
54#include "MRanForestGrow.h"
55
56#include "MData.h"
57#include "MParameters.h"
58
59
60ClassImp(MRFEnergyEst);
61
62using namespace std;
63
64static const TString gsDefName = "MRFEnergyEst";
65static const TString gsDefTitle = "RF for energy estimation";
66
67// --------------------------------------------------------------------------
68//
69//
70MRFEnergyEst::MRFEnergyEst(const char *name, const char *title):fNumTrees(-1)
71{
72 //
73 // set the name and title of this object
74 //
75 fName = name ? name : gsDefName.Data();
76 fTitle = title ? title : gsDefTitle.Data();
77}
78
79// --------------------------------------------------------------------------
80//
81// Delete the data chains
82//
83MRFEnergyEst::~MRFEnergyEst()
84{
85 // delete fData;
86}
87
88// --------------------------------------------------------------------------
89Int_t MRFEnergyEst::Train()
90{
91 if(!fMatrixTrain)
92 {
93 *fLog << err << dbginf << "fMatrixTrain not found... aborting." << endl;
94 return kFALSE;
95 }
96
97 if(!fMatrixTrain->GetColumns())
98 {
99 *fLog << err << dbginf << "fMatrixTrain does not contain rules... aborting." << endl;
100 return kFALSE;
101 }
102
103 const Int_t ncols = (fMatrixTrain->GetM()).GetNcols();
104 const Int_t nrows = (fMatrixTrain->GetM()).GetNrows();
105
106 cout<<"ncols="<<ncols<<" nrows="<<nrows<<endl<<flush;
107
108 if(ncols<=0 || nrows <=0)
109 {
110 *fLog << err << dbginf << "No. of columns or no. of rows of fMatrixTrain equal 0 ... aborting." << endl;
111 return kFALSE;
112 }
113
114 // rules (= combination of image par) to be used for energy estimation
115 MDataArray used_rules;
116 TString energy_rule;
117 for(Int_t i=0;i<ncols;i++)
118 {
119 MDataArray *rules=fMatrixTrain->GetColumns();
120 MData &data=(*rules)[i];
121
122 if(i<ncols-1)
123 used_rules.AddEntry(data.GetRule());
124 else
125 energy_rule=data.GetRule();
126 }
127
128 if(!energy_rule.Contains("fEnergy"))
129 {
130 *fLog << err << dbginf << "Can not accept "<<energy_rule<<" as true energy ... aborting." << endl;
131 return kFALSE;
132 }
133
134 TFile fileRF(fRFfileName,"recreate");
135 if(!fileRF.IsOpen())
136 {
137 *fLog << err << dbginf << "File to store RFs could not be opened... aborting." << endl;
138 return kFALSE;
139 }
140
141 TMatrix *mptr=(TMatrix*)&(fMatrixTrain->GetM());
142 const Int_t nbins = fEnergyGrid.GetSize()-1;
143 if(nbins<=0)
144 {
145 *fLog << err << dbginf << "Energy grid not vaild... aborting." << endl;
146 return kFALSE;
147 }
148
149 // training of RF for each energy bin
150 Int_t numbins=0;
151 for(Int_t ie=0;ie<nbins;ie++)
152 {
153 TMatrix mat1(nrows,ncols);
154 TMatrix mat0(nrows,ncols);
155
156 // prepare matrix for current energy bin
157 Int_t irow1=0;
158 Int_t irow0=0;
159
160 for(Int_t j=0;j<nrows;j++)
161 {
162 Double_t energy=(*mptr)(j,ncols-1);
163 if(energy>pow(10.,fEnergyGrid[ie]) && energy<=pow(10.,fEnergyGrid[ie+1]))
164 {
165 TMatrixRow(mat1, irow1) = TMatrixRow(*mptr,j);
166 irow1++;
167 }else{
168 TMatrixRow(mat0, irow0) = TMatrixRow(*mptr,j);
169 irow0++;
170 }
171 }
172 if(irow1*irow0==0)continue;
173
174 *fLog << inf << dbginf << "Training RF for energy bin "<<ie<< endl;
175
176 // last column must be removed (true energy col.)
177 mat1.ResizeTo(irow1,ncols-1);
178 mat0.ResizeTo(irow0,ncols-1);
179
180 // create MHMatrix as input for RF
181 MHMatrix matrix1(mat1,"MatrixHadrons");
182 MHMatrix matrix0(mat0,"MatrixGammas");
183
184 MDataArray *rules1=matrix1.GetColumns();
185 MDataArray *rules0=matrix0.GetColumns();
186 // rules of new matrices be re-set
187 if(rules1)delete rules1; rules1=new MDataArray();
188 if(rules0)delete rules0; rules0=new MDataArray();
189
190 for(Int_t i=0;i<ncols-1;i++)
191 {
192 //MDataArray *rules=fMatrixTrain->GetColumns();
193 //MData &data=(*rules)[i];
194 MData &data=used_rules[i];
195 rules1->AddEntry(data.GetRule());
196 rules0->AddEntry(data.GetRule());
197 }
198
199 // training of RF
200 MParList plist;
201 MTaskList tlist;
202 plist.AddToList(&tlist);
203 plist.AddToList(&matrix0);
204 plist.AddToList(&matrix1);
205
206 MRanForestGrow rfgrow;
207 rfgrow.SetNumTrees(fNumTrees); // number of trees
208 rfgrow.SetNumTry(fNumTry); // number of trials in random split selection
209 rfgrow.SetNdSize(fNdSize); // limit for nodesize
210
211 tlist.AddToList(&rfgrow);
212
213 MEvtLoop evtloop;
214 evtloop.SetParList(&plist);
215
216 //gLog.SetNullOutput(kTRUE);
217
218 if (!evtloop.Eventloop()) return kFALSE;
219
220 //gLog.SetNullOutput(kFALSE);
221
222 // save whole forest
223 MRanForest *forest=(MRanForest*)plist.FindObject("MRanForest");
224 forest->SetTitle(Form("%f",0.5*(fEnergyGrid[ie]+fEnergyGrid[ie+1])));
225 forest->Write(Form("%d",numbins));
226 numbins++;
227 }
228
229 // save rules
230 used_rules.Write("rules");
231
232 fileRF.Close();
233
234 return kTRUE;
235}
236
237Int_t MRFEnergyEst::Test()
238{
239 if(!fMatrixTest)
240 {
241 *fLog << err << dbginf << "fMatrixTest not found... aborting." << endl;
242 return kFALSE;
243 }
244
245 const Int_t ncols = (fMatrixTest->GetM()).GetNcols();
246 const Int_t nrows = (fMatrixTest->GetM()).GetNrows();
247
248 if(ncols<=0 || nrows <=0)
249 {
250 *fLog << err << dbginf << "No. of columns or no. of rows of fMatrixTrain equal 0 ... aborting." << endl;
251 return kFALSE;
252 }
253
254 TMatrix *mptr=(TMatrix*)&(fMatrixTest->GetM());
255
256 if(!ReadForests())
257 {
258 *fLog << err << dbginf << "Reading RFs failed... aborting." << endl;
259 return kFALSE;
260 }
261
262 const Int_t nbins=fEForests.GetSize();
263
264 Double_t e_low = 100;
265 Double_t e_up = 0;
266
267 for(Int_t j=0;j<nbins;j++)
268 {
269 e_low = TMath::Min(atof((fEForests[j])->GetTitle()),e_low);
270 e_up = TMath::Max(atof((fEForests[j])->GetTitle()),e_up);
271 }
272
273 TH1F hres("hres","",1000,-10,10);
274 TH2F henergy("henergy","",100,e_low,e_up,100,e_low,e_up);
275
276 for(Int_t i=0;i<nrows;i++)
277 {
278 Double_t e_true = (*mptr)(i,ncols-1);
279 Double_t e_est = 0;
280 //Double_t hmax = 0;
281 Double_t hsum = 0;
282
283 for(Int_t j=0;j<nbins;j++)
284 {
285 const TVector v=TMatrixRow(*mptr,i);
286 Double_t h = ((MRanForest*) (fEForests[j]))->CalcHadroness(v);
287 Double_t e = atof((fEForests[j])->GetTitle());
288 /*if(h>=hmax)
289 {
290 hmax=h;
291 e_est=pow(10.,e);
292 }*/
293 hsum+=h;
294 e_est+=h*e;
295 }
296 e_est/=hsum;
297 e_est=pow(10.,e_est);
298
299 if(e_true>80.) hres.Fill((e_est-e_true)/e_true);
300 henergy.Fill(log10(e_true),log10(e_est));
301 }
302
303 if(TestBit(kEnableGraphicalOutput))
304 {
305 if(gStyle) gStyle->SetOptFit(1);
306 // show results
307 TCanvas *c=MH::MakeDefCanvas("c","",300,900);
308
309 c->Divide(1,3);
310 c->cd(1);
311 henergy.SetTitle("Estimated vs true energy");
312 henergy.GetXaxis()->SetTitle("log10(E_{true}[GeV])");
313 henergy.GetYaxis()->SetTitle("log10(E_{est}[GeV])");
314 henergy.DrawCopy();
315
316 c->cd(2);
317
318 TH1F *hptr=(TH1F*)henergy.ProfileX("_px",1,100,"S");
319 hptr->SetTitle("Estimated vs true energy - profile");
320 hptr->GetXaxis()->SetTitle("log10(E_{true}[GeV])");
321 hptr->GetYaxis()->SetTitle("log10(E_{est}[GeV])");
322 hptr->DrawCopy();
323
324 c->cd(3);
325 hres.Fit("gaus");
326 hres.SetTitle("Energy resolution for E_{true}>80Gev");
327 hres.GetXaxis()->SetTitle("(E_{estimated}-E_{true})/E_{true})");
328 hres.GetYaxis()->SetTitle("counts");
329 hres.DrawCopy();
330
331
332 c->GetPad(1)->SetGrid();
333 c->GetPad(2)->SetGrid();
334 c->GetPad(3)->SetGrid();
335
336 }
337
338 return kTRUE;
339}
340
341Int_t MRFEnergyEst::ReadForests(MParList *plist)
342{
343 TFile fileRF(fRFfileName,"read");
344 if(!fileRF.IsOpen())
345 {
346 *fLog << err << dbginf << "File containing RFs could not be opened... aborting." << endl;
347 return kFALSE;
348 }
349
350 TList *list=(TList*)fileRF.GetListOfKeys();
351 const Int_t n=list->GetSize()-1;// subtract 1 since 1 key belongs to MDataArray
352
353 fEForests.Expand(n);
354
355 MRanForest forest;
356 for(Int_t i=0;i<n;i++)
357 {
358 forest.Read(Form("%d",i));
359 MRanForest *curforest=(MRanForest*)forest.Clone();
360 const char *energy=list->FindObject(Form("%d",i))->GetTitle();
361 const char *bin =list->FindObject(Form("%d",i))->GetName();
362
363 curforest->SetTitle(energy);
364 curforest->SetName(bin);
365
366 fEForests[i]=curforest;
367 }
368
369 if(plist)
370 {
371 fData=(MDataArray*)plist->FindCreateObj("MDataArray");
372 fData->Read("rules");
373 }
374 fileRF.Close();
375
376 return kTRUE;
377}
378
379Int_t MRFEnergyEst::PreProcess(MParList *plist)
380{
381 fEnergyEst = (MParameterD*)plist->FindCreateObj("MParameterD", "MEnergyEst");
382 if (!fEnergyEst)
383 {
384 *fLog << err << dbginf << "MEnergyEst [MParameterD] not found... aborting." << endl;
385 return kFALSE;
386 }
387
388 if(!ReadForests(plist))
389 {
390 *fLog << err << dbginf << "Reading RFs failed... aborting." << endl;
391 return kFALSE;
392 }
393
394 if (!fData)
395 {
396 *fLog << err << dbginf << "MDataArray not found... aborting." << endl;
397 return kFALSE;
398 }
399
400 if (!fData->PreProcess(plist))
401 {
402 *fLog << err << dbginf << "PreProcessing of the MDataArray failed for the columns failed... aborting." << endl;
403 return kFALSE;
404 }
405
406 return kTRUE;
407}
408
409// --------------------------------------------------------------------------
410//
411//
412Int_t MRFEnergyEst::Process()
413{
414 TVector event;
415 *fData >> event;
416
417 Double_t eest = 0;
418 //Double_t hmax = 0;
419 Double_t hsum = 0;
420
421 for(Int_t j=0;j<fEForests.GetSize();j++)
422 {
423 Double_t h = ((MRanForest*) (fEForests[j]))->CalcHadroness(event);
424 Double_t e = atof((fEForests[j])->GetTitle());
425 /*if(h>=hmax)
426 {
427 hmax=h;
428 e_est=pow(10.,e);
429 }*/
430 hsum+=h;
431 eest+=h*e;
432 }
433 eest/=hsum;
434 eest=pow(10.,eest);
435
436 fEnergyEst->SetVal(eest);
437 fEnergyEst->SetReadyToSave();
438
439 return kTRUE;
440}
Note: See TracBrowser for help on using the repository browser.