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

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