source: trunk/MagicSoft/Mars/mhist/MHRanForest.cc@ 1880

Last change on this file since 1880 was 1880, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 5.7 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@alwa02.physik.uni-siegen.de>
19!
20! Copyright: MAGIC Software Development, 2000-2003
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHRanForest
28//
29// This histogram shows the evolution of the standard deviation
30// of est. hadronness as the number of trees increases.
31// It helps you to find out how many trees are really needed in g/h-sep.
32//
33////////////////////////////////////////////////////////////////////////////
34#include "MHRanForest.h"
35
36#include <TPad.h>
37#include <TGraph.h>
38#include <TStyle.h>
39#include <TCanvas.h>
40#include <TMarker.h>
41
42#include "MParList.h"
43#include "MBinning.h"
44#include "MRanForest.h"
45
46#include "MLog.h"
47#include "MLogManip.h"
48
49#include "MMcEvt.hxx"
50
51ClassImp(MHRanForest);
52
53// --------------------------------------------------------------------------
54//
55// Setup histograms, nbins is the number of bins used for the evaluation.
56// The default is 100 bins.
57//
58MHRanForest::MHRanForest(Int_t nbins, const char *name, const char *title)
59{
60 //
61 // set the name and title of this object
62 //
63 fName = name ? name : "MHRanForest";
64 fTitle = title ? title : "Histogram showing Standard deviation of estimated hadronness";
65
66 fGraphSigma = new TGraph;
67 fGraphSigma->SetTitle("Evolution of Standard deviation of estimated hadronness in tree combination");
68 fGraphSigma->SetMaximum(1);
69}
70
71// --------------------------------------------------------------------------
72//
73// Delete the histograms.
74//
75MHRanForest::~MHRanForest()
76{
77 delete fGraphSigma;
78}
79
80// --------------------------------------------------------------------------
81//
82// Setup Filling of the histograms. It needs:
83// MMcEvt and MRanForest
84//
85Bool_t MHRanForest::SetupFill(const MParList *plist)
86{
87 fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
88 if (!fMcEvt)
89 {
90 *fLog << err << dbginf << "MMcEvt not found... aborting." << endl;
91 return kFALSE;
92 }
93
94 fRanForest = (MRanForest*)plist->FindObject("MRanForest");
95 if (!fRanForest)
96 {
97 *fLog << err << dbginf << "MRanForest not found... aborting." << endl;
98 return kFALSE;
99 }
100
101 fSigma.Set(fRanForest->GetNumTrees());
102 fNumEvent=0;
103 return kTRUE;
104}
105
106// --------------------------------------------------------------------------
107//
108//
109Bool_t MHRanForest::Fill(const MParContainer *par)
110{
111 fNumEvent++;
112 Double_t hest=0;
113 Double_t htrue=fMcEvt->GetPartId()==kGAMMA ? 0. : 1.;
114
115 Int_t ntrees=fRanForest->GetNumTrees();
116
117 for (Int_t i=0;i<ntrees;i++)
118 {
119 for(Int_t j=0;j<=i;j++)
120 hest+=fRanForest->GetTreeHad(j);
121
122 hest/=i+1;
123 fSigma[i]+=(htrue-hest)*(htrue-hest);
124 }
125
126 return kTRUE;
127}
128
129// --------------------------------------------------------------------------
130//
131// Finalize the histogram:
132// calculate standard deviation and set histogram max and min
133//
134Bool_t MHRanForest::Finalize()
135{
136 Int_t n = fSigma.GetSize();
137
138 fGraphSigma->Set(n);
139
140 Stat_t max=0.;
141 Stat_t min=0.;
142 for (Int_t i=0; i<n; i++)
143 {
144 Stat_t ip = i+1.;
145 fSigma[i] = TMath::Sqrt(fSigma[i]/Stat_t(fNumEvent));
146 Stat_t ig = fSigma[i];
147 max=TMath::Max(max,ig);
148 min=TMath::Min(min,ig);
149 fGraphSigma->SetPoint(i,ip,ig);
150 }
151 fGraphSigma->SetMaximum(1.05*max);
152 fGraphSigma->SetMinimum(0.95*min);
153
154 return kTRUE;
155}
156
157// --------------------------------------------------------------------------
158//
159// Draw clone of histogram
160//
161TObject *MHRanForest::DrawClone(Option_t *opt) const
162{
163 if (fGraphSigma->GetN()==0)
164 return NULL;
165
166 TCanvas &c = *MakeDefCanvas("RanForest", fTitle);
167 gROOT->SetSelectedPad(NULL);
168
169 gStyle->SetOptStat(10);
170 TGraph &g = (TGraph&)*fGraphSigma->DrawClone("AL");
171 g.SetBit(kCanDelete);
172 gPad->Modified();
173 gPad->Update();
174 if (g.GetHistogram())
175 {
176 g.GetXaxis()->SetRangeUser(0, fNumEvent);
177 g.GetXaxis()->SetTitle("Number of Trees");
178 g.GetYaxis()->SetTitle("Standard deviation of estimated hadronness");
179 g.SetMarkerStyle(kFullDotMedium);
180 gPad->Modified();
181 gPad->Update();
182 //g.Draw("P");
183 }
184 gPad->SetGrid();
185
186 return &c;
187}
188
189// --------------------------------------------------------------------------
190//
191// Draw histogram. (For the Meaning see class description)
192//
193void MHRanForest::Draw(Option_t *)
194{
195 if (fGraphSigma->GetN()==0)
196 return;
197
198 if (!gPad)
199 MakeDefCanvas("RanForest", fTitle);
200
201 gStyle->SetOptStat(10);
202 fGraphSigma->Draw("ALP");
203 gPad->Modified();
204 gPad->Update();
205 if (fGraphSigma->GetHistogram())
206 {
207 fGraphSigma->GetXaxis()->SetRangeUser(0, 1);
208 fGraphSigma->GetXaxis()->SetTitle("Number of Trees");
209 fGraphSigma->GetYaxis()->SetTitle("Standard deviation of estimated hadronness");
210
211 fGraphSigma->SetMarkerStyle(kFullDotSmall);
212 //fGraphSigma->Draw("P");
213 gPad->Modified();
214 gPad->Update();
215 }
216}
Note: See TracBrowser for help on using the repository browser.