source: trunk/MagicSoft/Mars/manalysis/MCalibrationCam.cc@ 2716

Last change on this file since 2716 was 2715, checked in by gaug, 22 years ago
*** empty log message ***
File size: 16.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// 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// Return TRUE if:
211//
212// 1) Pixel is inside the range of the TClonesArray
213// 2) Pixel has a fitted charge greater than 3*PedRMS
214// 3) Pixel has a charge sigma bigger than its Pedestal RMS
215// 4) Pixel has a fit error greater than 0.
216// 5) Pixel has a fit Probability greater than 0.0001
217// 6) If FitTimes is used,
218// the mean arrival time is at least 1.5 slices from the edge
219//
220Bool_t MCalibrationCam::IsPixelFitValid(Int_t idx) const
221{
222
223 if (!CheckBounds(idx))
224 return kFALSE;
225
226 return (*this)[idx].IsFitValid();
227}
228
229// --------------------------------------
230//
231void MCalibrationCam::Clear(Option_t *o)
232{
233 fPixels->ForEach(TObject, Clear)();
234}
235
236// --------------------------------------------------------------------------
237//
238// Sets the user ranges of all histograms such that
239// empty bins at the edges are not used. Additionally, it rebins the
240// histograms such that in total, 50 bins are used.
241//
242void MCalibrationCam::CutEdges()
243{
244
245 fBlindPixel->GetHist()->CutAllEdges();
246 fPINDiode->GetHist()->CutAllEdges();
247
248 TIter Next(fPixels);
249 MCalibrationPix *pix;
250 while ((pix=(MCalibrationPix*)Next()))
251 {
252 pix->GetHist()->CutAllEdges();
253 }
254
255 return;
256}
257
258// --------------------------------------------------------------------------
259//
260//
261//
262//
263//
264void MCalibrationCam::Print(Option_t *o) const
265{
266
267 *fLog << all << GetDescriptor() << ":" << endl;
268 int id = 0;
269
270 *fLog << "Succesfully calibrated pixels:" << endl;
271 *fLog << endl;
272
273 TIter Next(fPixels);
274 MCalibrationPix *pix;
275 while ((pix=(MCalibrationPix*)Next()))
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
287 *fLog << id << " succesful pixels :-))" << endl;
288 id = 0;
289
290 *fLog << endl;
291 *fLog << "Pixels with errors:" << endl;
292 *fLog << endl;
293
294 TIter Next2(fPixels);
295 while ((pix=(MCalibrationPix*)Next2()))
296 {
297
298 if (!pix->IsFitValid())
299 {
300 *fLog << pix->GetPixId() << " Pedestals: " << pix->GetPed() << " +- " << pix->GetPedRms()
301 << " Reduced Charge: " << pix->GetCharge() << " +- "
302 << pix->GetSigmaCharge() << " Reduced Sigma: " << TMath::Sqrt(pix->GetRSigmaSquare()) << endl;
303 id++;
304 }
305 }
306 *fLog << id << " pixels with errors :-((" << endl;
307
308}
309
310
311Bool_t MCalibrationCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
312{
313
314 if (idx > GetSize())
315 return kFALSE;
316
317 switch (type)
318 {
319 case 0:
320 val = (*this)[idx].GetCharge();
321 break;
322 case 1:
323 val = (*this)[idx].GetErrCharge();
324 break;
325 case 2:
326 val = (*this)[idx].GetSigmaCharge();
327 break;
328 case 3:
329 val = (*this)[idx].GetErrSigmaCharge();
330 break;
331 case 4:
332 val = (*this)[idx].GetChargeProb();
333 break;
334 case 5:
335 val = (*this)[idx].GetTime();
336 break;
337 case 6:
338 val = (*this)[idx].GetSigmaTime();
339 break;
340 case 7:
341 val = (*this)[idx].GetTimeChiSquare();
342 break;
343 case 8:
344 val = (*this)[idx].GetPed();
345 break;
346 case 9:
347 val = (*this)[idx].GetPedRms();
348 break;
349 case 10:
350 if ((*this)[idx].GetRSigmaSquare() > 0.)
351 val = TMath::Sqrt((*this)[idx].GetRSigmaSquare());
352 else
353 val = -1.;
354 break;
355 case 11:
356 val = (*this)[idx].GetPheFFactorMethod();
357 break;
358 case 12:
359 val = (*this)[idx].GetPheFFactorMethodError();
360 break;
361 case 13:
362 val = (*this)[idx].GetMeanConversionFFactorMethod();
363 break;
364 case 14:
365 val = (*this)[idx].GetErrorConversionFFactorMethod();
366 break;
367 case 15:
368 if (idx < 397)
369 val = (double)fMeanPhotInsidePlexiglass;
370 else
371 val = (double)fMeanPhotInsidePlexiglass*gkCalibrationOutervsInnerPixelArea;
372 break;
373 case 16:
374 if (idx < 397)
375 val = (*this)[idx].GetMeanConversionBlindPixelMethod();
376 else
377 val = (*this)[idx].GetMeanConversionBlindPixelMethod()*gkCalibrationOutervsInnerPixelArea;
378 break;
379 case 17:
380 if ((*this)[idx].GetCharge() != 0.)
381 val = ((*this)[idx].GetSigmaCharge()/(*this)[idx].GetCharge())*
382 ((*this)[idx].GetSigmaCharge()/(*this)[idx].GetCharge());
383 else
384 val = -1.;
385 break;
386 default:
387 return kFALSE;
388 }
389 return val>=0;
390}
391
392void MCalibrationCam::DrawPixelContent(Int_t idx) const
393{
394
395 (*this)[idx].Draw();
396
397}
398
399Bool_t MCalibrationCam::CalcNumPhotInsidePlexiglass()
400{
401
402 if (!fBlindPixel->IsFitOK())
403 return kFALSE;
404
405 const Float_t mean = fBlindPixel->GetLambda();
406 const Float_t merr = fBlindPixel->GetErrLambda();
407
408 switch (fColor)
409 {
410 case kECGreen:
411 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEGreen) // real photons
412 *TMath::Power(10,gkCalibrationBlindPixelAttGreen) // correct for absorption
413 * gkCalibrationInnerPixelArea; // correct for area
414 break;
415 case kECBlue:
416 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEBlue )
417 *TMath::Power(10,gkCalibrationBlindPixelAttBlue)
418 * gkCalibrationInnerPixelArea;
419 break;
420 case kECUV:
421 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQEUV )
422 *TMath::Power(10,gkCalibrationBlindPixelAttUV)
423 * gkCalibrationInnerPixelArea;
424 break;
425 case kECCT1:
426 default:
427 fMeanPhotInsidePlexiglass = (mean / gkCalibrationBlindPixelQECT1 )
428 *TMath::Power(10,gkCalibrationBlindPixelAttCT1)
429 * gkCalibrationInnerPixelArea;
430 break;
431 }
432
433 fNumPhotInsidePlexiglassAvailable = kTRUE;
434
435 *fLog << endl;
436 *fLog << mean << " Mean number of Photons for an Inner Pixel: " << fMeanPhotInsidePlexiglass << endl;
437 *fLog << endl;
438
439 TIter Next(fPixels);
440 MCalibrationPix *pix;
441 while ((pix=(MCalibrationPix*)Next()))
442 {
443 if((pix->GetCharge() > 0.) && (fMeanPhotInsidePlexiglass > 0.))
444 pix->SetConversionBlindPixelMethod(fMeanPhotInsidePlexiglass/pix->GetCharge(), 0., 0.);
445 }
446 return kTRUE;
447}
448
449
450Bool_t MCalibrationCam::CalcNumPhotOutsidePlexiglass()
451{
452
453 if (!fPINDiode->IsFitOK())
454 return kFALSE;
455
456 const Float_t mean = fPINDiode->GetMean();
457 const Float_t merr = fPINDiode->GetMeanError();
458
459 switch (fColor)
460 {
461 case kECGreen:
462 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEGreen) // real photons
463 * gkCalibrationInnerPixelvsPINDiodeArea; // correct for area
464 break;
465 case kECBlue:
466 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEBlue )
467 * gkCalibrationInnerPixelvsPINDiodeArea;
468 break;
469 case kECUV:
470 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQEUV )
471 * gkCalibrationInnerPixelvsPINDiodeArea;
472 break;
473 case kECCT1:
474 default:
475 fMeanPhotOutsidePlexiglass = (mean / gkCalibrationPINDiodeQECT1 )
476 * gkCalibrationInnerPixelvsPINDiodeArea;
477 break;
478 }
479
480 fNumPhotOutsidePlexiglassAvailable = kTRUE;
481
482 TIter Next(fPixels);
483 MCalibrationPix *pix;
484 while ((pix=(MCalibrationPix*)Next()))
485 {
486
487 if((pix->GetCharge() > 0.) && (fMeanPhotInsidePlexiglass > 0.))
488 pix->SetConversionPINDiodeMethod(fMeanPhotOutsidePlexiglass/pix->GetCharge(), 0., 0.);
489 }
490 return kTRUE;
491}
492
493
494
495Bool_t MCalibrationCam::GetConversionFactorBlindPixel(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
496{
497
498 if (ipx < 0 || !IsPixelFitted(ipx))
499 return kFALSE;
500
501 if (!fNumPhotInsidePlexiglassAvailable)
502 if (!CalcNumPhotInsidePlexiglass())
503 return kFALSE;
504
505 mean = (*this)[ipx].GetMeanConversionBlindPixelMethod();
506 err = (*this)[ipx].GetErrorConversionBlindPixelMethod();
507 sigma = (*this)[ipx].GetSigmaConversionBlindPixelMethod();
508
509 return kTRUE;
510}
511
512
513Bool_t MCalibrationCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
514{
515
516 if (ipx < 0 || !IsPixelFitted(ipx))
517 return kFALSE;
518
519 Float_t conv = (*this)[ipx].GetMeanConversionFFactorMethod();
520
521 if (conv < 0.)
522 return kFALSE;
523
524 mean = conv;
525 err = (*this)[ipx].GetErrorConversionFFactorMethod();
526 sigma = (*this)[ipx].GetSigmaConversionFFactorMethod();
527
528 return kTRUE;
529}
530
531
532//-----------------------------------------------------------------------------------
533//
534// Calculates the conversion factor between the integral of FADCs slices
535// (as defined in the signal extractor MExtractSignal.cc)
536// and the number of photons reaching the plexiglass for one Inner Pixel
537//
538// FIXME: The PINDiode is still not working and so is the code
539//
540Bool_t MCalibrationCam::GetConversionFactorPINDiode(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
541{
542
543 if (ipx < 0 || !IsPixelFitted(ipx))
544 return kFALSE;
545
546 return kFALSE;
547
548}
549
550//-----------------------------------------------------------------------------------
551//
552// Calculates the best combination of the three used methods possible
553// between the integral of FADCs slices
554// (as defined in the signal extractor MExtractSignal.cc)
555// and the number of photons reaching one Inner Pixel.
556// The procedure is not yet defined.
557//
558// FIXME: The PINDiode is still not working and so is the code
559//
560Bool_t MCalibrationCam::GetConversionFactorCombined(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
561{
562
563 if (ipx < 0 || !IsPixelFitted(ipx))
564 return kFALSE;
565
566 return kFALSE;
567
568}
569
570
571void MCalibrationCam::DrawHiLoFits()
572{
573
574 if (!fOffsets)
575 fOffsets = new TH1D("pp","Offsets of the HiGain LoGain Fit",100,-600.,400.);
576 if (!fSlopes)
577 fSlopes = new TH1D("mm","Slopes of the HiGain LoGain Fit",100,-2.,2.);
578 if (!fOffvsSlope)
579 fOffvsSlope = new TH2D("aa","Slopes vs Offsets of the HiGain LoGain Fit",100,-600.,400.,100,-2.,2.);
580
581 TIter Next(fPixels);
582 MCalibrationPix *pix;
583 MHCalibrationPixel *hist;
584 while ((pix=(MCalibrationPix*)Next()))
585 {
586 hist = pix->GetHist();
587 hist->FitHiGainvsLoGain();
588 fOffsets->Fill(hist->GetOffset(),1.);
589 fSlopes->Fill(hist->GetSlope(),1.);
590 fOffvsSlope->Fill(hist->GetOffset(),hist->GetSlope(),1.);
591 }
592
593 TCanvas *c1 = new TCanvas();
594
595 c1->Divide(1,3);
596 c1->cd(1);
597 fOffsets->Draw();
598 gPad->Modified();
599 gPad->Update();
600
601 c1->cd(2);
602 fSlopes->Draw();
603 gPad->Modified();
604 gPad->Update();
605
606 c1->cd(3);
607 fOffvsSlope->Draw("col1");
608 gPad->Modified();
609 gPad->Update();
610}
611
Note: See TracBrowser for help on using the repository browser.