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

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