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

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