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

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