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 ,04/2004 <mailto:markus@ifae.es>
|
---|
20 | !
|
---|
21 | ! Copyright: MAGIC Software Development, 2000-2004
|
---|
22 | !
|
---|
23 | !
|
---|
24 | \* ======================================================================== */
|
---|
25 |
|
---|
26 | /////////////////////////////////////////////////////////////////////////////
|
---|
27 | //
|
---|
28 | // MJPedestal
|
---|
29 | //
|
---|
30 | /////////////////////////////////////////////////////////////////////////////
|
---|
31 | #include "MJPedestal.h"
|
---|
32 |
|
---|
33 | #include <TF1.h>
|
---|
34 | #include <TEnv.h>
|
---|
35 | #include <TFile.h>
|
---|
36 | #include <TLine.h>
|
---|
37 | #include <TLatex.h>
|
---|
38 | #include <TString.h>
|
---|
39 | #include <TCanvas.h>
|
---|
40 | #include <TSystem.h>
|
---|
41 | #include <TLegend.h>
|
---|
42 |
|
---|
43 | #include "MLog.h"
|
---|
44 | #include "MLogManip.h"
|
---|
45 |
|
---|
46 | #include "MTaskEnv.h"
|
---|
47 | #include "MSequence.h"
|
---|
48 | #include "MRunIter.h"
|
---|
49 | #include "MParList.h"
|
---|
50 | #include "MTaskList.h"
|
---|
51 | #include "MEvtLoop.h"
|
---|
52 | #include "MExtractor.h"
|
---|
53 |
|
---|
54 | #include "MStatusDisplay.h"
|
---|
55 |
|
---|
56 | #include "MGeomCam.h"
|
---|
57 | #include "MHCamera.h"
|
---|
58 | #include "MPedestalCam.h"
|
---|
59 | #include "MBadPixelsCam.h"
|
---|
60 |
|
---|
61 | #include "MReadMarsFile.h"
|
---|
62 | #include "MRawFileRead.h"
|
---|
63 | #include "MGeomApply.h"
|
---|
64 | #include "MBadPixelsMerge.h"
|
---|
65 | #include "MPedCalcPedRun.h"
|
---|
66 | #include "MPedCalcFromLoGain.h"
|
---|
67 |
|
---|
68 | ClassImp(MJPedestal);
|
---|
69 |
|
---|
70 | using namespace std;
|
---|
71 |
|
---|
72 | const Double_t MJPedestal::fgPedestalMin = 4.;
|
---|
73 | const Double_t MJPedestal::fgPedestalMax = 16.;
|
---|
74 | const Double_t MJPedestal::fgPedRmsMin = 0.;
|
---|
75 | const Double_t MJPedestal::fgPedRmsMax = 20.;
|
---|
76 |
|
---|
77 | const Float_t MJPedestal::fgRefPedClosedLids = 9.635;
|
---|
78 | const Float_t MJPedestal::fgRefPedExtraGalactic = 9.93;
|
---|
79 | const Float_t MJPedestal::fgRefPedGalactic = 10.03;
|
---|
80 | const Float_t MJPedestal::fgRefPedRmsClosedLidsInner = 1.7;
|
---|
81 | const Float_t MJPedestal::fgRefPedRmsExtraGalacticInner = 5.6;
|
---|
82 | const Float_t MJPedestal::fgRefPedRmsGalacticInner = 6.92;
|
---|
83 | const Float_t MJPedestal::fgRefPedRmsClosedLidsOuter = 1.7;
|
---|
84 | const Float_t MJPedestal::fgRefPedRmsExtraGalacticOuter = 3.35;
|
---|
85 | const Float_t MJPedestal::fgRefPedRmsGalacticOuter = 4.2;
|
---|
86 | // --------------------------------------------------------------------------
|
---|
87 | //
|
---|
88 | // Default constructor.
|
---|
89 | //
|
---|
90 | // Sets fRuns to 0, fExtractor to NULL, fDataCheck to kFALSE
|
---|
91 | //
|
---|
92 | MJPedestal::MJPedestal(const char *name, const char *title)
|
---|
93 | : fRuns(0), fExtractor(NULL), fDisplayType(kNormalDisplay),
|
---|
94 | fDataCheck(kFALSE), fUseData(kFALSE)
|
---|
95 | {
|
---|
96 | fName = name ? name : "MJPedestal";
|
---|
97 | fTitle = title ? title : "Tool to create a pedestal file (MPedestalCam)";
|
---|
98 | }
|
---|
99 |
|
---|
100 | const char* MJPedestal::GetOutputFile() const
|
---|
101 | {
|
---|
102 | if (fSequence.IsValid())
|
---|
103 | return Form("%s/pedest%06d.root", (const char*)fPathOut, fSequence.GetSequence());
|
---|
104 |
|
---|
105 | if (!fRuns)
|
---|
106 | return "";
|
---|
107 |
|
---|
108 | return Form("%s/%s-F0.root", (const char*)fPathOut, (const char*)fRuns->GetRunsAsFileName());
|
---|
109 | }
|
---|
110 |
|
---|
111 | Bool_t MJPedestal::ReadPedestalCam()
|
---|
112 | {
|
---|
113 | const TString fname = GetOutputFile();
|
---|
114 |
|
---|
115 | if (gSystem->AccessPathName(fname, kFileExists))
|
---|
116 | {
|
---|
117 | *fLog << warn << "Input file " << fname << " doesn't exist, will create it." << endl;
|
---|
118 | return kFALSE;
|
---|
119 | }
|
---|
120 |
|
---|
121 | *fLog << inf << "Reading from file: " << fname << endl;
|
---|
122 |
|
---|
123 | TFile file(fname, "READ");
|
---|
124 | if (fPedestalCam.Read()<=0)
|
---|
125 | {
|
---|
126 | *fLog << err << "Unable to read MPedestalCam from " << fname << endl;
|
---|
127 | return kFALSE;
|
---|
128 | }
|
---|
129 |
|
---|
130 | if (file.FindKey("MBadPixelsCam"))
|
---|
131 | {
|
---|
132 | MBadPixelsCam bad;
|
---|
133 | if (bad.Read()<=0)
|
---|
134 | {
|
---|
135 | *fLog << err << "Unable to read MBadPixelsCam from " << fname << endl;
|
---|
136 | return kFALSE;
|
---|
137 | }
|
---|
138 | fBadPixels.Merge(bad);
|
---|
139 | }
|
---|
140 |
|
---|
141 | if (fDisplay && !fDisplay->GetCanvas("Pedestals"))
|
---|
142 | fDisplay->Read();
|
---|
143 |
|
---|
144 | return kTRUE;
|
---|
145 | }
|
---|
146 |
|
---|
147 | void MJPedestal::DisplayResult(MParList &plist)
|
---|
148 | {
|
---|
149 | if (!fDisplay)
|
---|
150 | return;
|
---|
151 |
|
---|
152 | //
|
---|
153 | // Update display
|
---|
154 | //
|
---|
155 | TString title = fDisplay->GetTitle();
|
---|
156 | title += "-- Pedestal ";
|
---|
157 | if (fSequence.IsValid())
|
---|
158 | title += fSequence.GetName();
|
---|
159 | else
|
---|
160 | if (fRuns) // FIXME: What to do if an environmentfile was used?
|
---|
161 | title += fRuns->GetRunsAsString();
|
---|
162 | title += " --";
|
---|
163 | fDisplay->SetTitle(title);
|
---|
164 |
|
---|
165 | //
|
---|
166 | // Get container from list
|
---|
167 | //
|
---|
168 | MGeomCam &geomcam = *(MGeomCam*)plist.FindObject("MGeomCam");
|
---|
169 |
|
---|
170 | //
|
---|
171 | // Create container to display
|
---|
172 | //
|
---|
173 | MHCamera disp0(geomcam, "MPedestalCam;ped", "Mean Pedestal");
|
---|
174 | MHCamera disp1(geomcam, "MPedestalCam;RMS", "Pedestal RMS");
|
---|
175 |
|
---|
176 | disp0.SetCamContent(fPedestalCam, 0);
|
---|
177 | disp0.SetCamError (fPedestalCam, 1);
|
---|
178 |
|
---|
179 | disp1.SetCamContent(fPedestalCam, 2);
|
---|
180 | disp1.SetCamError (fPedestalCam, 3);
|
---|
181 |
|
---|
182 | disp0.SetYTitle("Pedestal [counts/slice]");
|
---|
183 | disp1.SetYTitle("RMS [counts/slice]");
|
---|
184 |
|
---|
185 | //
|
---|
186 | // Display data
|
---|
187 | //
|
---|
188 | TCanvas &c3 = fDisplay->AddTab("Pedestals");
|
---|
189 | c3.Divide(2,3);
|
---|
190 |
|
---|
191 | if (fDisplayType != kDataCheckDisplay)
|
---|
192 | {
|
---|
193 | disp0.CamDraw(c3, 1, 2, 1);
|
---|
194 | disp1.CamDraw(c3, 2, 2, 6);
|
---|
195 | return;
|
---|
196 | }
|
---|
197 |
|
---|
198 | c3.cd(1);
|
---|
199 | gPad->SetBorderMode(0);
|
---|
200 | gPad->SetTicks();
|
---|
201 | MHCamera *obj1=(MHCamera*)disp0.DrawCopy("hist");
|
---|
202 | //
|
---|
203 | // for the datacheck, fix the ranges!!
|
---|
204 | //
|
---|
205 | obj1->SetMinimum(fgPedestalMin);
|
---|
206 | obj1->SetMaximum(fgPedestalMax);
|
---|
207 |
|
---|
208 | //
|
---|
209 | // Set the datacheck sizes:
|
---|
210 | //
|
---|
211 | FixDataCheckHist((TH1D*)obj1);
|
---|
212 | //
|
---|
213 | // set reference lines
|
---|
214 | //
|
---|
215 | DisplayReferenceLines(obj1,0);
|
---|
216 |
|
---|
217 | //
|
---|
218 | // end reference lines
|
---|
219 | //
|
---|
220 | c3.cd(3);
|
---|
221 | gPad->SetBorderMode(0);
|
---|
222 | obj1->SetPrettyPalette();
|
---|
223 | obj1->Draw();
|
---|
224 |
|
---|
225 | c3.cd(5);
|
---|
226 | gPad->SetBorderMode(0);
|
---|
227 | gPad->SetTicks();
|
---|
228 | TH1D *obj2 = (TH1D*)obj1->Projection(obj1->GetName());
|
---|
229 | obj2->Draw();
|
---|
230 | obj2->SetBit(kCanDelete);
|
---|
231 | obj2->Fit("gaus","Q");
|
---|
232 | obj2->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
233 | //
|
---|
234 | // Set the datacheck sizes:
|
---|
235 | //
|
---|
236 | FixDataCheckHist(obj2);
|
---|
237 | obj2->SetStats(1);
|
---|
238 |
|
---|
239 | c3.cd(2);
|
---|
240 | gPad->SetBorderMode(0);
|
---|
241 | gPad->SetTicks();
|
---|
242 | MHCamera *obj3=(MHCamera*)disp1.DrawCopy("hist");
|
---|
243 | //
|
---|
244 | // for the datacheck, fix the ranges!!
|
---|
245 | //
|
---|
246 | obj3->SetMinimum(fgPedRmsMin);
|
---|
247 | obj3->SetMaximum(fgPedRmsMax);
|
---|
248 | //
|
---|
249 | // Set the datacheck sizes:
|
---|
250 | //
|
---|
251 | FixDataCheckHist((TH1D*)obj3);
|
---|
252 | //
|
---|
253 | // set reference lines
|
---|
254 | //
|
---|
255 | DisplayReferenceLines(obj1,1);
|
---|
256 |
|
---|
257 | c3.cd(4);
|
---|
258 | gPad->SetBorderMode(0);
|
---|
259 | obj3->SetPrettyPalette();
|
---|
260 | obj3->Draw();
|
---|
261 |
|
---|
262 | c3.cd(6);
|
---|
263 | gPad->SetBorderMode(0);
|
---|
264 |
|
---|
265 | if (geomcam.InheritsFrom("MGeomCamMagic"))
|
---|
266 | {
|
---|
267 | TArrayI inner(1);
|
---|
268 | inner[0] = 0;
|
---|
269 |
|
---|
270 | TArrayI outer(1);
|
---|
271 | outer[0] = 1;
|
---|
272 |
|
---|
273 | TArrayI s0(6);
|
---|
274 | s0[0] = 6;
|
---|
275 | s0[1] = 1;
|
---|
276 | s0[2] = 2;
|
---|
277 | s0[3] = 3;
|
---|
278 | s0[4] = 4;
|
---|
279 | s0[5] = 5;
|
---|
280 |
|
---|
281 | TArrayI s1(3);
|
---|
282 | s1[0] = 6;
|
---|
283 | s1[1] = 1;
|
---|
284 | s1[2] = 2;
|
---|
285 |
|
---|
286 | TArrayI s2(3);
|
---|
287 | s2[0] = 3;
|
---|
288 | s2[1] = 4;
|
---|
289 | s2[2] = 5;
|
---|
290 |
|
---|
291 | TVirtualPad *pad = gPad;
|
---|
292 | pad->Divide(2,1);
|
---|
293 |
|
---|
294 | TH1D *inout[2];
|
---|
295 | inout[0] = disp1.ProjectionS(s0, inner, "Inner");
|
---|
296 | inout[1] = disp1.ProjectionS(s0, outer, "Outer");
|
---|
297 | FixDataCheckHist(inout[0]);
|
---|
298 | FixDataCheckHist(inout[1]);
|
---|
299 |
|
---|
300 | inout[0]->SetTitle(Form("%s %s",disp1.GetTitle(),"Inner"));
|
---|
301 | inout[1]->SetTitle(Form("%s %s",disp1.GetTitle(),"Outer"));
|
---|
302 |
|
---|
303 |
|
---|
304 | for (int i=0; i<2; i++)
|
---|
305 | {
|
---|
306 | pad->cd(i+1);
|
---|
307 | gPad->SetBorderMode(0);
|
---|
308 | gPad->SetTicks();
|
---|
309 |
|
---|
310 | inout[i]->SetDirectory(NULL);
|
---|
311 | inout[i]->SetLineColor(kRed+i);
|
---|
312 | inout[i]->SetBit(kCanDelete);
|
---|
313 | inout[i]->Draw();
|
---|
314 | inout[i]->Fit("gaus", "Q");
|
---|
315 |
|
---|
316 | TLegend *leg2 = new TLegend(0.6,0.2,0.9,0.55);
|
---|
317 | leg2->SetHeader(inout[i]->GetName());
|
---|
318 | leg2->AddEntry(inout[i], inout[i]->GetName(), "l");
|
---|
319 |
|
---|
320 | //
|
---|
321 | // Display the outliers as dead and noisy pixels
|
---|
322 | //
|
---|
323 | DisplayOutliers(inout[i]);
|
---|
324 |
|
---|
325 | //
|
---|
326 | // Display the two half of the camera separately
|
---|
327 | //
|
---|
328 | TH1D *half[2];
|
---|
329 | half[0] = disp1.ProjectionS(s1, i==0 ? inner : outer , "Sector 6-1-2");
|
---|
330 | half[1] = disp1.ProjectionS(s2, i==0 ? inner : outer , "Sector 3-4-5");
|
---|
331 |
|
---|
332 | for (int j=0; j<2; j++)
|
---|
333 | {
|
---|
334 | half[j]->SetLineColor(kRed+i+2*j+1);
|
---|
335 | half[j]->SetDirectory(NULL);
|
---|
336 | half[j]->SetBit(kCanDelete);
|
---|
337 | half[j]->Draw("same");
|
---|
338 | leg2->AddEntry(half[j], half[j]->GetName(), "l");
|
---|
339 | }
|
---|
340 | leg2->Draw();
|
---|
341 | }
|
---|
342 | }
|
---|
343 | }
|
---|
344 |
|
---|
345 | void MJPedestal::DisplayReferenceLines(MHCamera *cam, const Int_t what) const
|
---|
346 | {
|
---|
347 |
|
---|
348 | Double_t x = cam->GetNbinsX();
|
---|
349 |
|
---|
350 | const MGeomCam *geom = cam->GetGeometry();
|
---|
351 |
|
---|
352 | if (geom->InheritsFrom("MGeomCamMagic"))
|
---|
353 | x = what ? 397 : cam->GetNbinsX();
|
---|
354 |
|
---|
355 | TLine line;
|
---|
356 | line.SetLineStyle(kDashed);
|
---|
357 | line.SetLineWidth(3);
|
---|
358 |
|
---|
359 | line.SetLineColor(kBlue);
|
---|
360 | TLine *l1 = line.DrawLine(0, what ? fgRefPedRmsGalacticInner : fgRefPedGalactic,
|
---|
361 | x, what ? fgRefPedRmsGalacticInner : fgRefPedGalactic);
|
---|
362 |
|
---|
363 | line.SetLineColor(kYellow);
|
---|
364 | TLine *l2 = line.DrawLine(0, what ? fgRefPedRmsExtraGalacticInner : fgRefPedExtraGalactic,
|
---|
365 | x, what ? fgRefPedRmsExtraGalacticInner : fgRefPedExtraGalactic);
|
---|
366 |
|
---|
367 | line.SetLineColor(kMagenta);
|
---|
368 | TLine *l3 = line.DrawLine(0, what ? fgRefPedRmsClosedLidsInner : fgRefPedClosedLids,
|
---|
369 | x, what ? fgRefPedRmsClosedLidsInner : fgRefPedClosedLids);
|
---|
370 |
|
---|
371 | if (geom->InheritsFrom("MGeomCamMagic"))
|
---|
372 | if (what)
|
---|
373 | {
|
---|
374 | const Double_t x2 = cam->GetNbinsX();
|
---|
375 |
|
---|
376 | line.SetLineColor(kBlue);
|
---|
377 | line.DrawLine(398, fgRefPedRmsGalacticOuter,
|
---|
378 | x2, fgRefPedRmsGalacticOuter);
|
---|
379 |
|
---|
380 | line.SetLineColor(kYellow);
|
---|
381 | line.DrawLine(398, fgRefPedRmsExtraGalacticOuter,
|
---|
382 | x2, fgRefPedRmsExtraGalacticOuter);
|
---|
383 |
|
---|
384 | line.SetLineColor(kMagenta);
|
---|
385 | line.DrawLine(398, fgRefPedRmsClosedLidsOuter,
|
---|
386 | x2, fgRefPedRmsClosedLidsOuter);
|
---|
387 | }
|
---|
388 |
|
---|
389 |
|
---|
390 | TLegend *leg = new TLegend(0.4,0.75,0.7,0.99);
|
---|
391 | leg->SetBit(kCanDelete);
|
---|
392 | leg->AddEntry(l1, "Galactic Source","l");
|
---|
393 | leg->AddEntry(l2, "Extra-Galactic Source","l");
|
---|
394 | leg->AddEntry(l3, "Closed Lids","l");
|
---|
395 | leg->Draw();
|
---|
396 | }
|
---|
397 |
|
---|
398 | void MJPedestal::DisplayOutliers(TH1D *hist) const
|
---|
399 | {
|
---|
400 | const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
|
---|
401 | const Float_t lolim = mean - 3.5*hist->GetFunction("gaus")->GetParameter(2);
|
---|
402 | const Float_t uplim = mean + 3.5*hist->GetFunction("gaus")->GetParameter(2);
|
---|
403 | const Stat_t dead = hist->Integral(0,hist->FindBin(lolim)-1);
|
---|
404 | const Stat_t noisy = hist->Integral(hist->FindBin(uplim)+1,hist->GetNbinsX()+1);
|
---|
405 |
|
---|
406 | TLatex deadtex;
|
---|
407 | deadtex.SetTextSize(0.06);
|
---|
408 | deadtex.DrawLatex(0.1,hist->GetBinContent(hist->GetMaximumBin())/1.1,Form("%3i dead pixels",(Int_t)dead));
|
---|
409 |
|
---|
410 | TLatex noisytex;
|
---|
411 | noisytex.SetTextSize(0.06);
|
---|
412 | noisytex.DrawLatex(0.1,hist->GetBinContent(hist->GetMaximumBin())/1.2,Form("%3i noisy pixels",(Int_t)noisy));
|
---|
413 | }
|
---|
414 |
|
---|
415 | void MJPedestal::FixDataCheckHist(TH1D *hist) const
|
---|
416 | {
|
---|
417 |
|
---|
418 | hist->SetDirectory(NULL);
|
---|
419 | hist->SetStats(0);
|
---|
420 |
|
---|
421 | //
|
---|
422 | // set the labels bigger
|
---|
423 | //
|
---|
424 | TAxis *xaxe = hist->GetXaxis();
|
---|
425 | TAxis *yaxe = hist->GetYaxis();
|
---|
426 |
|
---|
427 | xaxe->CenterTitle();
|
---|
428 | yaxe->CenterTitle();
|
---|
429 | xaxe->SetTitleSize(0.06);
|
---|
430 | yaxe->SetTitleSize(0.06);
|
---|
431 | xaxe->SetTitleOffset(0.8);
|
---|
432 | yaxe->SetTitleOffset(0.5);
|
---|
433 | xaxe->SetLabelSize(0.05);
|
---|
434 | yaxe->SetLabelSize(0.05);
|
---|
435 |
|
---|
436 | }
|
---|
437 | /*
|
---|
438 | Bool_t MJPedestal::WriteEventloop(MEvtLoop &evtloop) const
|
---|
439 | {
|
---|
440 | if (fOutputPath.IsNull())
|
---|
441 | return kTRUE;
|
---|
442 |
|
---|
443 | const TString oname(GetOutputFile());
|
---|
444 |
|
---|
445 | *fLog << inf << "Writing to file: " << oname << endl;
|
---|
446 |
|
---|
447 | TFile file(oname, fOverwrite?"RECREATE":"NEW", "File created by MJPedestal", 9);
|
---|
448 | if (!file.IsOpen())
|
---|
449 | {
|
---|
450 | *fLog << err << "ERROR - Couldn't open file " << oname << " for writing..." << endl;
|
---|
451 | return kFALSE;
|
---|
452 | }
|
---|
453 |
|
---|
454 | if (evtloop.Write(fName)<=0)
|
---|
455 | {
|
---|
456 | *fLog << err << "Unable to write MEvtloop to " << oname << endl;
|
---|
457 | return kFALSE;
|
---|
458 | }
|
---|
459 |
|
---|
460 | return kTRUE;
|
---|
461 | }
|
---|
462 | */
|
---|
463 |
|
---|
464 | // --------------------------------------------------------------------------
|
---|
465 | //
|
---|
466 | // The following resource options are available:
|
---|
467 | // Prefix.DataCheckDisplay: Yes, No
|
---|
468 | // Prefix.DataCheck: Yes, No
|
---|
469 | // Prefix.UseData: Yes, No
|
---|
470 | //
|
---|
471 | Bool_t MJPedestal::CheckEnvLocal()
|
---|
472 | {
|
---|
473 | if (HasEnv("DataCheckDisplay"))
|
---|
474 | fDisplayType = GetEnv("DataCheckDisplay", kFALSE) ? kDataCheckDisplay : kNormalDisplay;
|
---|
475 |
|
---|
476 | SetDataCheck(GetEnv("DataCheck", fDataCheck));
|
---|
477 | SetUseData(GetEnv("UseData", fUseData));
|
---|
478 |
|
---|
479 | return kTRUE;
|
---|
480 | }
|
---|
481 |
|
---|
482 | Bool_t MJPedestal::WriteResult()
|
---|
483 | {
|
---|
484 | if (fPathOut.IsNull())
|
---|
485 | return kTRUE;
|
---|
486 |
|
---|
487 | const TString oname(GetOutputFile());
|
---|
488 |
|
---|
489 | *fLog << inf << "Writing to file: " << oname << endl;
|
---|
490 |
|
---|
491 | TFile file(oname, "UPDATE", "File created by MJPedestal", 9);
|
---|
492 | if (!file.IsOpen())
|
---|
493 | {
|
---|
494 | *fLog << err << "ERROR - Couldn't open file " << oname << " for writing..." << endl;
|
---|
495 | return kFALSE;
|
---|
496 | }
|
---|
497 |
|
---|
498 | if (fDisplay && fDisplay->Write()<=0)
|
---|
499 | {
|
---|
500 | *fLog << err << "Unable to write MStatusDisplay to " << oname << endl;
|
---|
501 | return kFALSE;
|
---|
502 | }
|
---|
503 |
|
---|
504 | if (fPedestalCam.Write()<=0)
|
---|
505 | {
|
---|
506 | *fLog << err << "Unable to write MPedestalCam to " << oname << endl;
|
---|
507 | return kFALSE;
|
---|
508 | }
|
---|
509 |
|
---|
510 | if (fBadPixels.Write()<=0)
|
---|
511 | {
|
---|
512 | *fLog << err << "Unable to write MBadPixelsCam to " << oname << endl;
|
---|
513 | return kFALSE;
|
---|
514 | }
|
---|
515 |
|
---|
516 | return kTRUE;
|
---|
517 | }
|
---|
518 |
|
---|
519 | Bool_t MJPedestal::Process()
|
---|
520 | {
|
---|
521 | if (!ReadPedestalCam())
|
---|
522 | return ProcessFile();
|
---|
523 |
|
---|
524 | return kTRUE;
|
---|
525 | }
|
---|
526 |
|
---|
527 | Bool_t MJPedestal::ProcessFile()
|
---|
528 | {
|
---|
529 | if (!fSequence.IsValid())
|
---|
530 | {
|
---|
531 | if (!fRuns)
|
---|
532 | {
|
---|
533 | *fLog << err << "Neither AddRuns nor SetSequence nor SetEnv was called... abort." << endl;
|
---|
534 | return kFALSE;
|
---|
535 | }
|
---|
536 | if (fRuns && fRuns->GetNumRuns() != fRuns->GetNumEntries())
|
---|
537 | {
|
---|
538 | *fLog << err << "Number of files found doesn't match number of runs... abort." << endl;
|
---|
539 | return kFALSE;
|
---|
540 | }
|
---|
541 | }
|
---|
542 |
|
---|
543 | //if (!CheckEnv())
|
---|
544 | // return kFALSE;
|
---|
545 |
|
---|
546 | CheckEnv();
|
---|
547 |
|
---|
548 | // --------------------------------------------------------------------------------
|
---|
549 |
|
---|
550 | const TString type = fUseData ? "data" : "pedestal";
|
---|
551 |
|
---|
552 | *fLog << inf;
|
---|
553 | fLog->Separator(GetDescriptor());
|
---|
554 | *fLog << "Calculate MPedestalCam from " << type << "-runs ";
|
---|
555 | if (fSequence.IsValid())
|
---|
556 | *fLog << fSequence.GetName() << endl;
|
---|
557 | else
|
---|
558 | if (fRuns)
|
---|
559 | *fLog << fRuns->GetRunsAsString() << endl;
|
---|
560 | else
|
---|
561 | *fLog << "in Resource File." << endl;
|
---|
562 | *fLog << endl;
|
---|
563 |
|
---|
564 | // --------------------------------------------------------------------------------
|
---|
565 |
|
---|
566 | MParList plist;
|
---|
567 | MTaskList tlist;
|
---|
568 | plist.AddToList(&tlist);
|
---|
569 | plist.AddToList(this); // take care of fDisplay!
|
---|
570 |
|
---|
571 | MReadMarsFile read("Events");
|
---|
572 | MRawFileRead rawread(NULL);
|
---|
573 |
|
---|
574 | MDirIter iter;
|
---|
575 | if (fSequence.IsValid())
|
---|
576 | {
|
---|
577 | const Int_t n0 = fUseData ? fSequence.SetupDatRuns(iter, fPathData, fDataCheck) : fSequence.SetupPedRuns(iter, fPathData, fDataCheck);
|
---|
578 | const Int_t n1 = fUseData ? fSequence.GetNumDatRuns() : fSequence.GetNumPedRuns();
|
---|
579 | if (n0==0)
|
---|
580 | {
|
---|
581 | *fLog << err << "ERROR - No " << type << " input files of sequence found in " << fPathData << endl;
|
---|
582 | return kFALSE;
|
---|
583 | }
|
---|
584 | if (n0!=n1)
|
---|
585 | {
|
---|
586 | *fLog << err << "ERROR - Number of " << type << " files found (" << n0 << ") in " << fPathData << " doesn't match number of files in sequence (" << n1 << ")" << endl;
|
---|
587 | return kFALSE;
|
---|
588 | }
|
---|
589 | }
|
---|
590 |
|
---|
591 | if (fDataCheck)
|
---|
592 | {
|
---|
593 | if (fRuns || fSequence.IsValid())
|
---|
594 | rawread.AddFiles(fSequence.IsValid() ? iter : *fRuns);
|
---|
595 | tlist.AddToList(&rawread);
|
---|
596 | }
|
---|
597 | else
|
---|
598 | {
|
---|
599 | read.DisableAutoScheme();
|
---|
600 | if (fRuns || fSequence.IsValid())
|
---|
601 | read.AddFiles(fSequence.IsValid() ? iter : *fRuns);
|
---|
602 | tlist.AddToList(&read);
|
---|
603 | }
|
---|
604 |
|
---|
605 | // Setup Tasklist
|
---|
606 | plist.AddToList(&fPedestalCam);
|
---|
607 | plist.AddToList(&fBadPixels);
|
---|
608 |
|
---|
609 | MGeomApply geomapl;
|
---|
610 | MBadPixelsMerge merge(&fBadPixels);
|
---|
611 |
|
---|
612 | MPedCalcPedRun pedcalc;
|
---|
613 | MPedCalcFromLoGain pedlogain;
|
---|
614 | pedlogain.SetPedestalUpdate(kFALSE);
|
---|
615 |
|
---|
616 | MTaskEnv taskenv("ExtractPedestal");
|
---|
617 | if (IsUseData())
|
---|
618 | taskenv.SetDefault(&pedlogain);
|
---|
619 | else
|
---|
620 | taskenv.SetDefault(&pedcalc);
|
---|
621 |
|
---|
622 | if (fExtractor)
|
---|
623 | {
|
---|
624 | if (IsUseData())
|
---|
625 | pedlogain.SetWindowSize(12,(Int_t)fExtractor->GetNumHiGainSamples());
|
---|
626 | else
|
---|
627 | {
|
---|
628 | pedcalc.SetWindowSize((Int_t)fExtractor->GetNumHiGainSamples());
|
---|
629 | pedcalc.SetRange(fExtractor->GetHiGainFirst(), fExtractor->GetHiGainLast());
|
---|
630 | }
|
---|
631 | }
|
---|
632 | /*
|
---|
633 | else
|
---|
634 | {
|
---|
635 | *fLog << warn << GetDescriptor();
|
---|
636 | *fLog << ": No extractor has been chosen, take default number of FADC samples " << endl;
|
---|
637 | }
|
---|
638 | */
|
---|
639 |
|
---|
640 | tlist.AddToList(&geomapl);
|
---|
641 | tlist.AddToList(&merge);
|
---|
642 | tlist.AddToList(&taskenv);
|
---|
643 |
|
---|
644 | //
|
---|
645 | // Execute the eventloop
|
---|
646 | //
|
---|
647 | MEvtLoop evtloop(fName);
|
---|
648 | evtloop.SetParList(&plist);
|
---|
649 | evtloop.SetDisplay(fDisplay);
|
---|
650 | evtloop.SetLogStream(fLog);
|
---|
651 | if (!SetupEnv(evtloop))
|
---|
652 | return kFALSE;
|
---|
653 |
|
---|
654 | // if (!WriteEventloop(evtloop))
|
---|
655 | // return kFALSE;
|
---|
656 |
|
---|
657 | // Execute first analysis
|
---|
658 | if (!evtloop.Eventloop(fMaxEvents))
|
---|
659 | {
|
---|
660 | *fLog << err << GetDescriptor() << ": Failed." << endl;
|
---|
661 | return kFALSE;
|
---|
662 | }
|
---|
663 |
|
---|
664 | tlist.PrintStatistics();
|
---|
665 |
|
---|
666 | DisplayResult(plist);
|
---|
667 |
|
---|
668 | if (!WriteResult())
|
---|
669 | return kFALSE;
|
---|
670 |
|
---|
671 | *fLog << all << GetDescriptor() << ": Done." << endl;
|
---|
672 | *fLog << endl << endl;
|
---|
673 |
|
---|
674 | return kTRUE;
|
---|
675 | }
|
---|