source: trunk/MagicSoft/Mars/mjobs/MJPedestal.cc@ 7023

Last change on this file since 7023 was 7013, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 38.6 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! Author(s): Markus Gaug, 4/2004 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2005
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MJPedestal
29//
30// Resource file entries are case sensitive!
31//
32/////////////////////////////////////////////////////////////////////////////
33#include "MJPedestal.h"
34
35// C/C++ includes
36#include <fstream>
37
38// root classes
39#include <TF1.h>
40#include <TEnv.h>
41#include <TFile.h>
42#include <TLine.h>
43#include <TLatex.h>
44#include <TString.h>
45#include <TCanvas.h>
46#include <TSystem.h>
47#include <TLegend.h>
48#include <TPad.h>
49#include <TEnv.h>
50#include <TH2F.h>
51
52// mars core
53#include "MLog.h"
54#include "MLogManip.h"
55
56#include "MTaskEnv.h"
57#include "MSequence.h"
58#include "MRunIter.h"
59#include "MParList.h"
60#include "MTaskList.h"
61#include "MEvtLoop.h"
62
63#include "MStatusDisplay.h"
64#include "MFEventSelector.h"
65
66// Other basic classes
67#include "MExtractTimeAndCharge.h"
68
69// parameter containers
70#include "MGeomCam.h"
71#include "MHCamera.h"
72#include "MPedestalPix.h"
73
74#include "MCalibrationPedCam.h"
75#include "MCalibrationPix.h"
76#include "MHPedestalPix.h"
77#include "MHCalibrationPulseTimeCam.h"
78#include "MCalibrationPulseTimeCam.h"
79
80// tasks
81#include "MReadMarsFile.h"
82#include "MRawFileRead.h"
83#include "MRawEvtData.h"
84#include "MGeomApply.h"
85#include "MTriggerPatternDecode.h"
86#include "MBadPixelsMerge.h"
87#include "MFillH.h"
88#include "MPedCalcPedRun.h"
89#include "MPedCalcFromLoGain.h"
90#include "MFTriggerPattern.h"
91
92ClassImp(MJPedestal);
93
94using namespace std;
95
96const TString MJPedestal::fgReferenceFile = "mjobs/pedestalref.rc";
97const TString MJPedestal::fgBadPixelsFile = "mjobs/badpixels_0_559.rc";
98const Float_t MJPedestal::fgExtractWinLeft = 2.5;
99const Float_t MJPedestal::fgExtractWinRight = 4.5;
100
101// --------------------------------------------------------------------------
102//
103// Default constructor.
104//
105// Sets:
106// - fExtractor to NULL,
107// - fExtractType to kUsePedRun
108// - fStorage to Normal Storage
109// - fExtractorResolution to kFALSE
110//
111MJPedestal::MJPedestal(const char *name, const char *title)
112 : fExtractor(NULL), fDisplayType(kDisplayDataCheck),
113 fExtractType(kUsePedRun), fExtractionType(kFundamental), fIsUseHists(kFALSE)
114{
115 fName = name ? name : "MJPedestal";
116 fTitle = title ? title : "Tool to create a pedestal file (MPedestalCam)";
117
118 SetUsePedRun();
119 SetPathIn("");
120 SetReferenceFile();
121 SetBadPixelsFile();
122
123 SetExtractWinLeft();
124 SetExtractWinRight();
125 //
126 // Default references for case that no value references file is there
127 // (should not occur)
128 //
129
130 fPedestalMin = 4.;
131 fPedestalMax = 16.;
132 fPedRmsMin = 0.;
133 fPedRmsMax = 20.;
134 fRefPedClosedLids = 9.635;
135 fRefPedExtraGalactic = 9.93;
136 fRefPedGalactic = 10.03;
137 fRefPedRmsClosedLidsInner = 1.7;
138 fRefPedRmsExtraGalacticInner = 5.6;
139 fRefPedRmsGalacticInner = 6.92;
140 fRefPedRmsClosedLidsOuter = 1.7;
141 fRefPedRmsExtraGalacticOuter = 3.35;
142 fRefPedRmsGalacticOuter = 4.2;
143}
144
145MJPedestal::~MJPedestal()
146{
147 if (fExtractor)
148 delete fExtractor;
149}
150
151const char* MJPedestal::GetOutputFileName() const
152{
153 return Form("pedest%08d.root", fSequence.GetSequence());
154}
155
156MExtractor *MJPedestal::ReadCalibration() const
157{
158 const TString fname = Form("%s/calib%08d.root", fPathIn.Data(), fSequence.GetSequence());
159
160 *fLog << inf << "Reading extractor from file: " << fname << endl;
161
162 TFile file(fname, "READ");
163 if (!file.IsOpen())
164 {
165 *fLog << err << dbginf << "ERROR - Could not open file " << fname << endl;
166 return NULL;
167 }
168
169 TObject *o = file.Get("ExtractSignal");
170 if (o && !o->InheritsFrom(MExtractor::Class()))
171 {
172 *fLog << err << dbginf << "ERROR - ExtractSignal read from " << fname << " doesn't inherit from MExtractor!" << endl;
173 return NULL;
174 }
175
176 return o ? (MExtractor*)o->Clone("ExtractSignal") : NULL;
177}
178
179//---------------------------------------------------------------------------------
180//
181// Display the results.
182// If Display type "kDataCheck" was chosen, also the reference lines are displayed.
183//
184void MJPedestal::DisplayResult(MParList &plist)
185{
186 if (!fDisplay)
187 return;
188
189 //
190 // Update display
191 //
192 TString title = fDisplay->GetTitle();
193 title += "-- Pedestal ";
194 title += fSequence.GetName();
195 title += " --";
196 fDisplay->SetTitle(title);
197
198 //
199 // Get container from list
200 //
201 MGeomCam &geomcam = *(MGeomCam*)plist.FindObject("MGeomCam");
202 // MCalibrationPedCam &calpedcam = *(MCalibrationPedCam*)plist.FindObject("MCalibrationPedCam");
203
204 //
205 // Create container to display
206 //
207 MHCamera disp0 (geomcam, "MPedestalCam;ped", "Mean Pedestal");
208 MHCamera disp1 (geomcam, "MPedestalCam;RMS", "Pedestal RMS");
209 MHCamera disp2 (geomcam, "MCalibPedCam;histmean", "Mean Pedestal (Hist.)");
210 MHCamera disp3 (geomcam, "MCalibPedCam;histsigma", "Pedestal RMS (Hist.)");
211 MHCamera disp4 (geomcam, "MCalibPedCam;ped", "Mean Pedestal");
212 MHCamera disp5 (geomcam, "MCalibPedCam;RMS", "Pedestal RMS");
213 MHCamera disp6 (geomcam, "MCalibDiffCam;ped", "Diff. Mean Pedestal (Hist.)");
214 MHCamera disp7 (geomcam, "MCalibDiffCam;RMS", "Diff. Pedestal RMS (Hist.)");
215 MHCamera disp8 (geomcam, "MCalibDiffCam;ped", "Diff. Mean Pedestal");
216 MHCamera disp9 (geomcam, "MCalibDiffCam;AbsRMS", "Diff. Abs. Pedestal RMS");
217 MHCamera disp10(geomcam, "MCalibDiffCam;RelRMS", "Diff. Rel. Pedestal RMS");
218
219 disp0.SetCamContent(fPedestalCamOut, 0);
220 disp0.SetCamError (fPedestalCamOut, 1);
221
222 disp1.SetCamContent(fPedestalCamOut, 2);
223 disp1.SetCamError (fPedestalCamOut, 3);
224
225 /*
226 if (fIsUseHists)
227 {
228 disp2.SetCamContent(calpedcam, 0);
229 disp2.SetCamError (calpedcam, 1);
230
231 disp3.SetCamContent(calpedcam, 2);
232 disp3.SetCamError (calpedcam, 3);
233
234 disp4.SetCamContent(calpedcam, 5);
235 disp4.SetCamError (calpedcam, 6);
236
237 disp5.SetCamContent(calpedcam, 7);
238 disp5.SetCamError (calpedcam, 8);
239
240 for (UInt_t i=0;i<geomcam.GetNumPixels();i++)
241 {
242
243 MPedestalPix &ped = fPedestalCamOut[i];
244 MCalibrationPix &hist = calpedcam [i];
245 MBadPixelsPix &bad = fBadPixels[i];
246
247 if (bad.IsUnsuitable())
248 continue;
249
250 disp6.Fill(i,ped.GetPedestal()-hist.GetHiGainMean());
251 disp6.SetUsed(i);
252
253 disp7.Fill(i,hist.GetHiGainSigma()-ped.GetPedestalRms());
254 if (TMath::Abs(ped.GetPedestalRms()-hist.GetHiGainSigma()) < 4.0)
255 disp7.SetUsed(i);
256
257 disp8.Fill(i,ped.GetPedestal()-hist.GetLoGainMean());
258 disp8.SetUsed(i);
259
260 disp9.Fill(i,hist.GetLoGainSigma()-ped.GetPedestalRms());
261 if (TMath::Abs(hist.GetLoGainSigma() - ped.GetPedestalRms()) < 4.0)
262 disp9.SetUsed(i);
263 }
264 }
265 */
266
267 if (fExtractionType!=kFundamental/*fExtractorResolution*/)
268 {
269 for (UInt_t i=0;i<geomcam.GetNumPixels();i++)
270 {
271
272 MPedestalPix &pedo = fPedestalCamOut[i];
273 MPedestalPix &pedi = fPedestalCamIn[i];
274 MBadPixelsPix &bad = fBadPixels[i];
275
276 if (bad.IsUnsuitable())
277 continue;
278
279 const Float_t diff = pedo.GetPedestalRms()-pedi.GetPedestalRms();
280 const Float_t sum = 0.5*(pedo.GetPedestalRms()+pedi.GetPedestalRms());
281
282 disp9.Fill(i,pedo.GetPedestalRms()-pedi.GetPedestalRms());
283 if (pedo.IsValid() && pedi.IsValid())
284 disp9.SetUsed(i);
285
286 disp10.Fill(i,sum == 0. ? 0. : diff/sum);
287 if (pedo.IsValid() && pedi.IsValid() && sum != 0.)
288 disp10.SetUsed(i);
289 }
290 }
291
292 disp0.SetYTitle("P [cts/slice]");
293 disp1.SetYTitle("P_{rms} [cts/slice]");
294 disp2.SetYTitle("Hist. Mean [cts/slice]");
295 disp3.SetYTitle("Hist. Sigma [cts/slice]");
296 disp4.SetYTitle("Calc. Mean [cts/slice]");
297 disp5.SetYTitle("Calc. RMS [cts/slice]");
298 disp6.SetYTitle("Diff. Mean [cts/slice]");
299 disp7.SetYTitle("Diff. RMS [cts/slice]");
300 disp8.SetYTitle("Diff. Mean [cts/slice]");
301 disp9.SetYTitle("Abs.Diff.RMS [cts/slice]");
302 disp10.SetYTitle("Rel.Diff.RMS [1]");
303
304 //
305 // Display data
306 //
307 if (fDisplayType != kDisplayDataCheck && fExtractionType==kFundamental/*fExtractorResolution*/)
308 {
309 TCanvas &c3 = fDisplay->AddTab("Pedestals");
310 c3.Divide(2,3);
311
312 disp0.CamDraw(c3, 1, 2, 1);
313 disp1.CamDraw(c3, 2, 2, 6);
314 return;
315 }
316
317#if 0
318 if (fIsUseHists)
319 {
320
321 TCanvas &c3 = fDisplay->AddTab("Extractor Hist.");
322 c3.Divide(2,3);
323
324 disp2.CamDraw(c3, 1, 2, 1);
325 disp3.CamDraw(c3, 2, 2, 5);
326
327 TCanvas &c4 = fDisplay->AddTab("Extractor Calc.");
328 c4.Divide(2,3);
329
330 disp4.CamDraw(c4, 1, 2, 1);
331 disp5.CamDraw(c4, 2, 2, 5);
332
333 /*
334 TCanvas &c5 = fDisplay->AddTab("Difference Hist.");
335 c5.Divide(2,3);
336
337 disp6.CamDraw(c5, 1, 2, 1);
338 disp7.CamDraw(c5, 2, 2, 5);
339 */
340
341 TCanvas &c6 = fDisplay->AddTab("Difference Calc.");
342 c6.Divide(2,3);
343
344 disp8.CamDraw(c6, 1, 2, 1);
345 disp9.CamDraw(c6, 2, 2, 5);
346 return;
347 }
348#endif
349 if (fDisplayType == kDisplayDataCheck)
350 {
351
352 TCanvas &c3 = fDisplay->AddTab(fExtractionType!=kFundamental/*fExtractorResolution*/ ? "PedExtrd" : "Ped");
353 c3.Divide(2,3);
354
355 c3.cd(1);
356 gPad->SetBorderMode(0);
357 gPad->SetTicks();
358 MHCamera *obj1=(MHCamera*)disp0.DrawCopy("hist");
359 //
360 // for the datacheck, fix the ranges!!
361 //
362 if (fExtractionType==kFundamental/*!fExtractorResolution*/)
363 {
364 obj1->SetMinimum(fPedestalMin);
365 obj1->SetMaximum(fPedestalMax);
366 }
367 //
368 // Set the datacheck sizes:
369 //
370 FixDataCheckHist((TH1D*)obj1);
371 //
372 // set reference lines
373 //
374 DisplayReferenceLines(obj1,0);
375 //
376 // end reference lines
377 //
378 c3.cd(3);
379 gPad->SetBorderMode(0);
380 obj1->SetPrettyPalette();
381 obj1->Draw();
382
383 c3.cd(5);
384 gPad->SetBorderMode(0);
385 gPad->SetTicks();
386 TH1D *obj2 = (TH1D*)obj1->Projection(obj1->GetName());
387 obj2->Draw();
388 obj2->SetBit(kCanDelete);
389 obj2->Fit("gaus","Q");
390 obj2->GetFunction("gaus")->SetLineColor(kYellow);
391 //
392 // Set the datacheck sizes:
393 //
394 FixDataCheckHist(obj2);
395 obj2->SetStats(1);
396
397 c3.cd(2);
398 gPad->SetBorderMode(0);
399 gPad->SetTicks();
400 MHCamera *obj3=(MHCamera*)disp1.DrawCopy("hist");
401 //
402 // for the datacheck, fix the ranges!!
403 //
404 obj3->SetMinimum(fPedRmsMin);
405 obj3->SetMaximum(fPedRmsMax);
406 //
407 // Set the datacheck sizes:
408 //
409 FixDataCheckHist((TH1D*)obj3);
410 //
411 // set reference lines
412 //
413 DisplayReferenceLines(obj3,1);
414
415 c3.cd(4);
416 gPad->SetBorderMode(0);
417 obj3->SetPrettyPalette();
418 obj3->Draw();
419
420 c3.cd(6);
421 gPad->SetBorderMode(0);
422
423 if (geomcam.InheritsFrom("MGeomCamMagic"))
424 {
425 TArrayI inner(1);
426 inner[0] = 0;
427
428 TArrayI outer(1);
429 outer[0] = 1;
430
431 TArrayI s0(6);
432 s0[0] = 6;
433 s0[1] = 1;
434 s0[2] = 2;
435 s0[3] = 3;
436 s0[4] = 4;
437 s0[5] = 5;
438
439 TArrayI s1(3);
440 s1[0] = 6;
441 s1[1] = 1;
442 s1[2] = 2;
443
444 TArrayI s2(3);
445 s2[0] = 3;
446 s2[1] = 4;
447 s2[2] = 5;
448
449 TVirtualPad *pad = gPad;
450 pad->Divide(2,1);
451
452 TH1D *inout[2];
453 inout[0] = disp1.ProjectionS(s0, inner, "Inner");
454 inout[1] = disp1.ProjectionS(s0, outer, "Outer");
455 FixDataCheckHist(inout[0]);
456 FixDataCheckHist(inout[1]);
457
458 inout[0]->SetTitle(Form("%s %s",disp1.GetTitle(),"Inner"));
459 inout[1]->SetTitle(Form("%s %s",disp1.GetTitle(),"Outer"));
460
461
462 for (int i=0; i<2; i++)
463 {
464 pad->cd(i+1);
465 gPad->SetBorderMode(0);
466 gPad->SetTicks();
467
468 inout[i]->SetDirectory(NULL);
469 inout[i]->SetLineColor(kRed+i);
470 inout[i]->SetBit(kCanDelete);
471 inout[i]->Draw();
472 inout[i]->Fit("gaus", "Q");
473
474 TLegend *leg2 = new TLegend(0.6,0.2,0.9,0.55);
475 leg2->SetHeader(inout[i]->GetName());
476 leg2->AddEntry(inout[i], inout[i]->GetName(), "l");
477
478 //
479 // Display the outliers as dead and noisy pixels
480 //
481 DisplayOutliers(inout[i]);
482
483 //
484 // Display the two half of the camera separately
485 //
486 TH1D *half[2];
487 half[0] = disp1.ProjectionS(s1, i==0 ? inner : outer , "Sector 6-1-2");
488 half[1] = disp1.ProjectionS(s2, i==0 ? inner : outer , "Sector 3-4-5");
489
490 for (int j=0; j<2; j++)
491 {
492 half[j]->SetLineColor(kRed+i+2*j+1);
493 half[j]->SetDirectory(NULL);
494 half[j]->SetBit(kCanDelete);
495 half[j]->Draw("same");
496 leg2->AddEntry(half[j], half[j]->GetName(), "l");
497 }
498 leg2->Draw();
499 delete leg2;
500 }
501 return;
502 }
503 }
504
505 if (fExtractionType!=kFundamental/*fExtractorResolution*/)
506 {
507
508 TCanvas &c3 = fDisplay->AddTab(fExtractionType==kWithExtractor?"PedExtrd":"PedRndm");
509 c3.Divide(2,3);
510
511 disp0.CamDraw(c3, 1, 2, 1);
512 disp1.CamDraw(c3, 2, 2, 6);
513
514 TCanvas &c13 = fDisplay->AddTab(fExtractionType==kWithExtractor?"DiffExtrd":"DiffRndm");
515 c13.Divide(2,3);
516
517 disp9.CamDraw(c13, 1, 2, 5);
518 disp10.CamDraw(c13, 2, 2, 5);
519 return;
520 }
521}
522
523
524void MJPedestal::DisplayReferenceLines(MHCamera *cam, const Int_t what) const
525{
526
527 Double_t x = cam->GetNbinsX();
528
529 const MGeomCam *geom = cam->GetGeometry();
530
531 if (geom->InheritsFrom("MGeomCamMagic"))
532 x = what ? 397 : cam->GetNbinsX();
533
534 TLine line;
535 line.SetLineStyle(kDashed);
536 line.SetLineWidth(3);
537 line.SetLineColor(kBlue);
538
539 TLegend *leg = new TLegend(0.6,0.75,0.9,0.99);
540 leg->SetBit(kCanDelete);
541
542 if (fExtractionType==kWithExtractorRndm && !(what))
543 {
544 TLine *l0 = line.DrawLine(0,0.,cam->GetNbinsX(),0.);
545 l0->SetBit(kCanDelete);
546 leg->AddEntry(l0, "Reference","l");
547 leg->Draw();
548 return;
549 }
550
551 line.SetLineColor(kBlue);
552 TLine *l1 = line.DrawLine(0, what ? fRefPedRmsGalacticInner : fRefPedGalactic,
553 x, what ? fRefPedRmsGalacticInner : fRefPedGalactic);
554 l1->SetBit(kCanDelete);
555 line.SetLineColor(kYellow);
556 TLine *l2 = line.DrawLine(0, what ? fRefPedRmsExtraGalacticInner : fRefPedExtraGalactic,
557 x, what ? fRefPedRmsExtraGalacticInner : fRefPedExtraGalactic);
558 l2->SetBit(kCanDelete);
559 line.SetLineColor(kMagenta);
560 TLine *l3 = line.DrawLine(0, what ? fRefPedRmsClosedLidsInner : fRefPedClosedLids,
561 x, what ? fRefPedRmsClosedLidsInner : fRefPedClosedLids);
562 l3->SetBit(kCanDelete);
563
564 if (geom->InheritsFrom("MGeomCamMagic"))
565 if (what)
566 {
567 const Double_t x2 = cam->GetNbinsX();
568
569 line.SetLineColor(kBlue);
570 line.DrawLine(398, fRefPedRmsGalacticOuter,
571 x2, fRefPedRmsGalacticOuter);
572
573 line.SetLineColor(kYellow);
574 line.DrawLine(398, fRefPedRmsExtraGalacticOuter,
575 x2, fRefPedRmsExtraGalacticOuter);
576
577 line.SetLineColor(kMagenta);
578 line.DrawLine(398, fRefPedRmsClosedLidsOuter,
579 x2, fRefPedRmsClosedLidsOuter);
580 }
581
582
583 leg->AddEntry(l1, "Galactic Source","l");
584 leg->AddEntry(l2, "Extra-Galactic Source","l");
585 leg->AddEntry(l3, "Closed Lids","l");
586 leg->Draw();
587}
588
589void MJPedestal::DisplayOutliers(TH1D *hist) const
590{
591 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
592 const Float_t lolim = mean - 3.5*hist->GetFunction("gaus")->GetParameter(2);
593 const Float_t uplim = mean + 3.5*hist->GetFunction("gaus")->GetParameter(2);
594 const Stat_t dead = hist->Integral(0,hist->FindBin(lolim)-1);
595 const Stat_t noisy = hist->Integral(hist->FindBin(uplim)+1,hist->GetNbinsX()+1);
596
597 TLatex deadtex;
598 deadtex.SetTextSize(0.06);
599 deadtex.DrawLatex(0.1,hist->GetBinContent(hist->GetMaximumBin())/1.1,Form("%3i dead pixels",(Int_t)dead));
600
601 TLatex noisytex;
602 noisytex.SetTextSize(0.06);
603 noisytex.DrawLatex(0.1,hist->GetBinContent(hist->GetMaximumBin())/1.2,Form("%3i noisy pixels",(Int_t)noisy));
604}
605
606void MJPedestal::FixDataCheckHist(TH1D *hist) const
607{
608 hist->SetDirectory(NULL);
609 hist->SetStats(0);
610
611 //
612 // set the labels bigger
613 //
614 TAxis *xaxe = hist->GetXaxis();
615 TAxis *yaxe = hist->GetYaxis();
616
617 xaxe->CenterTitle();
618 yaxe->CenterTitle();
619 xaxe->SetTitleSize(0.06);
620 yaxe->SetTitleSize(0.06);
621 xaxe->SetTitleOffset(0.8);
622 yaxe->SetTitleOffset(0.5);
623 xaxe->SetLabelSize(0.05);
624 yaxe->SetLabelSize(0.05);
625}
626
627/*
628Bool_t MJPedestal::WriteEventloop(MEvtLoop &evtloop) const
629{
630 if (fOutputPath.IsNull())
631 return kTRUE;
632
633 const TString oname(GetOutputFile());
634
635 *fLog << inf << "Writing to file: " << oname << endl;
636
637 TFile file(oname, fOverwrite?"RECREATE":"NEW", "File created by MJPedestal", 9);
638 if (!file.IsOpen())
639 {
640 *fLog << err << "ERROR - Couldn't open file " << oname << " for writing..." << endl;
641 return kFALSE;
642 }
643
644 if (evtloop.Write(fName)<=0)
645 {
646 *fLog << err << "Unable to write MEvtloop to " << oname << endl;
647 return kFALSE;
648 }
649
650 return kTRUE;
651}
652*/
653
654void MJPedestal::SetExtractor(MExtractor* ext)
655{
656 if (ext)
657 {
658 if (fExtractor)
659 delete fExtractor;
660 fExtractor = ext ? (MExtractor*)ext->Clone(ext->GetName()) : NULL;
661 }
662 else
663 fExtractor = 0;
664}
665
666// --------------------------------------------------------------------------
667//
668// Read the following values from resource file:
669//
670// PedestalMin
671// PedestalMax
672//
673// PedRmsMin
674// PedRmsMax
675//
676// RefPedClosedLids
677// RefPedExtraGalactic
678// RefPedGalactic
679//
680// RefPedRmsClosedLidsInner
681// RefPedRmsExtraGalacticInner
682// RefPedRmsGalacticInner
683// RefPedRmsClosedLidsOuter
684// RefPedRmsExtraGalacticOuter
685// RefPedRmsGalacticOuter
686//
687void MJPedestal::ReadReferenceFile()
688{
689 TEnv refenv(fReferenceFile);
690
691 fPedestalMin = refenv.GetValue("PedestalMin",fPedestalMin);
692 fPedestalMax = refenv.GetValue("PedestalMax",fPedestalMax);
693 fPedRmsMin = refenv.GetValue("PedRmsMin",fPedRmsMin);
694 fPedRmsMax = refenv.GetValue("PedRmsMax",fPedRmsMax);
695 fRefPedClosedLids = refenv.GetValue("RefPedClosedLids",fRefPedClosedLids);
696 fRefPedExtraGalactic = refenv.GetValue("RefPedExtraGalactic",fRefPedExtraGalactic);
697 fRefPedGalactic = refenv.GetValue("RefPedGalactic",fRefPedGalactic);
698 fRefPedRmsClosedLidsInner = refenv.GetValue("RefPedRmsClosedLidsInner",fRefPedRmsClosedLidsInner);
699 fRefPedRmsExtraGalacticInner = refenv.GetValue("RefPedRmsExtraGalacticInner",fRefPedRmsExtraGalacticInner);
700 fRefPedRmsGalacticInner = refenv.GetValue("RefPedRmsGalacticInner",fRefPedRmsGalacticInner);
701 fRefPedRmsClosedLidsOuter = refenv.GetValue("RefPedRmsClosedLidsOuter",fRefPedRmsClosedLidsOuter);
702 fRefPedRmsExtraGalacticOuter = refenv.GetValue("RefPedRmsExtraGalacticOuter",fRefPedRmsExtraGalacticOuter);
703 fRefPedRmsGalacticOuter = refenv.GetValue("RefPedRmsGalacticOuter",fRefPedRmsGalacticOuter);
704}
705
706// --------------------------------------------------------------------------
707//
708// The following resource options are available:
709//
710// Do a datacheck run (read raw-data and enable display)
711// Prefix.DataCheck: Yes, No <default>
712//
713// Setup display type
714// Prefix.Display: normal <default>, datacheck, none
715//
716// Use cosmic data instead of pedestal data (DatRuns)
717// Prefix.UseData: Yes, No <default>
718//
719// Write an output file with pedestals and status-display
720// Prefix.DisableOutput: Yes, No <default>
721//
722// Name of a file containing reference values (see ReadReferenceFile)
723// Prefix.ReferenceFile: filename
724// (see ReadReferenceFile)
725//
726Bool_t MJPedestal::CheckEnvLocal()
727{
728 if (HasEnv("Display"))
729 {
730 TString type = GetEnv("Display", "normal");
731 type.ToLower();
732 if (type==(TString)"normal")
733 fDisplayType = kDisplayNormal;
734 if (type==(TString)"datacheck")
735 fDisplayType = kDisplayDataCheck;
736 if (type==(TString)"none")
737 fDisplayType = kDisplayNone;
738 }
739
740
741 SetExtractWinLeft (GetEnv("ExtractWinLeft", fExtractWinLeft ));
742 SetExtractWinRight(GetEnv("ExtractWinRight", fExtractWinRight));
743
744 if (!MJCalib::CheckEnvLocal())
745 return kFALSE;
746
747 if (HasEnv("UseData"))
748 fExtractType = GetEnv("UseData",kFALSE) ? kUseData : kUsePedRun;
749
750 if (IsUseMC() && fExtractType==kUseData)
751 {
752 // The reason is, that the standard data files contains empty
753 // (untriggered) events. If we would loop over the default 500
754 // first events of the data file you would calculate the
755 // pedestal from only some single events...
756 *fLog << inf;
757 *fLog << "Sorry, you cannot extract the starting pedestal from the first" << endl;
758 *fLog << "events in your data files... using pedestal file instead. The" << endl;
759 *fLog << "result should not differ..." << endl;
760 fExtractType = kUsePedRun;
761 }
762
763
764 if (HasEnv("UseHists"))
765 if (GetEnv("UseHists",kFALSE))
766 fIsUseHists = kTRUE;
767
768 SetNoStorage(GetEnv("DisableOutput", IsNoStorage()));
769
770 MTaskEnv tenv("ExtractSignal");
771 tenv.SetDefault(fExtractor);
772
773 if (tenv.ReadEnv(*GetEnv(), GetEnvPrefix()+".ExtractSignal", GetEnvDebug()>2)==kERROR)
774 return kFALSE;
775
776 if (fExtractor==tenv.GetTask())
777 return kTRUE;
778
779 if (!tenv.GetTask()->InheritsFrom(MExtractor::Class()))
780 {
781 *fLog << err << "ERROR: ExtractSignal from resource file doesn't inherit from MExtractor.... abort." << endl;
782 return kFALSE;
783 }
784
785 SetExtractor((MExtractor*)tenv.GetTask());
786
787 fBadPixelsFile = GetEnv("BadPixelsFile",fBadPixelsFile.Data());
788 fReferenceFile = GetEnv("ReferenceFile",fReferenceFile.Data());
789 ReadReferenceFile();
790
791 return kTRUE;
792}
793
794//---------------------------------------------------------------------------------
795//
796// Try to write the created MPedestalCam in the output file.
797// If possible, also an MBadPixelsCam and the corresponding display is written.
798//
799// In case of Storage type "kNoStorage" or if any of the containers
800// cannot be written, return kFALSE, otherwise kTRUE.
801//
802Bool_t MJPedestal::WriteResult()
803{
804 if (IsNoStorage())
805 return kTRUE;
806
807 TObjArray cont;
808
809 cont.Add(&fPedestalCamOut);
810 cont.Add(&fBadPixels);
811
812 return WriteContainer(cont, GetOutputFileName(), fOverwrite?"RECREATE":"NEW");
813}
814
815Bool_t MJPedestal::Process()
816{
817 if (!fSequence.IsValid())
818 {
819 *fLog << err << "ERROR - Sequence invalid..." << endl;
820 return kFALSE;
821 }
822
823 *fLog << inf;
824 fLog->Separator(GetDescriptor());
825 *fLog << "Calculate pedestal from Sequence #";
826 *fLog << fSequence.GetSequence() << endl << endl;
827
828 if (!CheckEnv())
829 return kFALSE;
830
831 // --------------------------------------------------------------------------------
832
833 const TString type = IsUseData() ? "data" : "pedestal";
834
835 *fLog << inf;
836 fLog->Separator(GetDescriptor());
837 *fLog << "Calculate MPedestalCam from " << type << "-runs ";
838 *fLog << fSequence.GetName() << endl;
839 *fLog << endl;
840
841 // --------------------------------------------------------------------------------
842
843 MParList plist;
844 MTaskList tlist;
845 plist.AddToList(&tlist);
846 plist.AddToList(this); // take care of fDisplay!
847
848 MReadMarsFile read("Events");
849 MRawFileRead rawread(NULL);
850
851 MDirIter iter;
852 if (fSequence.IsValid())
853 {
854 const Int_t n0 = IsUseData()
855 ? fSequence.SetupDatRuns(iter, fPathData, IsUseRawData())
856 : fSequence.SetupPedRuns(iter, fPathData, IsUseRawData());
857 const Int_t n1 = IsUseData()
858 ? fSequence.GetNumDatRuns()
859 : fSequence.GetNumPedRuns();
860 if (n0==0)
861 {
862 *fLog << err << "ERROR - No " << type << " input files of sequence found in " << (fPathData.IsNull()?"<default>":fPathData.Data()) << endl;
863 return kFALSE;
864 }
865 if (n0!=n1)
866 {
867 *fLog << err << "ERROR - Number of " << type << " files found (" << n0 << ") in " << (fPathData.IsNull()?"<default>":fPathData.Data()) << " doesn't match number of files in sequence (" << n1 << ")" << endl;
868 if (fLog->GetDebugLevel()>4)
869 {
870 *fLog << dbg << "Files which are searched:" << endl;
871 iter.Print();
872 }
873 return kFALSE;
874 }
875 }
876
877 if (IsUseRawData())
878 {
879 rawread.AddFiles(iter);
880 tlist.AddToList(&rawread);
881 }
882 else
883 {
884 read.DisableAutoScheme();
885 read.AddFiles(iter);
886 tlist.AddToList(&read);
887 }
888
889 // Setup Tasklist
890 plist.AddToList(&fPedestalCamOut);
891 plist.AddToList(&fBadPixels);
892
893 //
894 // Read bad pixels from outside
895 //
896 if (!fBadPixelsFile.IsNull())
897 {
898 *fLog << inf << "Excluding: " << fBadPixelsFile << endl;
899 ifstream fin(fBadPixelsFile.Data());
900 fBadPixels.AsciiRead((istream&)fin);
901 }
902
903 MGeomApply geomapl;
904 MBadPixelsMerge merge(&fBadPixels);
905
906 MPedCalcPedRun pedcalc;
907 pedcalc.SetPedestalUpdate(kFALSE);
908
909 MPedCalcFromLoGain pedlogain;
910 pedlogain.SetPedestalUpdate(kFALSE);
911
912 MHPedestalCam hpedcam;
913 if (fExtractionType != kFundamental)
914 hpedcam.SetRenorm(kTRUE);
915 // fPedestalHist.SetRenorm(kTRUE);
916 // fPedestalHist.SetPedestalsOut(&fPedestalCamOut);
917 hpedcam.SetPedestalsOut(&fPedestalCamOut);
918
919 MPedestalCam pedinter;
920 pedinter.SetName("MPedestalCamIntermediate");
921
922 MFillH fillped(&hpedcam, "MPedestalCamIntermediate", "FillPedCam");
923 // MFillH fillped(&fPedestalHist, "MPedestalCamIntermediate", "FillPedCam");
924 MFillH fillpul("MHCalibrationPulseTimeCam", "MRawEvtData", "FillPulseTime");
925 fillped.SetBit(MFillH::kDoNotDisplay);
926 fillpul.SetBit(MFillH::kDoNotDisplay);
927
928 tlist.AddToList(&geomapl);
929 tlist.AddToList(&merge);
930
931 if (!fPathIn.IsNull())
932 {
933 fExtractor = ReadCalibration();
934 if (!fExtractor)
935 return kFALSE;
936
937 *fLog << all;
938 *fLog << underline << "Signal Extractor found in calibration file:" << endl;
939 fExtractor->Print();
940 *fLog << endl;
941 }
942
943 // This will make that for data with version less than 5, where
944 // trigger patterns were not yet correct, all the events in the real
945 // data file will be processed. In any case there are no interleaved
946 // calibration events in such data, so this is fine.
947 MTriggerPatternDecode decode;
948 MFTriggerPattern fcalib("CalibFilter");
949 fcalib.SetDefault(kFALSE);
950 fcalib.RequireCalibration();
951 fcalib.SetInverted();
952
953 if (fIsPulsePosCheck)
954 {
955 fillpul.SetFilter(&fcalib);
956 tlist.AddToList(&decode);
957 tlist.AddToList(&fcalib);
958 tlist.AddToList(&fillpul);
959 }
960
961 // ----------------------------------------------------------------------
962 // Now we make sure, that in all cases the ranges are setup correctly
963 // ----------------------------------------------------------------------
964 MTaskEnv taskenv("ExtractPedestal");
965 switch (fExtractType)
966 {
967 case kUsePedRun:
968 // In case other than 'fundamental' second argument is obsolete
969 // pedcalc.SetExtractWindow(0,14); // kUsePedRun (take default from class)
970 taskenv.SetDefault(&pedcalc);
971 tlist.AddToList(&taskenv);
972 break;
973
974 case kUseData:
975 // In case other than 'fundamental' second argument is obsolete
976 // pedlogain.SetExtractWindow(15,14); // kUseData (take default from class)
977 taskenv.SetDefault(&pedlogain);
978 tlist.AddToList(&taskenv);
979 break;
980 }
981
982 if (fIsUseHists)
983 {
984 if (fExtractor->InheritsFrom("MExtractTimeAndCharge"))
985 hpedcam.SetFitStart(-10.*((MExtractTimeAndCharge*)fExtractor)->GetWindowSizeHiGain());
986
987 pedcalc.SetIntermediateStorage();
988 pedlogain.SetIntermediateStorage();
989 plist.AddToList(&pedinter);
990 plist.AddToList(&hpedcam);
991 // plist.AddToList(&fPedestalHist);
992 tlist.AddToList(&fillped);
993 }
994
995 pedcalc.SetPedestalsIn(&fPedestalCamIn);
996 pedlogain.SetPedestalsIn(&fPedestalCamIn);
997
998 pedcalc.SetPedestalsInter(&pedinter);
999 pedlogain.SetPedestalsInter(&pedinter);
1000 pedcalc.SetPedestalsOut(&fPedestalCamOut);
1001 pedlogain.SetPedestalsOut(&fPedestalCamOut);
1002
1003 // kFundamental
1004 if (fExtractor)
1005 {
1006
1007 if (fExtractionType!=kFundamental)
1008 {
1009 pedcalc.SetRandomCalculation(fExtractionType==kWithExtractorRndm);
1010 pedlogain.SetRandomCalculation(fExtractionType==kWithExtractorRndm);
1011
1012 pedcalc.SetExtractor((MExtractTimeAndCharge*)fExtractor);
1013 pedlogain.SetExtractor((MExtractTimeAndCharge*)fExtractor);
1014 }
1015
1016 if (fExtractor->InheritsFrom("MExtractTimeAndCharge"))
1017 {
1018
1019 const Float_t f = 0.1+fExtractor->GetHiGainFirst();
1020 const Int_t win = ((MExtractTimeAndCharge*)fExtractor)->GetWindowSizeHiGain();
1021 pedcalc.SetExtractWindow((Int_t)f, win);
1022 pedlogain.SetExtractWindow((Int_t)(15+f), win);
1023
1024 }
1025 else
1026 {
1027 const Float_t f = 0.1+fExtractor->GetHiGainFirst();
1028 const Float_t n = 0.1+fExtractor->GetNumHiGainSamples();
1029 pedcalc.SetExtractWindow((Int_t)f, (Int_t)n);
1030 pedlogain.SetExtractWindow((Int_t)(15+f), (Int_t)n);
1031
1032 if (fExtractionType!=kFundamental)
1033 {
1034 *fLog << inf;
1035 *fLog << "Signal extractor doesn't inherit from MExtractTimeAndCharge..." << endl;
1036 *fLog << " --> falling back to fundamental pedestal extraction." << endl;
1037 fExtractionType=kFundamental;
1038 }
1039 }
1040
1041
1042 }
1043 else
1044 {
1045 *fLog << warn << GetDescriptor() << ": WARNING - No extractor has been handed over! " << endl;
1046 *fLog << "Taking default window for pedestal extraction. The calculated pedestal RMS" << endl;
1047 *fLog << "will probably not match with future pedestal RMS' from different extraction" << endl;
1048 *fLog << "windows." << endl;
1049 }
1050
1051 /*
1052 switch (fExtractType)
1053 {
1054 case kUseData:
1055 *fLog << all << "TYPE: USEDATA " << fExtractor << endl;
1056 taskenv.SetDefault(&pedlogain);
1057 tlist.AddToList(&taskenv);
1058
1059 if (!SetupExtractor(plist, pedlogain))
1060 {
1061 *fLog << all << "SETTING TO: " << fExtractor << " " << fExtractor->GetNumHiGainSamples() << endl;
1062 fExtractor->Print();
1063 pedlogain.SetExtractWindow(15, (Int_t)TMath::Nint(fExtractor->GetNumHiGainSamples()));
1064 }
1065 break;
1066
1067 case kUsePedRun:
1068 *fLog << all << "TYPE: USEPEDRUN " << fExtractor << endl;
1069 taskenv.SetDefault(&pedcalc);
1070 tlist.AddToList(&taskenv);
1071
1072 if (!SetupExtractor(plist, pedcalc))
1073 pedcalc.SetExtractWindow(fExtractor->GetHiGainFirst(), TMath::Nint(fExtractor->GetNumHiGainSamples()));
1074 break;
1075
1076 case kUseHists:
1077 if (!fExtractor)
1078 {
1079 *fLog << err << GetDescriptor() << " - ERROR: ";
1080 *fLog << "Extraction Type is kUseHists, but no extractor was set" << endl;
1081 return kFALSE;
1082 }
1083
1084 tlist.AddToList(fExtractor);
1085 tlist.AddToList(&fillped);
1086 break;
1087 } */
1088
1089 /*
1090 if (!fPathIn.IsNull())
1091 {
1092 delete fExtractor;
1093 fExtractor = 0;
1094 }
1095 */
1096
1097
1098 //
1099 // Execute the eventloop
1100 //
1101 MEvtLoop evtloop(fName);
1102 evtloop.SetParList(&plist);
1103 evtloop.SetDisplay(fDisplay);
1104 evtloop.SetLogStream(fLog);
1105 if (!SetupEnv(evtloop))
1106 return kFALSE;
1107
1108 // if (!WriteEventloop(evtloop))
1109 // return kFALSE;
1110
1111 // Execute first analysis
1112 if (!evtloop.Eventloop(fMaxEvents))
1113 {
1114 *fLog << err << GetDescriptor() << ": Failed." << endl;
1115 return kFALSE;
1116 }
1117
1118 tlist.PrintStatistics();
1119
1120 if (fIsPixelCheck)
1121 {
1122 MHPedestalCam *hcam = (MHPedestalCam*)plist.FindObject("MHPedestalCam");
1123 if (hcam)
1124 {
1125 MHPedestalPix &pix1 = (MHPedestalPix&)(*hcam)[fCheckedPixId];
1126 pix1.DrawClone("");
1127 }
1128 }
1129
1130 if (fDisplayType!=kDisplayNone)
1131 DisplayResult(plist);
1132
1133 if (!WriteResult())
1134 return kFALSE;
1135
1136 if (fIsPulsePosCheck)
1137 {
1138
1139 Int_t numhigainsamples = 0;
1140 Int_t numlogainsamples = 0;
1141 Float_t meanpulsetime = 0.;
1142 Float_t rmspulsetime = 0.;
1143
1144 if (IsUseMC())
1145 {
1146 //
1147 // FIXME:
1148 // The MC cannot run over the first 2000 pedestal events since almost all
1149 // events are empty, therefore a pulse pos. check is not possible, either.
1150 // For the moment, have to fix the problem hardcoded...
1151 //
1152 // MMcEvt *evt = (MMcEvt*)plist.FindObject("MMcEvt");
1153 // const Float_t meanpulsetime = evt->GetFadcTimeJitter();
1154 meanpulsetime = 4.5;
1155 rmspulsetime = 1.0;
1156
1157 *fLog << all << "Mean pulse time (MC): " << meanpulsetime << "+-" << rmspulsetime << endl;
1158
1159 numhigainsamples = 15;
1160 numlogainsamples = 15;
1161
1162 }
1163 else
1164 {
1165 MHCalibrationPulseTimeCam *hcam = (MHCalibrationPulseTimeCam*)plist.FindObject("MHCalibrationPulseTimeCam");
1166 if (fIsPixelCheck)
1167 {
1168 hcam->DrawClone();
1169 gPad->SaveAs(Form("%s/PulsePosTest_all.root",fPathOut.Data()));
1170 MHCalibrationPix &pix = (*hcam)[fCheckedPixId];
1171 pix.DrawClone();
1172 gPad->SaveAs(Form("%s/PulsePosTest_Pixel%04d.root",fPathOut.Data(),fCheckedPixId));
1173 }
1174
1175 MCalibrationPulseTimeCam *cam = (MCalibrationPulseTimeCam*)plist.FindObject("MCalibrationPulseTimeCam");
1176 if (!cam)
1177 {
1178 *fLog << err << "Could not determine mean pulse position, abort... " << endl;
1179 return kFALSE;
1180 }
1181
1182 meanpulsetime = cam->GetAverageArea(0).GetHiGainMean();
1183 rmspulsetime = cam->GetAverageArea(0).GetHiGainRms();
1184
1185 *fLog << all << "Mean pulse time (cosmics): " << meanpulsetime << "+-" << rmspulsetime << endl;
1186
1187 MRawEvtData *data = (MRawEvtData*)plist.FindObject("MRawEvtData");
1188
1189 numhigainsamples = data->GetNumHiGainSamples();
1190 numlogainsamples = data->GetNumLoGainSamples();
1191 }
1192
1193 //
1194 // Get the ranges for the new extractor setting
1195 //
1196 const Int_t newfirst = TMath::Nint(meanpulsetime-fExtractWinLeft);
1197 Int_t wshigain = fExtractor->InheritsFrom("MExtractTimeAndCharge")
1198 ? ((MExtractTimeAndCharge*)fExtractor)->GetWindowSizeHiGain()
1199 : 6;
1200
1201 if (wshigain > 6)
1202 wshigain = 6;
1203
1204 Int_t wslogain = fExtractor->InheritsFrom("MExtractTimeAndCharge")
1205 ? ((MExtractTimeAndCharge*)fExtractor)->GetWindowSizeLoGain()
1206 : 6;
1207 if (wslogain > 4)
1208 wslogain = 4;
1209
1210 const Int_t newlast = TMath::Nint(meanpulsetime+fExtractWinRight);
1211
1212 *fLog << all << underline
1213 << "Try to set new range limits: ("<<newfirst<<","<<newlast<<"+"<<wshigain
1214 <<","<<newfirst-1<<","<<newlast<<"+"<<wslogain<<")"<<endl;
1215 //
1216 // Check the ranges for the new extractor setting
1217 //
1218 if (newfirst < 0)
1219 {
1220 *fLog << err << "Pulse is too much to the left, cannot go below 0! " << endl;
1221 return kFALSE;
1222 }
1223 if (newlast+wshigain > numhigainsamples+numlogainsamples-1)
1224 {
1225 *fLog << err << "Pulse is too much to the right, cannot go beyond limits: "
1226 << numhigainsamples << "+" << numlogainsamples << "-1" << endl;
1227 *fLog << " Cannot extract at all! ... " << endl;
1228 return kFALSE;
1229 }
1230 if (newlast+wslogain > numlogainsamples)
1231 {
1232 *fLog << err << "Pulse is too much to the right, cannot go beyond logain limits! " << endl;
1233 *fLog << endl;
1234 *fLog << "Try to use a different extractor (e.g. with a window size of only 4 sl.) or:" << endl;
1235 *fLog << "Set the limit to a lower value (callisto.rc: line 329): " << endl;
1236 *fLog << " MJPedestalY2:ExtractWinRight: 4.5 " << endl;
1237 *fLog << "(ATTENTION, you will lose late cosmics pulses!)" << endl;
1238 *fLog << endl;
1239 return kFALSE;
1240 }
1241 //
1242 // Set and store the new ranges
1243 //
1244 fExtractor->SetRange(newfirst,newlast+wshigain,
1245 newfirst>0?newfirst-1:newfirst,numlogainsamples-1);
1246 }
1247
1248
1249 *fLog << all << GetDescriptor() << ": Done." << endl;
1250 *fLog << endl << endl;
1251
1252 return kTRUE;
1253}
Note: See TracBrowser for help on using the repository browser.