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

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