source: trunk/MagicSoft/Mars/mjobs/MGCamDisplays.cc@ 3875

Last change on this file since 3875 was 3875, checked in by gaug, 21 years ago
*** empty log message ***
File size: 11.3 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 Bretz, 1/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Markus Gaug, 02/2004 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25/////////////////////////////////////////////////////////////////////////////
26//
27// MGCamDisplays
28//
29// Graphical interfaces to display the camera with fits and projections
30//
31/////////////////////////////////////////////////////////////////////////////
32#include "MGCamDisplays.h"
33
34#include <TStyle.h>
35#include <TCanvas.h>
36
37#include "MH.h"
38#include "MHCamera.h"
39#include "MGeomCam.h"
40#include "TVirtualPad.h"
41#include "TProfile.h"
42#include "TF1.h"
43
44#include "MLog.h"
45#include "MLogManip.h"
46
47#include "MStatusDisplay.h"
48
49ClassImp(MGCamDisplays);
50
51using namespace std;
52
53// --------------------------------------------------------------------------
54//
55// Default constructor.
56//
57MGCamDisplays::MGCamDisplays()
58{
59}
60
61// --------------------------------------------------------------------------
62//
63// Draw the MHCamera into the MStatusDisplay:
64//
65// 1) Draw it as histogram (MHCamera::DrawCopy("hist")
66// 2) Draw it as a camera, with MHCamera::SetPrettyPalette() set.
67// 3) If "rad" is not zero, draw its values vs. the radius from the camera center.
68// (DrawRadialProfile())
69// 4) Depending on the variable "fit", draw the values projection on the y-axis
70// (DrawProjection()):
71// 0: don't draw
72// 1: Draw fit to Single Gauss (for distributions flat-fielded over the whole camera)
73// 2: Draw and fit to Double Gauss (for distributions different for inner and outer pixels)
74// 3: Draw and fit to Triple Gauss (for distributions with inner, outer pixels and outliers)
75// 4: Draw and fit to Polynomial grade 0: (for the probability distributions)
76// >4: Draw and don;t fit.
77//
78void MGCamDisplays::CamDraw(TCanvas &c, const Int_t x, const Int_t y, const MHCamera &cam1,
79 const Int_t fit, const Int_t rad, TObject *notify)
80{
81
82 c.cd(x);
83 gPad->SetBorderMode(0);
84 gPad->SetTicks();
85 MHCamera *obj1=(MHCamera*)cam1.DrawCopy("hist");
86 obj1->SetDirectory(NULL);
87
88 if (notify)
89 obj1->AddNotify(notify);
90
91 c.cd(x+y);
92 gPad->SetBorderMode(0);
93 obj1->SetPrettyPalette();
94 obj1->Draw();
95
96 if (rad)
97 {
98 c.cd(x+2*y);
99 gPad->SetBorderMode(0);
100 gPad->SetTicks();
101 DrawRadialProfile(obj1);
102 }
103
104
105 if (!fit)
106 return;
107
108 c.cd(rad ? x+3*y : x+2*y);
109 gPad->SetBorderMode(0);
110 gPad->SetTicks();
111 DrawProjection(obj1, fit);
112}
113
114// --------------------------------------------------------------------------
115//
116// Draw a projection of MHCamera onto the y-axis values. Depending on the
117// variable fit, the following fits are performed:
118//
119// 1: Single Gauss (for distributions flat-fielded over the whole camera)
120// 2: Double Gauss (for distributions different for inner and outer pixels)
121// 3: Triple Gauss (for distributions with inner, outer pixels and outliers)
122// 4: flat (for the probability distributions)
123// 5: Fit Inner and Outer pixels separately by a single Gaussian
124// (only for MAGIC cameras)
125//
126// Moreover, sectors 6,1 and 2 of the camera and sectors 3,4 and 5 are
127// drawn separately, for inner and outer pixels.
128//
129void MGCamDisplays::DrawProjection(MHCamera *obj, Int_t fit) const
130{
131
132 TArrayI inner(1);
133 inner[0] = 0;
134
135 TArrayI outer(1);
136 outer[0] = 1;
137
138 if (fit==5)
139 {
140
141 if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
142 {
143 TArrayI s0(6);
144 s0[0] = 6;
145 s0[1] = 1;
146 s0[2] = 2;
147 s0[3] = 3;
148 s0[4] = 4;
149 s0[5] = 5;
150
151 gPad->Clear();
152 TVirtualPad *pad = gPad;
153 pad->Divide(2,1);
154
155 TH1D *half[2];
156 half[0] = obj->ProjectionS(s0, inner, "Inner");
157 half[1] = obj->ProjectionS(s0, outer, "Outer");
158
159 half[0]->SetDirectory(NULL);
160 half[1]->SetDirectory(NULL);
161
162 for (int i=0; i<2; i++)
163 {
164 pad->cd(i+1);
165 half[i]->SetLineColor(kRed+i);
166 half[i]->SetBit(kCanDelete);
167 half[i]->Draw();
168 half[i]->Fit("gaus","Q");
169 }
170
171 gLog << all << obj->GetName()
172 << Form("%s%5.3f%s%3.2f"," Mean: Inner Pixels: ",
173 half[0]->GetFunction("gaus")->GetParameter(1),"+-",
174 half[0]->GetFunction("gaus")->GetParError(1));
175 gLog << Form("%s%5.3f%s%3.2f"," Outer Pixels: ",
176 half[1]->GetFunction("gaus")->GetParameter(1),"+-",
177 half[1]->GetFunction("gaus")->GetParError(1)) << endl;
178 gLog << all << obj->GetName()
179 << Form("%s%5.3f%s%3.2f"," Sigma: Inner Pixels: ",
180 half[0]->GetFunction("gaus")->GetParameter(2),"+-",
181 half[0]->GetFunction("gaus")->GetParError(2));
182 gLog << Form("%s%5.3f%s%3.2f"," Outer Pixels: ",
183 half[1]->GetFunction("gaus")->GetParameter(2),"+-",
184 half[1]->GetFunction("gaus")->GetParError(2)) << endl;
185
186 }
187 return;
188 }
189
190 TH1D *obj2 = (TH1D*)obj->Projection(obj->GetName());
191 obj2->SetDirectory(0);
192 obj2->Draw();
193 obj2->SetBit(kCanDelete);
194
195
196 if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
197 {
198 TArrayI s0(3);
199 s0[0] = 6;
200 s0[1] = 1;
201 s0[2] = 2;
202
203 TArrayI s1(3);
204 s1[0] = 3;
205 s1[1] = 4;
206 s1[2] = 5;
207
208
209 TH1D *halfInOut[4];
210
211 // Just to get the right (maximum) binning
212 halfInOut[0] = obj->ProjectionS(s0, inner, "Sector 6-1-2 Inner");
213 halfInOut[1] = obj->ProjectionS(s1, inner, "Sector 3-4-5 Inner");
214 halfInOut[2] = obj->ProjectionS(s0, outer, "Sector 6-1-2 Outer");
215 halfInOut[3] = obj->ProjectionS(s1, outer, "Sector 3-4-5 Outer");
216
217 for (int i=0; i<4; i++)
218 {
219 halfInOut[i]->SetLineColor(kRed+i);
220 halfInOut[i]->SetDirectory(0);
221 halfInOut[i]->SetBit(kCanDelete);
222 halfInOut[i]->Draw("same");
223 }
224 }
225
226 const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
227 const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
228 const Double_t integ = obj2->Integral("width")/2.5;
229 const Double_t mean = obj2->GetMean();
230 const Double_t rms = obj2->GetRMS();
231 const Double_t width = max-min;
232
233 const TString dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
234 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
235
236 const TString tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
237 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"
238 "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])";
239 TF1 *f=0;
240 switch (fit)
241 {
242 case 1:
243 f = new TF1("sgaus", "gaus(0)", min, max);
244 f->SetLineColor(kYellow);
245 f->SetBit(kCanDelete);
246 f->SetParNames("Area", "#mu", "#sigma");
247 f->SetParameters(integ/rms, mean, rms);
248 f->SetParLimits(0, 0, integ);
249 f->SetParLimits(1, min, max);
250 f->SetParLimits(2, 0, width/1.5);
251
252 obj2->Fit(f, "QLR");
253 break;
254
255 case 2:
256 f = new TF1("dgaus",dgausformula.Data(),min,max);
257 f->SetLineColor(kYellow);
258 f->SetBit(kCanDelete);
259 f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2");
260 f->SetParameters(integ,(min+mean)/2.,width/4.,
261 integ/width/2.,(max+mean)/2.,width/4.);
262 // The left-sided Gauss
263 f->SetParLimits(0,integ-1.5 , integ+1.5);
264 f->SetParLimits(1,min+(width/10.), mean);
265 f->SetParLimits(2,0 , width/2.);
266 // The right-sided Gauss
267 f->SetParLimits(3,0 , integ);
268 f->SetParLimits(4,mean, max-(width/10.));
269 f->SetParLimits(5,0 , width/2.);
270 obj2->Fit(f,"QLRM");
271 break;
272
273 case 3:
274 f = new TF1("tgaus",tgausformula.Data(),min,max);
275 f->SetLineColor(kYellow);
276 f->SetBit(kCanDelete);
277 f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}",
278 "A_{2}","#mu_{2}","#sigma_{2}",
279 "A_{3}","#mu_{3}","#sigma_{3}");
280 f->SetParameters(integ,(min+mean)/2,width/4.,
281 integ/width/3.,(max+mean)/2.,width/4.,
282 integ/width/3.,mean,width/2.);
283 // The left-sided Gauss
284 f->SetParLimits(0,integ-1.5,integ+1.5);
285 f->SetParLimits(1,min+(width/10.),mean);
286 f->SetParLimits(2,width/15.,width/2.);
287 // The right-sided Gauss
288 f->SetParLimits(3,0.,integ);
289 f->SetParLimits(4,mean,max-(width/10.));
290 f->SetParLimits(5,width/15.,width/2.);
291 // The Gauss describing the outliers
292 f->SetParLimits(6,0.,integ);
293 f->SetParLimits(7,min,max);
294 f->SetParLimits(8,width/4.,width/1.5);
295 obj2->Fit(f,"QLRM");
296 break;
297
298 case 4:
299 obj2->Fit("pol0", "Q");
300 obj2->GetFunction("pol0")->SetLineColor(kYellow);
301 break;
302
303 case 9:
304 break;
305
306 default:
307 obj2->Fit("gaus", "Q");
308 obj2->GetFunction("gaus")->SetLineColor(kYellow);
309 break;
310 }
311}
312
313// --------------------------------------------------------------------------
314//
315// Draw a projection of MHCamera vs. the radius from the central pixel.
316//
317// The inner and outer pixels are drawn separately, both fitted by a polynomial
318// of grade 1.
319//
320void MGCamDisplays::DrawRadialProfile(MHCamera *obj) const
321{
322
323 TProfile *obj2 = (TProfile*)obj->RadialProfile(obj->GetName());
324 obj2->SetDirectory(0);
325 obj2->Draw();
326 obj2->SetBit(kCanDelete);
327
328 if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
329 {
330
331 TArrayI s0(6);
332 s0[0] = 1;
333 s0[1] = 2;
334 s0[2] = 3;
335 s0[3] = 4;
336 s0[4] = 5;
337 s0[5] = 6;
338
339 TArrayI inner(1);
340 inner[0] = 0;
341
342 TArrayI outer(1);
343 outer[0] = 1;
344
345 // Just to get the right (maximum) binning
346 TProfile *half[2];
347 half[0] = obj->RadialProfileS(s0, inner,Form("%s%s",obj->GetName(),"Inner"));
348 half[1] = obj->RadialProfileS(s0, outer,Form("%s%s",obj->GetName(),"Outer"));
349
350 for (Int_t i=0; i<2; i++)
351 {
352 Double_t min = obj->GetGeomCam().GetMinRadius(i);
353 Double_t max = obj->GetGeomCam().GetMaxRadius(i);
354
355 half[i]->SetLineColor(kRed+i);
356 half[i]->SetDirectory(0);
357 half[i]->SetBit(kCanDelete);
358 half[i]->Draw("same");
359 half[i]->Fit("pol1","Q","",min,max);
360 half[i]->GetFunction("pol1")->SetLineColor(kRed+i);
361 half[i]->GetFunction("pol1")->SetLineWidth(1);
362 }
363 }
364}
365
Note: See TracBrowser for help on using the repository browser.