source: trunk/MagicSoft/Mars/mcalib/MCalibrationChargeCam.cc@ 3698

Last change on this file since 3698 was 3683, checked in by gaug, 21 years ago
*** empty log message ***
File size: 22.7 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! Author(s): Markus Gaug 11/2003 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24/////////////////////////////////////////////////////////////////////////////
25//
26// MCalibrationChargeCam
27//
28// Storage container for charge calibration results of the whole camera.
29//
30// Individual pixels have to be cast when retrieved e.g.:
31// MCalibrationChargePix &avpix = (MCalibrationChargePix&)(*fChargeCam)[i]
32//
33// The following calibration constants can be retrieved from each pixel:
34//
35// - GetConversionFactorFFactor ( Int_t idx, Float_t &mean, Float_t &err, Float_t &sigma );
36// - GetConversionFactorBlindPixel ( Int_t idx, Float_t &mean, Float_t &err, Float_t &sigma );
37// - GetConversionFactorPINDiode ( Int_t idx, Float_t &mean, Float_t &err, Float_t &sigma );
38// - GetConversionFactorCombined ( Int_t idx, Float_t &mean, Float_t &err, Float_t &sigma );
39//
40// where:
41// - idx is the pixel software ID
42// - "mean" is the mean conversion constant, to be multiplied with the retrieved signal
43// in order to get a calibrated number of PHOTONS.
44// - "err" is the pure statistical uncertainty about the MEAN
45// - "sigma", if mulitplied with the square root of signal, gives the approximate sigma of the
46// retrieved mean number of incident Cherenkov photons.
47//
48// Note, Conversion is ALWAYS (included the F-Factor method!) from summed FADC slices to PHOTONS.
49//
50// Averaged values over one whole area index (e.g. inner or outer pixels for
51// the MAGIC camera), can be retrieved via:
52// MCalibrationChargePix &avpix = (MCalibrationChargePix&)fChargeCam->GetAverageArea(i)
53//
54// Averaged values over one whole camera sector can be retrieved via:
55// MCalibrationChargePix &avpix = (MCalibrationChargePix&)fChargeCam->GetAverageSector(i)
56//
57// Note the averageing has been done on an event-by-event basis. Resulting
58// Sigma's of the Gauss fit have been multiplied with the square root of the number
59// of involved pixels in order to make a direct comparison possible with the mean of
60// sigmas.
61//
62// See also: MCalibrationChargePix, MCalibrationChargeCalc
63// MHCalibrationChargePix, MHCalibrationChargeCam
64//
65// The types of GetPixelContent() are as follows:
66//
67// Fitted values:
68// ==============
69//
70// 0: Fitted Charge (see MCalibrationPix::GetMean())
71// 1: Error of fitted Charge (see MCalibrationPix::GetMeanErr())
72// 2: Sigma of fitted Charge (see MCalibrationPix::GetSigma())
73// 3: Error of Sigma of fitted Charge (see MCalibrationPix::GetSigmaErr())
74//
75// Useful variables derived from the fit results:
76// =============================================
77//
78// 4: Probability Gauss fit Charge distribution (see MCalibrationPix::GetProb())
79// 5: Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigma())
80// 6: Error Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigmaErr())
81// 7: Reduced Sigma per Charge
82// 8: Error of Reduced Sigma per Charge
83//
84// Results of the different calibration methods:
85// =============================================
86//
87// 9: Nr. Photo-electrons from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethod())
88// 10: Error Nr. Photo-el. from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethodErr())
89// 11: Conversion factor from F-Factor Method (see MCalibrationChargePix::GetMeanConversionFFactorMethod())
90// 12: Error conv. factor from F-Factor Method (see MCalibrationChargePix::GetConversionFFactorMethodErr())
91// 13: Overall F-Factor from F-Factor Method (see MCalibrationChargePix::GetTotalFFactorFFactorMethod())
92// 14: Error F-Factor from F-Factor Method (see MCalibrationChargePix::GetTotalFFactorFFactorMethodErr())
93// 15: Pixels valid calibration F-Factor-Method (see MCalibrationChargePix::IsFFactorMethodValid())
94// 16: Conversion factor from BlindPixel Method (see MCalibrationChargePix::GetMeanConversionBlindPixelMethod())
95// 17: Error conv. factor from BlindPixel Method (see MCalibrationChargePix::GetConversionBlindPixelMethodErr())
96// 18: Overall F-Factor from BlindPixel Method (see MCalibrationChargePix::GetTotalFFactorBlindPixelMethod())
97// 19: Error F-Factor from BlindPixel Method (see MCalibrationChargePix::GetTotalFFactorBlindPixelMethodErr())
98// 20: Pixels valid calibration BlindPixel-Method (see MCalibrationChargePix::IsBlindPixelMethodValid())
99// 21: Conversion factor from PINDiode Method (see MCalibrationChargePix::GetMeanConversionPINDiodeMethod())
100// 22: Error conv. factor from PINDiode Method (see MCalibrationChargePix::GetConversionPINDiodeMethodErr())
101// 23: Overall F-Factor from PINDiode Method (see MCalibrationChargePix::GetTotalFFactorPINDiodeMethod())
102// 24: Error F-Factor from PINDiode Method (see MCalibrationChargePix::GetTotalFFactorPINDiodeMethodErr())
103// 25: Pixels valid calibration PINDiode-Method (see MCalibrationChargePix::IsPINDiodeMethodValid())
104//
105// Localized defects:
106// ==================
107//
108// 26: Excluded Pixels
109// 27: Number of pickup events in the Hi Gain (see MCalibrationPix::GetHiGainNumPickup())
110// 28: Number of pickup events in the Lo Gain (see MCalibrationPix::GetLoGainNumPickup())
111//
112// Other classifications of pixels:
113// ================================
114//
115// 29: Pixels with saturated Hi-Gain (see MCalibrationPix::IsHighGainSaturation())
116//
117// Used Pedestals:
118// ===============
119//
120// 30: Pedestal for entire signal extr. range (see MCalibrationChargePix::Ped())
121// 31: Error Pedestal entire signal extr. range (see MCalibrationChargePix::PedErr())
122// 32: Ped. RMS entire signal extraction range (see MCalibrationChargePix::PedRms())
123// 33: Error Ped. RMS entire signal extr. range (see MCalibrationChargePix::PedRmsErr())
124//
125// Calculated absolute arrival times (very low precision!):
126// ========================================================
127//
128// 34: Absolute Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeMean())
129// 35: RMS Ab. Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeRms())
130//
131/////////////////////////////////////////////////////////////////////////////
132#include "MCalibrationChargeCam.h"
133#include "MCalibrationCam.h"
134
135#include <TH2.h>
136#include <TCanvas.h>
137#include <TClonesArray.h>
138
139#include "MLog.h"
140#include "MLogManip.h"
141
142#include "MGeomCam.h"
143#include "MGeomPix.h"
144
145#include "MBadPixelsCam.h"
146#include "MBadPixelsPix.h"
147
148#include "MCalibrationChargePix.h"
149#include "MCalibrationChargeBlindPix.h"
150#include "MCalibrationChargePINDiode.h"
151
152ClassImp(MCalibrationChargeCam);
153
154using namespace std;
155// --------------------------------------------------------------------------
156//
157// Default constructor.
158//
159// Sets all pointers to 0
160//
161// Creates a TClonesArray of MCalibrationChargePix containers, initialized to 1 entry, destinated
162// to hold one container per pixel. Later, a call to MCalibrationChargeCam::InitSize()
163// has to be performed (in MGeomApply).
164//
165// Creates a TClonesArray of MCalibrationChargePix containers, initialized to 1 entry, destinated
166// to hold one container per pixel AREA. Later, a call to MCalibrationChargeCam::InitAreas()
167// has to be performed (in MGeomApply).
168//
169// Creates a TClonesArray of MCalibrationChargePix containers, initialized to 1 entry, destinated
170// to hold one container per camera SECTOR. Later, a call to MCalibrationChargeCam::InitSectors()
171// has to be performed (in MGeomApply).
172//
173// Calls:
174// - Clear()
175//
176MCalibrationChargeCam::MCalibrationChargeCam(const char *name, const char *title)
177 : fOffsets(NULL),
178 fSlopes(NULL),
179 fOffvsSlope(NULL)
180{
181 fName = name ? name : "MCalibrationChargeCam";
182 fTitle = title ? title : "Storage container for the Calibration Information in the camera";
183
184 fPixels = new TClonesArray("MCalibrationChargePix",1);
185 fAverageAreas = new TClonesArray("MCalibrationChargePix",1);
186 fAverageSectors = new TClonesArray("MCalibrationChargePix",1);
187
188 Clear();
189}
190
191// --------------------------------------------------------------------------
192//
193// Deletes the histograms if they exist
194//
195MCalibrationChargeCam::~MCalibrationChargeCam()
196{
197
198 if (fOffsets)
199 delete fOffsets;
200 if (fSlopes)
201 delete fSlopes;
202 if (fOffvsSlope)
203 delete fOffvsSlope;
204}
205
206
207// --------------------------------------
208//
209// Sets all variable to 0.
210// Sets all flags to kFALSE
211// Calls MCalibrationCam::Clear()
212//
213void MCalibrationChargeCam::Clear(Option_t *o)
214{
215
216 SetFFactorMethodValid ( kFALSE );
217
218 MCalibrationCam::Clear();
219
220 return;
221}
222
223
224// -----------------------------------------------
225//
226// Sets the kFFactorMethodValid bit from outside
227//
228void MCalibrationChargeCam::SetFFactorMethodValid(const Bool_t b)
229{
230 b ? SETBIT(fFlags, kFFactorMethodValid) : CLRBIT(fFlags, kFFactorMethodValid);
231}
232
233
234// --------------------------------------------------------------------------
235//
236// Test bit kFFactorMethodValid
237//
238Bool_t MCalibrationChargeCam::IsFFactorMethodValid() const
239{
240 return TESTBIT(fFlags,kFFactorMethodValid);
241}
242
243// --------------------------------------------------------------------------
244//
245// Print first the well fitted pixels
246// and then the ones which are not FitValid
247//
248void MCalibrationChargeCam::Print(Option_t *o) const
249{
250
251 *fLog << all << GetDescriptor() << ":" << endl;
252 int id = 0;
253
254 *fLog << all << "Calibrated pixels:" << endl;
255 *fLog << all << endl;
256
257 TIter Next(fPixels);
258 MCalibrationChargePix *pix;
259 while ((pix=(MCalibrationChargePix*)Next()))
260 {
261
262 if (!pix->IsExcluded())
263 {
264
265 *fLog << all << "Pix " << pix->GetPixId()
266 << ": Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
267 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetSigma()
268 << " Reduced Sigma: " << pix->GetRSigma()
269 << " Nr Phe's: " << pix->GetPheFFactorMethod()
270 << " Saturated? :" << pix->IsHiGainSaturation()
271 << endl;
272 id++;
273 }
274 }
275
276 *fLog << all << id << " pixels" << endl;
277 id = 0;
278
279
280 *fLog << all << endl;
281 *fLog << all << "Excluded pixels:" << endl;
282 *fLog << all << endl;
283
284 id = 0;
285
286 TIter Next4(fPixels);
287 while ((pix=(MCalibrationChargePix*)Next4()))
288 {
289 if (pix->IsExcluded())
290 {
291 *fLog << all << pix->GetPixId() << endl;
292 id++;
293 }
294 }
295 *fLog << all << id << " Excluded pixels " << endl;
296 *fLog << endl;
297
298 TIter Next5(fAverageAreas);
299 while ((pix=(MCalibrationChargePix*)Next5()))
300 {
301 *fLog << all << "Average Area:"
302 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
303 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
304 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
305 << " Reduced Sigma: " << pix->GetRSigma()
306 << " Nr Phe's: " << pix->GetPheFFactorMethod()
307 << endl;
308 }
309
310 TIter Next6(fAverageSectors);
311 while ((pix=(MCalibrationChargePix*)Next5()))
312 {
313 *fLog << all << "Average Sector:"
314 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
315 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
316 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
317 << " Reduced Sigma: " << pix->GetRSigma()
318 << " Nr Phe's: " << pix->GetPheFFactorMethod()
319 << endl;
320 }
321
322}
323
324
325// --------------------------------------------------------------------------
326//
327// The types are as follows:
328//
329// Fitted values:
330// ==============
331//
332// 0: Fitted Charge (see MCalibrationPix::GetMean())
333// 1: Error of fitted Charge (see MCalibrationPix::GetMeanErr())
334// 2: Sigma of fitted Charge (see MCalibrationPix::GetSigma())
335// 3: Error of Sigma of fitted Charge (see MCalibrationPix::GetSigmaErr())
336//
337// Useful variables derived from the fit results:
338// =============================================
339//
340// 4: Probability Gauss fit Charge distribution (see MCalibrationPix::GetProb())
341// 5: Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigma())
342// 6: Error Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigmaErr())
343// 7: Reduced Sigma per Charge
344// 8: Error of Reduced Sigma per Charge
345//
346// Results of the different calibration methods:
347// =============================================
348//
349// 9: Nr. Photo-electrons from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethod())
350// 10: Error Nr. Photo-el. from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethodErr())
351// 11: Conversion factor from F-Factor Method (see MCalibrationChargePix::GetMeanConversionFFactorMethod())
352// 12: Error conv. factor from F-Factor Method (see MCalibrationChargePix::GetConversionFFactorMethodErr())
353// 13: Overall F-Factor from F-Factor Method (see MCalibrationChargePix::GetTotalFFactorFFactorMethod())
354// 14: Error F-Factor from F-Factor Method (see MCalibrationChargePix::GetTotalFFactorFFactorMethodErr())
355// 15: Pixels valid calibration F-Factor-Method (see MCalibrationChargePix::IsFFactorMethodValid())
356// 16: Conversion factor from BlindPixel Method (see MCalibrationChargePix::GetMeanConversionBlindPixelMethod())
357// 17: Error conv. factor from BlindPixel Method (see MCalibrationChargePix::GetConversionBlindPixelMethodErr())
358// 18: Overall F-Factor from BlindPixel Method (see MCalibrationChargePix::GetTotalFFactorBlindPixelMethod())
359// 19: Error F-Factor from BlindPixel Method (see MCalibrationChargePix::GetTotalFFactorBlindPixelMethodErr())
360// 20: Pixels valid calibration BlindPixel-Method (see MCalibrationChargePix::IsBlindPixelMethodValid())
361// 21: Conversion factor from PINDiode Method (see MCalibrationChargePix::GetMeanConversionPINDiodeMethod())
362// 22: Error conv. factor from PINDiode Method (see MCalibrationChargePix::GetConversionPINDiodeMethodErr())
363// 23: Overall F-Factor from PINDiode Method (see MCalibrationChargePix::GetTotalFFactorPINDiodeMethod())
364// 24: Error F-Factor from PINDiode Method (see MCalibrationChargePix::GetTotalFFactorPINDiodeMethodErr())
365// 25: Pixels valid calibration PINDiode-Method (see MCalibrationChargePix::IsPINDiodeMethodValid())
366//
367// Localized defects:
368// ==================
369//
370// 26: Excluded Pixels
371// 27: Number of pickup events in the Hi Gain (see MCalibrationPix::GetHiGainNumPickup())
372// 28: Number of pickup events in the Lo Gain (see MCalibrationPix::GetLoGainNumPickup())
373//
374// Other classifications of pixels:
375// ================================
376//
377// 29: Pixels with saturated Hi-Gain (see MCalibrationPix::IsHighGainSaturation())
378//
379// Used Pedestals:
380// ===============
381//
382// 30: Pedestal for entire signal extr. range (see MCalibrationChargePix::Ped())
383// 31: Error Pedestal entire signal extr. range (see MCalibrationChargePix::PedErr())
384// 32: Ped. RMS entire signal extraction range (see MCalibrationChargePix::PedRms())
385// 33: Error Ped. RMS entire signal extr. range (see MCalibrationChargePix::PedRmsErr())
386//
387// Calculated absolute arrival times (very low precision!):
388// ========================================================
389//
390// 34: Absolute Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeMean())
391// 35: RMS Ab. Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeRms())
392//
393Bool_t MCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
394{
395
396 if (idx > GetSize())
397 return kFALSE;
398
399 Float_t area = cam[idx].GetA();
400
401 if (area == 0)
402 return kFALSE;
403
404 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[idx];
405
406 switch (type)
407 {
408 case 0:
409 if (pix.IsExcluded())
410 return kFALSE;
411 val = pix.GetMean();
412 break;
413 case 1:
414 if (pix.IsExcluded())
415 return kFALSE;
416 val = pix.GetMeanErr();
417 break;
418 case 2:
419 if (pix.IsExcluded())
420 return kFALSE;
421 val = pix.GetSigma();
422 break;
423 case 3:
424 if (pix.IsExcluded())
425 return kFALSE;
426 val = pix.GetSigmaErr();
427 break;
428 case 4:
429 if (pix.IsExcluded())
430 return kFALSE;
431 val = pix.GetProb();
432 break;
433 case 5:
434 if (pix.IsExcluded())
435 return kFALSE;
436 if (pix.GetRSigma() == -1.)
437 return kFALSE;
438 val = pix.GetRSigma();
439 break;
440 case 6:
441 if (pix.IsExcluded())
442 return kFALSE;
443 if (pix.GetRSigma() == -1.)
444 return kFALSE;
445 val = pix.GetRSigmaErr();
446 break;
447 case 7:
448 if (pix.IsExcluded())
449 return kFALSE;
450 if (pix.GetRSigma() == -1.)
451 return kFALSE;
452 if (pix.GetMean() == 0.)
453 return kFALSE;
454 val = pix.GetRSigma() / pix.GetMean();
455 break;
456 case 8:
457 if (pix.IsExcluded())
458 return kFALSE;
459 if (pix.GetRSigmaRelVar() <= 0.)
460 return kFALSE;
461 if (pix.GetMeanRelVar() <= 0.)
462 return kFALSE;
463 // relative error Rsigma square
464 val = pix.GetRSigmaRelVar() + pix.GetMeanRelVar();
465 //a calculate relative error out of squares
466 if (val < 0)
467 val = -1.;
468 else
469 val = TMath::Sqrt(val) * pix.GetRSigma() / pix.GetMean();
470 break;
471 case 9:
472 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
473 return kFALSE;
474 val = pix.GetPheFFactorMethod();
475 break;
476 case 10:
477 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
478 return kFALSE;
479 val = pix.GetPheFFactorMethodErr();
480 break;
481 case 11:
482 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
483 return kFALSE;
484 val = pix.GetMeanConvFADC2Phe();
485 break;
486 case 12:
487 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
488 return kFALSE;
489 val = pix.GetMeanConvFADC2PheErr();
490 break;
491 case 13:
492 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
493 return kFALSE;
494 val = pix.GetMeanFFactorFADC2Phot();
495 break;
496 case 14:
497 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
498 return kFALSE;
499 val = pix.GetMeanFFactorFADC2PhotErr();
500 break;
501 case 15:
502 if (pix.IsExcluded())
503 return kFALSE;
504 if (pix.IsFFactorMethodValid())
505 val = 1;
506 else
507 return kFALSE;
508 break;
509 case 16:
510 return kFALSE;
511 break;
512 case 17:
513 return kFALSE;
514 break;
515 case 18:
516 return kFALSE;
517 break;
518 case 19:
519 return kFALSE;
520 break;
521 case 20:
522 return kFALSE;
523 break;
524 case 21:
525 return kFALSE;
526 break;
527 case 22:
528 return kFALSE;
529 break;
530 case 23:
531 return kFALSE;
532 break;
533 case 24:
534 return kFALSE;
535 break;
536 case 25:
537 return kFALSE;
538 break;
539 case 26:
540 if (pix.IsExcluded())
541 val = 1.;
542 else
543 return kFALSE;
544 break;
545 case 27:
546 if (pix.IsExcluded())
547 return kFALSE;
548 val = pix.GetHiGainNumPickup();
549 break;
550 case 28:
551 if (pix.IsExcluded())
552 return kFALSE;
553 val = pix.GetLoGainNumPickup();
554 break;
555 case 29:
556 if (pix.IsExcluded())
557 return kFALSE;
558 val = pix.IsHiGainSaturation();
559 break;
560 case 30:
561 if (pix.IsExcluded())
562 return kFALSE;
563 val = pix.GetPed();
564 break;
565 case 31:
566 if (pix.IsExcluded())
567 return kFALSE;
568 val = pix.GetPedErr();
569 break;
570 case 32:
571 if (pix.IsExcluded())
572 return kFALSE;
573 val = pix.GetPedRms();
574 break;
575 case 33:
576 if (pix.IsExcluded())
577 return kFALSE;
578 val = pix.GetPedErr()/2.;
579 break;
580 case 34:
581 if (pix.IsExcluded())
582 return kFALSE;
583 val = pix.GetAbsTimeMean();
584 break;
585 case 35:
586 if (pix.IsExcluded())
587 return kFALSE;
588 val = pix.GetAbsTimeRms();
589 break;
590 default:
591 return kFALSE;
592 }
593
594 return val!=-1.;
595
596}
597
598// --------------------------------------------------------------------------
599//
600// Calls MCalibrationChargePix::DrawClone()
601//
602void MCalibrationChargeCam::DrawPixelContent(Int_t idx) const
603{
604 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[idx];
605 pix.DrawClone();
606}
607
608
609
610Bool_t MCalibrationChargeCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &ffactor)
611{
612
613 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[ipx];
614
615 Float_t conv = pix.GetMeanConvFADC2Phe();
616
617 if (conv < 0.)
618 return kFALSE;
619
620 mean = conv;
621 err = pix.GetMeanConvFADC2PheErr();
622 ffactor = pix.GetMeanFFactorFADC2Phot();
623
624 return kTRUE;
625}
626
627
628/*
629void MCalibrationChargeCam::DrawHiLoFits()
630{
631
632 if (!fOffsets)
633 fOffsets = new TH1D("pp","Offsets of the HiGain LoGain Fit",100,-600.,400.);
634 if (!fSlopes)
635 fSlopes = new TH1D("mm","Slopes of the HiGain LoGain Fit",100,-2.,2.);
636 if (!fOffvsSlope)
637 fOffvsSlope = new TH2D("aa","Slopes vs Offsets of the HiGain LoGain Fit",100,-600.,400.,100,-2.,2.);
638
639 TIter Next(fPixels);
640 MCalibrationPix *pix;
641 MHCalibrationPixel *hist;
642 while ((pix=(MCalibrationPix*)Next()))
643 {
644 hist = pix->GetHist();
645 hist->FitHiGainvsLoGain();
646 fOffsets->Fill(hist->GetOffset(),1.);
647 fSlopes->Fill(hist->GetSlope(),1.);
648 fOffvsSlope->Fill(hist->GetOffset(),hist->GetSlope(),1.);
649 }
650
651 TCanvas *c1 = new TCanvas();
652
653 c1->Divide(1,3);
654 c1->cd(1);
655 fOffsets->Draw();
656 gPad->Modified();
657 gPad->Update();
658
659 c1->cd(2);
660 fSlopes->Draw();
661 gPad->Modified();
662 gPad->Update();
663
664 c1->cd(3);
665 fOffvsSlope->Draw("col1");
666 gPad->Modified();
667 gPad->Update();
668}
669
670*/
Note: See TracBrowser for help on using the repository browser.