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

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