source: trunk/MagicSoft/Mars/mhistmc/MHMcTriggerLvl2.cc@ 2173

Last change on this file since 2173 was 2173, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 9.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): Nicola Galante, 2003 <mailto:nicola.galante@pi.infn.it>
19!
20! Copyright: MAGIC Software Development, 2000-2003
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHMcTriggerLvl2
28//
29// This class contains histograms for the Trigger Level2 image parameters
30//
31/////////////////////////////////////////////////////////////////////////////
32
33#include "MHMcTriggerLvl2.h"
34
35#include <math.h>
36
37#include <TH1.h>
38#include <TF1.h>
39#include <TPad.h>
40#include <TStyle.h>
41#include <TCanvas.h>
42#include <TPaveLabel.h>
43
44#include "MLog.h"
45#include "MLogManip.h"
46
47#include "MParList.h"
48
49#include "MMcTriggerLvl2.h"
50#include "MGeomCam.h"
51#include "MBinning.h"
52
53using namespace std;
54
55/*
56 Please, DON'T USE IFDEFS IN SUCH A CONTEXT, Thomas.
57 --------------------------------------------------
58#ifndef COLOR_LINELPS
59#define COLOR_LINELPS Int_t colorlps = 1
60COLOR_LINELPS;
61#endif
62
63#ifndef COLOR_LINESBC
64#define COLOR_LINESBC Int_t colorsbc = 1
65COLOR_LINESBC;
66#endif
67
68#ifndef COLOR_LINEPS
69#define COLOR_LINEPS Int_t colorps = 1
70COLOR_LINEPS;
71#endif
72*/
73
74/* Use this insteadif you want to have some value which is the same for all
75 your instances of MHMcTriggerLvl2, if not define the values in the constructor
76 and remove the 'static'
77 */
78
79Int_t MHMcTriggerLvl2::fColorLps = 1;
80Int_t MHMcTriggerLvl2::fColorSbc = 1;
81Int_t MHMcTriggerLvl2::fColorPs = 1;
82
83ClassImp(MHMcTriggerLvl2);
84
85// --------------------------------------------------------------------------
86//
87// Setup three histograms for fLutPseudoSize, fPseudoSize, fSizeBiggerCell
88//
89MHMcTriggerLvl2::MHMcTriggerLvl2(const char *name, const char *title)
90{
91 //
92 // set the name and title of this object
93 //
94 fName = name ? name : "MHMcTriggerLvl2";
95 fTitle = title ? title : "Trigger L2 image parameters";
96
97 fHistLutPseudoSize = new TH1F("fHistLutPseudoSize", "number of compact pixels in one lut", 13, 0, 12);
98 fHistPseudoSize = new TH1F("fHistPseudoSize", "Multiplicity of the cluster identified by the L2T", 41, 0, 40);
99 fHistSizeBiggerCell = new TH1F("fHistSizeBiggerCell", "Number of fired pixel in bigger cell", 37, 0, 36);
100
101 fHistLutPseudoSizeNorm = new TH1F("fHistLutPseudoSizeNorm", "Normalized Number of compact pixels in one lut", 13, 0, 12);
102 fHistPseudoSizeNorm = new TH1F("fHistPseudoSizeNorm", "Normalized Multiplicity of the cluster identified by the L2T", 41, 0, 40);
103 fHistSizeBiggerCellNorm = new TH1F("fHistSizeBiggerCellNorm", "Normalized Number of fired pixel in bigger cell", 37, 0, 36);
104
105 fHistLutPseudoSize->SetDirectory(NULL);
106 fHistLutPseudoSizeNorm->SetDirectory(NULL);
107 fHistPseudoSize->SetDirectory(NULL);
108 fHistPseudoSizeNorm->SetDirectory(NULL);
109 fHistSizeBiggerCell->SetDirectory(NULL);
110 fHistSizeBiggerCellNorm->SetDirectory(NULL);
111
112 fHistLutPseudoSize->SetXTitle("Number of Pixels");
113 fHistLutPseudoSizeNorm->SetXTitle("Number of Pixels");
114 fHistPseudoSize->SetXTitle("Number of Pixels");
115 fHistPseudoSizeNorm->SetXTitle("Number of Pixels");
116 fHistSizeBiggerCell->SetXTitle("Number of Pixels");
117 fHistSizeBiggerCellNorm->SetXTitle("Number of Pixels");
118
119 fHistLutPseudoSize->SetYTitle("Counts");
120 fHistLutPseudoSizeNorm->SetYTitle("Counts/Total Counts");
121 fHistPseudoSize->SetYTitle("Counts");
122 fHistPseudoSizeNorm->SetYTitle("Counts/Total Counts");
123 fHistSizeBiggerCell->SetYTitle("Counts");
124 fHistSizeBiggerCellNorm->SetYTitle("Counts/Total Counts");
125
126 fFNorm = new TF1("FNorm", "1", -1, 40);
127}
128
129// --------------------------------------------------------------------------
130//
131// Delete the histograms
132//
133MHMcTriggerLvl2::~MHMcTriggerLvl2()
134{
135 delete fHistLutPseudoSize;
136 delete fHistPseudoSize;
137 delete fHistSizeBiggerCell;
138 delete fHistLutPseudoSizeNorm;
139 delete fHistPseudoSizeNorm;
140 delete fHistSizeBiggerCellNorm;
141
142 delete fFNorm;
143}
144
145// --------------------------------------------------------------------------
146//
147// Fill the histograms with data from a MMcTriggerLvl2-Container.
148// Be careful: Only call this with an object of type MMcTriggerLvl2
149//
150Bool_t MHMcTriggerLvl2::Fill(const MParContainer *par, const Stat_t w)
151{
152 const MMcTriggerLvl2 &h = *(MMcTriggerLvl2 *)par;
153
154 fHistLutPseudoSize->Fill(h.GetLutPseudoSize());
155 fHistPseudoSize->Fill(h.GetPseudoSize());
156 fHistSizeBiggerCell->Fill(h.GetSizeBiggerCell());
157
158 return kTRUE;
159}
160
161
162// --------------------------------------------------------------------------
163//
164// This is the private function member which draw a clone of a histogram.
165// This method is called by the DrawClone method.
166//
167TObject *MHMcTriggerLvl2::DrawHist(TH1 &hist, TH1 &histNorm, const TString &canvasname, Int_t &col) const
168{
169 col++;
170
171 TCanvas *c = (TCanvas*)gROOT->FindObject(canvasname);
172
173 Bool_t same = kTRUE;
174 if (!c)
175 {
176 c = MakeDefCanvas(canvasname,canvasname, 800, 610);
177 c->Divide(2,1);
178 same = kFALSE;
179 }
180
181 c->cd(1);
182 hist.SetLineColor(col);
183 hist.DrawCopy(same?"same":"");
184
185 c->cd(2);
186 histNorm.SetLineColor(col);
187 histNorm.DrawCopy(same?"same":"");
188
189 return c;
190}
191
192
193// --------------------------------------------------------------------------
194//
195// Draw a clone of a data member histogram. So that the object can be deleted
196// and the histogram is still visible in the canvas.
197// The cloned object are deleted together with the canvas if the canvas is
198// destroyed. If you want to handle dostroying the canvas you can get a
199// pointer to it from this function
200// Possible options are:
201// "lps" for the fLutPseudoSize histogram;
202// "sbc" for the fSizeBiggerCell histogram;
203// "ps" for the fPseudoSize histogram;
204// default: "ps"
205//
206TObject *MHMcTriggerLvl2::DrawClone(Option_t *opt) const
207{
208 TString str(opt);
209
210 if (str.IsNull())
211 str = "ps";
212
213 if (!str.Contains("lps", TString::kIgnoreCase) &&
214 !str.Contains("sbc", TString::kIgnoreCase) &&
215 !str.Contains("ps", TString::kIgnoreCase))
216 {
217 *fLog << "ARGH!@! Possible options are \"lps\", \"sbc\", \"ps\" or NULL!" <<endl;
218 return NULL;
219 }
220
221 TH1 *hist=NormalizeHist(fHistLutPseudoSizeNorm, fHistLutPseudoSize);
222
223 if (!hist)
224 return NULL;
225
226 if (str.Contains("lps",TString::kIgnoreCase))
227 return DrawHist(*fHistLutPseudoSize, *hist, "CanvasLPS", fColorLps);
228
229 if (str.Contains("sbc",TString::kIgnoreCase))
230 return DrawHist(*fHistSizeBiggerCell, *hist, "CanvasSBC", fColorSbc);
231
232 if (str.Contains("ps",TString::kIgnoreCase))
233 return DrawHist(*fHistPseudoSize, *hist, "CanvasPS", fColorPs);
234
235 return NULL;
236}
237
238
239
240// --------------------------------------------------------------------------
241//
242// Creates a new canvas and draws the three histograms into it.
243// Be careful: The histograms belongs to this object and won't get deleted
244// together with the canvas.
245//
246void MHMcTriggerLvl2::Draw(Option_t *)
247{
248 MakeDefCanvas("c","L2T Parameters", 720, 810);
249
250 TPad* pad1 = new TPad("pad1","Pad with fLutPseudoSize", 0.003, 0.7, 0.4, 0.997);
251 TPad* pad2 = new TPad("pad2","Pad with fSizeBiggerCell", 0.403, 0.7, 0.997, 0.997);
252 TPad* pad3 = new TPad("pad3","Pad with fPseudoSize", 0.003, 0.003, 0.997, 0.697);
253 pad1->Draw();
254 pad2->Draw();
255 pad3->Draw();
256
257 pad1->cd();
258 fHistLutPseudoSize->Draw();
259
260 pad2->cd();
261 fHistSizeBiggerCell->Draw();
262
263 pad3->cd();
264 fHistPseudoSize->Draw();
265}
266
267
268
269// --------------------------------------------------------------------------
270//
271// Return the histogram by its name
272//
273TH1 *MHMcTriggerLvl2::GetHistByName(const TString name)
274{
275 if (name.Contains("fHistLutPseudoSize", TString::kIgnoreCase))
276 return fHistLutPseudoSize;
277 if (name.Contains("fHistSizeBiggerCell", TString::kIgnoreCase))
278 return fHistSizeBiggerCell;
279 if (name.Contains("fHistPseudoSize", TString::kIgnoreCase))
280 return fHistPseudoSize;
281
282 if (name.Contains("fHistLutPseudoSizeNorm", TString::kIgnoreCase))
283 return fHistLutPseudoSizeNorm;
284 if (name.Contains("fHistSizeBiggerCellNorm", TString::kIgnoreCase))
285 return fHistSizeBiggerCellNorm;
286 if (name.Contains("fHistPseudoSizeNorm", TString::kIgnoreCase))
287 return fHistPseudoSizeNorm;
288
289 return NULL;
290}
291
292
293
294// --------------------------------------------------------------------------
295//
296// Normalize the histogram on its integral, i.e. normalize the
297// values of the histogram on the statistics. This method creates
298// a copy (histNorm) of the histogram to normalize (hist) and normalize
299// the copy (histNorm) on its integral. It returns histNorm.
300//
301TH1 *MHMcTriggerLvl2::NormalizeHist(TH1 *histNorm, TH1 *hist) const
302{
303 if (histNorm == hist){
304 *fLog << "ARGH!@! You cannot pass the same histogram into each argument!" << endl;
305 return NULL;
306 }
307
308 if ((hist == NULL) || (hist->Integral() == (Stat_t)0)){
309 *fLog << "ARGH!@! You are trying to normalize the histogram to 0!" << endl;
310 return NULL;
311 }
312
313 histNorm->Reset("ICE");
314 histNorm->Add(hist, 1);
315 histNorm->Divide(fFNorm, (Double_t)(hist->Integral()));
316
317 return histNorm;
318}
Note: See TracBrowser for help on using the repository browser.