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

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