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

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