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

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