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

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