source: trunk/MagicSoft/Mars/mtemp/mpisa/macros/production.C@ 5138

Last change on this file since 5138 was 4055, checked in by stamerra, 21 years ago
*** empty log message ***
File size: 9.0 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
26///////////////////////////////////////////////////////////////////////////
27//
28// production.C
29// ============
30//
31// This macro takes calibration runs and pedestal runs
32// to set pedestals and calibration factors
33// Then raw data files are read and calibrated.
34// The output file contains the calibrated data, pedestal subtracted,
35// with the signal (in photons) for each pixels
36//
37// input: calibration file(s)
38// pedestal file(s)
39// raw data files
40//
41// output: root files with following containers (branches)
42// - MCerPhotEvt
43// - MPedestalCam
44// - MCalibrationRelTimeCam
45// - MCerPhotEvt
46// - MRawEvtHeader
47// - MBadPixelsCam
48// - MPedPhotCam
49//
50///////////////////////////////////////////////////////////////////////////
51
52//
53// Set path for files and run numbers
54//
55const TString defpath = "/data1/earth/magic/data/Period014/rootdata/2004_02_15/";
56const TString defrout = "output_Mrk421_20040215.root";
57
58const Int_t defpedr [] = {17284};
59const Int_t defcalr [] = {17285};
60const Int_t defdatar[] = {17286,17287,17288,17289,17290,17291,17292,17293,
61 17294,17295,17296,17297,17298,17299,17300,17301,
62 17302,17303,17304,17305,17306,17307,17308,17309,
63 17310,17311,17312,17313,17314,17315,17316,17317,
64 17318,17319,17320,17321,17322,17323,17324,17325,
65 17326,17327,17328,17329,17330,17331,17332,17333,
66 17334,17335,17336,17337,17338,17339,17340,17341,
67 17342,17343,17344,17345,17346,17347,17348,17349,
68 17350,17351,17352,17353,17354,17355,17356,17357,
69 17358,17359,17360,17361,17362,17363,17364,17365,
70 17366,17367,17368,17369,17370,17371,17372,17373,
71 17374,17375,17376,17377,17378,17379,17380,17381,
72 17382,17383,17384,17385,17386,17387,17388,17389,
73 17390,17391,17392,17393,17394,17395,17396,17397,
74 17398,17399,17400,17401,17402,17403,17404,17405,
75 17406,17407,17408,17409,17410,17411,17412,17413,
76 17414,17415,17416,17417,17418,17419,17420,17421,
77 17422,17423,17424,17425,17426,17427,17428,17429,
78 17430,17431,17432,17433,17434,17435,17436,17437,
79 17438,17439,17440,17441,17442,17443,17444,17445,
80 17446,17447,17448,17449,17450,17451,17452,17453,
81 17454,17455,17456};
82
83
84void production(const TString inpath=defpath,
85 const Int_t psize=1, const Int_t pedruns[]=defpedr,
86 const Int_t csize=1, const Int_t calruns[]=defcalr,
87 const Int_t dsize=171, const Int_t dataruns[]=defdatar,
88 const TString resname=defrout)
89
90{
91
92 MExtractSlidingWindow extractor;
93
94 MRunIter pruns;
95 MRunIter cruns;
96 MRunIter druns;
97
98 for (Int_t i=0;i<psize;i++) {
99 cout << "Adding pedestal run: " << pedruns[i] << endl;
100 pruns.AddRun(pedruns[i],inpath);
101 }
102 for (Int_t i=0;i<csize;i++) {
103 cout << "Adding calibration run: " << calruns[i] << endl;
104 cruns.AddRun(calruns[i],inpath);
105 }
106 for (Int_t i=0;i<dsize;i++) {
107 cout << "Adding data run: " << dataruns[i] << endl;
108 druns.AddRun(dataruns[i],inpath);
109 }
110
111
112 MStatusDisplay *display = new MStatusDisplay;
113 display->SetUpdateTime(3000);
114 display->Resize(850,700);
115
116 gStyle->SetOptStat(1111);
117 gStyle->SetOptFit();
118
119 /************************************/
120 /* FIRST LOOP: PEDESTAL COMPUTATION */
121 /************************************/
122
123 MParList plist1;
124 MTaskList tlist1;
125 plist1.AddToList(&tlist1);
126
127 // containers
128 MPedestalCam pedcam;
129 MBadPixelsCam badcam;
130 //
131 // for excluding pixels from the beginning:
132 //
133 // badcam.AsciiRead("badpixels.dat");
134
135
136 plist1.AddToList(&pedcam);
137 plist1.AddToList(&badcam);
138
139 //tasks
140 MReadMarsFile read("Events");
141 read.DisableAutoScheme();
142 static_cast<MRead&>(read).AddFiles(pruns);
143
144 MGeomApply geomapl;
145 MPedCalcPedRun pedcalc;
146 MGeomCamMagic geomcam;
147
148 tlist1.AddToList(&read);
149 tlist1.AddToList(&geomapl);
150 tlist1.AddToList(&pedcalc);
151
152 // Create and execute the event looper
153 MEvtLoop pedloop;
154 pedloop.SetParList(&plist1);
155 pedloop.SetDisplay(display);
156
157 cout << "*************************" << endl;
158 cout << "** COMPUTING PEDESTALS **" << endl;
159 cout << "*************************" << endl;
160
161 if (!pedloop.Eventloop())
162 return;
163
164 tlist1.PrintStatistics();
165
166 //
167 // Now the short version:
168 //
169 //
170 // Now setup the new tasks for the calibration:
171 // ---------------------------------------------------
172 //
173 MJCalibration calloop;
174 calloop.SetInput(&cruns);
175 // calloop.SetFullDisplay();
176 //
177 calloop.SetExtractor(&extractor);
178 //
179 // Set the corr. cams:
180 //
181 calloop.SetBadPixels(badcam);
182 //
183 // The next two commands are for the display:
184 //
185 calloop.SetDisplay(display);
186
187 //
188 // Apply rel. time calibration:
189 //
190 calloop.SetRelTimeCalibration();
191 //
192 // Use as arrival time extractor MArrivalTimeCalc2:
193 //
194 calloop.SetArrivalTimeLevel(2);
195
196 //
197 // Do the event-loop:
198 //
199 cout << "***************************" << endl;
200 cout << "** COMPUTING CALIBRATION **" << endl;
201 cout << "***************************" << endl;
202
203 if (!calloop.Process(pedcam))
204 return;
205
206 badcam.Print();
207
208 MBadPixelsCam &badbad = calloop.GetBadPixels();
209 MCalibrationChargeCam &calcam = calloop.GetCalibrationCam();
210 MCalibrationRelTimeCam &timecam = calloop.GetRelTimeCam();
211 MCalibrationQECam &qecam = calloop.GetQECam();
212
213 badbad.Print();
214
215 /************************************************************************/
216 /* THIRD LOOP: DATA CALIBRATION INTO PHOTONS */
217 /************************************************************************/
218
219 // Create an empty Parameter List and an empty Task List
220 MParList plist3;
221 MTaskList tlist3;
222 plist3.AddToList(&tlist3);
223
224 // containers
225 MCerPhotEvt photevt;
226 MPedPhotCam pedphotcam;
227 MSrcPosCam srccam;
228 MRawRunHeader runhead;
229 MExtractedSignalCam sigcam;
230
231 plist3.AddToList(&geomcam );
232 plist3.AddToList(&pedcam );
233 plist3.AddToList(&calcam );
234 plist3.AddToList(&qecam );
235 plist3.AddToList(&badbad );
236 plist3.AddToList(&timecam );
237 plist3.AddToList(&sigcam );
238 plist3.AddToList(&photevt);
239 plist3.AddToList(&pedphotcam);
240 plist3.AddToList(&srccam);
241 plist3.AddToList(&runhead);
242
243 //tasks
244 MReadMarsFile read3("Events");
245 read3.DisableAutoScheme();
246 static_cast<MRead&>(read3).AddFiles(druns);
247
248 MArrivalTimeCalc2 timecalc;
249 MCalibrateData photcalc;
250 photcalc.SetCalibrationMode(MCalibrateData::kFfactor); // !!! was only MCalibrate
251 // MPedPhotCalc pedphotcalc; // already done by MCalibrate Data
252 // MCerPhotCalc cerphotcalc; // already done by MCalibrate Data
253
254 tlist3.AddToList(&read3);
255 tlist3.AddToList(&geomapl);
256 tlist3.AddToList(&extractor);
257 tlist3.AddToList(&timecalc);
258 // tlist3.AddToList(&cerphotcalc); // already done by MCalibrate Data
259 tlist3.AddToList(&photcalc);
260 // tlist3.AddToList(&pedphotcalc); // already done by MCalibrate Data
261
262 MWriteRootFile write(resname);
263
264 write.AddContainer("MGeomCam" , "RunHeaders");
265 write.AddContainer("MRawRunHeader" , "RunHeaders");
266 write.AddContainer("MSrcPosCam" , "RunHeaders");
267 write.AddContainer("MCalibrationChargeCam" , "RunHeaders");
268 write.AddContainer("MCalibrationQECam" , "RunHeaders");
269 // write.AddContainer("MPedPhotCam","RunHeaders"); // Attention, was in Events - Tree!!
270 write.AddContainer("MPedestalCam" , "RunHeaders");
271 write.AddContainer("MCalibrationRelTimeCam", "RunHeaders");
272
273 write.AddContainer("MCerPhotEvt" , "Events");
274 write.AddContainer("MRawEvtHeader" , "Events");
275 write.AddContainer("MBadPixelsCam" , "Events");
276 write.AddContainer("MPedPhotCam" , "Events");
277
278 tlist3.AddToList(&write);
279
280 // Create and execute eventloop
281 MEvtLoop evtloop3;
282 evtloop3.SetParList(&plist3);
283
284 cout << "*************************************************************" << endl;
285 cout << "*** COMPUTING DATA USING EXTRACTED SIGNAL (IN PHOTONS) ***" << endl;
286 cout << "*************************************************************" << endl;
287
288 if (!evtloop3.Eventloop())
289 return;
290 tlist3.PrintStatistics();
291
292}
293
294
295
Note: See TracBrowser for help on using the repository browser.