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

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