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

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