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

Last change on this file since 8785 was 8704, checked in by tbretz, 17 years ago
*** empty log message ***
File size: 17.1 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-2007
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 7
50// ------------------
51// + UShort_t fNumBytesPerSample;
52//
53// Version 6
54// ------------------
55// - The data can now be stoe in a single array keeping he full
56// compatibility
57//
58// Version 5 (0.8.5):
59// ------------------
60// - Changed type of ABFlags from TArrayC to MArrayB
61//
62// Version 4 (0.8.4):
63// ------------------
64// - Changed derivement from MCamEvent to MParContainer and MCamEvent
65//
66// Version 3 (0.8.4):
67// ------------------
68// - Added fABFlags
69//
70// Version 2:
71// ----------
72// - Derives from MCamEvent now
73//
74// Version 1:
75// ----------
76// - First implementation
77//
78/////////////////////////////////////////////////////////////////////////////
79#include "MRawEvtData.h"
80
81#include <fstream>
82
83#include <TH1.h>
84#include <TGraph.h>
85#include <TArrayC.h>
86#include <TVirtualPad.h>
87
88#include "MLog.h"
89#include "MLogManip.h"
90
91#include "MArrayS.h"
92#include "MArrayB.h"
93#include "MGeomCam.h"
94
95#include "MRawCrateArray.h"
96#include "MRawCrateData.h"
97
98#include "MRawRunHeader.h"
99#include "MRawEvtPixelIter.h"
100
101ClassImp(MRawEvtData);
102
103using namespace std;
104
105// --------------------------------------------------------------------------
106//
107// Default constructor. It initializes all arrays with zero size.
108//
109MRawEvtData::MRawEvtData(const char *name, const char *title)
110 : fRunHeader(0), fNumBytesPerSample(1)
111{
112 fName = name ? name : "MRawEvtData";
113 fTitle = title ? title : "Raw Event Data Information";
114
115 InitArrays();
116}
117
118// --------------------------------------------------------------------------
119//
120// Destructor. Deletes all the arrays.
121//
122MRawEvtData::~MRawEvtData()
123{
124 DeleteArrays();
125}
126
127// --------------------------------------------------------------------------
128//
129// reset all arrays
130//
131void MRawEvtData::Clear(Option_t *)
132{
133 /*
134 FIXME:
135 Is Reset (set all entries to zero) what you want to do
136 or Set(0) (delete the array)
137 */
138 fHiGainPixId->Reset();
139 fLoGainPixId->Reset();
140 fHiGainFadcSamples->Reset();
141 fLoGainFadcSamples->Reset();
142 fABFlags->Reset();
143}
144
145// --------------------------------------------------------------------------
146//
147// return the number of hi gain samples per pixel
148//
149UShort_t MRawEvtData::GetNumHiGainSamples() const
150{
151 // If value<0 we are reading old MC files which don't have the
152 // value set so far. So we use the old methid to determin the
153 // numbers and calculate them from the length of the arrays.
154 return fHiGainPixId->GetSize() ? fHiGainFadcSamples->GetSize()/(fHiGainPixId->GetSize()*fNumBytesPerSample) : 0;
155}
156
157// --------------------------------------------------------------------------
158//
159// return the number of lo gain samples per pixel
160//
161UShort_t MRawEvtData::GetNumLoGainSamples() const
162{
163 // If value<0 we are reading old MC files which don't have the
164 // value set so far. So we use the old methid to determin the
165 // numbers and calculate them from the length of the arrays.
166 return fLoGainPixId->GetSize() ? fLoGainFadcSamples->GetSize()/(fLoGainPixId->GetSize()*fNumBytesPerSample) : 0;
167}
168
169// --------------------------------------------------------------------------
170//
171// return the number of stored pixel
172//
173UShort_t MRawEvtData::GetNumPixels() const
174{
175 return fHiGainPixId ? fHiGainPixId->GetSize() : 0;
176}
177
178// --------------------------------------------------------------------------
179//
180// Print out the onformation to *fLog.
181// Options:
182// "hex" Prints the time slices hexadecimal (default)
183// "dec" Prints the time slices decimal
184//
185void MRawEvtData::Print(Option_t *opt) const
186{
187 //
188 // print fadc inforation to screen
189 // Possible Options:
190 // - DEC: Print values decimal instead of hexadecimal (default)
191 //
192 const UShort_t nHiSamp = GetNumHiGainSamples();
193 const UShort_t nLoSamp = GetNumLoGainSamples();
194
195 const UShort_t bps = GetNumBytesPerSample();
196
197 fLog->unsetf(ios::showbase);
198
199 *fLog << dec << all;
200 *fLog << GetDescriptor() << ": " << endl;
201 *fLog << GetNumPixels() << " Pixels with " << (Int_t)nHiSamp << "/" << (Int_t)nLoSamp << " samples" << endl;
202
203 TString str(opt);
204 Int_t manip = str.Contains("dec", TString::kIgnoreCase);
205
206 MRawEvtPixelIter pixel(const_cast<MRawEvtData*>(this));
207 Int_t i = 0;
208 while (pixel.Next())
209 {
210 const UShort_t idx = pixel.GetPixelId();
211
212 *fLog << endl << dec << setfill(' ');
213 *fLog << " " << setw(3) << i++ << " - " << setw(3) << idx << " ";
214
215 if (pixel.IsABFlagValid())
216 *fLog << "<" << (pixel.HasABFlag()?"B":"A") << "> ";
217
218 *fLog << (manip?dec:hex) << setfill(manip?' ':'0');
219
220 const Byte_t *hi = (Byte_t*)pixel.GetHiGainSamples();
221 const Byte_t *lo = (Byte_t*)pixel.GetLoGainSamples();
222
223 for (int j=0; j<nHiSamp*bps; j++)
224 {
225 *fLog << setw(manip?3:2);
226 *fLog << (hi[j]&0xff);
227 if (manip)
228 *fLog << ' ';
229 }
230
231 for (int j=0; j<nLoSamp*bps; j++)
232 {
233 *fLog << setw(manip?3:2);
234 *fLog << ((UShort_t)lo[j]&0xff);
235 if (manip)
236 *fLog << ' ';
237 }
238 }
239 *fLog << dec << endl;
240}
241
242// --------------------------------------------------------------------------
243//
244// Draw a pixel. A Histogram or Graph is created and it's draw function is
245// called.
246// Options:
247// "GRAPH" A graph is drawn
248// "HIST" A histogram is drawn
249// <index> The pixel with the given index is drawn
250//
251void MRawEvtData::Draw(Option_t *opt)
252{
253 if (GetNumPixels()==0)
254 {
255 *fLog << warn << "Sorry, no pixel to draw!" << endl;
256 return;
257 }
258
259 TString str(opt);
260 str.ToLower();
261
262 UInt_t id = 0;
263
264 if (str.BeginsWith("graph"))
265 if (str.Length()>5)
266 sscanf(&str[5], "%d", &id);
267 if (str.BeginsWith("hist"))
268 if (str.Length()>4)
269 sscanf(&str[4], "%d", &id);
270
271 MRawEvtPixelIter pix(this);
272 if (!pix.Jump(id))
273 {
274 *fLog << warn << dec << "Pixel Idx #" << id << " doesn't exist!" << endl;
275 return;
276 }
277
278 const void *higains = pix.GetHiGainSamples();
279 const void *logains = pix.GetLoGainSamples();
280
281 const Int_t nh = GetNumHiGainSamples();
282 const Int_t nl = GetNumLoGainSamples();
283
284 TString name = "Pixel Idx.";
285 name += pix.GetPixelId();
286
287 Bool_t same = str.Contains("same");
288
289 if (str.BeginsWith("graph"))
290 {
291 *fLog << inf << "Drawing Graph: Pixel Idx #" << dec << pix.GetPixelId();
292 *fLog << " of " << (int)GetNumPixels() << "Pixels" << endl;
293
294 TGraph *graphhi = new TGraph;
295
296 for (int i=0; i<nh; i++)
297 graphhi->SetPoint(graphhi->GetN(), i, GetSample(higains, i));
298 for (int i=0; i<nl; i++)
299 graphhi->SetPoint(graphhi->GetN(), i+nh, GetSample(logains, i));
300
301 graphhi->SetMaximum(256);
302 graphhi->SetMinimum(0);
303
304 graphhi->SetBit(kCanDelete);
305 graphhi->Draw(same ? "C*" : "AC*");
306
307 TH1F *histhi = graphhi->GetHistogram();
308 histhi->SetMinimum(0);
309 histhi->SetMaximum(255);
310
311 histhi->SetXTitle("Time/FADC Slices");
312 histhi->SetYTitle("Signal/FADC Units");
313
314 return;
315 }
316
317 if (str.BeginsWith("hist"))
318 {
319 // FIXME: Add Legend
320 *fLog << inf << "Drawing Histogram of Pixel with Idx #" << dec << pix.GetPixelId() << " to " << gPad << endl;
321
322 TH1F *histh = new TH1F(name, "FADC Samples", nh, -0.5, nh-.5);
323 histh->SetMinimum(0);
324 histh->SetMaximum(255);
325 histh->SetXTitle("Time [FADC Slices]");
326 histh->SetYTitle("Signal [FADC Units]");
327 histh->SetDirectory(NULL);
328 for (int i=0; i<nh; i++)
329 histh->Fill(i, GetSample(higains, i));
330 for (int i=0; i<nl; i++)
331 histh->Fill(i, GetSample(logains, i));
332 histh->SetBit(kCanDelete);
333 histh->Draw(same ? "same" : "");
334
335 return;
336 }
337
338 *fLog << warn << dbginf << "WARNING - You must specify either 'GRAPH' or 'HIST'" << endl;
339}
340
341// --------------------------------------------------------------------------
342//
343// Deletes all the arrays
344// The flag is for future usage.
345//
346void MRawEvtData::InitArrays(UShort_t numconnected, UShort_t maxid)
347{
348 // fRunHeader should not be set only in the constructor!
349 const Int_t numhi = fRunHeader ? fRunHeader->GetNumSamplesHiGain() : 0;
350 const Int_t numlo = fRunHeader ? fRunHeader->GetNumSamplesLoGain() : 0;
351
352 fNumBytesPerSample = fRunHeader ? fRunHeader->GetNumBytesPerSample() : 1;
353
354 fHiGainPixId = new MArrayS(numconnected);
355 fLoGainPixId = new MArrayS(0);
356 fHiGainFadcSamples = new MArrayB(numconnected*(numhi+numlo)*fNumBytesPerSample);
357 fLoGainFadcSamples = new MArrayB(0);
358
359 fABFlags = new MArrayB(maxid==0 ? 0 : maxid/8+1);
360
361 fConnectedPixels = 0;
362}
363
364// --------------------------------------------------------------------------
365//
366// Deletes all the arrays
367//
368void MRawEvtData::DeleteArrays()
369{
370 delete fHiGainPixId;
371 delete fLoGainPixId;
372 delete fHiGainFadcSamples;
373 delete fLoGainFadcSamples;
374 delete fABFlags;
375}
376
377// --------------------------------------------------------------------------
378//
379// Deletes all arrays describing the pixel Id and Samples in pixels.
380// The flag is for future usage.
381//
382void MRawEvtData::ResetPixels(UShort_t numconnected, UShort_t maxid)
383{
384 //const int npix = fRunHeader->GetNumCrates()*fRunHeader->GetNumPixInCrate();
385 if (fHiGainPixId && fHiGainPixId->GetSize()==numconnected && (UShort_t)fABFlags->GetSize()==(maxid/8+1))
386 {
387 fConnectedPixels = 0;
388 return;
389 }
390
391 DeleteArrays();
392 InitArrays(numconnected, maxid);
393}
394
395// --------------------------------------------------------------------------
396//
397// This is to fill the data of one pixel to the MRawEvtHeader Class.
398// The parameters are the pixelnumber and the FADC_SLICES values of ADCs
399//
400void MRawEvtData::AddPixel(UShort_t nOfPixel, TArrayC *data)
401{
402 const Int_t n = fRunHeader->GetNumSamples();
403 if (data->GetSize()!=n)
404 {
405 *fLog << err << "RawEvtData::AddPixel: Error, number of samples in ";
406 *fLog << "TArrayC " << data->GetSize() << " doesn't match current number " << n << endl;
407 return;
408 }
409
410 //
411 // enhance pixel array by one
412 //
413 fHiGainPixId->Set(fHiGainPixId->GetSize()+1);
414
415 //
416 // add the number of the new pixel to the array as last entry
417 //
418 fHiGainPixId->AddAt(nOfPixel, fHiGainPixId->GetSize()-1);
419
420 //
421 // enhance the array by the number of new samples
422 //
423 fHiGainFadcSamples->Set(fHiGainFadcSamples->GetSize()+n);
424
425 //
426 // add the new slices as last entries to array
427 //
428 fHiGainFadcSamples->AddAt((Byte_t*)data->GetArray(), fHiGainFadcSamples->GetSize()-n, n);
429}
430/*
431void MRawEvtData::AddPixel(UShort_t nOfPixel, TArrayC *data, Bool_t lflag)
432{
433 MArrayS *arrpix = lflag ? fLoGainPixId : fHiGainPixId;
434 MArrayB *arrsam = lflag ? fLoGainFadcSamples : fHiGainFadcSamples;
435
436 //
437 // check whether we got the right number of new samples
438 // if there are no samples already stored: this is the new number of samples
439 //
440 const UShort_t ns = data->GetSize();
441 const UShort_t nSamp = lflag ? GetNumLoGainSamples() : GetNumHiGainSamples();
442 if (nSamp && ns!=nSamp)
443 {
444 *fLog << err << "RawEvtData::AddPixel: Error, number of samples in ";
445 *fLog << "TArrayC " << ns << " doesn't match current number " << nSamp << endl;
446 return;
447 }
448
449 //
450 // enhance pixel array by one
451 //
452 arrpix->Set(arrpix->GetSize()+1);
453
454 //
455 // add the number of the new pixel to the array as last entry
456 //
457 arrpix->AddAt(nOfPixel, arrpix->GetSize()-1);
458
459 //
460 // enhance the array by the number of new samples
461 //
462 arrsam->Set(arrsam->GetSize()+ns);
463
464 //
465 // add the new slices as last entries to array
466 //
467 arrsam->AddAt((Byte_t*)data->GetArray(), arrsam->GetSize()-ns, ns);
468}
469*/
470
471void MRawEvtData::ReadPixel(istream &fin, Int_t npix)
472{
473 // number of samples
474 const UInt_t ns = fRunHeader->GetNumSamples();
475
476 // bytes per sample
477 const Int_t bps = fRunHeader->GetNumBytesPerSample();
478
479 // Number of bytes
480 const Int_t nb = ns*bps;
481
482 // position in higain array
483 Byte_t *pos = fHiGainFadcSamples->GetArray() + fConnectedPixels*nb;
484
485 // Set pixel index
486 fHiGainPixId->AddAt(npix, fConnectedPixels++);
487
488 // Read data for one pixel
489 fin.read((char*)pos, nb);
490}
491
492void MRawEvtData::SetABFlag(Int_t npix, Bool_t ab)
493{
494 if (ab)
495 SETBIT((*fABFlags)[npix/8], npix%8);
496 else
497 CLRBIT((*fABFlags)[npix/8], npix%8);
498}
499
500// --------------------------------------------------------------------------
501//
502// Fills members with information from a magic binary file.
503// WARNING: you have to use Init() before you can do this
504//
505/*
506void MRawEvtData::ReadEvt(istream &fin, Int_t posinarray)
507{
508
509 const UShort_t npic = fRunHeader->GetNumPixInCrate();
510
511 const UShort_t npos = npic*posinarray;
512
513 //const Byte_t ab = fCrateArray->GetEntry(posinarray)->GetABFlags();
514
515 for (int i=npos; i<npic+npos; i++)
516 {
517 // calc the spiral hardware pixel number
518 const UShort_t ipos = i;
519
520 // Get Hardware Id
521 const Short_t hwid = fRunHeader->GetPixAssignment(ipos);
522
523 // Check whether the pixel is connected or not
524 if (hwid==0)
525 {
526 const UShort_t n = fRunHeader->GetNumSamplesLoGain()+fRunHeader->GetNumSamplesHiGain();
527 fin.seekg(n, ios::cur);
528 return;
529 }
530
531 // -1 converts the hardware pixel Id into the software pixel index
532 const Int_t npix = (Int_t)hwid-1;
533
534 const Byte_t ab = fCrateArray->GetEntry(posinarray)->GetABFlags();
535 AddPixel(fin, npix, TESTBIT(ab, i-npos));
536 }
537}
538*/
539
540// --------------------------------------------------------------------------
541//
542// Return the size in bytes of one event data block
543//
544Int_t MRawEvtData::GetNumBytes() const
545{
546 const UShort_t nlo = fRunHeader->GetNumSamplesLoGain();
547 const UShort_t nhi = fRunHeader->GetNumSamplesHiGain();
548 const UShort_t npic = fRunHeader->GetNumPixInCrate();
549 const UShort_t nbps = fRunHeader->GetNumBytesPerSample();
550
551 return (nhi+nlo)*npic*nbps;
552}
553
554// --------------------------------------------------------------------------
555//
556// Make sure, that you skip the whole event. This function only skips a part
557// of the event - see MRawRead::SkipEvent
558//
559void MRawEvtData::SkipEvt(istream &fin)
560{
561 fin.seekg(GetNumBytes(), ios::cur);
562}
563
564Bool_t MRawEvtData::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
565{
566 *fLog << warn << "WARNING - The use of MRawEvtData::GetPixelContent is deprecated!" << endl;
567
568 /*
569 MRawEvtPixelIter Next(const_cast<MRawEvtData*>(this));
570 if (!Next.Jump(idx))
571 return kFALSE;
572
573 switch (type)
574 {
575 case 0:
576 val = Next.GetSumHiGainSamples()-(float)GetNumHiGainSamples()*fHiGainFadcSamples->GetArray()[0];
577 val*= cam.GetPixRatio(idx);
578 break;
579 case 1:
580 val = Next.GetMaxHiGainSample();
581 break;
582 case 2:
583 val = Next.GetMaxLoGainSample();
584 break;
585 case 3:
586 val = Next.GetIdxMaxHiGainSample();
587 break;
588 case 4:
589 val = Next.GetIdxMaxLoGainSample();
590 break;
591 case 5:
592 val = Next.GetIdxMaxHiLoGainSample();
593 return val >= 0;
594 }
595*/
596 return kTRUE;
597}
598
599void MRawEvtData::Copy(TObject &named)
600#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
601const
602#endif
603{
604 MRawEvtData &evt = (MRawEvtData &)named;
605
606 *evt.fHiGainPixId = *fHiGainPixId;
607 *evt.fLoGainPixId = *fLoGainPixId;
608
609 *evt.fHiGainFadcSamples = *fHiGainFadcSamples;
610 *evt.fLoGainFadcSamples = *fLoGainFadcSamples;
611
612 *evt.fABFlags = *fABFlags;
613
614 evt.fConnectedPixels = fConnectedPixels;
615}
Note: See TracBrowser for help on using the repository browser.