source: trunk/MagicSoft/Mars/mraw/MRawEvtData.cc@ 839

Last change on this file since 839 was 761, checked in by tbretz, 24 years ago
*** empty log message ***
File size: 11.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! Author(s): Thomas Bretz 12/2000 (tbretz@uni-sw.gwdg.de)
19!
20! Copyright: MAGIC Software Development, 2000-2001
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MRawEvtData
28//
29// Storage container to store the raw FADC values.
30//
31// MArrayS fHiGainPixId
32// ---------------------
33// Array of Pixel Numbers for their high voltage channel in the order the
34// FADC values are stored in fHiGainFadcSamples
35//
36// MArrayB fHiGainFadcSaples
37// -------------------------
38// FADC samples (hi gain) of all pixels
39//
40// MArrayS fLoGainPixId
41// --------------------
42// see fHiGainPixId
43//
44// MArrayB fLoGainFadcSamples
45// --------------------------
46// see fHiGainFadcSaples
47//
48/////////////////////////////////////////////////////////////////////////////
49
50#include "MRawEvtData.h"
51
52#include <iomanip.h>
53
54#include <fstream.h>
55
56#include <TH1.h>
57#include <TGraph.h>
58#include <TArrayC.h>
59
60#include "MLog.h"
61#include "MArrayS.h"
62#include "MArrayB.h"
63#include "MRawRunHeader.h"
64
65ClassImp(MRawEvtData)
66
67// --------------------------------------------------------------------------
68//
69// Default constructor. It initializes all arrays with zero size.
70//
71MRawEvtData::MRawEvtData(const char *name, const char *title)
72{
73 *fName = name ? name : "MRawEvtData";
74 *fTitle = title ? title : "Raw Event Data Information";
75
76 InitArrays();
77}
78
79// --------------------------------------------------------------------------
80//
81// Destructor. Deletes all the arrays.
82//
83MRawEvtData::~MRawEvtData()
84{
85 DeleteArrays();
86}
87
88// --------------------------------------------------------------------------
89//
90// reset all arrays
91//
92void MRawEvtData::Clear(Option_t *)
93{
94 /*
95 FIXME:
96 Is Reset (set all entries to zero) what you want to do
97 or Set(0) (delete the array)
98 */
99 fHiGainPixId->Reset();
100 fLoGainPixId->Reset();
101 fHiGainFadcSamples->Reset();
102 fLoGainFadcSamples->Reset();
103}
104
105// --------------------------------------------------------------------------
106//
107// return the number of hi gain samples per pixel
108//
109Byte_t MRawEvtData::GetNumHiGainSamples() const
110{
111 return fHiGainPixId->GetSize() ? fHiGainFadcSamples->GetSize()/fHiGainPixId->GetSize() : 0;
112}
113
114// --------------------------------------------------------------------------
115//
116// return the number of lo gain samples per pixel
117//
118Byte_t MRawEvtData::GetNumLoGainSamples() const
119{
120 return fLoGainPixId->GetSize() ? fLoGainFadcSamples->GetSize()/fLoGainPixId->GetSize() : 0;
121}
122
123// --------------------------------------------------------------------------
124//
125// return the number of stored pixel
126//
127UShort_t MRawEvtData::GetNumPixels() const
128{
129 return fHiGainPixId->GetSize();
130}
131
132
133// --------------------------------------------------------------------------
134//
135// Print out the onformation to *fLog.
136// Options:
137// "hex" Prints the time slices hexadecimal (default)
138// "dec" Prints the time slices decimal
139//
140void MRawEvtData::Print(Option_t *opt)
141{
142 //
143 // print fadc inforation to screen
144 // Possible Options:
145 // - DEC: Print values decimal instead of hexadecimal (default)
146 //
147 const Byte_t nHiSamp = GetNumHiGainSamples();
148 const Byte_t nLoSamp = GetNumLoGainSamples();
149
150 const UShort_t nHiPix = fHiGainPixId->GetSize();;
151 const UShort_t nLoPix = fLoGainPixId->GetSize();;
152
153 fLog->unsetf(ios::showbase);
154
155 *fLog << dec;
156 *fLog << "HiGain: " << nHiPix << " Pixels with " << (Int_t)nHiSamp << " Samples" << endl;
157 *fLog << "LoGain: " << nLoPix << " Pixels with " << (Int_t)nLoSamp << " Samples";;
158
159 TString str(opt);
160 Int_t manip = str.Contains("dec", TString::kIgnoreCase);
161
162 Int_t l=0;
163 for (int i=0; i<nHiPix; i++)
164 {
165 *fLog << endl;
166 *fLog << " " << setfill(' ') << setw(3) << dec << i << ": ";
167 *fLog << (manip?dec:hex) << flush;
168
169 if (!manip)
170 *fLog << setfill('0');
171
172 for (int j=0; j<nHiSamp; j++)
173 {
174 *fLog << setw(manip?3:2);
175 *fLog << ((UShort_t)(*fHiGainFadcSamples)[j+i*nHiSamp]&0xff);
176 if (manip)
177 *fLog << ' ';
178 *fLog << flush;
179 }
180
181 if (!(l<nLoPix && (*fLoGainPixId)[l]==(*fHiGainPixId)[i]))
182 continue;
183
184 if (manip)
185 *fLog << "/ ";
186
187 for (int j=0; j<nLoSamp; j++)
188 {
189 *fLog << setw(manip?3:2);
190 *fLog << ((UShort_t)(*fLoGainFadcSamples)[j+i*nLoSamp]&0xff);
191 if (manip)
192 *fLog << ' ';
193 *fLog << flush;
194 }
195 l++;
196 }
197 *fLog << endl;
198}
199
200// --------------------------------------------------------------------------
201//
202// Draw a pixel. A Histogram or Graph is created and it's draw function is
203// called.
204// Options:
205// "GRAPH" A graph is drawn
206// "HIST" A histogram is drawn
207// number The pixel with the given number is drawn
208//
209void MRawEvtData::Draw(Option_t *opt)
210{
211 // ----- AppendPad(opt);
212
213 //
214 // FIXME: BIG MEMORY LEAK! (( How and when are the objects deleted?)
215 //
216
217 TString str(opt);
218
219 UInt_t num = 0;
220
221 if (str.BeginsWith("GRAPH", TString::kIgnoreCase))
222 {
223 if (str.Length()>5)
224 sscanf(&str[5], "%d", &num);
225
226 if (num>=GetNumPixels())
227 num= GetNumPixels();
228
229 *fLog << "Drawing Graph: Pixel #" << num << " of " << (int)GetNumPixels() << endl;
230
231 const Int_t n = GetNumHiGainSamples();
232
233 Float_t *x = new Float_t[n];
234 Float_t *y = new Float_t[n];
235
236 for (int i=0; i<n; i++)
237 {
238 x[i] = i;
239 y[i] = (*fHiGainFadcSamples)[i + num*GetNumHiGainSamples()];
240 }
241
242 TGraph *graph = new TGraph(n, x, y);
243 graph->SetMaximum (256) ;
244 graph->SetMinimum (0) ;
245
246 graph->Draw("AC*");
247
248 return;
249 }
250
251 if (str.BeginsWith("HIST", TString::kIgnoreCase))
252 {
253 *fLog << "Length: " << str.Length() << endl;
254
255 if (str.Length()>4)
256 sscanf(&str[4], "%d", &num);
257
258 if (num>=GetNumPixels())
259 num= GetNumPixels();
260
261 *fLog << "Drawing Histogram of Pixel " << num << endl;
262
263 const Int_t n = GetNumHiGainSamples();
264
265 char *name = new char[16];
266
267 sprintf(name, "Pixel No.%d", (*fHiGainPixId)[num]);
268
269 TH1F *hist = new TH1F(name, "Hi Gain Samples FADC", n, 0, n);
270
271 for (int i=0; i<n; i++)
272 hist->Fill(0.5+i, (*fHiGainFadcSamples)[i + num*GetNumHiGainSamples()]);
273
274 hist->Draw();
275
276 return;
277 }
278
279 *fLog << "MRawEvtData::Draw: Warning: You must specify either 'GRAPH' or 'HIST'" << endl;
280}
281
282// --------------------------------------------------------------------------
283//
284// Deletes all arrays describing the pixel Id and Samples in pixels.
285// The flag is for future usage.
286//
287void MRawEvtData::DeletePixels(Bool_t flag)
288{
289 DeleteArrays();
290 InitArrays(flag);
291}
292
293// --------------------------------------------------------------------------
294//
295// Deletes all the arrays
296//
297void MRawEvtData::DeleteArrays()
298{
299 delete fHiGainPixId;
300 delete fLoGainPixId;
301 delete fHiGainFadcSamples;
302 delete fLoGainFadcSamples;
303}
304
305// --------------------------------------------------------------------------
306//
307// Deletes all the arrays
308// The flag is for future usage.
309//
310void MRawEvtData::InitArrays(Bool_t flag)
311{
312 // const int npix = !flag ? 0 : fRunHeader->GetNumCrates()*fRunHeader->GetNumPixInCrate();
313
314 fHiGainPixId = new MArrayS(0);//npix);
315 fLoGainPixId = new MArrayS(0);//npix);
316 fHiGainFadcSamples = new MArrayB(0);//npix*fRunHeader->GetNumSamplesHiGain());
317 fLoGainFadcSamples = new MArrayB(0);//npix*fRunHeader->GetNumSamplesLoGain());
318}
319
320// --------------------------------------------------------------------------
321//
322// This is to fill the data of one pixel to the MRawEvtHeader Class.
323// The parameters are the pixelnumber and the FADC_SLICES values of ADCs
324// Add to lo gains if lflag = 1
325//
326void MRawEvtData::AddPixel(UShort_t nOfPixel, TArrayC *data, Bool_t lflag)
327{
328 MArrayS *arrpix = lflag ? fLoGainPixId : fHiGainPixId;
329 MArrayB *arrsam = lflag ? fLoGainFadcSamples : fHiGainFadcSamples;
330
331 //
332 // check whether we got the right number of new samples
333 // if there are no samples already stored: this is the new number of samples
334 //
335 const Byte_t ns = data->GetSize();
336 const Byte_t nSamp = lflag ? GetNumLoGainSamples() : GetNumHiGainSamples();
337 if (nSamp && ns!=nSamp)
338 {
339 *fLog << "RawEvtData::FillPixel: Error, number of samples in ";
340 *fLog << "TArrayC doesn't match actual number" << endl;
341 return;
342 }
343
344 //
345 // enhance pixel array by one
346 //
347 arrpix->Set(arrpix->GetSize()+1);
348
349 //
350 // add the number of the new pixel to the array as last entry
351 //
352 arrpix->AddAt(nOfPixel, arrpix->GetSize()-1);
353
354 //
355 // enhance the array by the number of new samples
356 //
357 arrsam->Set(arrsam->GetSize()+ns);
358
359 //
360 // add the new slices as last entries to array
361 //
362 arrsam->AddAt((Byte_t*)data->GetArray(), arrsam->GetSize()-ns, ns);
363}
364
365// --------------------------------------------------------------------------
366//
367// Fills members with information from a magic binary file.
368// WARNING: you have to use Init() before you can do this
369//
370void MRawEvtData::ReadEvt(istream &fin)
371{
372 const UShort_t nlo = fRunHeader->GetNumSamplesLoGain();
373 const UShort_t nhi = fRunHeader->GetNumSamplesHiGain();
374
375 const UShort_t npic = fRunHeader->GetNumPixInCrate();
376
377 //
378 // Enhance array by the size which we'll read now
379 //
380 const int npixhi = fHiGainPixId->GetSize();
381 const int npixlo = fLoGainPixId->GetSize();
382
383 fHiGainPixId->Set(npixhi+npic);
384 fLoGainPixId->Set(npixlo+npic);
385
386 const int nsamhi = fHiGainFadcSamples->GetSize();
387 const int nsamlo = fLoGainFadcSamples->GetSize();
388
389 fHiGainFadcSamples->Set(nsamhi+nhi*npic);
390 fLoGainFadcSamples->Set(nsamlo+nlo*npic);
391
392 Byte_t *higainsam = fHiGainFadcSamples->GetArray()+nsamhi;
393 Byte_t *logainsam = fLoGainFadcSamples->GetArray()+nsamlo;
394
395 for (int i=0; i<npic; i++)
396 {
397 //
398 // get the spiral pixel number from the run header
399 //
400 const UShort_t npix = fRunHeader->GetPixAssignment(i);
401
402 //
403 // This is to fill the data of one pixel to the MRawEvtHeader Class.
404 // The parameters are the pixelnumber and the FADC_SLICES values of ADCs
405 // Add to lo gains if lflag = 1
406 //
407 fHiGainPixId->AddAt(npix, npixhi+i);
408 fin.read(higainsam, nhi);
409 higainsam += nhi;
410
411 // FIXME: Not implemented in the raw files yet
412 //if (IsLoGainOn(i, j))
413 //{
414 fLoGainPixId->AddAt(npix, npixlo+i);
415 fin.read(logainsam, nlo);
416 logainsam += nlo;
417 //}
418 }
419}
420
Note: See TracBrowser for help on using the repository browser.