source: branches/Mars_McMismatchStudy/mcalib/MCalibrationIntensityConstCam.cc@ 18066

Last change on this file since 18066 was 7804, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 11.2 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 06/2005 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2005
21!
22!
23\* ======================================================================== */
24/////////////////////////////////////////////////////////////////////////////
25//
26// MCalibrationIntensityConstCam
27//
28// Storage container for intensity charge calibration results.
29//
30// Individual MCalibrationConstCam's can be retrieved with:
31// - GetCam() yielding the current cam.
32// - GetCam("name") yielding the current camera with name "name".
33// - GetCam(i) yielding the i-th camera.
34//
35// See also: MCalibrationIntensityCam, MCalibConstCam,
36// MCalibrationChargePix, MCalibrationChargeCalc, MCalibrationQECam
37// MHCalibrationPix, MHCalibrationCam
38//
39/////////////////////////////////////////////////////////////////////////////
40#include "MCalibrationIntensityConstCam.h"
41#include "MCalibConstCam.h"
42#include "MCalibConstPix.h"
43
44#include "MGeomCam.h"
45#include "MGeomPix.h"
46
47#include "MBadPixelsCam.h"
48#include "MBadPixelsPix.h"
49
50#include "MLogManip.h"
51#include "MHCamera.h"
52
53#include <TH1D.h>
54#include <TF1.h>
55#include <TOrdCollection.h>
56#include <TGraphErrors.h>
57
58ClassImp(MCalibrationIntensityConstCam);
59
60using namespace std;
61// --------------------------------------------------------------------------
62//
63// Default constructor.
64//
65MCalibrationIntensityConstCam::MCalibrationIntensityConstCam(const char *name, const char *title)
66{
67
68 fName = name ? name : "MCalibrationIntensityConstCam";
69 fTitle = title ? title : "Results of the Intensity Calibration";
70
71 fCams = new TOrdCollection;
72 fCams->SetOwner();
73
74 InitSize(1);
75}
76
77// --------------------------------------------------------------------------
78//
79// Deletes the histograms if they exist
80//
81MCalibrationIntensityConstCam::~MCalibrationIntensityConstCam()
82{
83 if (fCams)
84 delete fCams;
85}
86
87// -------------------------------------------------------------------
88//
89// Add MCalibrationConstCam's in the ranges from - to.
90//
91void MCalibrationIntensityConstCam::Add(const UInt_t from, const UInt_t to)
92{
93 for (UInt_t i=from; i<to; i++)
94 fCams->AddAt(new MCalibConstCam,i);
95}
96
97// --------------------------------------------------------------------------
98//
99// Add a new MCalibConstCam to fCams, give it the name "name" and initialize
100// it with geom.
101//
102void MCalibrationIntensityConstCam::AddToList( const char* name, const MGeomCam &geom)
103{
104 InitSize(GetSize()+1);
105 GetCam()->SetName(name);
106 GetCam()->Init(geom);
107}
108
109// -----------------------------------------------------
110//
111// Calls Clear() for all entries fCams
112//
113void MCalibrationIntensityConstCam::Clear(Option_t *o)
114{
115 fCams->R__FOR_EACH(MCalibConstCam, Clear)();
116}
117
118// -----------------------------------------------------
119//
120// Calls Print(o) for all entries fCams
121//
122void MCalibrationIntensityConstCam::Print(Option_t *o) const
123{
124 fCams->R__FOR_EACH(MCalibConstCam, Print)(o);
125}
126
127
128// -------------------------------------------------------------------
129//
130// Initialize the objects inside the TOrdCollection using the
131// virtual function Add().
132//
133// InitSize can only increase the size, but not shrink.
134//
135// It can be called more than one time. New Containers are
136// added only from the current size to the argument i.
137//
138void MCalibrationIntensityConstCam::InitSize(const UInt_t i)
139{
140
141 const UInt_t save = GetSize();
142
143 if (i==save)
144 return;
145
146 if (i>save)
147 Add(save,i);
148}
149
150// -------------------------------------------------------------------
151//
152// If size is still 0, Intialize a first Cam.
153// Calls Init(geom) for all fCams
154//
155void MCalibrationIntensityConstCam::Init(const MGeomCam &geom)
156{
157 if (GetSize() == 0)
158 InitSize(1);
159
160 fCams->R__FOR_EACH(MCalibConstCam,Init)(geom);
161}
162
163
164
165// --------------------------------------------------------------------------
166//
167// Returns the current size of the TOrdCollection fCams
168// independently if the MCalibrationCam is filled with values or not.
169//
170const Int_t MCalibrationIntensityConstCam::GetSize() const
171{
172 return fCams->GetSize();
173}
174
175// --------------------------------------------------------------------------
176//
177// Get i-th pixel from current camera
178//
179MCalibConstPix &MCalibrationIntensityConstCam::operator[](UInt_t i)
180{
181 return (*GetCam())[i];
182}
183
184// --------------------------------------------------------------------------
185//
186// Get i-th pixel from current camera
187//
188const MCalibConstPix &MCalibrationIntensityConstCam::operator[](UInt_t i) const
189{
190 return (*GetCam())[i];
191}
192
193// --------------------------------------------------------------------------
194//
195// Get i-th camera
196//
197MCalibConstCam *MCalibrationIntensityConstCam::GetCam(Int_t i)
198{
199 return static_cast<MCalibConstCam*>(i==-1 ? fCams->Last() : fCams->At(i));
200}
201
202// --------------------------------------------------------------------------
203//
204// Get i-th camera
205//
206const MCalibConstCam *MCalibrationIntensityConstCam::GetCam(Int_t i) const
207{
208 return static_cast<MCalibConstCam*>(i==-1 ? fCams->Last() : fCams->At(i));
209}
210
211Bool_t MCalibrationIntensityConstCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
212{
213 if ((*fBadPixels)[idx].IsUnsuitable())
214 return kFALSE;
215
216 return GetCam()->GetPixelContent(val,idx,cam,type);
217}
218
219TGraphErrors *MCalibrationIntensityConstCam::GetConvFactorPerAreaVsTime( const Int_t aidx, const MGeomCam &geom)
220{
221
222 const Int_t size = GetSize();
223
224 if (size == 0)
225 return NULL;
226
227 TArrayF convarea(size);
228 TArrayF convareaerr(size);
229 TArrayF time(size);
230 TArrayF timeerr(size);
231
232 TH1D *h = 0;
233 MCalibConstCam *cam = 0;
234
235 for (Int_t i=0;i<GetSize();i++)
236 {
237 //
238 // Get the calibration cam from the intensity cam
239 //
240 cam = GetCam(i);
241
242 //
243 // Get the calibration pix from the calibration cam
244 //
245 Double_t variab = 0.;
246 Double_t variab2 = 0.;
247 Double_t variance = 0.;
248 Int_t num = 0;
249 Float_t pvar = 0.;
250
251 MHCamera camconv(geom,"CamConv","Variable;;channels");
252 //
253 // Get the area calibration pix from the calibration cam
254 //
255 for (Int_t j=0; j<cam->GetSize(); j++)
256 {
257 const MCalibConstPix &pix = (*cam)[j];
258 //
259 // Use only pixels of a same area index
260 //
261 if (aidx != geom[j].GetAidx())
262 continue;
263
264 pvar = pix.GetCalibConst();
265
266 if (pvar < 0)
267 continue;
268
269 variab += pvar;
270 variab2 += pvar*pvar;
271 num++;
272
273 camconv.Fill(j,pvar);
274 camconv.SetUsed(j);
275 }
276
277 if (num > 1)
278 {
279 variab /= num;
280 variance = (variab2 - variab*variab*num) / (num-1);
281
282 convarea[i] = variab;
283 if (variance > 0.)
284 convareaerr[i] = TMath::Sqrt(variance);
285 else
286 convareaerr[i] = 999999999.;
287
288 //
289 // Make also a Gauss-fit to the distributions. The RMS can be determined by
290 // outlier, thus we look at the sigma and the RMS and take the smaller one, afterwards.
291 //
292 h = camconv.ProjectionS(TArrayI(),TArrayI(1,&aidx),"_py",500);
293 h->SetDirectory(NULL);
294 h->Fit("gaus","QL");
295 TF1 *fit = h->GetFunction("gaus");
296
297 Float_t ci2 = fit->GetChisquare();
298 Float_t sigma = fit->GetParameter(2);
299
300 if (ci2 > 500. || sigma > convareaerr[i])
301 {
302 h->Fit("gaus","QLM");
303 fit = h->GetFunction("gaus");
304
305 ci2 = fit->GetChisquare();
306 sigma = fit->GetParameter(2);
307 }
308
309 const Float_t mean = fit->GetParameter(1);
310 const Float_t meanerr = fit->GetParError(1);
311 const Float_t ndf = fit->GetNDF();
312
313 *fLog << inf << "Camera Nr: " << i << endl;
314 *fLog << inf << " area idx: " << aidx << " Results: " << endl;
315 *fLog << inf << "Mean: " << Form("%4.3f",mean)
316 << "+-" << Form("%4.3f",meanerr)
317 << " Sigma: " << Form("%4.3f",sigma) << "+-" << Form("%4.3f",fit->GetParError(2))
318 << " Chisquare: " << Form("%4.3f",ci2) << " NDF : " << ndf << endl;
319 delete h;
320 gROOT->GetListOfFunctions()->Remove(fit);
321
322 if (meanerr < convareaerr[i] && ndf > 2 && ci2/ndf < 10.)
323 {
324 convarea [i] = mean;
325 convareaerr[i] = meanerr;
326 }
327 else
328 {
329 convareaerr[i] /= TMath::Sqrt((Float_t)num);
330 }
331 }
332 else
333 {
334 convarea[i] = -1.;
335 convareaerr[i] = 0.;
336 }
337
338 time[i] = i;
339 timeerr[i] = 0.;
340 }
341
342 TGraphErrors *gr = new TGraphErrors(size,
343 time.GetArray(),convarea.GetArray(),
344 timeerr.GetArray(),convareaerr.GetArray());
345 gr->SetTitle(Form("Conversion Factors Area %3i Average",aidx));
346 gr->GetXaxis()->SetTitle("Camera Nr.");
347 gr->GetYaxis()->SetTitle("<N_{phe}>/<Q> [FADC cnt^{-1}]");
348 return gr;
349}
350
351TGraph *MCalibrationIntensityConstCam::GetConvFactorVsTime ( const Int_t pixid )
352{
353
354 const Int_t size = GetSize();
355
356 if (size == 0)
357 return NULL;
358
359 TArrayF time(size);
360 TArrayF conv(size);
361
362 MCalibConstCam *cam = 0;
363
364 for (Int_t i=0;i<size;i++)
365 {
366 //
367 // Get the calibration cam from the intensity cam
368 //
369 cam = GetCam(i);
370
371 if (cam->GetSize()-1 < pixid)
372 {
373 *fLog << err << "Pixel " << pixid << " out of range " << endl;
374 return NULL;
375 }
376
377 //
378 // Get the calibration pix from the calibration cam
379 //
380 MCalibConstPix &pix = (*cam)[pixid];
381 //
382 time[i] = i;
383 conv[i] = pix.GetCalibConst();
384 //
385 }
386
387
388 TGraph *gr = new TGraph(size,time.GetArray(),conv.GetArray());
389 gr->SetTitle(Form("%s%3i","Pixel ",pixid));
390 gr->GetXaxis()->SetTitle("Camera Nr.");
391 gr->GetYaxis()->SetTitle("<N_{phe}>/<Q> [FADC cnts^{-1}]");
392 return gr;
393
394
395}
396
397void MCalibrationIntensityConstCam::DrawConvFactorPerAreaVsTime( const Int_t aidx )
398{
399
400 TGraphErrors *gr = GetConvFactorPerAreaVsTime(aidx,MGeomCamMagic());
401
402 if (gr)
403 {
404 gr->SetBit(kCanDelete);
405 gr->Draw("A*");
406 }
407}
408
409void MCalibrationIntensityConstCam::DrawConvFactorVsTime( const Int_t idx )
410{
411 TGraph *gr = GetConvFactorVsTime(idx);
412 if (gr)
413 {
414 gr->SetBit(kCanDelete);
415 gr->Draw("A*");
416 }
417}
418
Note: See TracBrowser for help on using the repository browser.