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

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