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

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