source: trunk/MagicSoft/Mars/manalysis/MPedCalcPedRun.cc@ 2827

Last change on this file since 2827 was 2821, checked in by rico, 21 years ago
*** empty log message ***
File size: 5.6 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): Josep Flix 04/2001 <mailto:jflix@ifae.es>
19! Author(s): Thomas Bretz 05/2001 <mailto:tbretz@uni-sw.gwdg.de>
20! Author(s): Sebastian Commichau 12/2003
21! Author(s): Javier Rico 01/2004
22!
23! Copyright: MAGIC Software Development, 2000-2004
24!
25!
26\* ======================================================================== */
27
28/////////////////////////////////////////////////////////////////////////////
29// //
30// MPedCalcPedRun //
31// //
32// This task takes a pedestal run file and fills MPedestalCam during //
33// the Process() with the pedestal and rms computed in an event basis. //
34// In the PostProcess() MPedestalCam is finally filled with the pedestal //
35// mean and rms computed in a run basis. //
36// More than one run (file) can be merged //
37// //
38// Input Containers: //
39// MRawEvtData //
40// //
41// Output Containers: //
42// MPedestalCam //
43// //
44// //
45/////////////////////////////////////////////////////////////////////////////
46
47#include "MPedCalcPedRun.h"
48
49#include "MParList.h"
50
51#include "MLog.h"
52#include "MLogManip.h"
53
54#include "MRawRunHeader.h"
55#include "MRawEvtPixelIter.h"
56#include "MRawEvtData.h"
57
58#include "MPedestalPix.h"
59#include "MPedestalCam.h"
60
61#include "MGeomCamMagic.h"
62
63ClassImp(MPedCalcPedRun);
64
65using namespace std;
66
67MPedCalcPedRun::MPedCalcPedRun(const char *name, const char *title)
68{
69 fName = name ? name : "MPedCalcPedRun";
70 fTitle = title ? title : "Task to calculate pedestals from pedestal runs raw data";
71
72 AddToBranchList("fHiGainPixId");
73 AddToBranchList("fHiGainFadcSamples");
74}
75
76Int_t MPedCalcPedRun::PreProcess( MParList *pList )
77{
78 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
79 if (!fRawEvt)
80 {
81 *fLog << err << "MRawEvtData not found... aborting." << endl;
82 return kFALSE;
83 }
84
85 fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam");
86 if (!fPedestals)
87 return kFALSE;
88
89 fNumSamplesTot=0;
90
91 return kTRUE;
92}
93
94Bool_t MPedCalcPedRun::ReInit(MParList *pList)
95{
96 const MRawRunHeader *runheader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
97 if (!runheader)
98 {
99 *fLog << warn << dbginf;
100 *fLog << "Warning - cannot check file type, MRawRunHeader not found." << endl;
101 }
102 else
103 if (runheader->GetRunType() == kRTMonteCarlo)
104 return kTRUE;
105
106 fNumPixels = fPedestals->GetSize();
107
108 if(fSumx.GetSize()==0)
109 {
110 fSumx.Set(fNumPixels);
111 fSumx2.Set(fNumPixels);
112
113 memset(fSumx.GetArray(), 0, sizeof(Float_t)*fNumPixels);
114 memset(fSumx2.GetArray(), 0, sizeof(Float_t)*fNumPixels);
115 }
116
117 // Calculate an even number for the hi gain samples to avoid
118 // biases due to the fluctuation in pedestal from one slice to
119 // the other one
120 fNumHiGainSamples = runheader->GetNumSamplesHiGain() & ~1;
121
122 return kTRUE;
123}
124
125Int_t MPedCalcPedRun::Process()
126{
127 MRawEvtPixelIter pixel(fRawEvt);
128
129 while (pixel.Next())
130 {
131 Byte_t *ptr = pixel.GetHiGainSamples();
132 const Byte_t *end = ptr + fNumHiGainSamples;
133
134 UInt_t sum = 0;
135 UInt_t sqr = 0;
136
137 do
138 {
139 sum += *ptr;
140 sqr += *ptr * *ptr;
141 }
142 while (++ptr != end);
143
144 const Float_t msum = (Float_t)sum;
145 const Float_t msqr = (Float_t)sqr;
146
147 const Float_t higainped = msum/fNumHiGainSamples;
148 const Float_t higainrms = TMath::Sqrt((msqr-msum*msum/fNumHiGainSamples)/(fNumHiGainSamples-1.));
149
150 const UInt_t idx = pixel.GetPixelId();
151 (*fPedestals)[idx].Set(higainped, higainrms);
152
153 fSumx[idx] += msum;
154 fSumx2[idx] += sqr;
155 }
156
157 fPedestals->SetReadyToSave();
158 fNumSamplesTot+=fNumHiGainSamples;
159
160 return kTRUE;
161}
162
163Int_t MPedCalcPedRun::PostProcess()
164{
165 // Compute pedestals and rms from the whole run
166 const Int_t n = fNumSamplesTot;
167
168 MRawEvtPixelIter pixel(fRawEvt);
169
170 while (pixel.Next())
171 {
172 const UInt_t pixid = pixel.GetPixelId();
173
174 const Float_t sum = fSumx.At(pixid);
175 const Float_t sum2 = fSumx2.At(pixid);
176
177 const Float_t higainped = sum/n;
178 const Float_t higainrms = TMath::Sqrt((sum2-sum*sum/n)/(n-1.));
179
180 (*fPedestals)[pixid].Set(higainped, higainrms);
181 }
182
183 return kTRUE;
184}
185
186
Note: See TracBrowser for help on using the repository browser.