source: trunk/Mars/fact/analysis/mc/callisto.C@ 17667

Last change on this file since 17667 was 17641, checked in by ftemme, 11 years ago
Added the calculation of Disp and ThetaSq in the MC-Star macro and commented the lightpulser calibration in the MC-Callisto macro out
File size: 21.4 KB
Line 
1using namespace std;
2
3#include <sstream>
4#include <iostream>
5
6#include "TH1F.h"
7#include "TFile.h"
8#include "TStyle.h"
9#include "TGraph.h"
10#include "TLine.h"
11
12#include "MDrsCalibration.h"
13//#include "MLogManip.h"
14#include "MExtralgoSpline.h"
15#include "MSequence.h"
16#include "MStatusArray.h"
17#include "MHCamera.h"
18#include "MJob.h"
19#include "MWriteRootFile.h"
20#include "MHCamera.h"
21#include "MBadPixelsCam.h"
22#include "MBadPixelsPix.h"
23#include "MDirIter.h"
24#include "MTaskList.h"
25#include "MFDataPhrase.h"
26#include "MArrayF.h"
27#include "MBadPixelsTreat.h"
28#include "MCalibrateDrsTimes.h"
29#include "MHSectorVsTime.h"
30#include "MHCamEvent.h"
31#include "MExtractTimeAndChargeSpline.h"
32#include "MFillH.h"
33#include "MDrsCalibApply.h"
34#include "MGeomApply.h"
35#include "MContinue.h"
36#include "MRawFitsRead.h"
37#include "MEvtLoop.h"
38#include "MParList.h"
39#include "MStatusDisplay.h"
40#include "MDrsCalibrationTime.h"
41#include "MH3.h"
42#include "MGeomCamFACT.h"
43#include "MCalibrateFact.h"
44#include "MParameters.h"
45#include "MWriteAsciiFile.h"
46
47/* Maybe you wanna use this macro like this:
48 *
49 * 0.) ---- call root ----
50 * root -b
51 *
52 * 1.) ---- compile the stuff ----
53 * .L fact/analysis/callisto_buildable_no_sequence_file.C++
54 * <read a lot of warnings>
55 *
56 * 2.) ---- you can call it then ----
57 * Therefore you need to specify all the paths ... see below.
58 *
59 * When you wanna call the stuff directly from the bash make sure to
60 * escape the bracets and quotes correctly.
61 *
62 * your can do:
63 * root -b -q callisto_buildable_no_sequence_file.C++'("path1","path2",...)'
64 * or:
65 * root -b -q callisto_buildable_no_sequence_file.C++(\"path1\",\"$HOME\",...)
66 * using bash enviroment variables like $HOME is not possible in the upper variant.
67 */
68
69//using namespace std;
70
71int callisto(
72 const char *drs_calib_300_path="/fhgfs/groups/app/fact/mc_test/testingdrsfile/test300samples.drs.fits",
73 const char *pedestal_file_path="/fhgfs/groups/app/fact/mc_test/ceresfitstest/mcFilesForTests/mcNsbPedestal/00000001.001_P_MonteCarlo000_Events.fits",
74 const char *data_file_path="/fhgfs/groups/app/fact/mc_test/ceresfitstest/mcFilesForTests/mcProton/00000003.387_D_MonteCarlo010_Events.fits",
75
76 const char *root_file_output_path = "/fhgfs/groups/app/callisto_star_test/callisto_for_mc_test_output/callisto.root",
77 const char *status_display_output_path = "/fhgfs/groups/app/callisto_star_test/callisto_for_mc_test_output/callisto_status_display.root",
78 const char *status_display_title = "callisto_status_display")
79{
80
81 // ======================================================
82
83 // true: Display correctly mapped pixels in the camera displays
84 // but the value-vs-index plot is in software/spiral indices
85 // false: Display pixels in hardware/linear indices,
86 // but the order is the camera display is distorted
87 bool usemap = true;
88
89 // map file to use (get that from La Palma!)
90 const char *map = usemap ? "/fhgfs/groups/app/fact/resources/monte_carlo_FACTmap.txt" : NULL;
91
92 Bool_t maximum = kTRUE;
93
94 //const char *lp_template = maximum ?
95 // "/cm/shared/apps/fact/Mars_svn_LP/template-lp-extractor-maximum.root" :
96 // "/cm/shared/apps/fact/Mars_svn_LP/template-lp-extractor-leading-edge.root";
97
98 const char *lp_template = "/cm/shared/apps/fact/Mars_Trunk/template-lp.root";
99
100 //const char *pulse_template = "/cm/shared/apps/fact/Mars_svn_LP/template-pulse.root";
101
102 const char *pulse_template = "/cm/shared/apps/fact/Mars_Trunk/template-pulse.root";
103
104 // ------------------------------------------------------
105
106 // Calib: 51 / 90 / 197 (20% TH)
107 // Data: 52 / 64 / 104 (20% TH)
108
109 // Extraction range in slices. It will always(!) contain the full range
110 // of integration
111 const int first_slice = 20; // 10ns
112 const int last_slice = 250; // 125ns
113
114 // Note that rise and fall time mean different things whether you use IntegralFixed or IntegraRel:
115 //
116 // IntegralFixed:
117 // * fRiseTime: Number of slices left from arrival time
118 // * fFallTime: Number of slices right from arrival time
119 // IntegralRel:
120 // * fRiseTime: Number of slices left from maximum time
121 // * fFallTime: Number of slices right from maximum time
122 //
123 const int rise_time_cal = maximum ? 40 : 10; // was 13; 5ns
124 const int fall_time_cal = maximum ? 120 : 160; // was 23; 80ns
125
126 const int rise_time_dat = maximum ? 10 : 2; // was 13; was 10; 1ns
127 const int fall_time_dat = maximum ? 40 : 48; // was 23; was 40; 24ns
128
129 // Extraction type: Extract integral and half leading edge
130
131 const MExtralgoSpline::ExtractionType_t type = maximum ? (MExtralgoSpline::kIntegralRel) : (MExtralgoSpline::kIntegralFixed);
132 //const int type = MExtralgoSpline::kIntegralFixed;
133
134
135 const double heighttm = 0.5; // IntegralAbs { 1.5pe * 9.6mV/pe } / IntegralRel { 0.5 }
136
137 Long_t max = 0; // All
138 Long_t max3 = max; // Pedestal Rndm
139 Long_t max4 = max; // Pedestal Ext
140
141 // ======================================================
142
143 if (map && gSystem->AccessPathName(map, kFileExists))
144 {
145 gLog << err << "ERROR - Cannot access mapping file '" << map << "'" << endl;
146 return 1;
147 }
148
149 TString datfile = TString(data_file_path);
150 TString drsfile = TString(drs_calib_300_path);
151 TString pedfile = TString(pedestal_file_path);
152
153 gLog.Separator("Callisto");
154 gLog << all;
155 gLog << "Data File : " << datfile << "\n";
156 gLog << "DRS calib 300: " << drsfile << '\n';
157
158 MDrsCalibration drscalib300;
159 if (!drscalib300.ReadFits(drsfile.Data())) {
160 gLog << err << "ERROR - Cannot access drscallib300 file '" << drsfile << "'" << endl;
161 return 5;
162 }
163 gLog << all;
164 gLog << "Pedestal file: " << pedfile << '\n';
165
166 gLog << "root_file_output_path: " << root_file_output_path << endl;
167 gLog << "status_display_output_path: " << status_display_output_path << endl;
168 gLog << "status_display_title: " << status_display_title << endl;
169
170 // ------------------------------------------------------
171 MStatusArray arrt, arrp;
172
173 TFile ft(lp_template);
174 if (arrt.Read()<=0)
175 {
176 gLog << err << "ERROR - Reading LP template from " << lp_template << endl;
177 return 100;
178 }
179
180 MHCamera *lpref = (MHCamera*)arrt.FindObjectInCanvas("ExtCalSig;avg", "MHCamera", "Cam");
181 if (!lpref)
182 {
183 gLog << err << "ERROR - LP Template not found in " << lp_template << endl;
184 return 101;
185 }
186 lpref->SetDirectory(0);
187
188 // MHCamera *gain = (MHCamera*)arrt.FindObjectInCanvas("gain", "MHCamera", "Gain");
189 // if (!gain)
190 // {
191 // gLog << err << "ERROR - Gain not found in " << lp_template << endl;
192 // return 101;
193 // }
194 // gain->SetDirectory(0);
195
196 TFile fp(pulse_template);
197 if (arrp.Read()<=0)
198 {
199 gLog << err << "ERROR - Reading Pulse template from " << pulse_template << endl;
200 return 102;
201 }
202
203 TH1F *hpulse = (TH1F*)arrp.FindObjectInCanvas("hPixelEdgeMean0_0", "TH1F", "cgpPixelPulses0");
204 if (!hpulse)
205 {
206 gLog << err << "ERROR - Pulse Template not found in " << pulse_template << endl;
207 return 103;
208 }
209 hpulse->SetDirectory(0);
210 // ======================================================
211
212 MStatusDisplay *d = new MStatusDisplay;
213
214 MBadPixelsCam badpixels;
215 badpixels.InitSize(1440);
216 /*
217 badpixels[ 424].SetUnsuitable(MBadPixelsPix::kUnsuitable);
218 badpixels[ 583].SetUnsuitable(MBadPixelsPix::kUnsuitable);
219 badpixels[ 830].SetUnsuitable(MBadPixelsPix::kUnsuitable);
220 badpixels[ 923].SetUnsuitable(MBadPixelsPix::kUnsuitable);
221 badpixels[1208].SetUnsuitable(MBadPixelsPix::kUnsuitable);
222 badpixels[1399].SetUnsuitable(MBadPixelsPix::kUnsuitable);
223 */
224 // Twin pixel
225 // 113
226 // 115
227 // 354
228 // 423
229 // 1195
230 // 1393
231
232 MDrsCalibrationTime timecam;
233
234 // Plot the trigger pattern rates vs. run-number
235 MH3 hrate("MRawRunHeader.GetFileID", "MRawEvtHeader.GetTriggerID&0xff00");
236 hrate.SetWeight("1./TMath::Max(MRawRunHeader.GetRunLength,1)");
237 hrate.SetName("Rate");
238 hrate.SetTitle("Event rate [Hz];File Id;Trigger Type;");
239 hrate.InitLabels(MH3::kLabelsXY);
240 hrate.DefineLabelY( 0, "Data"); // What if TriggerID==0 already???
241 hrate.DefineLabelY(0x100, "Cal");
242 hrate.DefineLabelY(0x400, "Ped");
243 // hrate.DefaultLabelY("ERROR");
244 gStyle->SetOptFit(kTRUE);
245
246
247 // ========================= Result ==================================
248
249 //~ Double_t avgS = evt1f.GetHist()->GetMean();
250 //~ Double_t medS = evt1f.GetHist()->GetMedian();
251 //~ Double_t rmsS = evt1f.GetHist()->GetRMS();
252 //~ Double_t maxS = evt1f.GetHist()->GetMaximum();
253
254 MArrayF der1(hpulse->GetNbinsX());
255 MArrayF der2(hpulse->GetNbinsX());
256
257 MExtralgoSpline spline(hpulse->GetArray()+1, hpulse->GetNbinsX(),
258 der1.GetArray(), der2.GetArray());
259 spline.SetRiseFallTime(rise_time_dat, fall_time_dat);
260 spline.SetExtractionType(type);
261 spline.SetHeightTm(heighttm);
262
263 spline.Extract(hpulse->GetMaximumBin()-1);
264
265 // The pulser signal is most probably around 400mV/9.5mV
266 // IntegraFixed 2/48 corresponds to roughly 215mV*50slices
267 Double_t scale = 1./spline.GetSignal();
268
269 MArrayD calib(1440);
270 for (int i=0; i<1440; i++)
271 calib[i] =1.;
272
273 gROOT->SetSelectedPad(0);
274 d->AddTab("PulseTemp");
275 gPad->SetGrid();
276 hpulse->SetNameTitle("Pulse", "Single p.e. pulse template");
277 hpulse->SetDirectory(0);
278 hpulse->SetLineColor(kBlack);
279 hpulse->DrawCopy();
280
281 TAxis *ax = hpulse->GetXaxis();
282
283 Double_t w = hpulse->GetBinWidth(1);
284 Double_t T = w*(spline.GetTime()+0.5) +ax->GetXmin();
285 //~ Double_t H = w*(hpulse->GetMaximumBin()+0.5)+ax->GetXmin();
286
287 TLine line;
288 line.SetLineColor(kRed);
289 line.DrawLine(T-rise_time_dat*w, spline.GetHeight(),
290 T+fall_time_dat*w, spline.GetHeight());
291 line.DrawLine(T, spline.GetHeight()/4, T, 3*spline.GetHeight()/4);
292 line.DrawLine(T-rise_time_dat*w, 0,
293 T-rise_time_dat*w, spline.GetHeight());
294 line.DrawLine(T+fall_time_dat*w, 0,
295 T+fall_time_dat*w, spline.GetHeight());
296
297 TGraph gg;
298 for (int ix=1; ix<=hpulse->GetNbinsX(); ix++)
299 for (int i=0; i<10; i++)
300 {
301 Double_t x = hpulse->GetBinLowEdge(ix)+i*hpulse->GetBinWidth(ix)/10.;
302 gg.SetPoint(gg.GetN(), x+w/2, spline.EvalAt(ix-1+i/10.));
303 }
304
305 gg.SetLineColor(kBlue);
306 gg.SetMarkerColor(kBlue);
307 gg.SetMarkerStyle(kFullDotMedium);
308 gg.DrawClone("L");
309
310 gROOT->SetSelectedPad(0);
311 d->AddTab("CalConst");
312 MGeomCamFACT fact;
313 MHCamera hcalco(fact);
314 hcalco.SetName("CalConst");
315 hcalco.SetTitle(Form("Relative calibration constant [%.0f/pe]", 1./scale));
316 hcalco.SetCamContent(calib);
317 hcalco.SetAllUsed();
318 //hcalco.Scale(scale);
319 hcalco.DrawCopy();
320
321 // ======================================================
322
323 gLog << endl;
324 gLog.Separator("Extracting random pedestal");
325
326 MTaskList tlist3;
327
328 MParList plist3;
329 plist3.AddToList(&tlist3);
330 plist3.AddToList(&drscalib300);
331 plist3.AddToList(&badpixels);
332 plist3.AddToList(&timecam);
333
334 MEvtLoop loop3("DetermineRndmPed");
335 loop3.SetDisplay(d);
336 loop3.SetParList(&plist3);
337
338 // ------------------ Setup the tasks ---------------
339
340 MRawFitsRead read3;
341 read3.LoadMap(map);
342 read3.AddFile(pedfile);
343
344 MFillH fill3a(&hrate);
345
346 MContinue cont3("(MRawEvtHeader.GetTriggerID&0xff00)!=0x400", "SelectPed");
347
348 MGeomApply apply3;
349
350 MDrsCalibApply drsapply3;
351
352 //---
353
354 MExtractTimeAndChargeSpline extractor3;
355 extractor3.SetRange(first_slice, last_slice);
356 extractor3.SetRiseTimeHiGain(rise_time_dat);
357 extractor3.SetFallTimeHiGain(fall_time_dat);
358 extractor3.SetHeightTm(heighttm);
359 extractor3.SetChargeType(type);
360 extractor3.SetSaturationLimit(600000);
361 extractor3.SetNoiseCalculation(kTRUE);
362
363// MHCamEvent evt2a(0, "PedRdm", "Extracted Pedestal Signal;;S");
364
365// MFillH fill2a(&evt2a, "MExtractedSignalCam", "FillPedRndm");
366
367 // Use this for data, but not for calibration events
368// evt2a.SetErrorSpread(kFALSE);
369
370 /*
371 MCalibrateData conv3;
372 conv3.SetCalibrationMode(MCalibrateData::kNone);
373 conv3.SetPedestalFlag(MCalibrateData::kNo);
374 conv3.SetCalibConvMinLimit(0);
375 conv3.SetCalibConvMaxLimit(10000);
376 conv3.SetScaleFactor(scale);
377 */
378
379 MCalibrateFact conv3;
380 conv3.SetScale(scale);
381 conv3.SetCalibConst(calib);
382
383 MBadPixelsTreat treat3;
384 treat3.SetProcessPedestalRun(kFALSE);
385 treat3.SetProcessPedestalEvt(kFALSE);
386 treat3.SetProcessTimes(kFALSE);
387
388 MHCamEvent evt3b(0, "PedRdm","Interpolated random pedestal;;Signal [~phe]");
389 //evt2b.SetErrorSpread(kFALSE);
390
391 MFillH fill3b(&evt3b, "MSignalCam", "FillPedRdm");
392 fill3b.SetDrawOption("gaus");
393
394 // ------------------ Setup eventloop and run analysis ---------------
395
396 tlist3.AddToList(&read3);
397 tlist3.AddToList(&apply3);
398 tlist3.AddToList(&drsapply3);
399 tlist3.AddToList(&cont3);
400 tlist3.AddToList(&extractor3);
401// tlist3.AddToList(&fill3a);
402 tlist3.AddToList(&conv3);
403 tlist3.AddToList(&treat3);
404 tlist3.AddToList(&fill3b);
405
406 if (!loop3.Eventloop(max3))
407 return 14;
408
409 if (!loop3.GetDisplay())
410 return 15;
411
412 // ======================================================
413
414 gLog << endl;
415 gLog.Separator("Extracting pedestal");
416
417 MTaskList tlist4;
418
419 MParList plist4;
420 plist4.AddToList(&tlist4);
421 plist4.AddToList(&drscalib300);
422 plist4.AddToList(&badpixels);
423 plist4.AddToList(&timecam);
424
425 MEvtLoop loop4("DetermineExtractedPed");
426 loop4.SetDisplay(d);
427 loop4.SetParList(&plist4);
428
429 // ------------------ Setup the tasks ---------------
430
431 MRawFitsRead read4;
432 read4.LoadMap(map);
433 read4.AddFile(pedfile);
434
435 MContinue cont4("(MRawEvtHeader.GetTriggerID&0xff00)!=0x400", "SelectPed");
436
437 MGeomApply apply4;
438
439 MDrsCalibApply drsapply4;
440
441 MExtractTimeAndChargeSpline extractor4;
442 extractor4.SetRange(first_slice, last_slice);
443 extractor4.SetRiseTimeHiGain(rise_time_dat);
444 extractor4.SetFallTimeHiGain(fall_time_dat);
445 extractor4.SetHeightTm(heighttm);
446 extractor4.SetChargeType(type);
447 extractor4.SetSaturationLimit(600000);
448 extractor4.SetNoiseCalculation(kFALSE);
449
450 // MHCamEvent evt3a(0, "PedExt", "Extracted Pedestal Signal;;S");
451
452 // MFillH fill3a(&evt3a, "MExtractedSignalCam", "FillPedExt");
453
454 // Use this for data, but not for calibration events
455// evt3a.SetErrorSpread(kFALSE);
456/*
457 MCalibrateData conv4;
458 conv4.SetCalibrationMode(MCalibrateData::kNone);
459 conv4.SetPedestalFlag(MCalibrateData::kNo);
460 conv4.SetCalibConvMinLimit(0);
461 conv4.SetCalibConvMaxLimit(10000);
462 conv4.SetScaleFactor(scale);
463*/
464 MCalibrateFact conv4;
465 conv4.SetScale(scale);
466 conv4.SetCalibConst(calib);
467
468 MBadPixelsTreat treat4;
469 treat4.SetProcessPedestalRun(kFALSE);
470 treat4.SetProcessPedestalEvt(kFALSE);
471
472 MHCamEvent evt4b(0, "PedExt","Interpolated extracted pedestal;;Signal [~phe]");
473 //evt4b.SetErrorSpread(kFALSE);
474
475 MFillH fill4b(&evt4b, "MSignalCam", "FillPedExt");
476 fill4b.SetDrawOption("gaus");
477
478 // ------------------ Setup eventloop and run analysis ---------------
479
480 tlist4.AddToList(&read4);
481 tlist4.AddToList(&apply4);
482 tlist4.AddToList(&drsapply4);
483 tlist4.AddToList(&cont4);
484 tlist4.AddToList(&extractor4);
485// tlist4.AddToList(&fill4a);
486 tlist4.AddToList(&conv4);
487 tlist4.AddToList(&treat4);
488 tlist4.AddToList(&fill4b);
489
490 if (!loop4.Eventloop(max4))
491 return 15;
492
493 if (!loop4.GetDisplay())
494 return 16;
495
496 // ===================================================================
497
498 gLog << endl;
499 gLog.Separator("Extracting and calibration data");
500
501 MTaskList tlist5;
502
503 MParList plist5;
504 plist5.AddToList(&tlist5);
505 plist5.AddToList(&drscalib300);
506 plist5.AddToList(&badpixels);
507 plist5.AddToList(&timecam);
508
509 MEvtLoop loop5("CalibratingData");
510 loop5.SetDisplay(d);
511 loop5.SetParList(&plist5);
512
513 // ------------------ Setup the tasks ---------------
514
515 MRawFitsRead read5;
516 read5.LoadMap(map);
517 read5.AddFile(datfile);
518
519 MFillH fill5a(&hrate);
520
521 MGeomApply apply5;
522
523 MDrsCalibApply drsapply5;
524
525 MFDataPhrase filterdat("(MRawEvtHeader.GetTriggerID&0xff00)==0", "SelectDat");
526 MFDataPhrase filtercal("(MRawEvtHeader.GetTriggerID&0xff00)==0x100", "SelectCal");
527 MFDataPhrase filterped("(MRawEvtHeader.GetTriggerID&0xff00)==0x400", "SelectPed");
528 MFDataPhrase filterncl("(MRawEvtHeader.GetTriggerID&0xff00)!=0x100", "SelectNonCal");
529
530 //MContinue cont4("MRawEvtHeader.GetTriggerID!=4", "SelectData");
531
532 // ---
533
534 MExtractTimeAndChargeSpline extractor5dat;
535 extractor5dat.SetRange(first_slice, last_slice);
536 extractor5dat.SetRiseTimeHiGain(rise_time_dat);
537 extractor5dat.SetFallTimeHiGain(fall_time_dat);
538 extractor5dat.SetHeightTm(heighttm);
539 extractor5dat.SetChargeType(type);
540 extractor5dat.SetSaturationLimit(600000);
541 extractor5dat.SetNoiseCalculation(kFALSE);
542
543 MExtractTimeAndChargeSpline extractor5cal;
544 extractor5cal.SetRange(first_slice, last_slice);
545 extractor5cal.SetRiseTimeHiGain(rise_time_cal);
546 extractor5cal.SetFallTimeHiGain(fall_time_cal);
547 extractor5cal.SetHeightTm(heighttm);
548 extractor5cal.SetChargeType(type);
549 extractor5cal.SetSaturationLimit(600000);
550 extractor5cal.SetNoiseCalculation(kFALSE);
551
552 MExtractTimeAndChargeSpline extractor5tm("ExtractTM");
553 extractor5tm.SetRange(last_slice, 294);
554 extractor5tm.SetRiseTimeHiGain(1);
555 extractor5tm.SetFallTimeHiGain(1);
556 extractor5tm.SetHeightTm(0.5);
557 extractor5tm.SetChargeType(MExtralgoSpline::kAmplitudeRel);
558 extractor5tm.SetSaturationLimit(600000);
559 extractor5tm.SetNoiseCalculation(kFALSE);
560 extractor5tm.SetNameSignalCam("TimeMarkerAmplitude");
561 extractor5tm.SetNameTimeCam("TimeMarkerTime");
562
563 extractor5dat.SetFilter(&filterncl);
564 extractor5cal.SetFilter(&filtercal);
565 //extractor4tm.SetFilter(&filtercal);
566
567 // ---
568/*
569 MCalibrateData conv5;
570 conv5.SetCalibrationMode(MCalibrateData::kNone);
571 conv5.SetPedestalFlag(MCalibrateData::kNo);
572 conv5.SetCalibConvMinLimit(0);
573 conv5.SetCalibConvMaxLimit(10000);
574 conv5.SetScaleFactor(scale);
575*/
576 MCalibrateFact conv5;
577 conv5.SetScale(scale);
578 conv5.SetCalibConst(calib);
579
580 MCalibrateDrsTimes calctm5;
581 calctm5.SetNameUncalibrated("UncalibratedTimes");
582
583 MCalibrateDrsTimes calctm5tm("CalibrateTimeMarker");
584 calctm5tm.SetNameArrivalTime("TimeMarkerTime");
585 calctm5tm.SetNameUncalibrated("UncalTimeMarker");
586 calctm5tm.SetNameCalibrated("TimeMarker");
587 calctm5tm.SetTimeMarker();
588 //calctm4tm.SetFilter(&filtercal);
589
590 MBadPixelsTreat treat5;
591 treat5.SetProcessPedestalRun(kFALSE);
592 treat5.SetProcessPedestalEvt(kFALSE);
593
594 MHCamEvent evt5b(0, "ExtSig", "Extracted signal;;S [mV·sl]");
595 MHCamEvent evt5c(0, "CalSig", "Calibrated and interpolated signal;;S [~phe]");
596 MHCamEvent evt5d(4, "ExtSigTm", "Extracted time;;T [sl]");
597 MHCamEvent evt5e(6, "CalSigTm", "Calibrated and interpolated time;;T [ns]");
598
599 MFillH fill5b(&evt5b, "MExtractedSignalCam", "FillExtSig");
600 MFillH fill5c(&evt5c, "MSignalCam", "FillCalSig");
601 MFillH fill5d(&evt5d, "MArrivalTimeCam", "FillExtTm");
602 MFillH fill5e(&evt5e, "MSignalCam", "FillCalTm");
603
604 fill5c.SetDrawOption("gaus");
605 fill5d.SetDrawOption("gaus");
606 fill5e.SetDrawOption("gaus");
607
608 /*
609 fill4b.SetFilter(&filterdat);
610 fill4c.SetFilter(&filterdat);
611 fill4d.SetFilter(&filterdat);
612 fill4e.SetFilter(&filterdat);
613 */
614
615 //MFSoftwareTrigger swtrig;
616 //MContinue contsw(&swtrig, "FilterSwTrigger", "Software trigger");
617 //contsw.SetInverted();
618
619 // The second rule is for the case reading raw-files!
620
621 MWriteRootFile write5(root_file_output_path, "RECREATE", "Calibrated Data", 2);
622 write5.AddContainer("MRawRunHeader", "RunHeaders");
623 write5.AddContainer("MGeomCam", "RunHeaders");
624 write5.AddContainer("MSignalCam", "Events");
625 write5.AddContainer("MTime", "Events");
626 write5.AddContainer("MRawEvtHeader", "Events");
627 //write.AddContainer("MTriggerPattern", "Events");
628
629
630 // ------------------ Setup histograms and fill tasks ----------------
631
632 MContinue test;
633 test.SetFilter(&filterncl);
634
635 MTaskList tlist5tm;
636 tlist5tm.AddToList(&extractor5tm);
637 tlist5tm.AddToList(&calctm5tm);
638 tlist5tm.SetFilter(&filtercal);
639
640 MTaskList tlist5dat;
641 tlist5dat.AddToList(&fill5b);
642 tlist5dat.AddToList(&fill5c);
643 tlist5dat.AddToList(&fill5d);
644 tlist5dat.AddToList(&fill5e);
645 tlist5dat.SetFilter(&filterdat);
646
647 tlist5.AddToList(&read5);
648 tlist5.AddToList(&apply5);
649 tlist5.AddToList(&drsapply5);
650 tlist5.AddToList(&filterncl);
651 //tlist5.AddToList(&test);
652 tlist5.AddToList(&filterdat);
653 tlist5.AddToList(&filtercal);
654 tlist5.AddToList(&filterped);
655 tlist5.AddToList(&fill5a);
656 tlist5.AddToList(&extractor5dat);
657 tlist5.AddToList(&extractor5cal);
658 tlist5.AddToList(&calctm5);
659 tlist5.AddToList(&tlist5tm);
660 tlist5.AddToList(&conv5);
661 tlist5.AddToList(&treat5);
662 tlist5.AddToList(&tlist5dat);
663 tlist5.AddToList(&write5);
664
665 if (!loop5.Eventloop(max4))
666 return 18;
667
668 if (!loop5.GetDisplay())
669 return 19;
670
671 d->SetTitle(status_display_title, kFALSE);
672 d->SaveAs(status_display_output_path);
673
674 return 0;
675}
Note: See TracBrowser for help on using the repository browser.