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