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

Last change on this file since 3873 was 3871, checked in by gaug, 21 years ago
*** empty log message ***
File size: 11.2 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)
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 c.cd(x+y);
89 gPad->SetBorderMode(0);
90 obj1->SetPrettyPalette();
91 obj1->Draw();
92
93 if (rad)
94 {
95 c.cd(x+2*y);
96 gPad->SetBorderMode(0);
97 gPad->SetTicks();
98 DrawRadialProfile(obj1);
99 }
100
101
102 if (!fit)
103 return;
104
105 c.cd(rad ? x+3*y : x+2*y);
106 gPad->SetBorderMode(0);
107 gPad->SetTicks();
108 DrawProjection(obj1, fit);
109}
110
111// --------------------------------------------------------------------------
112//
113// Draw a projection of MHCamera onto the y-axis values. Depending on the
114// variable fit, the following fits are performed:
115//
116// 1: Single Gauss (for distributions flat-fielded over the whole camera)
117// 2: Double Gauss (for distributions different for inner and outer pixels)
118// 3: Triple Gauss (for distributions with inner, outer pixels and outliers)
119// 4: flat (for the probability distributions)
120// 5: Fit Inner and Outer pixels separately by a single Gaussian
121// (only for MAGIC cameras)
122//
123// Moreover, sectors 6,1 and 2 of the camera and sectors 3,4 and 5 are
124// drawn separately, for inner and outer pixels.
125//
126void MGCamDisplays::DrawProjection(MHCamera *obj, Int_t fit) const
127{
128
129 TArrayI inner(1);
130 inner[0] = 0;
131
132 TArrayI outer(1);
133 outer[0] = 1;
134
135 if (fit==5)
136 {
137
138 if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
139 {
140 TArrayI s0(6);
141 s0[0] = 6;
142 s0[1] = 1;
143 s0[2] = 2;
144 s0[3] = 3;
145 s0[4] = 4;
146 s0[5] = 5;
147
148 gPad->Clear();
149 TVirtualPad *pad = gPad;
150 pad->Divide(2,1);
151
152 TH1D *half[2];
153 half[0] = obj->ProjectionS(s0, inner, "Inner");
154 half[1] = obj->ProjectionS(s0, outer, "Outer");
155
156 half[0]->SetDirectory(NULL);
157 half[1]->SetDirectory(NULL);
158
159 for (int i=0; i<2; i++)
160 {
161 pad->cd(i+1);
162 half[i]->SetLineColor(kRed+i);
163 half[i]->SetBit(kCanDelete);
164 half[i]->Draw();
165 half[i]->Fit("gaus","Q");
166 }
167
168 gLog << all << obj->GetName()
169 << Form("%s%5.3f%s%3.2f"," Mean: Inner Pixels: ",
170 half[0]->GetFunction("gaus")->GetParameter(1),"+-",
171 half[0]->GetFunction("gaus")->GetParError(1));
172 gLog << Form("%s%5.3f%s%3.2f"," Outer Pixels: ",
173 half[1]->GetFunction("gaus")->GetParameter(1),"+-",
174 half[1]->GetFunction("gaus")->GetParError(1)) << endl;
175 gLog << all << obj->GetName()
176 << Form("%s%5.3f%s%3.2f"," Sigma: Inner Pixels: ",
177 half[0]->GetFunction("gaus")->GetParameter(2),"+-",
178 half[0]->GetFunction("gaus")->GetParError(2));
179 gLog << Form("%s%5.3f%s%3.2f"," Outer Pixels: ",
180 half[1]->GetFunction("gaus")->GetParameter(2),"+-",
181 half[1]->GetFunction("gaus")->GetParError(2)) << endl;
182
183 }
184 return;
185 }
186
187 TH1D *obj2 = (TH1D*)obj->Projection(obj->GetName());
188 obj2->SetDirectory(0);
189 obj2->Draw();
190 obj2->SetBit(kCanDelete);
191
192
193 if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
194 {
195 TArrayI s0(3);
196 s0[0] = 6;
197 s0[1] = 1;
198 s0[2] = 2;
199
200 TArrayI s1(3);
201 s1[0] = 3;
202 s1[1] = 4;
203 s1[2] = 5;
204
205
206 TH1D *halfInOut[4];
207
208 // Just to get the right (maximum) binning
209 halfInOut[0] = obj->ProjectionS(s0, inner, "Sector 6-1-2 Inner");
210 halfInOut[1] = obj->ProjectionS(s1, inner, "Sector 3-4-5 Inner");
211 halfInOut[2] = obj->ProjectionS(s0, outer, "Sector 6-1-2 Outer");
212 halfInOut[3] = obj->ProjectionS(s1, outer, "Sector 3-4-5 Outer");
213
214 for (int i=0; i<4; i++)
215 {
216 halfInOut[i]->SetLineColor(kRed+i);
217 halfInOut[i]->SetDirectory(0);
218 halfInOut[i]->SetBit(kCanDelete);
219 halfInOut[i]->Draw("same");
220 }
221 }
222
223 const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
224 const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
225 const Double_t integ = obj2->Integral("width")/2.5;
226 const Double_t mean = obj2->GetMean();
227 const Double_t rms = obj2->GetRMS();
228 const Double_t width = max-min;
229
230 const TString dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
231 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
232
233 const TString tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
234 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"
235 "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])";
236 TF1 *f=0;
237 switch (fit)
238 {
239 case 1:
240 f = new TF1("sgaus", "gaus(0)", min, max);
241 f->SetLineColor(kYellow);
242 f->SetBit(kCanDelete);
243 f->SetParNames("Area", "#mu", "#sigma");
244 f->SetParameters(integ/rms, mean, rms);
245 f->SetParLimits(0, 0, integ);
246 f->SetParLimits(1, min, max);
247 f->SetParLimits(2, 0, width/1.5);
248
249 obj2->Fit(f, "QLR");
250 break;
251
252 case 2:
253 f = new TF1("dgaus",dgausformula.Data(),min,max);
254 f->SetLineColor(kYellow);
255 f->SetBit(kCanDelete);
256 f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2");
257 f->SetParameters(integ,(min+mean)/2.,width/4.,
258 integ/width/2.,(max+mean)/2.,width/4.);
259 // The left-sided Gauss
260 f->SetParLimits(0,integ-1.5 , integ+1.5);
261 f->SetParLimits(1,min+(width/10.), mean);
262 f->SetParLimits(2,0 , width/2.);
263 // The right-sided Gauss
264 f->SetParLimits(3,0 , integ);
265 f->SetParLimits(4,mean, max-(width/10.));
266 f->SetParLimits(5,0 , width/2.);
267 obj2->Fit(f,"QLRM");
268 break;
269
270 case 3:
271 f = new TF1("tgaus",tgausformula.Data(),min,max);
272 f->SetLineColor(kYellow);
273 f->SetBit(kCanDelete);
274 f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}",
275 "A_{2}","#mu_{2}","#sigma_{2}",
276 "A_{3}","#mu_{3}","#sigma_{3}");
277 f->SetParameters(integ,(min+mean)/2,width/4.,
278 integ/width/3.,(max+mean)/2.,width/4.,
279 integ/width/3.,mean,width/2.);
280 // The left-sided Gauss
281 f->SetParLimits(0,integ-1.5,integ+1.5);
282 f->SetParLimits(1,min+(width/10.),mean);
283 f->SetParLimits(2,width/15.,width/2.);
284 // The right-sided Gauss
285 f->SetParLimits(3,0.,integ);
286 f->SetParLimits(4,mean,max-(width/10.));
287 f->SetParLimits(5,width/15.,width/2.);
288 // The Gauss describing the outliers
289 f->SetParLimits(6,0.,integ);
290 f->SetParLimits(7,min,max);
291 f->SetParLimits(8,width/4.,width/1.5);
292 obj2->Fit(f,"QLRM");
293 break;
294
295 case 4:
296 obj2->Fit("pol0", "Q");
297 obj2->GetFunction("pol0")->SetLineColor(kYellow);
298 break;
299
300 case 9:
301 break;
302
303 default:
304 obj2->Fit("gaus", "Q");
305 obj2->GetFunction("gaus")->SetLineColor(kYellow);
306 break;
307 }
308}
309
310// --------------------------------------------------------------------------
311//
312// Draw a projection of MHCamera vs. the radius from the central pixel.
313//
314// The inner and outer pixels are drawn separately, both fitted by a polynomial
315// of grade 1.
316//
317void MGCamDisplays::DrawRadialProfile(MHCamera *obj) const
318{
319
320 TProfile *obj2 = (TProfile*)obj->RadialProfile(obj->GetName());
321 obj2->SetDirectory(0);
322 obj2->Draw();
323 obj2->SetBit(kCanDelete);
324
325 if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
326 {
327
328 TArrayI s0(6);
329 s0[0] = 1;
330 s0[1] = 2;
331 s0[2] = 3;
332 s0[3] = 4;
333 s0[4] = 5;
334 s0[5] = 6;
335
336 TArrayI inner(1);
337 inner[0] = 0;
338
339 TArrayI outer(1);
340 outer[0] = 1;
341
342 // Just to get the right (maximum) binning
343 TProfile *half[2];
344 half[0] = obj->RadialProfileS(s0, inner,Form("%s%s",obj->GetName(),"Inner"));
345 half[1] = obj->RadialProfileS(s0, outer,Form("%s%s",obj->GetName(),"Outer"));
346
347 for (Int_t i=0; i<2; i++)
348 {
349 Double_t min = obj->GetGeomCam().GetMinRadius(i);
350 Double_t max = obj->GetGeomCam().GetMaxRadius(i);
351
352 half[i]->SetLineColor(kRed+i);
353 half[i]->SetDirectory(0);
354 half[i]->SetBit(kCanDelete);
355 half[i]->Draw("same");
356 half[i]->Fit("pol1","Q","",min,max);
357 half[i]->GetFunction("pol1")->SetLineColor(kRed+i);
358 half[i]->GetFunction("pol1")->SetLineWidth(1);
359 }
360 }
361}
362
Note: See TracBrowser for help on using the repository browser.