source: trunk/MagicSoft/Mars/macros/calibrate_data.C@ 5293

Last change on this file since 5293 was 4120, checked in by moralejo, 20 years ago
*** empty log message ***
File size: 11.9 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): Hendrik Bartko, 03/2004 <mailto:hbartko@mppmu.mpg.de>
19! Markus Gaug, 03/2004 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25
26const TString defpath = "./";
27const TString defrout = "output_test.root";
28
29const Int_t defpedr [] = {20122};
30const Int_t defcalr [] = {20125};
31const Int_t defdatar[] = {20122};
32
33void calibrate_data(const TString inpath=defpath,
34 const Int_t psize=1, const Int_t pedruns[]=defpedr,
35 const Int_t csize=1, const Int_t calruns[]=defcalr,
36 const Int_t dsize=1, const Int_t dataruns[]=defdatar,
37 const TString resname=defrout)
38
39{
40
41// MExtractSlidingWindow extractor;
42
43 MExtractFixedWindowPeakSearch extractor;
44 MExtractTimeFastSpline timeext;
45
46 MRunIter pruns;
47 MRunIter cruns;
48 MRunIter druns;
49
50 for (Int_t i=0;i<psize;i++) {
51 cout << "Adding pedestal run: " << pedruns[i] << endl;
52 pruns.AddRun(pedruns[i],inpath);
53 }
54 for (Int_t i=0;i<csize;i++) {
55 cout << "Adding calibration run: " << calruns[i] << endl;
56 cruns.AddRun(calruns[i],inpath);
57 }
58 for (Int_t i=0;i<dsize;i++) {
59 cout << "Adding data run: " << dataruns[i] << endl;
60 druns.AddRun(dataruns[i],inpath);
61 }
62
63
64 MStatusDisplay *display = new MStatusDisplay;
65 display->SetUpdateTime(3000);
66 display->Resize(850,700);
67
68 gStyle->SetOptStat(1111);
69 gStyle->SetOptFit();
70
71 /************************************/
72 /* FIRST LOOP: PEDESTAL COMPUTATION */
73 /************************************/
74
75 MParList plist1;
76 MTaskList tlist1;
77 plist1.AddToList(&tlist1);
78
79 // containers
80 MPedestalCam pedcam;
81 MBadPixelsCam badcam;
82 //
83 // for excluding pixels from the beginning:
84 //
85 // badcam.AsciiRead("badpixels.dat");
86
87
88 plist1.AddToList(&pedcam);
89 plist1.AddToList(&badcam);
90
91 //tasks
92 MReadMarsFile read("Events");
93 read.DisableAutoScheme();
94 static_cast<MRead&>(read).AddFiles(pruns);
95
96 MGeomApply geomapl;
97 MPedCalcPedRun pedcalc;
98 MGeomCamMagic geomcam;
99
100 tlist1.AddToList(&read);
101 tlist1.AddToList(&geomapl);
102 tlist1.AddToList(&pedcalc);
103
104 // Create and execute the event looper
105 MEvtLoop pedloop;
106 pedloop.SetParList(&plist1);
107 pedloop.SetDisplay(display);
108
109 cout << "*************************" << endl;
110 cout << "** COMPUTING PEDESTALS **" << endl;
111 cout << "*************************" << endl;
112
113 if (!pedloop.Eventloop())
114 return;
115
116 tlist1.PrintStatistics();
117
118 //
119 // Now the short version:
120 //
121 //
122 // Now setup the new tasks for the calibration:
123 // ---------------------------------------------------
124 //
125 MJCalibration calloop;
126 calloop.SetInput(&cruns);
127 // calloop.SetFullDisplay();
128 //
129 calloop.SetExtractor(&extractor);
130 //
131 // Set the corr. cams:
132 //
133 calloop.SetBadPixels(badcam);
134 //
135 // The next two commands are for the display:
136 //
137 calloop.SetDisplay(display);
138
139 //
140 // Apply rel. time calibration:
141 //
142 calloop.SetRelTimeCalibration();
143 calloop.SetTimeExtractor(&timeext);
144 //
145 // Do the event-loop:
146 //
147 cout << "***************************" << endl;
148 cout << "** COMPUTING CALIBRATION **" << endl;
149 cout << "***************************" << endl;
150
151 if (!calloop.Process(pedcam))
152 return;
153
154 badcam.Print();
155
156 MBadPixelsCam &badbad = calloop.GetBadPixels();
157 MCalibrationChargeCam &calcam = calloop.GetCalibrationCam();
158 MCalibrationRelTimeCam &timecam = calloop.GetRelTimeCam();
159 MCalibrationQECam &qecam = calloop.GetQECam();
160
161 badbad.Print();
162
163 /************************************************************************/
164 /* THIRD LOOP: DATA CALIBRATION INTO PHOTONS */
165 /************************************************************************/
166
167 // Create an empty Parameter List and an empty Task List
168 MParList plist3;
169 MTaskList tlist3;
170 plist3.AddToList(&tlist3);
171
172 // containers
173 MCerPhotEvt photevt;
174 MPedPhotCam pedphotcam;
175 MSrcPosCam srccam;
176 MRawRunHeader runhead;
177 MExtractedSignalCam sigcam;
178
179 plist3.AddToList(&geomcam );
180 plist3.AddToList(&pedcam );
181 plist3.AddToList(&calcam );
182 plist3.AddToList(&qecam );
183 plist3.AddToList(&badbad );
184 plist3.AddToList(&timecam );
185 plist3.AddToList(&sigcam );
186 plist3.AddToList(&photevt);
187 plist3.AddToList(&pedphotcam);
188 plist3.AddToList(&srccam);
189 plist3.AddToList(&runhead);
190
191 //tasks
192 MReadMarsFile read3("Events");
193 read3.DisableAutoScheme();
194 static_cast<MRead&>(read3).AddFiles(druns);
195
196 MArrivalTimeCalc2 timecalc;
197 MCalibrateData photcalc;
198 photcalc.SetCalibrationMode(MCalibrateData::kFfactor); // !!! was only MCalibrate
199 // MPedPhotCalc pedphotcalc; // already done by MCalibrate Data
200 // MCerPhotCalc cerphotcalc; // already done by MCalibrate Data
201
202 tlist3.AddToList(&read3);
203 tlist3.AddToList(&geomapl);
204 tlist3.AddToList(&extractor);
205 tlist3.AddToList(&timecalc);
206 // tlist3.AddToList(&cerphotcalc); // already done by MCalibrate Data
207 tlist3.AddToList(&photcalc);
208 // tlist3.AddToList(&pedphotcalc); // already done by MCalibrate Data
209
210 MWriteRootFile write(resname);
211
212 write.AddContainer("MGeomCam" , "RunHeaders");
213 write.AddContainer("MRawRunHeader" , "RunHeaders");
214 write.AddContainer("MSrcPosCam" , "RunHeaders");
215 write.AddContainer("MCalibrationChargeCam" , "RunHeaders");
216 write.AddContainer("MCalibrationQECam" , "RunHeaders");
217 // write.AddContainer("MPedPhotCam","RunHeaders"); // Attention, was in Events - Tree!!
218 write.AddContainer("MPedestalCam" , "RunHeaders");
219 write.AddContainer("MCalibrationRelTimeCam", "RunHeaders");
220
221 write.AddContainer("MCerPhotEvt" , "Events");
222 write.AddContainer("MRawEvtHeader" , "Events");
223 write.AddContainer("MBadPixelsCam" , "Events");
224 write.AddContainer("MPedPhotCam" , "Events");
225
226 tlist3.AddToList(&write);
227
228 // Create and execute eventloop
229 MEvtLoop evtloop3;
230 evtloop3.SetParList(&plist3);
231
232 cout << "*************************************************************" << endl;
233 cout << "*** COMPUTING DATA USING EXTRACTED SIGNAL (IN PHOTONS) ***" << endl;
234 cout << "*************************************************************" << endl;
235
236 if (!evtloop3.Eventloop())
237 return;
238 tlist3.PrintStatistics();
239
240}
241
242void CamDraw(TCanvas &c, MHCamera &cam, MCamEvent &evt, Int_t i, Int_t j, Int_t fit)
243{
244
245 c.cd(i);
246 gPad->SetBorderMode(0);
247 MHCamera *obj1=(MHCamera*)cam.DrawCopy("hist");
248 // obj1->AddNotify(evt);
249
250 c.cd(i+j);
251 gPad->SetBorderMode(0);
252 obj1->Draw();
253 ((MHCamera*)obj1)->SetPrettyPalette();
254
255 if (fit != 0)
256 {
257 c.cd(i+2*j);
258 gPad->SetBorderMode(0);
259 TH1D *obj2 = (TH1D*)obj1->Projection(obj1.GetName());
260
261// obj2->Sumw2();
262 obj2->Draw();
263 obj2->SetBit(kCanDelete);
264
265 const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
266 const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
267 const Double_t integ = obj2->Integral("width")/2.5066283;
268 const Double_t mean = obj2->GetMean();
269 const Double_t rms = obj2->GetRMS();
270 const Double_t width = max-min;
271
272 if (rms == 0. || width == 0. )
273 return;
274
275 switch (fit)
276 {
277 case 1:
278 TF1 *sgaus = new TF1("sgaus","gaus(0)",min,max);
279 sgaus->SetBit(kCanDelete);
280 sgaus->SetParNames("Area","#mu","#sigma");
281 sgaus->SetParameters(integ/rms,mean,rms);
282 sgaus->SetParLimits(0,0.,integ);
283 sgaus->SetParLimits(1,min,max);
284 sgaus->SetParLimits(2,0,width/1.5);
285 obj2->Fit("sgaus","QLR");
286 obj2->GetFunction("sgaus")->SetLineColor(kYellow);
287 break;
288
289 case 2:
290 TString dgausform = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])";
291 dgausform += "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
292 TF1 *dgaus = new TF1("dgaus",dgausform.Data(),min,max);
293 dgaus->SetBit(kCanDelete);
294 dgaus->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}","A_{2}","#mu_{2}","#sigma_{2}");
295 dgaus->SetParameters(integ,(min+mean)/2.,width/4.,
296 integ/width/2.,(max+mean)/2.,width/4.);
297 // The left-sided Gauss
298 dgaus->SetParLimits(0,integ-1.5,integ+1.5);
299 dgaus->SetParLimits(1,min+(width/10.),mean);
300 dgaus->SetParLimits(2,0,width/2.);
301 // The right-sided Gauss
302 dgaus->SetParLimits(3,0,integ);
303 dgaus->SetParLimits(4,mean,max-(width/10.));
304 dgaus->SetParLimits(5,0,width/2.);
305 obj2->Fit("dgaus","QLRM");
306 obj2->GetFunction("dgaus")->SetLineColor(kYellow);
307 break;
308
309 case 3:
310 TString tgausform = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])";
311 tgausform += "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
312 tgausform += "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])";
313 TF1 *tgaus = new TF1("tgaus",tgausform.Data(),min,max);
314 tgaus->SetBit(kCanDelete);
315 tgaus->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}",
316 "A_{2}","#mu_{2}","#sigma_{2}",
317 "A_{3}","#mu_{3}","#sigma_{3}");
318 tgaus->SetParameters(integ,(min+mean)/2,width/4.,
319 integ/width/3.,(max+mean)/2.,width/4.,
320 integ/width/3.,mean,width/2.);
321 // The left-sided Gauss
322 tgaus->SetParLimits(0,integ-1.5,integ+1.5);
323 tgaus->SetParLimits(1,min+(width/10.),mean);
324 tgaus->SetParLimits(2,width/15.,width/2.);
325 // The right-sided Gauss
326 tgaus->SetParLimits(3,0.,integ);
327 tgaus->SetParLimits(4,mean,max-(width/10.));
328 tgaus->SetParLimits(5,width/15.,width/2.);
329 // The Gauss describing the outliers
330 tgaus->SetParLimits(6,0.,integ);
331 tgaus->SetParLimits(7,min,max);
332 tgaus->SetParLimits(8,width/4.,width/1.5);
333 obj2->Fit("tgaus","QLRM");
334 obj2->GetFunction("tgaus")->SetLineColor(kYellow);
335 break;
336 case 4:
337 obj2->Fit("pol0","Q");
338 obj2->GetFunction("pol0")->SetLineColor(kYellow);
339 break;
340 case 9:
341 break;
342 default:
343 obj2->Fit("gaus","Q");
344 obj2->GetFunction("gaus")->SetLineColor(kYellow);
345 break;
346 }
347
348 TArrayI s0(3);
349 s0[0] = 6;
350 s0[1] = 1;
351 s0[2] = 2;
352
353 TArrayI s1(3);
354 s1[0] = 3;
355 s1[1] = 4;
356 s1[2] = 5;
357
358 TArrayI inner(1);
359 inner[0] = 0;
360
361 TArrayI outer(1);
362 outer[0] = 1;
363
364 // Just to get the right (maximum) binning
365 TH1D *half[4];
366 half[0] = obj1->ProjectionS(s0, inner, "Sector 6-1-2 Inner");
367 half[1] = obj1->ProjectionS(s1, inner, "Sector 3-4-5 Inner");
368 half[2] = obj1->ProjectionS(s0, outer, "Sector 6-1-2 Outer");
369 half[3] = obj1->ProjectionS(s1, outer, "Sector 3-4-5 Outer");
370
371 for (int i=0; i<4; i++)
372 {
373 half[i]->SetLineColor(kRed+i);
374 half[i]->SetDirectory(0);
375 half[i]->SetBit(kCanDelete);
376 half[i]->Draw("same");
377 }
378
379 gPad->Modified();
380 gPad->Update();
381
382 }
383}
384
385
Note: See TracBrowser for help on using the repository browser.