source: trunk/MagicSoft/Mars/mjobs/MJCalibration.cc@ 3718

Last change on this file since 3718 was 3709, checked in by gaug, 20 years ago
*** empty log message ***
File size: 18.8 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! Markus Gaug, 02/2004 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MJCalibration
29//
30//Begin_Html
31/*
32<img src="images/CalibClasses.gif">
33*/
34//End_Html
35//
36/////////////////////////////////////////////////////////////////////////////
37#include "MJCalibration.h"
38
39#include <TF1.h>
40#include <TFile.h>
41#include <TStyle.h>
42#include <TCanvas.h>
43#include <TSystem.h>
44#include <TProfile.h>
45
46#include "MLog.h"
47#include "MLogManip.h"
48
49#include "MRunIter.h"
50#include "MParList.h"
51#include "MTaskList.h"
52#include "MEvtLoop.h"
53
54#include "MHCamera.h"
55#include "MGeomCam.h"
56
57#include "MPedestalCam.h"
58#include "MCalibrationChargeCam.h"
59#include "MCalibrationChargePINDiode.h"
60#include "MCalibrationChargeBlindPix.h"
61#include "MCalibrationChargeCalc.h"
62
63#include "MReadMarsFile.h"
64#include "MGeomApply.h"
65#include "MBadPixelsMerge.h"
66#include "MBadPixelsCam.h"
67#include "MExtractSignal.h"
68#include "MExtractPINDiode.h"
69#include "MExtractBlindPixel.h"
70#include "MExtractSignal2.h"
71#include "MFCosmics.h"
72#include "MContinue.h"
73#include "MFillH.h"
74
75#include "MJCalibration.h"
76#include "MStatusDisplay.h"
77
78ClassImp(MJCalibration);
79
80using namespace std;
81
82MJCalibration::MJCalibration(const char *name, const char *title) : fRuns(0)
83{
84 fName = name ? name : "MJCalibration";
85 fTitle = title ? title : "Tool to create a pedestal file (MPedestalCam)";
86}
87
88void MJCalibration::DrawProjection(MHCamera *obj, Int_t fit) const
89{
90
91 TH1D *obj2 = (TH1D*)obj->Projection(obj->GetName());
92 obj2->SetDirectory(0);
93 obj2->Draw();
94 obj2->SetBit(kCanDelete);
95
96 if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
97 {
98 TArrayI s0(3);
99 s0[0] = 6;
100 s0[1] = 1;
101 s0[2] = 2;
102
103 TArrayI s1(3);
104 s1[0] = 3;
105 s1[1] = 4;
106 s1[2] = 5;
107
108 TArrayI inner(1);
109 inner[0] = 0;
110
111 TArrayI outer(1);
112 outer[0] = 1;
113
114 // Just to get the right (maximum) binning
115 TH1D *half[4];
116 half[0] = obj->ProjectionS(s0, inner, "Sector 6-1-2 Inner");
117 half[1] = obj->ProjectionS(s1, inner, "Sector 3-4-5 Inner");
118 half[2] = obj->ProjectionS(s0, outer, "Sector 6-1-2 Outer");
119 half[3] = obj->ProjectionS(s1, outer, "Sector 3-4-5 Outer");
120
121 for (int i=0; i<4; i++)
122 {
123 half[i]->SetLineColor(kRed+i);
124 half[i]->SetDirectory(0);
125 half[i]->SetBit(kCanDelete);
126 half[i]->Draw("same");
127 }
128 }
129
130 const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
131 const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
132 const Double_t integ = obj2->Integral("width")/2.5;
133 const Double_t mean = obj2->GetMean();
134 const Double_t rms = obj2->GetRMS();
135 const Double_t width = max-min;
136
137 const TString dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
138 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
139
140 const TString tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
141 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"
142 "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])";
143 TF1 *f=0;
144 switch (fit)
145 {
146 case 1:
147 f = new TF1("sgaus", "gaus(0)", min, max);
148 f->SetLineColor(kYellow);
149 f->SetBit(kCanDelete);
150 f->SetParNames("Area", "#mu", "#sigma");
151 f->SetParameters(integ/rms, mean, rms);
152 f->SetParLimits(0, 0, integ);
153 f->SetParLimits(1, min, max);
154 f->SetParLimits(2, 0, width/1.5);
155
156 obj2->Fit(f, "QLR");
157 break;
158
159 case 2:
160 f = new TF1("dgaus",dgausformula.Data(),min,max);
161 f->SetLineColor(kYellow);
162 f->SetBit(kCanDelete);
163 f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2");
164 f->SetParameters(integ,(min+mean)/2.,width/4.,
165 integ/width/2.,(max+mean)/2.,width/4.);
166 // The left-sided Gauss
167 f->SetParLimits(0,integ-1.5 , integ+1.5);
168 f->SetParLimits(1,min+(width/10.), mean);
169 f->SetParLimits(2,0 , width/2.);
170 // The right-sided Gauss
171 f->SetParLimits(3,0 , integ);
172 f->SetParLimits(4,mean, max-(width/10.));
173 f->SetParLimits(5,0 , width/2.);
174 obj2->Fit(f,"QLRM");
175 break;
176
177 case 3:
178 f = new TF1("tgaus",tgausformula.Data(),min,max);
179 f->SetLineColor(kYellow);
180 f->SetBit(kCanDelete);
181 f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}",
182 "A_{2}","#mu_{2}","#sigma_{2}",
183 "A_{3}","#mu_{3}","#sigma_{3}");
184 f->SetParameters(integ,(min+mean)/2,width/4.,
185 integ/width/3.,(max+mean)/2.,width/4.,
186 integ/width/3.,mean,width/2.);
187 // The left-sided Gauss
188 f->SetParLimits(0,integ-1.5,integ+1.5);
189 f->SetParLimits(1,min+(width/10.),mean);
190 f->SetParLimits(2,width/15.,width/2.);
191 // The right-sided Gauss
192 f->SetParLimits(3,0.,integ);
193 f->SetParLimits(4,mean,max-(width/10.));
194 f->SetParLimits(5,width/15.,width/2.);
195 // The Gauss describing the outliers
196 f->SetParLimits(6,0.,integ);
197 f->SetParLimits(7,min,max);
198 f->SetParLimits(8,width/4.,width/1.5);
199 obj2->Fit(f,"QLRM");
200 break;
201
202 case 4:
203 obj2->Fit("pol0", "Q");
204 obj2->GetFunction("pol0")->SetLineColor(kYellow);
205 break;
206
207 case 9:
208 break;
209
210 default:
211 obj2->Fit("gaus", "Q");
212 obj2->GetFunction("gaus")->SetLineColor(kYellow);
213 break;
214 }
215}
216
217void MJCalibration::DrawRadialProfile(MHCamera *obj) const
218{
219
220 TProfile *obj2 = (TProfile*)obj->RadialProfile(obj->GetName());
221 obj2->SetDirectory(0);
222 obj2->Draw();
223 obj2->SetBit(kCanDelete);
224
225 if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
226 {
227
228 TArrayI s0(6);
229 s0[0] = 1;
230 s0[1] = 2;
231 s0[2] = 3;
232 s0[3] = 4;
233 s0[4] = 5;
234 s0[5] = 6;
235
236 TArrayI inner(1);
237 inner[0] = 0;
238
239 TArrayI outer(1);
240 outer[0] = 1;
241
242 // Just to get the right (maximum) binning
243 TProfile *half[2];
244 half[0] = obj->RadialProfileS(s0, inner,Form("%s%s",obj->GetName(),"Inner"));
245 half[1] = obj->RadialProfileS(s0, outer,Form("%s%s",obj->GetName(),"Outer"));
246
247 for (Int_t i=0; i<2; i++)
248 {
249 Double_t min = obj->GetGeomCam().GetMinRadius(i);
250 Double_t max = obj->GetGeomCam().GetMaxRadius(i);
251
252 half[i]->SetLineColor(kRed+i);
253 half[i]->SetDirectory(0);
254 half[i]->SetBit(kCanDelete);
255 half[i]->Draw("same");
256 half[i]->Fit("pol1","Q","",min,max);
257 half[i]->GetFunction("pol1")->SetLineColor(kRed+i);
258 half[i]->GetFunction("pol1")->SetLineWidth(1);
259 }
260 }
261}
262
263
264void MJCalibration::CamDraw(TCanvas &c, const Int_t x, const Int_t y, const MHCamera &cam1,
265 const Int_t fit, const Int_t rad)
266{
267 c.cd(x);
268 gPad->SetBorderMode(0);
269 gPad->SetTicks();
270 MHCamera *obj1=(MHCamera*)cam1.DrawCopy("hist");
271 obj1->SetDirectory(NULL);
272
273 c.cd(x+y);
274 gPad->SetBorderMode(0);
275 obj1->SetPrettyPalette();
276 obj1->AddNotify(&fCalibrationCam);
277 obj1->Draw();
278
279 if (rad)
280 {
281 c.cd(x+2*y);
282 gPad->SetBorderMode(0);
283 gPad->SetTicks();
284 DrawRadialProfile(obj1);
285 }
286
287
288 if (!fit)
289 return;
290
291 c.cd(rad ? x+3*y : x+2*y);
292 gPad->SetBorderMode(0);
293 gPad->SetTicks();
294 DrawProjection(obj1, fit);
295}
296
297
298void MJCalibration::DisplayResult(MParList &plist)
299{
300 if (!fDisplay)
301 return;
302
303 //
304 // Update display
305 //
306 TString title = fDisplay->GetTitle();
307 title += "-- Calibration ";
308 title += fRuns->GetRunsAsString();
309 title += " --";
310 fDisplay->SetTitle(title);
311
312 //
313 // Get container from list
314 //
315 MGeomCam &geomcam = *(MGeomCam*)plist.FindObject("MGeomCam");
316
317 // Create histograms to display
318 MHCamera disp1 (geomcam, "Cal;Charge", "Fitted Mean Charges");
319 MHCamera disp2 (geomcam, "Cal;SigmaCharge", "Sigma of Fitted Charges");
320 MHCamera disp3 (geomcam, "Cal;RSigma", "Reduced Sigmas");
321 MHCamera disp4 (geomcam, "Cal;RSigma/Charge", "Reduced Sigma per Charge");
322 MHCamera disp5 (geomcam, "Cal;FFactorPhe", "Nr. of Phe's (F-Factor Method)");
323 MHCamera disp6 (geomcam, "Cal;FFactorConv", "Conversion Factor (F-Factor Method)");
324 MHCamera disp7 (geomcam, "Cal;BlindPixConv", "Conversion Factor (Blind Pixel Method)");
325 MHCamera disp8 (geomcam, "Cal;PINDiodeConv", "Conversion Factor (PIN Diode Method)");
326 MHCamera disp9 (geomcam, "Cal;FFactorValid", "Pixels with valid F-Factor calibration");
327 MHCamera disp10(geomcam, "Cal;BlindPixelValid","Pixels with valid BlindPixel calibration");
328 MHCamera disp11(geomcam, "Cal;PINdiodeValid", "Pixels with valid PINDiode calibration");
329 MHCamera disp12(geomcam, "Cal;Excluded", "Pixels previously excluded");
330 MHCamera disp13(geomcam, "Cal;Saturation", "Pixels with saturated Hi Gain");
331 MHCamera disp14(geomcam, "Cal;HiGainPickup", "Number Pickup events Hi Gain");
332 MHCamera disp15(geomcam, "Cal;LoGainPickup", "Number Pickup events Lo Gain");
333 MHCamera disp16(geomcam, "Bad;UnSuitable", "Pixels not suited for further analysis");
334 MHCamera disp17(geomcam, "Bad;UnReliable", "Pixels not reliable for further analysis");
335 MHCamera disp18(geomcam, "Bad;Oscillation", "Oscillating Pixels");
336 MHCamera disp19(geomcam, "Cal;AbsTimeMean", "Abs. Arrival Times");
337 MHCamera disp20(geomcam, "Cal;AbsTimeRms", "RMS of Arrival Times");
338
339 // Fitted charge means and sigmas
340 disp1.SetCamContent(fCalibrationCam, 0);
341 disp1.SetCamError( fCalibrationCam, 1);
342 disp2.SetCamContent(fCalibrationCam, 2);
343 disp2.SetCamError( fCalibrationCam, 3);
344
345 // Reduced Sigmas and reduced sigmas per charge
346 disp3.SetCamContent(fCalibrationCam, 5);
347 disp3.SetCamError( fCalibrationCam, 6);
348 disp4.SetCamContent(fCalibrationCam, 7);
349 disp4.SetCamError( fCalibrationCam, 8);
350
351 // Conversion Factors, Methods
352 disp5.SetCamContent(fCalibrationCam, 9);
353 disp5.SetCamError( fCalibrationCam, 10);
354 disp6.SetCamContent(fCalibrationCam, 11);
355 disp6.SetCamError( fCalibrationCam, 12);
356 disp7.SetCamContent(fCalibrationCam, 13);
357 disp7.SetCamError( fCalibrationCam, 14);
358 disp8.SetCamContent(fCalibrationCam, 16);
359 disp8.SetCamError( fCalibrationCam, 17);
360
361 // Valid flags
362 disp9.SetCamContent (fCalibrationCam, 15);
363 disp10.SetCamContent(fCalibrationCam, 20);
364 disp11.SetCamContent(fCalibrationCam, 25);
365
366 // Pixels behavior
367 disp12.SetCamContent(fCalibrationCam, 26);
368 disp13.SetCamContent(fCalibrationCam, 29);
369 disp14.SetCamContent(fCalibrationCam, 27);
370 disp15.SetCamContent(fCalibrationCam, 28);
371
372 // Pixels with defects
373 disp16.SetCamContent(fBadPixels, 1);
374 disp17.SetCamContent(fBadPixels, 3);
375 disp18.SetCamContent(fBadPixels, 8);
376
377 // Arrival Times
378 disp19.SetCamContent(fCalibrationCam, 34);
379 disp19.SetCamError( fCalibrationCam, 35);
380 disp20.SetCamContent(fCalibrationCam, 35);
381
382 disp1.SetYTitle("Q [FADC units]");
383 disp2.SetYTitle("\\sigma_{Q} [FADC units]");
384
385 disp3.SetYTitle("\\sqrt{\\sigma^{2}_{Q} - RMS^{2}_{Ped}} [FADC Counts]");
386 disp4.SetYTitle("Red.Sigma/<Q> [1]");
387
388 disp5.SetYTitle("PhE [#]");
389 disp6.SetYTitle("Conv.Factor [PhE/FADC units]");
390
391 disp7.SetYTitle("Conv.Factor [Phot/FADC Count]");
392 disp8.SetYTitle("Conv.Factor [Phot/FADC Count]");
393
394 disp9.SetYTitle("[1]");
395 disp10.SetYTitle("[1]");
396 disp11.SetYTitle("[1]");
397 disp12.SetYTitle("[1]");
398 disp13.SetYTitle("[1]");
399 disp14.SetYTitle("[1]");
400 disp15.SetYTitle("[1]");
401 disp16.SetYTitle("[1]");
402 disp17.SetYTitle("[1]");
403 disp18.SetYTitle("[1]");
404 disp19.SetYTitle("Mean Abs. Time [FADC slice]");
405 disp20.SetYTitle("RMS Abs. Time [FADC slices]");
406
407 gStyle->SetOptStat(1111);
408 gStyle->SetOptFit();
409
410 // Charges
411 TCanvas &c1 = fDisplay->AddTab("Fit.Charge");
412 c1.Divide(2, 4);
413
414 CamDraw(c1, 1, 2, disp1, 2, 1);
415 CamDraw(c1, 2, 2, disp2, 2, 1);
416
417 // Reduced Sigmas
418 TCanvas &c3 = fDisplay->AddTab("Red.Sigma");
419 c3.Divide(2,4);
420
421 CamDraw(c3, 1, 2, disp3, 2, 1);
422 CamDraw(c3, 2, 2, disp4, 2, 1);
423
424 // Methods
425 TCanvas &c4 = fDisplay->AddTab("Methods");
426 c4.Divide(4,4);
427
428 CamDraw(c4, 1, 4, disp5, 2, 1);
429 CamDraw(c4, 2, 4, disp6, 2, 1);
430 CamDraw(c4, 3, 4, disp7, 2, 1);
431 CamDraw(c4, 4, 4, disp8, 2, 1);
432
433 // Validity
434 TCanvas &c5 = fDisplay->AddTab("Validity");
435 c5.Divide(3,2);
436
437 CamDraw(c5, 1, 3, disp9 , 0);
438 CamDraw(c5, 2, 3, disp10, 0);
439 CamDraw(c5, 3, 3, disp11, 0);
440
441 // Other info
442 TCanvas &c6 = fDisplay->AddTab("Behavior");
443 c6.Divide(4,2);
444
445 CamDraw(c6, 1, 4, disp12, 0);
446 CamDraw(c6, 2, 4, disp13, 0);
447 CamDraw(c6, 3, 4, disp14, 0);
448 CamDraw(c6, 4, 4, disp15, 0);
449
450 // Defects
451 TCanvas &c7 = fDisplay->AddTab("Defects");
452 c7.Divide(3,2);
453
454 CamDraw(c7, 1, 3, disp16, 0);
455 CamDraw(c7, 2, 3, disp17, 0);
456 CamDraw(c7, 3, 3, disp18, 0);
457
458 // Abs. Times
459 TCanvas &c8 = fDisplay->AddTab("Abs. Times");
460 c8.Divide(2,3);
461
462 CamDraw(c8, 1, 2, disp19, 2);
463 CamDraw(c8, 2, 2, disp20, 1);
464
465}
466
467Bool_t MJCalibration::WriteResult()
468{
469 if (fOutputPath.IsNull())
470 return kTRUE;
471
472 const TString oname(GetOutputFile());
473
474 *fLog << inf << "Writing to file: " << oname << endl;
475
476 TFile file(oname, "UPDATE");
477
478 if (fDisplay && fDisplay->Write()<=0)
479 {
480 *fLog << err << "Unable to write MStatusDisplay to " << oname << endl;
481 return kFALSE;
482 }
483
484 if (fCalibrationCam.Write()<=0)
485 {
486 *fLog << err << "Unable to write MCalibrationCam to " << oname << endl;
487 return kFALSE;
488 }
489
490 if (fBadPixels.Write()<=0)
491 {
492 *fLog << err << "Unable to write MBadPixelsCam to " << oname << endl;
493 return kFALSE;
494 }
495
496 return kTRUE;
497
498}
499
500void MJCalibration::SetOutputPath(const char *path)
501{
502 fOutputPath = path;
503 if (fOutputPath.EndsWith("/"))
504 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
505}
506
507Bool_t MJCalibration::Process(MPedestalCam &pedcam)
508{
509 if (!ReadCalibrationCam())
510 return ProcessFile(pedcam);
511
512 return kTRUE;
513}
514
515TString MJCalibration::GetOutputFile() const
516{
517 if (!fRuns)
518 return "";
519
520 return Form("%s/%s-F1.root", (const char*)fOutputPath, (const char*)fRuns->GetRunsAsFileName());
521}
522
523Bool_t MJCalibration::ReadCalibrationCam()
524{
525 const TString fname = GetOutputFile();
526
527 if (gSystem->AccessPathName(fname, kFileExists))
528 {
529 *fLog << err << "Input file " << fname << " doesn't exist." << endl;
530 return kFALSE;
531 }
532
533 *fLog << inf << "Reading from file: " << fname << endl;
534
535 TFile file(fname, "READ");
536 if (fCalibrationCam.Read()<=0)
537 {
538 *fLog << err << "Unable to read MCalibrationCam from " << fname << endl;
539 return kFALSE;
540 }
541
542 if (file.FindKey("MBadPixelsCam"))
543 {
544 MBadPixelsCam bad;
545 if (bad.Read()<=0)
546 {
547 *fLog << err << "Unable to read MBadPixelsCam from " << fname << endl;
548 return kFALSE;
549 }
550 fBadPixels.Merge(bad);
551 }
552
553 if (fDisplay /*&& !fDisplay->GetCanvas("Pedestals")*/) // FIXME!
554 fDisplay->Read();
555
556 return kTRUE;
557}
558
559
560Bool_t MJCalibration::ProcessFile(MPedestalCam &pedcam)
561{
562 if (!fRuns)
563 {
564 *fLog << err << "No Runs choosen... abort." << endl;
565 return kFALSE;
566 }
567 if (fRuns->GetNumRuns() != fRuns->GetNumEntries())
568 {
569 *fLog << err << "Number of files found doesn't match number of runs... abort." << endl;
570 return kFALSE;
571 }
572
573 *fLog << inf;
574 fLog->Separator(GetDescriptor());
575 *fLog << "Calculate MCalibrationCam from Runs " << fRuns->GetRunsAsString() << endl;
576 *fLog << endl;
577
578 MReadMarsFile read("Events");
579 read.DisableAutoScheme();
580 static_cast<MRead&>(read).AddFiles(*fRuns);
581
582 //
583 // As long, as we don't have digital modules,
584 // we have to set the color of the pulser LED by hand
585 //
586 MCalibrationChargePINDiode pindiode;
587 MCalibrationChargeBlindPix blindpix;
588
589 // Setup Tasklist
590 MParList plist;
591 plist.AddToList(&pedcam);
592 plist.AddToList(&fCalibrationCam);
593 plist.AddToList(&fBadPixels);
594 plist.AddToList(&pindiode);
595 plist.AddToList(&blindpix);
596
597 MTaskList tlist;
598 plist.AddToList(&tlist);
599
600 MGeomApply apply;
601 MBadPixelsMerge merge(&fBadPixels);
602 MExtractPINDiode pinext;
603 MExtractBlindPixel blindext;
604 // MExtractSignal extract; // Do not use this at the moment...
605 MExtractSignal2 extract;
606 MCalibrationChargeCalc calcalc;
607 calcalc.SetPulserColor(MCalibrationCam::kCT1);
608
609 MFillH fillpin("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode");
610 MFillH fillbnd("MHCalibrationChargeBlindPix", "MExtractedSignalBlindPixel");
611 MFillH fillcam("MHCalibrationChargeCam", "MExtractedSignalCam");
612 fillpin.SetNameTab("PINDiode");
613 fillbnd.SetNameTab("BlindPix");
614 fillcam.SetNameTab("Charge");
615
616 //
617 // Apply a filter against cosmics
618 // (will have to be needed in the future
619 // when the calibration hardware-trigger is working)
620 //
621 MFCosmics cosmics;
622 MContinue cont(&cosmics);
623
624 //calcalc.SkipBlindPixelFit();
625
626 tlist.AddToList(&read);
627 tlist.AddToList(&merge);
628 tlist.AddToList(&apply);
629 tlist.AddToList(&extract);
630 tlist.AddToList(&pinext);
631 tlist.AddToList(&blindext);
632 tlist.AddToList(&cont);
633 tlist.AddToList(&fillcam);
634 tlist.AddToList(&fillpin);
635 tlist.AddToList(&fillbnd);
636 tlist.AddToList(&calcalc);
637
638 // Create and setup the eventloop
639 MEvtLoop evtloop(fName);
640 evtloop.SetParList(&plist);
641 evtloop.SetDisplay(fDisplay);
642 evtloop.SetLogStream(fLog);
643
644 // Execute first analysis
645 if (!evtloop.Eventloop())
646 {
647 *fLog << err << GetDescriptor() << ": Failed." << endl;
648 return kFALSE;
649 }
650
651 tlist.PrintStatistics();
652
653 DisplayResult(plist);
654
655 if (!WriteResult())
656 return kFALSE;
657
658 *fLog << inf << GetDescriptor() << ": Done." << endl;
659
660 return kTRUE;
661}
Note: See TracBrowser for help on using the repository browser.