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

Last change on this file since 9029 was 8938, checked in by tbretz, 16 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-2008
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 sscanf(str.Data()+5, "%d", &id);
266 if (str.BeginsWith("hist"))
267 sscanf(str.Data()+4, "%d", &id);
268
269 MRawEvtPixelIter pix(this);
270 if (!pix.Jump(id))
271 {
272 *fLog << warn << dec << "Pixel Idx #" << id << " doesn't exist!" << endl;
273 return;
274 }
275
276 const void *higains = pix.GetHiGainSamples();
277 const void *logains = pix.GetLoGainSamples();
278
279 const Int_t nh = GetNumHiGainSamples();
280 const Int_t nl = GetNumLoGainSamples();
281
282 TString name = "Pixel Idx.";
283 name += pix.GetPixelId();
284
285 Bool_t same = str.Contains("same");
286
287 if (str.BeginsWith("graph"))
288 {
289 *fLog << inf << "Drawing Graph: Pixel Idx #" << dec << pix.GetPixelId();
290 *fLog << " of " << (int)GetNumPixels() << "Pixels" << endl;
291
292 TGraph *graphhi = new TGraph;
293
294 for (int i=0; i<nh; i++)
295 graphhi->SetPoint(graphhi->GetN(), i, GetSample(higains, i));
296 for (int i=0; i<nl; i++)
297 graphhi->SetPoint(graphhi->GetN(), i+nh, GetSample(logains, i));
298
299 graphhi->SetMaximum(GetMax()+0.5);
300 graphhi->SetMinimum(0);
301
302 graphhi->SetBit(kCanDelete);
303 graphhi->Draw(same ? "C*" : "AC*");
304
305 TH1F *histhi = graphhi->GetHistogram();
306 histhi->SetMinimum(0);
307 histhi->SetMaximum(GetMax()+0.5);
308
309 histhi->SetXTitle("Time/FADC Slices");
310 histhi->SetYTitle("Signal/FADC Units");
311
312 return;
313 }
314
315 if (str.BeginsWith("hist"))
316 {
317 // FIXME: Add Legend
318 *fLog << inf << "Drawing Histogram of Pixel with Idx #" << dec << pix.GetPixelId() << " to " << gPad << endl;
319
320 TH1F *histh = new TH1F(name, "FADC Samples", nh+nl, -0.5, nh+nl-.5);
321 histh->SetMinimum(0);
322 histh->SetMaximum(GetMax()+0.5);
323 histh->SetXTitle("Time [FADC Slices]");
324 histh->SetYTitle("Signal [FADC Units]");
325 histh->SetDirectory(NULL);
326 for (int i=0; i<nh; i++)
327 histh->Fill(i, GetSample(higains, i));
328 for (int i=0; i<nl; i++)
329 histh->Fill(i+nl, GetSample(logains, i));
330 histh->SetBit(kCanDelete);
331 histh->Draw(same ? "same" : "");
332
333 return;
334 }
335
336 *fLog << warn << dbginf << "WARNING - You must specify either 'GRAPH' or 'HIST'" << endl;
337}
338
339// --------------------------------------------------------------------------
340//
341// Deletes all the arrays
342// The flag is for future usage.
343//
344void MRawEvtData::InitArrays(UShort_t numconnected, UShort_t maxid)
345{
346 // fRunHeader should not be set only in the constructor!
347 const Int_t numhi = fRunHeader ? fRunHeader->GetNumSamplesHiGain() : 0;
348 const Int_t numlo = fRunHeader ? fRunHeader->GetNumSamplesLoGain() : 0;
349
350 fNumBytesPerSample = fRunHeader ? fRunHeader->GetNumBytesPerSample() : 1;
351
352 fHiGainPixId = new MArrayS(numconnected);
353 fLoGainPixId = new MArrayS(0);
354 fHiGainFadcSamples = new MArrayB(numconnected*(numhi+numlo)*fNumBytesPerSample);
355 fLoGainFadcSamples = new MArrayB(0);
356
357 fABFlags = new MArrayB(maxid==0 ? 0 : maxid/8+1);
358
359 fConnectedPixels = 0;
360}
361
362// --------------------------------------------------------------------------
363//
364// Deletes all the arrays
365//
366void MRawEvtData::DeleteArrays()
367{
368 delete fHiGainPixId;
369 delete fLoGainPixId;
370 delete fHiGainFadcSamples;
371 delete fLoGainFadcSamples;
372 delete fABFlags;
373}
374
375// --------------------------------------------------------------------------
376//
377// Deletes all arrays describing the pixel Id and Samples in pixels.
378// The flag is for future usage.
379//
380void MRawEvtData::ResetPixels(UShort_t numconnected, UShort_t maxid)
381{
382 //const int npix = fRunHeader->GetNumCrates()*fRunHeader->GetNumPixInCrate();
383 if (fHiGainPixId && fHiGainPixId->GetSize()==numconnected && (UShort_t)fABFlags->GetSize()==(maxid/8+1))
384 {
385 fConnectedPixels = 0;
386 return;
387 }
388
389 DeleteArrays();
390 InitArrays(numconnected, maxid);
391}
392
393// --------------------------------------------------------------------------
394//
395// This is to fill the data of one pixel to the MRawEvtHeader Class.
396// The parameters are the pixelnumber and the FADC_SLICES values of ADCs
397//
398void MRawEvtData::AddPixel(UShort_t nOfPixel, TArrayC *data)
399{
400 const Int_t n = fRunHeader->GetNumSamples();
401 if (data->GetSize()!=n)
402 {
403 *fLog << err << "RawEvtData::AddPixel: Error, number of samples in ";
404 *fLog << "TArrayC " << data->GetSize() << " doesn't match current number " << n << endl;
405 return;
406 }
407
408 //
409 // enhance pixel array by one
410 //
411 fHiGainPixId->Set(fHiGainPixId->GetSize()+1);
412
413 //
414 // add the number of the new pixel to the array as last entry
415 //
416 fHiGainPixId->AddAt(nOfPixel, fHiGainPixId->GetSize()-1);
417
418 //
419 // enhance the array by the number of new samples
420 //
421 fHiGainFadcSamples->Set(fHiGainFadcSamples->GetSize()+n);
422
423 //
424 // add the new slices as last entries to array
425 //
426 fHiGainFadcSamples->AddAt((Byte_t*)data->GetArray(), fHiGainFadcSamples->GetSize()-n, n);
427}
428/*
429void MRawEvtData::AddPixel(UShort_t nOfPixel, TArrayC *data, Bool_t lflag)
430{
431 MArrayS *arrpix = lflag ? fLoGainPixId : fHiGainPixId;
432 MArrayB *arrsam = lflag ? fLoGainFadcSamples : fHiGainFadcSamples;
433
434 //
435 // check whether we got the right number of new samples
436 // if there are no samples already stored: this is the new number of samples
437 //
438 const UShort_t ns = data->GetSize();
439 const UShort_t nSamp = lflag ? GetNumLoGainSamples() : GetNumHiGainSamples();
440 if (nSamp && ns!=nSamp)
441 {
442 *fLog << err << "RawEvtData::AddPixel: Error, number of samples in ";
443 *fLog << "TArrayC " << ns << " doesn't match current number " << nSamp << endl;
444 return;
445 }
446
447 //
448 // enhance pixel array by one
449 //
450 arrpix->Set(arrpix->GetSize()+1);
451
452 //
453 // add the number of the new pixel to the array as last entry
454 //
455 arrpix->AddAt(nOfPixel, arrpix->GetSize()-1);
456
457 //
458 // enhance the array by the number of new samples
459 //
460 arrsam->Set(arrsam->GetSize()+ns);
461
462 //
463 // add the new slices as last entries to array
464 //
465 arrsam->AddAt((Byte_t*)data->GetArray(), arrsam->GetSize()-ns, ns);
466}
467*/
468
469void MRawEvtData::ReadPixel(istream &fin, Int_t npix)
470{
471 // number of samples
472 const UInt_t ns = fRunHeader->GetNumSamples();
473
474 // bytes per sample
475 const Int_t bps = fRunHeader->GetNumBytesPerSample();
476
477 // Number of bytes
478 const Int_t nb = ns*bps;
479
480 // position in higain array
481 Byte_t *pos = fHiGainFadcSamples->GetArray() + fConnectedPixels*nb;
482
483 // Set pixel index
484 fHiGainPixId->AddAt(npix, fConnectedPixels++);
485
486 // Read data for one pixel
487 fin.read((char*)pos, nb);
488}
489
490void MRawEvtData::SetABFlag(Int_t npix, Bool_t ab)
491{
492 if (ab)
493 SETBIT((*fABFlags)[npix/8], npix%8);
494 else
495 CLRBIT((*fABFlags)[npix/8], npix%8);
496}
497
498// --------------------------------------------------------------------------
499//
500// Fills members with information from a magic binary file.
501// WARNING: you have to use Init() before you can do this
502//
503/*
504void MRawEvtData::ReadEvt(istream &fin, Int_t posinarray)
505{
506
507 const UShort_t npic = fRunHeader->GetNumPixInCrate();
508
509 const UShort_t npos = npic*posinarray;
510
511 //const Byte_t ab = fCrateArray->GetEntry(posinarray)->GetABFlags();
512
513 for (int i=npos; i<npic+npos; i++)
514 {
515 // calc the spiral hardware pixel number
516 const UShort_t ipos = i;
517
518 // Get Hardware Id
519 const Short_t hwid = fRunHeader->GetPixAssignment(ipos);
520
521 // Check whether the pixel is connected or not
522 if (hwid==0)
523 {
524 const UShort_t n = fRunHeader->GetNumSamplesLoGain()+fRunHeader->GetNumSamplesHiGain();
525 fin.seekg(n, ios::cur);
526 return;
527 }
528
529 // -1 converts the hardware pixel Id into the software pixel index
530 const Int_t npix = (Int_t)hwid-1;
531
532 const Byte_t ab = fCrateArray->GetEntry(posinarray)->GetABFlags();
533 AddPixel(fin, npix, TESTBIT(ab, i-npos));
534 }
535}
536*/
537
538// --------------------------------------------------------------------------
539//
540// Return the size in bytes of one event data block
541//
542Int_t MRawEvtData::GetNumBytes() const
543{
544 const UShort_t nlo = fRunHeader->GetNumSamplesLoGain();
545 const UShort_t nhi = fRunHeader->GetNumSamplesHiGain();
546 const UShort_t npic = fRunHeader->GetNumPixInCrate();
547 const UShort_t nbps = fRunHeader->GetNumBytesPerSample();
548
549 return (nhi+nlo)*npic*nbps;
550}
551
552// --------------------------------------------------------------------------
553//
554// Make sure, that you skip the whole event. This function only skips a part
555// of the event - see MRawRead::SkipEvent
556//
557void MRawEvtData::SkipEvt(istream &fin)
558{
559 fin.seekg(GetNumBytes(), ios::cur);
560}
561
562Bool_t MRawEvtData::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
563{
564 *fLog << warn << "WARNING - The use of MRawEvtData::GetPixelContent is deprecated!" << endl;
565
566 /*
567 MRawEvtPixelIter Next(const_cast<MRawEvtData*>(this));
568 if (!Next.Jump(idx))
569 return kFALSE;
570
571 switch (type)
572 {
573 case 0:
574 val = Next.GetSumHiGainSamples()-(float)GetNumHiGainSamples()*fHiGainFadcSamples->GetArray()[0];
575 val*= cam.GetPixRatio(idx);
576 break;
577 case 1:
578 val = Next.GetMaxHiGainSample();
579 break;
580 case 2:
581 val = Next.GetMaxLoGainSample();
582 break;
583 case 3:
584 val = Next.GetIdxMaxHiGainSample();
585 break;
586 case 4:
587 val = Next.GetIdxMaxLoGainSample();
588 break;
589 case 5:
590 val = Next.GetIdxMaxHiLoGainSample();
591 return val >= 0;
592 }
593*/
594 return kTRUE;
595}
596
597void MRawEvtData::Copy(TObject &named)
598#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
599const
600#endif
601{
602 MRawEvtData &evt = (MRawEvtData &)named;
603
604 *evt.fHiGainPixId = *fHiGainPixId;
605 *evt.fLoGainPixId = *fLoGainPixId;
606
607 *evt.fHiGainFadcSamples = *fHiGainFadcSamples;
608 *evt.fLoGainFadcSamples = *fLoGainFadcSamples;
609
610 *evt.fABFlags = *fABFlags;
611
612 evt.fConnectedPixels = fConnectedPixels;
613
614 evt.fNumBytesPerSample = fNumBytesPerSample;
615}
Note: See TracBrowser for help on using the repository browser.