source: trunk/MagicSoft/Mars/mcalib/MCalibrationTestCalc.cc@ 8149

Last change on this file since 8149 was 8049, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 20.3 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): Markus Gaug 08/2004 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MCalibrationTestCalc
28//
29// PreProcess(): Initialize pointers to MHCalibrationTestCam
30//
31// ReInit(): MCalibrationCam::InitSize(NumPixels) is called from MGeomApply (which allocates
32// memory in a TClonesArray of type MCalibrationChargePix)
33// Initializes pointer to MBadPixelsCam
34//
35// Process(): Nothing to be done, histograms getting filled by MHCalibrationTestCam
36//
37// PostProcess(): Print out interpolation results to file
38//
39// Input Containers:
40// MHCalibrationTestCam
41// MBadPixelsCam
42// MGeomCam
43//
44// Output Containers:
45// none
46//
47//////////////////////////////////////////////////////////////////////////////
48#include "MCalibrationTestCalc.h"
49
50#include <TSystem.h>
51#include <TH1.h>
52#include <TF1.h>
53
54#include "MLog.h"
55#include "MLogManip.h"
56
57#include "MParList.h"
58
59#include "MGeomCam.h"
60#include "MGeomPix.h"
61#include "MHCamera.h"
62
63#include "MHCalibrationTestCam.h"
64#include "MHCalibrationPix.h"
65
66#include "MCalibrationIntensityTestCam.h"
67#include "MCalibrationTestCam.h"
68#include "MCalibrationTestPix.h"
69
70#include "MBadPixelsIntensityCam.h"
71#include "MBadPixelsCam.h"
72#include "MBadPixelsPix.h"
73
74ClassImp(MCalibrationTestCalc);
75
76using namespace std;
77
78const Float_t MCalibrationTestCalc::fgPhotErrLimit = 4.5;
79// --------------------------------------------------------------------------
80//
81// Default constructor.
82//
83// Sets the pointer to fTestCam and fGeom to NULL
84// Sets outputpath to "."
85// Sets outputfile to "TestCalibStat.txt"
86// Sets fPhotErrLimit to fgPhotErrLimit
87//
88// Calls:
89// - Clear()
90//
91MCalibrationTestCalc::MCalibrationTestCalc(const char *name, const char *title)
92 : fIntensBad(NULL), fBadPixels(NULL),
93 fTestCam(NULL), fIntensCam(NULL), fCam(NULL),
94 fGeom(NULL)
95{
96
97 fName = name ? name : "MCalibrationTestCalc";
98 fTitle = title ? title : "Task to output the results of MHCalibrationTestCam ";
99
100 SetPhotErrLimit();
101
102 SetOutputPath();
103 SetOutputFile();
104
105}
106
107
108// --------------------------------------------------------------------------
109//
110// Search for the following input containers and abort if not existing:
111// - MGeomCam
112// - MHCalibrationTestCam
113// - MCalibrationIntensityTestCam or MCalibrationTestCam
114// - MBadPixelsIntensityCam or MBadPixelsCam
115//
116Bool_t MCalibrationTestCalc::ReInit(MParList *pList )
117{
118
119 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
120 if (!fGeom)
121 {
122 *fLog << err << "No MGeomCam found... aborting." << endl;
123 return kFALSE;
124 }
125
126 fTestCam = (MHCalibrationTestCam*)pList->FindObject("MHCalibrationTestCam");
127 if (!fTestCam)
128 {
129 *fLog << err << "Cannot find MHCalibrationTestCam... aborting" << endl;
130 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationTestCam before..." << endl;
131 return kFALSE;
132 }
133
134 fIntensCam = (MCalibrationIntensityTestCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityTestCam"));
135 if (fIntensCam)
136 *fLog << inf << "Found MCalibrationIntensityTestCam ... " << endl;
137 else
138 {
139 fCam = (MCalibrationTestCam*)pList->FindObject(AddSerialNumber("MCalibrationTestCam"));
140 if (!fCam)
141 {
142 *fLog << err << "Cannot find MCalibrationTestCam ... abort." << endl;
143 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationTestCam before..." << endl;
144 return kFALSE;
145 }
146 }
147
148 fIntensBad = (MBadPixelsIntensityCam*)pList->FindObject(AddSerialNumber("MBadPixelsIntensityCam"));
149 if (fIntensBad)
150 *fLog << inf << "Found MBadPixelsIntensityCam ... " << endl;
151 else
152 {
153 fBadPixels = (MBadPixelsCam*)pList->FindObject(AddSerialNumber("MBadPixelsCam"));
154 if (!fBadPixels)
155 {
156 *fLog << err << "Cannot find MBadPixelsCam ... abort." << endl;
157 return kFALSE;
158 }
159 }
160
161
162 return kTRUE;
163}
164
165// ----------------------------------------------------------------------------------
166//
167// Nothing to be done in Process, but have a look at MHCalibrationTestCam, instead
168//
169Int_t MCalibrationTestCalc::Process()
170{
171 return kTRUE;
172}
173
174// -----------------------------------------------------------------------
175//
176// Return if number of executions is null.
177//
178// Print out some statistics
179//
180Int_t MCalibrationTestCalc::PostProcess()
181{
182
183 if (GetNumExecutions()==0)
184 return kTRUE;
185
186 MCalibrationTestCam *testcam = fIntensCam
187 ? (MCalibrationTestCam*)fIntensCam->GetCam() : fCam;
188
189 //
190 // Re-direct the output to an ascii-file from now on:
191 //
192 MLog asciilog;
193 asciilog.SetOutputFile(GetOutputFile(),kTRUE);
194 SetLogStream(&asciilog);
195 //
196 // Finalize calibration statistics
197 //
198 FinalizeCalibratedPhotons();
199 FinalizeNotInterpolated();
200 const Int_t maxbad = CalcMaxNumBadPixelsCluster();
201
202
203 *fLog << inf << endl;
204 *fLog << GetDescriptor() << ": Pixel Interpolation status:" << endl;
205
206 if (fGeom->InheritsFrom("MGeomCamMagic"))
207 {
208 *fLog << " Not interpolateable Pixels:";
209 *fLog << " Inner: " << Form("%3i", testcam->GetNumUninterpolated(0));
210 *fLog << " Outer: " << Form("%3i", testcam->GetNumUninterpolated(1)) << endl;
211 *fLog << " Biggest not-interpolateable cluster: " << maxbad << endl;
212 }
213
214 testcam->SetNumUninterpolatedInMaxCluster(maxbad);
215
216 *fLog << endl;
217 SetLogStream(&gLog);
218
219 return kTRUE;
220}
221
222
223// ------------------------------------------------------------------------
224//
225//
226// First loop: Calculate a mean and mean RMS of calibrated photons per area index
227// Include only MHCalibrationTestPix's which are not empty (not interpolated)
228//
229// Second loop: Get weighted mean number of calibrated photons and its RMS
230// excluding those deviating by more than fPhotErrLimit
231// sigmas from the mean (obtained in first loop). Set
232// MBadPixelsPix::kDeviatingNumPhots if excluded.
233//
234void MCalibrationTestCalc::FinalizeCalibratedPhotons() const
235{
236
237 const UInt_t npixels = fGeom->GetNumPixels();
238 const UInt_t nareas = fGeom->GetNumAreas();
239 const UInt_t nsectors = fGeom->GetNumSectors();
240
241 MCalibrationTestCam *testcam = fIntensCam
242 ? (MCalibrationTestCam*)fIntensCam->GetCam() : fCam;
243 MBadPixelsCam *badcam = fIntensBad
244 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
245
246 TArrayD lowlim (nareas);
247 TArrayD upplim (nareas);
248 TArrayD areaphotons (nareas);
249 TArrayD sectorphotons(nsectors);
250 TArrayD areavars (nareas);
251 TArrayD sectorvars (nsectors);
252 TArrayD fittedmean (nareas);
253 TArrayD fittedsigma (nareas);
254 TArrayI numareavalid(nareas);
255 TArrayI numsectorvalid(nsectors);
256
257 //
258 // First loop: Get mean number of calibrated photons and the RMS
259 // The loop is only to recognize later pixels with very deviating numbers
260 //
261 MHCamera camphotons(*fGeom,"Camphotons","Photons in Camera");
262
263 for (UInt_t i=0; i<npixels; i++)
264 {
265
266 MHCalibrationPix &hist = (*fTestCam)[i];
267 MCalibrationTestPix &pix = (MCalibrationTestPix&)(*testcam)[i];
268 //
269 // We assume that the pixels have been interpolated so far.
270 // The MBadPixelsCam does not give any more information
271 //
272 if (hist.IsEmpty())
273 {
274 pix.SetExcluded();
275 continue;
276 }
277
278 const Double_t nphot = hist.GetMean();
279 const Double_t nphoterr = hist.GetMeanErr();
280 const Int_t aidx = (*fGeom)[i].GetAidx();
281
282 camphotons.Fill(i,nphot);
283 camphotons.SetUsed(i);
284
285 pix.SetNumPhotons ( nphot );
286 pix.SetNumPhotonsErr( nphoterr );
287
288 areaphotons [aidx] += nphot;
289 areavars [aidx] += nphot*nphot;
290 numareavalid[aidx] ++;
291 }
292
293 for (UInt_t aidx=0; aidx<nareas; aidx++)
294 {
295 if (numareavalid[aidx] == 0)
296 {
297 *fLog << warn << GetDescriptor() << ": No pixels with valid number of calibrated photons found "
298 << "in area index: " << aidx << endl;
299 continue;
300 }
301
302 if (numareavalid[aidx] == 1)
303 areavars[aidx] = 0.;
304 else if (numareavalid[aidx] == 0)
305 {
306 areaphotons[aidx] = -1.;
307 areavars[aidx] = -1.;
308 }
309 else
310 {
311 areavars[aidx] = (areavars[aidx] - areaphotons[aidx]*areaphotons[aidx]/numareavalid[aidx])
312 / (numareavalid[aidx]-1.);
313 areaphotons[aidx] = areaphotons[aidx] / numareavalid[aidx];
314 }
315
316 if (areavars[aidx] < 0.)
317 {
318 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of calibrated photons found "
319 << "in area index: " << aidx << endl;
320 continue;
321 }
322
323 const Float_t areamean = areaphotons[aidx];
324 const Float_t areaerr = TMath::Sqrt(areavars[aidx]);
325
326 MCalibrationTestPix &avpix = (MCalibrationTestPix&)testcam->GetAverageArea(aidx);
327 avpix.SetNumPhotons (areamean);
328 avpix.SetNumPhotonsErr(areaerr );
329
330 lowlim [aidx] = areamean - fPhotErrLimit*areaerr;
331 upplim [aidx] = areamean + fPhotErrLimit*areaerr;
332
333 TArrayI area(1);
334 area[0] = aidx;
335
336 TH1D *hist = camphotons.ProjectionS(TArrayI(),area,"_py",100);
337 hist->Fit("gaus","Q");
338 const Double_t mean = hist->GetFunction("gaus")->GetParameter(1);
339 const Double_t sigma = hist->GetFunction("gaus")->GetParameter(2);
340 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
341
342 if (ndf < 2)
343 {
344 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of calibrated photons "
345 << "in the camera with area index: " << aidx << endl;
346 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
347 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
348 delete hist;
349 continue;
350 }
351
352 const Double_t prob = hist->GetFunction("gaus")->GetProb();
353
354 if (prob < 0.001)
355 {
356 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of calibrated photons "
357 << "in the camera with area index: " << aidx << endl;
358 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
359 << " is smaller than 0.001 " << endl;
360 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
361 delete hist;
362 continue;
363 }
364
365 fittedmean [aidx] = mean;
366 fittedsigma[aidx] = sigma;
367
368 avpix.SetNumPhotons (mean );
369 avpix.SetNumPhotonsErr(sigma);
370
371 lowlim [aidx] = mean - fPhotErrLimit*sigma;
372 upplim [aidx] = mean + fPhotErrLimit*sigma;
373
374 *fLog << inf << GetDescriptor()
375 << ": Fitted number of calib. equiv. Cher. photons in area index " << aidx
376 << ": " << Form("%7.2f +- %6.2f",fittedmean[aidx],fittedsigma[aidx]) << endl;
377
378 delete hist;
379 }
380
381 *fLog << endl;
382
383 numareavalid.Reset();
384 areaphotons .Reset();
385 areavars .Reset();
386
387 //
388 // Second loop: Get mean number of calibrated photons and its RMS excluding
389 // pixels deviating by more than fPhotErrLimit sigma.
390 //
391 for (UInt_t i=0; i<npixels; i++)
392 {
393
394 MHCalibrationPix &hist = (*fTestCam)[i];
395 MCalibrationTestPix &pix = (MCalibrationTestPix&) (*testcam)[i];
396
397 const Int_t aidx = (*fGeom)[i].GetAidx();
398 const Int_t sector = (*fGeom)[i].GetSector();
399 const Double_t nphot = hist.GetMean();
400 const Double_t nphotpera = nphot / (*fGeom)[i].GetA();
401 const Double_t nphotperaerr = hist.GetMeanErr()/ (*fGeom)[i].GetA();
402
403 pix.SetNumPhotonsPerArea ( nphotpera );
404 pix.SetNumPhotonsPerAreaErr( nphotperaerr );
405
406 if ( nphot < lowlim[aidx] || nphot > upplim[aidx] )
407 {
408 *fLog << warn << GetDescriptor() << ": Number of photons: "
409 << Form("%8.2f out of %3.1f sigma limit: ",nphot,fPhotErrLimit)
410 << Form("[%8.2f,%8.2f] pixel%4i",lowlim[aidx],upplim[aidx],i) << endl;
411 MBadPixelsPix &bad = (*badcam)[i];
412 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhots );
413 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
414 continue;
415 }
416
417 areavars [aidx] += nphotpera*nphotpera;
418 areaphotons [aidx] += nphotpera;
419 numareavalid [aidx] ++;
420
421 sectorvars [sector] += nphotpera*nphotpera;
422 sectorphotons [sector] += nphotpera;
423 numsectorvalid[sector] ++;
424 }
425
426 *fLog << endl;
427
428 for (UInt_t aidx=0; aidx<nareas; aidx++)
429 {
430
431 if (numareavalid[aidx] == 1)
432 areavars[aidx] = 0.;
433 else if (numareavalid[aidx] == 0)
434 {
435 areaphotons[aidx] = -1.;
436 areavars[aidx] = -1.;
437 }
438 else
439 {
440 areavars[aidx] = (areavars[aidx] - areaphotons[aidx]*areaphotons[aidx]/numareavalid[aidx])
441 / (numareavalid[aidx]-1.);
442 areaphotons[aidx] /= numareavalid[aidx];
443 }
444
445
446 MCalibrationTestPix &avpix = (MCalibrationTestPix&)testcam->GetAverageArea(aidx);
447
448 if (areavars[aidx] < 0. || areaphotons[aidx] <= 0.)
449 {
450 *fLog << warn << GetDescriptor()
451 << ": Mean number of photons per area in area index "
452 << aidx << " could not be calculated! Mean: " << areaphotons[aidx]
453 << " Variance: " << areavars[aidx] << endl;
454 avpix.SetExcluded();
455 continue;
456 }
457
458 const Float_t areaerr = TMath::Sqrt(areavars[aidx]);
459
460 avpix.SetNumPhotonsPerArea (areaphotons[aidx]);
461 avpix.SetNumPhotonsPerAreaErr(areaerr );
462
463 *fLog << inf << GetDescriptor() << ": Mean number of equiv. Cher. photons "
464 << "per area in area idx " << aidx << ": "
465 << Form("%5.3f+-%5.4f [ph/mm^2]",areaphotons[aidx],areaerr) << endl;
466 }
467
468 *fLog << endl;
469
470 for (UInt_t sector=0; sector<nsectors; sector++)
471 {
472
473 if (numsectorvalid[sector] == 1)
474 sectorvars[sector] = 0.;
475 else if (numsectorvalid[sector] == 0)
476 {
477 sectorphotons[sector] = -1.;
478 sectorvars[sector] = -1.;
479 }
480 else
481 {
482 sectorvars[sector] = (sectorvars[sector]
483 - sectorphotons[sector]*sectorphotons[sector]/numsectorvalid[sector]
484 )
485 / (numsectorvalid[sector]-1.);
486 sectorphotons[sector] /= numsectorvalid[sector];
487 }
488
489 MCalibrationTestPix &avpix = (MCalibrationTestPix&)testcam->GetAverageSector(sector);
490
491 if (sectorvars[sector] < 0. || sectorphotons[sector] <= 0.)
492 {
493 *fLog << warn << GetDescriptor()
494 << ": Mean number of calibrated photons per area in sector "
495 << sector << " could not be calculated! Mean: " << sectorphotons[sector]
496 << " Variance: " << sectorvars[sector] << endl;
497 avpix.SetExcluded();
498 continue;
499 }
500
501
502 const Float_t sectorerr = TMath::Sqrt(sectorvars[sector]);
503
504 avpix.SetNumPhotonsPerArea (sectorphotons[sector]);
505 avpix.SetNumPhotonsPerAreaErr(sectorerr );
506
507 *fLog << inf << GetDescriptor() << ": Mean number of equiv. Cher. photons "
508 << "per area in sector " << sector << ": "
509 << Form("%5.3f+-%5.4f [ph/mm^2]",sectorphotons[sector],sectorerr) << endl;
510 }
511
512 return;
513}
514
515
516// -----------------------------------------------------------------------------------------------
517//
518// Print out statistics about not interpolated pixels
519//
520void MCalibrationTestCalc::FinalizeNotInterpolated()
521{
522
523 MCalibrationTestCam *testcam = fIntensCam
524 ? (MCalibrationTestCam*)fIntensCam->GetCam() : fCam;
525
526 const Int_t areas = fGeom->GetNumAreas();
527 TArrayI *newarr[areas];
528
529 for (Int_t aidx=0; aidx<areas; aidx++)
530 newarr[aidx] = new TArrayI(0);
531
532 for (Int_t i=0; i<testcam->GetSize(); i++)
533 {
534 const Int_t aidx = (*fGeom)[i].GetAidx();
535 MCalibrationTestPix &pix = (MCalibrationTestPix&)(*testcam)[i];
536 if (pix.IsExcluded())
537 {
538 const Int_t size = newarr[aidx]->GetSize();
539 newarr[aidx]->Set(size+1);
540 newarr[aidx]->AddAt(i,size);
541 }
542 }
543
544 Int_t num = 0;
545
546 for (Int_t aidx = 0; aidx<areas; aidx++)
547 {
548 *fLog << endl;
549 *fLog << " " << setw(7)
550 << "Not interpolated pixels by in area idx " << aidx << ": ";
551 for (Int_t i=0; i<newarr[aidx]->GetSize(); i++)
552 {
553 *fLog << newarr[aidx]->At(i) << " ";
554 num++;
555 }
556 testcam->SetNumUninterpolated(newarr[aidx]->GetSize(),aidx);
557 *fLog << endl;
558 }
559
560 *fLog << " " << setw(7) << num << " total not interpolateable pixels " << endl;
561
562}
563
564Int_t MCalibrationTestCalc::CalcMaxNumBadPixelsCluster()
565{
566
567 MCalibrationTestCam *testcam = fIntensCam
568 ? (MCalibrationTestCam*)fIntensCam->GetCam() : fCam;
569
570 TArrayI arr(0);
571
572 for (Int_t i=0; i<testcam->GetSize(); i++)
573 {
574 MCalibrationTestPix &pix = (MCalibrationTestPix&)(*testcam)[i];
575 if (pix.IsExcluded())
576 {
577 const Int_t s = arr.GetSize();
578 arr.Set(s+1);
579 arr[s] = i;
580 }
581 }
582
583 const Int_t size = arr.GetSize();
584
585 if (size == 0)
586 return 0;
587
588 if (size == 1)
589 return 1;
590
591 TArrayI knownpixels(0);
592 Int_t clustersize = 1;
593 Int_t oldclustersize = 0;
594 //
595 // Loop over the not-interpolateable pixels:
596 //
597 for (Int_t i=0; i<size; i++)
598 {
599
600 const Int_t id = arr[i];
601 const Int_t knownsize = knownpixels.GetSize();
602 knownpixels.Set(knownsize+1);
603 knownpixels[knownsize] = id;
604 LoopNeighbours(arr, knownpixels, clustersize, id);
605 if (clustersize > oldclustersize)
606 oldclustersize = clustersize;
607 clustersize = 1;
608 }
609
610 return oldclustersize;
611
612}
613
614
615void MCalibrationTestCalc::LoopNeighbours( const TArrayI &arr, TArrayI &knownpixels, Int_t &clustersize, const Int_t idx )
616{
617
618 const MGeomPix &pix = (*fGeom)[idx];
619 const Byte_t neighbours = pix.GetNumNeighbors();
620
621 //
622 // Loop over the next neighbours:
623 // - Check if they are already in the list of known pixels
624 // - If not, call loopneighbours for the rest
625 // - grow clustersize for those
626 //
627 for (Int_t i=0;i<neighbours;i++)
628 {
629 const Int_t newid = pix.GetNeighbor(i);
630 Bool_t known = kFALSE;
631
632 for (Int_t j=knownpixels.GetSize()-1;j>=0;j--)
633 if (newid == knownpixels.At(j))
634 {
635 known = kTRUE;
636 break;
637 }
638 if (known)
639 continue;
640
641 for (Int_t k=0;k<arr.GetSize();k++)
642 if (newid == arr.At(k))
643 {
644 // This is an unknown, new pixel in the cluster!!
645 clustersize++;
646 const Int_t knownsize = knownpixels.GetSize();
647 knownpixels.Set(knownsize+1);
648 knownpixels[knownsize] = newid;
649 LoopNeighbours(arr, knownpixels, clustersize, newid);
650 }
651 }
652}
653
654// --------------------------------------------------------------------------
655//
656// Set the path for output file
657//
658void MCalibrationTestCalc::SetOutputPath(TString path)
659{
660 fOutputPath = path;
661 if (fOutputPath.EndsWith("/"))
662 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
663}
664
665void MCalibrationTestCalc::SetOutputFile(TString file)
666{
667 fOutputFile = file;
668}
669
670// --------------------------------------------------------------------------
671//
672// Get the output file
673//
674const char* MCalibrationTestCalc::GetOutputFile()
675{
676 return Form("%s/%s", (const char*)fOutputPath, (const char*)fOutputFile);
677}
678
Note: See TracBrowser for help on using the repository browser.