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

Last change on this file since 2765 was 2765, checked in by gaug, 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//
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 *fLog << pix->GetPixId() << " Pedestals: " << pix->GetPed() << " +- " << pix->GetPedRms()
261 << " Reduced Charge: " << pix->GetCharge() << " +- "
262 << pix->GetSigmaCharge() << " Reduced Sigma: " << TMath::Sqrt(pix->GetRSigmaSquare()) << endl;
263 id++;
264 }
265 }
266
267 *fLog << id << " succesful pixels :-))" << endl;
268 id = 0;
269
270 *fLog << endl;
271 *fLog << "Pixels with errors:" << endl;
272 *fLog << endl;
273
274 TIter Next2(fPixels);
275 while ((pix=(MCalibrationPix*)Next2()))
276 {
277
278 if (!pix->IsFitValid())
279 {
280 *fLog << pix->GetPixId() << " Pedestals: " << pix->GetPed() << " +- " << pix->GetPedRms()
281 << " Reduced Charge: " << pix->GetCharge() << " +- "
282 << pix->GetSigmaCharge() << " Reduced Sigma: " << TMath::Sqrt(pix->GetRSigmaSquare()) << endl;
283 id++;
284 }
285 }
286 *fLog << id << " pixels with errors :-((" << endl;
287
288}
289
290// The types are as follows:
291//
292// 0: Fitted Charge
293// 1: Error of fitted Charge
294// 2: Sigma of fitted Charge
295// 3: Error of Sigma of fitted Charge
296// 4: Returned probability of Gauss fit to Charge distribution
297// 5: Mean arrival time
298// 6: Sigma of the arrival time
299// 7: Chi-square of the Gauss fit to the arrival times
300// 8: Pedestal
301// 9: Pedestal RMS
302// 10: Reduced Sigma Square
303// 11: Number of Photo-electrons after the F-Factor method
304// 12: Error on the Number of Photo-electrons after the F-Factor method
305// 13: Mean conversion factor after the F-Factor method
306// 14: Error on the conversion factor after the F-Factor method
307// 15: Number of Photons after the Blind Pixel method
308// 16: Mean conversion factor after the Blind Pixel method
309//
310Bool_t MCalibrationCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
311{
312
313 if (idx > GetSize())
314 return kFALSE;
315
316 switch (type)
317 {
318 case 0:
319 val = (*this)[idx].GetCharge();
320 break;
321 case 1:
322 val = (*this)[idx].GetErrCharge();
323 break;
324 case 2:
325 val = (*this)[idx].GetSigmaCharge();
326 break;
327 case 3:
328 val = (*this)[idx].GetErrSigmaCharge();
329 break;
330 case 4:
331 val = (*this)[idx].GetChargeProb();
332 break;
333 case 5:
334 val = (*this)[idx].GetTime();
335 break;
336 case 6:
337 val = (*this)[idx].GetSigmaTime();
338 break;
339 case 7:
340 val = (*this)[idx].GetTimeProb();
341 break;
342 case 8:
343 val = (*this)[idx].GetPed();
344 break;
345 case 9:
346 val = (*this)[idx].GetPedRms();
347 break;
348 case 10:
349 if ((*this)[idx].GetRSigmaSquare() > 0.)
350 val = TMath::Sqrt((*this)[idx].GetRSigmaSquare());
351 else
352 val = -1.;
353 break;
354 case 11:
355 val = (*this)[idx].GetPheFFactorMethod();
356 break;
357 case 12:
358 val = (*this)[idx].GetPheFFactorMethodError();
359 break;
360 case 13:
361 val = (*this)[idx].GetMeanConversionFFactorMethod();
362 break;
363 case 14:
364 val = (*this)[idx].GetErrorConversionFFactorMethod();
365 break;
366 case 15:
367 if (idx < 397)
368 val = (double)fMeanPhotInsidePlexiglass;
369 else
370 val = (double)fMeanPhotInsidePlexiglass*gkCalibrationOutervsInnerPixelArea;
371 break;
372 case 16:
373 if (idx < 397)
374 val = (*this)[idx].GetMeanConversionBlindPixelMethod();
375 else
376 val = (*this)[idx].GetMeanConversionBlindPixelMethod()*gkCalibrationOutervsInnerPixelArea;
377 break;
378 case 17:
379 if ((*this)[idx].GetCharge() != 0.)
380 val = ((*this)[idx].GetSigmaCharge()/(*this)[idx].GetCharge())*
381 ((*this)[idx].GetSigmaCharge()/(*this)[idx].GetCharge());
382 else
383 val = -1.;
384 break;
385 default:
386 return kFALSE;
387 }
388 return val>=0;
389}
390
391// --------------------------------------------------------------------------
392//
393// What MHCamera needs in order to draw an individual pixel in the camera
394//
395void MCalibrationCam::DrawPixelContent(Int_t idx) const
396{
397 (*this)[idx].Draw();
398}
399
400
401// --------------------------------------------------------------------------
402//
403
404//
405Bool_t MCalibrationCam::CalcNumPhotInsidePlexiglass()
406{
407
408 if (!fBlindPixel->IsFitOK())
409 return kFALSE;
410
411 const Float_t mean = fBlindPixel->GetLambda();
412 const Float_t merr = fBlindPixel->GetErrLambda();
413
414 switch (fColor)
415 {
416 case kECGreen:
417 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEGreen) // real photons
418 *TMath::Power(10,gkCalibrationBlindPixelAttGreen) // correct for absorption
419 * gkCalibrationInnerPixelArea; // correct for area
420 break;
421 case kECBlue:
422 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEBlue )
423 *TMath::Power(10,gkCalibrationBlindPixelAttBlue)
424 * gkCalibrationInnerPixelArea;
425 break;
426 case kECUV:
427 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEUV )
428 *TMath::Power(10,gkCalibrationBlindPixelAttUV)
429 * gkCalibrationInnerPixelArea;
430 break;
431 case kECCT1:
432 default:
433 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQECT1 )
434 *TMath::Power(10,gkCalibrationBlindPixelAttCT1)
435 * gkCalibrationInnerPixelArea;
436 break;
437 }
438
439 fNumPhotInsidePlexiglassAvailable = kTRUE;
440
441 *fLog << endl;
442 *fLog << mean << " Mean number of Photons for an Inner Pixel: " << fMeanPhotInsidePlexiglass << endl;
443 *fLog << endl;
444
445 TIter Next(fPixels);
446 MCalibrationPix *pix;
447 while ((pix=(MCalibrationPix*)Next()))
448 {
449 if((pix->GetCharge() > 0.) && (fMeanPhotInsidePlexiglass > 0.))
450 {
451
452 Float_t conversion = fMeanPhotInsidePlexiglass/pix->GetCharge();
453 Float_t conversionerr = 0.;
454 Float_t conversionsigma = 0.;
455 pix->SetConversionBlindPixelMethod(conversion, conversionerr, conversionsigma);
456
457 if (conversionerr/conversion < 0.1)
458 pix->SetBlindPixelMethodValid();
459 }
460 }
461 return kTRUE;
462}
463
464
465Bool_t MCalibrationCam::CalcNumPhotOutsidePlexiglass()
466{
467
468 if (!fPINDiode->IsFitOK())
469 return kFALSE;
470
471 const Float_t mean = fPINDiode->GetMean();
472 const Float_t merr = fPINDiode->GetMeanError();
473
474 switch (fColor)
475 {
476 case kECGreen:
477 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEGreen) // real photons
478 * gkCalibrationInnerPixelvsPINDiodeArea; // correct for area
479 break;
480 case kECBlue:
481 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEBlue )
482 * gkCalibrationInnerPixelvsPINDiodeArea;
483 break;
484 case kECUV:
485 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEUV )
486 * gkCalibrationInnerPixelvsPINDiodeArea;
487 break;
488 case kECCT1:
489 default:
490 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQECT1 )
491 * gkCalibrationInnerPixelvsPINDiodeArea;
492 break;
493 }
494
495 fNumPhotOutsidePlexiglassAvailable = kTRUE;
496
497 TIter Next(fPixels);
498 MCalibrationPix *pix;
499 while ((pix=(MCalibrationPix*)Next()))
500 {
501
502 if((pix->GetCharge() > 0.) && (fMeanPhotInsidePlexiglass > 0.))
503 pix->SetConversionPINDiodeMethod(fMeanPhotOutsidePlexiglass/pix->GetCharge(), 0., 0.);
504 }
505 return kTRUE;
506}
507
508
509
510Bool_t MCalibrationCam::GetConversionFactorBlindPixel(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
511{
512
513 if (ipx < 0 || !IsPixelFitted(ipx))
514 return kFALSE;
515
516 if (!fNumPhotInsidePlexiglassAvailable)
517 if (!CalcNumPhotInsidePlexiglass())
518 return kFALSE;
519
520 mean = (*this)[ipx].GetMeanConversionBlindPixelMethod();
521 err = (*this)[ipx].GetErrorConversionBlindPixelMethod();
522 sigma = (*this)[ipx].GetSigmaConversionBlindPixelMethod();
523
524 return kTRUE;
525}
526
527
528Bool_t MCalibrationCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
529{
530
531 if (ipx < 0 || !IsPixelFitted(ipx))
532 return kFALSE;
533
534 Float_t conv = (*this)[ipx].GetMeanConversionFFactorMethod();
535
536 if (conv < 0.)
537 return kFALSE;
538
539 mean = conv;
540 err = (*this)[ipx].GetErrorConversionFFactorMethod();
541 sigma = (*this)[ipx].GetSigmaConversionFFactorMethod();
542
543 return kTRUE;
544}
545
546
547//-----------------------------------------------------------------------------------
548//
549// Calculates the conversion factor between the integral of FADCs slices
550// (as defined in the signal extractor MExtractSignal.cc)
551// and the number of photons reaching the plexiglass for one Inner Pixel
552//
553// FIXME: The PINDiode is still not working and so is the code
554//
555Bool_t MCalibrationCam::GetConversionFactorPINDiode(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 return kFALSE;
562
563}
564
565//-----------------------------------------------------------------------------------
566//
567// Calculates the best combination of the three used methods possible
568// between the integral of FADCs slices
569// (as defined in the signal extractor MExtractSignal.cc)
570// and the number of photons reaching one Inner Pixel.
571// The procedure is not yet defined.
572//
573// FIXME: The PINDiode is still not working and so is the code
574//
575Bool_t MCalibrationCam::GetConversionFactorCombined(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
576{
577
578 if (ipx < 0 || !IsPixelFitted(ipx))
579 return kFALSE;
580
581 return kFALSE;
582
583}
584
585
586void MCalibrationCam::DrawHiLoFits()
587{
588
589 if (!fOffsets)
590 fOffsets = new TH1D("pp","Offsets of the HiGain LoGain Fit",100,-600.,400.);
591 if (!fSlopes)
592 fSlopes = new TH1D("mm","Slopes of the HiGain LoGain Fit",100,-2.,2.);
593 if (!fOffvsSlope)
594 fOffvsSlope = new TH2D("aa","Slopes vs Offsets of the HiGain LoGain Fit",100,-600.,400.,100,-2.,2.);
595
596 TIter Next(fPixels);
597 MCalibrationPix *pix;
598 MHCalibrationPixel *hist;
599 while ((pix=(MCalibrationPix*)Next()))
600 {
601 hist = pix->GetHist();
602 hist->FitHiGainvsLoGain();
603 fOffsets->Fill(hist->GetOffset(),1.);
604 fSlopes->Fill(hist->GetSlope(),1.);
605 fOffvsSlope->Fill(hist->GetOffset(),hist->GetSlope(),1.);
606 }
607
608 TCanvas *c1 = new TCanvas();
609
610 c1->Divide(1,3);
611 c1->cd(1);
612 fOffsets->Draw();
613 gPad->Modified();
614 gPad->Update();
615
616 c1->cd(2);
617 fSlopes->Draw();
618 gPad->Modified();
619 gPad->Update();
620
621 c1->cd(3);
622 fOffvsSlope->Draw("col1");
623 gPad->Modified();
624 gPad->Update();
625}
626
Note: See TracBrowser for help on using the repository browser.