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

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