source: trunk/MagicSoft/Mars/mcalib/MCalibrationCam.cc@ 2969

Last change on this file since 2969 was 2956, checked in by gaug, 21 years ago
*** empty log message ***
File size: 19.4 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!
19! Author(s): Markus Gaug 11/2003 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2001
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MCalibrationCam
29//
30// Hold the whole Calibration results of the camera:
31//
32// 1) MCalibrationCam initializes a TClonesArray whose elements are
33// pointers to MCalibrationPix Containers
34// 2) It initializes a pointer to an MCalibrationBlindPix container
35// 3) It initializes a pointer to an MCalibrationPINDiode container
36//
37// 4)
38//
39/////////////////////////////////////////////////////////////////////////////
40#include "MCalibrationCam.h"
41
42#include <TH2.h>
43#include <TCanvas.h>
44#include <TClonesArray.h>
45
46#include "MLog.h"
47#include "MLogManip.h"
48
49#include "MGeomCam.h"
50
51#include "MCalibrationPix.h"
52#include "MCalibrationConfig.h"
53#include "MCalibrationBlindPix.h"
54#include "MCalibrationPINDiode.h"
55
56#include "MHCalibrationPixel.h"
57
58ClassImp(MCalibrationCam);
59
60using namespace std;
61// --------------------------------------------------------------------------
62//
63// Default constructor.
64//
65// Creates a TClonesArray of MCalibrationPix containers, initialized to 1 entry
66// Later, a call to MCalibrationCam::InitSize(Int_t size) has to be performed
67//
68// Creates an MCalibrationBlindPix container
69// Creates an MCalibrationPINDiode container
70//
71MCalibrationCam::MCalibrationCam(const char *name, const char *title)
72 : fOffsets(NULL),
73 fSlopes(NULL),
74 fOffvsSlope(NULL),
75 fBlindPixelId(559),
76 fPINDiodeId(9999)
77{
78 fName = name ? name : "MCalibrationCam";
79 fTitle = title ? title : "Storage container for the Calibration Information in the camera";
80
81 fPixels = new TClonesArray("MCalibrationPix",1);
82 fBlindPixel = new MCalibrationBlindPix();
83 fPINDiode = new MCalibrationPINDiode();
84
85 Clear();
86}
87
88// --------------------------------------------------------------------------
89//
90// Delete the TClonesArray of MCalibrationPix containers
91// Delete the MCalibrationPINDiode and the MCalibrationBlindPix
92//
93// Delete the histograms if they exist
94//
95MCalibrationCam::~MCalibrationCam()
96{
97
98 //
99 // delete fPixels should delete all Objects stored inside
100 //
101 delete fPixels;
102 delete fBlindPixel;
103 delete fPINDiode;
104
105 if (fOffsets)
106 delete fOffsets;
107 if (fSlopes)
108 delete fSlopes;
109 if (fOffvsSlope)
110 delete fOffvsSlope;
111
112}
113
114// -------------------------------------------------------------------
115//
116// This function simply allocates memory via the ROOT command:
117// (TObject**) TStorage::ReAlloc(fCont, newSize * sizeof(TObject*),
118// fSize * sizeof(TObject*));
119// newSize corresponds to size in our case
120// fSize is the old size (in most cases: 1)
121//
122void MCalibrationCam::InitSize(const UInt_t i)
123{
124
125 //
126 // check if we have already initialized to size
127 //
128 if (CheckBounds(i))
129 return;
130
131 fPixels->ExpandCreate(i);
132
133}
134
135// --------------------------------------------------------------------------
136//
137// This function returns the current size of the TClonesArray
138// independently if the MCalibrationPix is filled with values or not.
139//
140// It is the size of the array fPixels.
141//
142Int_t MCalibrationCam::GetSize() const
143{
144 return fPixels->GetEntriesFast();
145}
146
147// --------------------------------------------------------------------------
148//
149// Check if position i is inside the current bounds of the TClonesArray
150//
151Bool_t MCalibrationCam::CheckBounds(Int_t i) const
152{
153 return i < GetSize();
154}
155
156
157// --------------------------------------------------------------------------
158//
159// Get i-th pixel (pixel number)
160//
161MCalibrationPix &MCalibrationCam::operator[](Int_t i)
162{
163 return *static_cast<MCalibrationPix*>(fPixels->UncheckedAt(i));
164}
165
166// --------------------------------------------------------------------------
167//
168// Get i-th pixel (pixel number)
169//
170MCalibrationPix &MCalibrationCam::operator[](Int_t i) const
171{
172 return *static_cast<MCalibrationPix*>(fPixels->UncheckedAt(i));
173}
174
175
176// --------------------------------------
177//
178void MCalibrationCam::Clear(Option_t *o)
179{
180
181 fPixels->ForEach(TObject, Clear)();
182 fBlindPixel->Clear();
183 fPINDiode->Clear();
184
185 fMeanPhotInsidePlexiglass = -1.;
186 fMeanPhotErrInsidePlexiglass = -1.;
187 fMeanPhotOutsidePlexiglass = -1.;
188 fMeanPhotErrOutsidePlexiglass = -1.;
189
190 fNumExcludedPixels = 0;
191
192 CLRBIT(fFlags,kBlindPixelMethodValid);
193 CLRBIT(fFlags,kPINDiodeMethodValid);
194 CLRBIT(fFlags,kNumPhotInsidePlexiglassAvailable);
195 CLRBIT(fFlags,kNumPhotOutsidePlexiglassAvailable);
196
197 return;
198}
199
200void MCalibrationCam::SetBlindPixelMethodValid(const Bool_t b)
201{
202
203 if (b)
204 SETBIT(fFlags, kBlindPixelMethodValid);
205 else
206 CLRBIT(fFlags, kBlindPixelMethodValid);
207
208}
209
210void MCalibrationCam::SetPINDiodeMethodValid(const Bool_t b)
211{
212
213 if (b)
214 SETBIT(fFlags, kPINDiodeMethodValid);
215 else
216 CLRBIT(fFlags, kPINDiodeMethodValid);
217
218
219}
220
221Bool_t MCalibrationCam::IsBlindPixelMethodValid() const
222{
223 return TESTBIT(fFlags,kBlindPixelMethodValid);
224}
225
226Bool_t MCalibrationCam::IsPINDiodeMethodValid() const
227{
228 return TESTBIT(fFlags,kPINDiodeMethodValid);
229}
230
231
232Bool_t MCalibrationCam::IsNumPhotInsidePlexiglassAvailable() const
233{
234 return TESTBIT(fFlags,kNumPhotInsidePlexiglassAvailable);
235}
236
237Bool_t MCalibrationCam::IsNumPhotOutsidePlexiglassAvailable() const
238{
239 return TESTBIT(fFlags,kNumPhotOutsidePlexiglassAvailable);
240}
241
242
243
244// --------------------------------------------------------------------------
245//
246// Print first the well fitted pixels
247// and then the ones which are not FitValid
248//
249void MCalibrationCam::Print(Option_t *o) const
250{
251
252 *fLog << all << GetDescriptor() << ":" << endl;
253 int id = 0;
254
255 *fLog << all << "Succesfully calibrated pixels:" << endl;
256 *fLog << all << endl;
257
258 TIter Next(fPixels);
259 MCalibrationPix *pix;
260 while ((pix=(MCalibrationPix*)Next()))
261 {
262
263 if (pix->IsFitValid() && !pix->IsExcluded())
264 {
265
266 Float_t rsigma = pix->GetRSigmaSquare();
267 if (rsigma > 0.)
268 rsigma = TMath::Sqrt(rsigma);
269
270 *fLog << all << pix->GetPixId() << " Pedestals: " << pix->GetPed() << " +- "
271 << pix->GetPedRms() << " Reduced Charge: " << pix->GetCharge() << " +- "
272 << pix->GetSigmaCharge() << " Reduced Sigma: " << rsigma
273 << " Nr Phe's: " << pix->GetPheFFactorMethod() << endl;
274 id++;
275 }
276 }
277
278 *fLog << all << id << " succesful pixels :-))" << endl;
279 id = 0;
280
281 *fLog << all << endl;
282 *fLog << all << "Pixels with errors:" << endl;
283 *fLog << all << endl;
284
285 TIter Next2(fPixels);
286 while ((pix=(MCalibrationPix*)Next2()))
287 {
288
289 if (!pix->IsFitValid() && !pix->IsExcluded())
290 {
291
292 Float_t rsigma = pix->GetRSigmaSquare();
293 if (rsigma > 0.)
294 rsigma = TMath::Sqrt(rsigma);
295
296 *fLog << all << pix->GetPixId() << " Pedestals: " << pix->GetPed() << " +- "
297 << pix->GetPedRms() << " Reduced Charge: " << pix->GetCharge() << " +- "
298 << pix->GetSigmaCharge() << " Reduced Sigma: " << rsigma << endl;
299 id++;
300 }
301 }
302 *fLog << all << id << " pixels with errors :-((" << endl;
303
304 *fLog << all << endl;
305 *fLog << all << "Excluded pixels:" << endl;
306 *fLog << all << endl;
307
308 TIter Next3(fPixels);
309 while ((pix=(MCalibrationPix*)Next3()))
310 if (pix->IsExcluded())
311 *fLog << all << pix->GetPixId() << endl;
312
313 *fLog << all << fNumExcludedPixels << " excluded pixels " << endl;
314}
315
316// --------------------------------------------------------------------------
317//
318// Return true if pixel is inside bounds of the TClonesArray fPixels
319//
320Bool_t MCalibrationCam::IsPixelUsed(Int_t idx) const
321{
322 if (!CheckBounds(idx))
323 return kFALSE;
324
325 return kTRUE;
326}
327
328// --------------------------------------------------------------------------
329//
330// Return true if pixel has already been fitted once (independent of the result)
331//
332Bool_t MCalibrationCam::IsPixelFitted(Int_t idx) const
333{
334
335 if (!CheckBounds(idx))
336 return kFALSE;
337
338 return (*this)[idx].IsFitted();
339}
340
341// --------------------------------------------------------------------------
342//
343// Sets the user ranges of all histograms such that
344// empty bins at the edges are not used. Additionally, it rebins the
345// histograms such that in total, 50 bins are used.
346//
347void MCalibrationCam::CutEdges()
348{
349
350 fBlindPixel->GetHist()->CutAllEdges();
351 fPINDiode->GetHist()->CutAllEdges();
352
353 TIter Next(fPixels);
354 MCalibrationPix *pix;
355 while ((pix=(MCalibrationPix*)Next()))
356 {
357 pix->GetHist()->CutAllEdges();
358 }
359
360 return;
361}
362
363
364// The types are as follows:
365//
366// 0: Fitted Charge
367// 1: Error of fitted Charge
368// 2: Sigma of fitted Charge
369// 3: Error of Sigma of fitted Charge
370// 4: Returned probability of Gauss fit to Charge distribution
371// 5: Mean arrival time
372// 6: Sigma of the arrival time
373// 7: Chi-square of the Gauss fit to the arrival times
374// 8: Pedestal
375// 9: Pedestal RMS
376// 10: Reduced Sigma Square
377// 11: Number of Photo-electrons after the F-Factor method
378// 12: Error on the Number of Photo-electrons after the F-Factor method
379// 13: Mean conversion factor after the F-Factor method
380// 14: Error on the conversion factor after the F-Factor method
381// 15: Number of Photons after the Blind Pixel method
382// 16: Mean conversion factor after the Blind Pixel method
383//
384Bool_t MCalibrationCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
385{
386
387 if (idx > GetSize())
388 return kFALSE;
389
390 if ( (!(*this)[idx].IsFitValid()) || (*this)[idx].IsExcluded())
391 return kFALSE;
392
393 if (idx == fBlindPixelId)
394 return kFALSE;
395
396 if (idx == fPINDiodeId)
397 return kFALSE;
398
399 switch (type)
400 {
401 case 0:
402 val = (*this)[idx].GetCharge();
403 break;
404 case 1:
405 val = (*this)[idx].GetErrCharge();
406 break;
407 case 2:
408 val = (*this)[idx].GetSigmaCharge();
409 break;
410 case 3:
411 val = (*this)[idx].GetErrSigmaCharge();
412 break;
413 case 4:
414 val = (*this)[idx].GetChargeProb();
415 break;
416 case 5:
417 val = (*this)[idx].GetTime();
418 break;
419 case 6:
420 val = (*this)[idx].GetSigmaTime();
421 break;
422 case 7:
423 val = (*this)[idx].GetTimeChiSquare();
424 break;
425 case 8:
426 val = (*this)[idx].GetPed();
427 break;
428 case 9:
429 val = (*this)[idx].GetPedRms();
430 break;
431 case 10:
432 if ((*this)[idx].GetRSigmaSquare() > 0.)
433 val = TMath::Sqrt((*this)[idx].GetRSigmaSquare());
434 else
435 val = -1.;
436 break;
437 case 11:
438 val = (*this)[idx].GetPheFFactorMethod();
439 break;
440 case 12:
441 val = (*this)[idx].GetPheFFactorMethodError();
442 break;
443 case 13:
444 val = (*this)[idx].GetMeanConversionFFactorMethod();
445 break;
446 case 14:
447 val = (*this)[idx].GetErrorConversionFFactorMethod();
448 break;
449 case 15:
450 if (idx < 397)
451 val = (double)fMeanPhotInsidePlexiglass;
452 else
453 val = (double)fMeanPhotInsidePlexiglass*gkCalibrationOutervsInnerPixelArea;
454 break;
455 case 16:
456 if (idx < 397)
457 val = (*this)[idx].GetMeanConversionBlindPixelMethod();
458 else
459 val = (*this)[idx].GetMeanConversionBlindPixelMethod()*gkCalibrationOutervsInnerPixelArea;
460 break;
461 case 17:
462 if ( (*this)[idx].GetRSigmaSquare() > 0. && (*this)[idx].GetCharge() > 0. )
463 val = TMath::Sqrt((*this)[idx].GetRSigmaSquare()) / (*this)[idx].GetCharge();
464 else
465 val = -1.;
466 break;
467 default:
468 return kFALSE;
469 }
470 return val!=-1.;
471}
472
473// --------------------------------------------------------------------------
474//
475// What MHCamera needs in order to draw an individual pixel in the camera
476//
477void MCalibrationCam::DrawPixelContent(Int_t idx) const
478{
479 (*this)[idx].Draw();
480}
481
482
483// --------------------------------------------------------------------------
484//
485//
486//
487Bool_t MCalibrationCam::CalcNumPhotInsidePlexiglass()
488{
489
490 if (!fBlindPixel->IsFitOK())
491 return kFALSE;
492
493 const Float_t mean = fBlindPixel->GetLambda();
494 const Float_t merr = fBlindPixel->GetErrLambda();
495
496 switch (fColor)
497 {
498 case kECGreen:
499 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEGreen) // real photons
500 *TMath::Power(10,gkCalibrationBlindPixelAttGreen) // correct for absorption
501 * gkCalibrationInnerPixelArea; // correct for area
502
503
504 break;
505 case kECBlue:
506 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEBlue )
507 *TMath::Power(10,gkCalibrationBlindPixelAttBlue)
508 * gkCalibrationInnerPixelArea;
509 break;
510 case kECUV:
511 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEUV )
512 *TMath::Power(10,gkCalibrationBlindPixelAttUV)
513 * gkCalibrationInnerPixelArea;
514 break;
515 case kECCT1:
516 default:
517 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQECT1 )
518 *TMath::Power(10,gkCalibrationBlindPixelAttCT1)
519 * gkCalibrationInnerPixelArea;
520 break;
521 }
522
523 SETBIT(fFlags,kNumPhotInsidePlexiglassAvailable);
524
525 *fLog << inf << endl;
526 *fLog << inf << "Mean number of Photons for an Inner Pixel (inside Plexiglass): "
527 << fMeanPhotInsidePlexiglass << endl;
528
529 TIter Next(fPixels);
530 MCalibrationPix *pix;
531 while ((pix=(MCalibrationPix*)Next()))
532 {
533 if((pix->GetCharge() > 0.) && (fMeanPhotInsidePlexiglass > 0.))
534 {
535
536 Float_t conversion = fMeanPhotInsidePlexiglass/pix->GetCharge();
537 Float_t conversionerr = 0.;
538 Float_t conversionsigma = 0.;
539 pix->SetConversionBlindPixelMethod(conversion, conversionerr, conversionsigma);
540
541 if (conversionerr/conversion < 0.1)
542 pix->SetBlindPixelMethodValid();
543 }
544 }
545 return kTRUE;
546}
547
548
549Bool_t MCalibrationCam::CalcNumPhotOutsidePlexiglass()
550{
551
552 if (!fPINDiode->IsFitValid())
553 return kFALSE;
554
555 const Float_t mean = fPINDiode->GetCharge();
556 const Float_t merr = fPINDiode->GetErrCharge();
557
558 switch (fColor)
559 {
560 case kECGreen:
561 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEGreen) // real photons
562 * gkCalibrationInnerPixelvsPINDiodeArea; // correct for area
563 break;
564 case kECBlue:
565 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEBlue )
566 * gkCalibrationInnerPixelvsPINDiodeArea;
567 break;
568 case kECUV:
569 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEUV )
570 * gkCalibrationInnerPixelvsPINDiodeArea;
571 break;
572 case kECCT1:
573 default:
574 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQECT1 )
575 * gkCalibrationInnerPixelvsPINDiodeArea;
576 break;
577 }
578
579 SETBIT(fFlags,kNumPhotOutsidePlexiglassAvailable);
580
581 *fLog << inf << endl;
582 *fLog << inf << mean << " Mean number of Photons for an Inner Pixel (outside Plexiglass): "
583 << fMeanPhotOutsidePlexiglass << endl;
584 *fLog << inf << endl;
585
586 TIter Next(fPixels);
587 MCalibrationPix *pix;
588 while ((pix=(MCalibrationPix*)Next()))
589 {
590
591 if((pix->GetCharge() > 0.) && (fMeanPhotInsidePlexiglass > 0.))
592 pix->SetConversionPINDiodeMethod(fMeanPhotOutsidePlexiglass/pix->GetCharge(), 0., 0.);
593 }
594 return kTRUE;
595}
596
597
598
599Bool_t MCalibrationCam::GetConversionFactorBlindPixel(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
600{
601
602 if (ipx < 0 || !IsPixelFitted(ipx))
603 return kFALSE;
604
605 if (!IsNumPhotInsidePlexiglassAvailable())
606 if (!CalcNumPhotInsidePlexiglass())
607 return kFALSE;
608
609 mean = (*this)[ipx].GetMeanConversionBlindPixelMethod();
610 err = (*this)[ipx].GetErrorConversionBlindPixelMethod();
611 sigma = (*this)[ipx].GetSigmaConversionBlindPixelMethod();
612
613 return kTRUE;
614}
615
616
617Bool_t MCalibrationCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
618{
619
620 if (ipx < 0 || !IsPixelFitted(ipx))
621 return kFALSE;
622
623 Float_t conv = (*this)[ipx].GetMeanConversionFFactorMethod();
624
625 if (conv < 0.)
626 return kFALSE;
627
628 mean = conv;
629 err = (*this)[ipx].GetErrorConversionFFactorMethod();
630 sigma = (*this)[ipx].GetSigmaConversionFFactorMethod();
631
632 return kTRUE;
633}
634
635
636//-----------------------------------------------------------------------------------
637//
638// Calculates the conversion factor between the integral of FADCs slices
639// (as defined in the signal extractor MExtractSignal.cc)
640// and the number of photons reaching the plexiglass for one Inner Pixel
641//
642// FIXME: The PINDiode is still not working and so is the code
643//
644Bool_t MCalibrationCam::GetConversionFactorPINDiode(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
645{
646
647 if (ipx < 0 || !IsPixelFitted(ipx))
648 return kFALSE;
649
650 if (!IsNumPhotOutsidePlexiglassAvailable())
651 if (!CalcNumPhotOutsidePlexiglass())
652 return kFALSE;
653
654 mean = (*this)[ipx].GetMeanConversionPINDiodeMethod();
655 err = (*this)[ipx].GetErrorConversionPINDiodeMethod();
656 sigma = (*this)[ipx].GetSigmaConversionPINDiodeMethod();
657
658 return kFALSE;
659
660}
661
662//-----------------------------------------------------------------------------------
663//
664// Calculates the best combination of the three used methods possible
665// between the integral of FADCs slices
666// (as defined in the signal extractor MExtractSignal.cc)
667// and the number of photons reaching one Inner Pixel.
668// The procedure is not yet defined.
669//
670// FIXME: The PINDiode is still not working and so is the code
671//
672Bool_t MCalibrationCam::GetConversionFactorCombined(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
673{
674
675 if (ipx < 0 || !IsPixelFitted(ipx))
676 return kFALSE;
677
678 return kFALSE;
679
680}
681
682
683void MCalibrationCam::DrawHiLoFits()
684{
685
686 if (!fOffsets)
687 fOffsets = new TH1D("pp","Offsets of the HiGain LoGain Fit",100,-600.,400.);
688 if (!fSlopes)
689 fSlopes = new TH1D("mm","Slopes of the HiGain LoGain Fit",100,-2.,2.);
690 if (!fOffvsSlope)
691 fOffvsSlope = new TH2D("aa","Slopes vs Offsets of the HiGain LoGain Fit",100,-600.,400.,100,-2.,2.);
692
693 TIter Next(fPixels);
694 MCalibrationPix *pix;
695 MHCalibrationPixel *hist;
696 while ((pix=(MCalibrationPix*)Next()))
697 {
698 hist = pix->GetHist();
699 hist->FitHiGainvsLoGain();
700 fOffsets->Fill(hist->GetOffset(),1.);
701 fSlopes->Fill(hist->GetSlope(),1.);
702 fOffvsSlope->Fill(hist->GetOffset(),hist->GetSlope(),1.);
703 }
704
705 TCanvas *c1 = new TCanvas();
706
707 c1->Divide(1,3);
708 c1->cd(1);
709 fOffsets->Draw();
710 gPad->Modified();
711 gPad->Update();
712
713 c1->cd(2);
714 fSlopes->Draw();
715 gPad->Modified();
716 gPad->Update();
717
718 c1->cd(3);
719 fOffvsSlope->Draw("col1");
720 gPad->Modified();
721 gPad->Update();
722}
723
Note: See TracBrowser for help on using the repository browser.