source: trunk/Mars/mhist/MHFadcPix.cc@ 19427

Last change on this file since 19427 was 2333, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 4.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! Author(s): Harald Kornmayer 1/2001
20!
21! Copyright: MAGIC Software Development, 2000-2002
22!
23!
24\* ======================================================================== */
25
26///////////////////////////////////////////////////////////////////////
27//
28// MHFadcPix
29//
30// This container stores a histogram to display an Fadc Spektrum.
31// The spektrum of all the values measured by the Fadcs.
32//
33///////////////////////////////////////////////////////////////////////
34#include "MHFadcPix.h"
35
36#include <TPad.h>
37
38#include "MH.h"
39#include "MBinning.h"
40
41#include "MRawEvtData.h"
42#include "MRawEvtPixelIter.h"
43
44ClassImp(MHFadcPix);
45
46// --------------------------------------------------------------------------
47//
48// Creates the histograms for lo and hi gain of one pixel
49//
50MHFadcPix::MHFadcPix(Int_t pixid, Type_t t)
51 : fPixId(pixid), fType(t)
52{
53 fHistHi.SetName(pixid>=0 ? Form("HiGain%03d", pixid) : "HiGain");
54 fHistHi.SetTitle(pixid>=0 ? Form("Hi Gain Pixel #%d", pixid) : "Hi Gain Samples");
55 fHistHi.SetDirectory(NULL);
56 fHistHi.UseCurrentStyle();
57
58 fHistLo.SetName(pixid>=0 ? Form("LoGain%03d", pixid) : "LoGain");
59 fHistLo.SetTitle(pixid>=0 ? Form("Lo Gain Pixel #%d", pixid) : "Lo Gain Samples");
60 fHistLo.SetDirectory(NULL);
61 fHistLo.UseCurrentStyle();
62
63 if (fType==kValue)
64 {
65 fHistHi.SetXTitle("Signal/FADC Units");
66 fHistLo.SetXTitle("Signal/FADC Units");
67 fHistHi.SetYTitle("Count");
68 fHistLo.SetYTitle("Count");
69
70 MBinning bins;
71 bins.SetEdges(255, -.5, 255.5);
72
73 bins.Apply(fHistHi);
74 bins.Apply(fHistLo);
75 }
76 else
77 {
78 fHistHi.SetXTitle("Time/FADC Slices");
79 fHistLo.SetXTitle("Time/FADC Slices");
80 fHistLo.SetYTitle("Sum Signal/FADC Units");
81 fHistHi.SetYTitle("Sum Signal/FADC Units");
82 }
83}
84
85void MHFadcPix::Init(Byte_t nhi, Byte_t nlo)
86{
87 MBinning bins;
88
89 bins.SetEdges(nhi, -.5, -.5+nhi);
90 bins.Apply(fHistHi);
91
92 bins.SetEdges(nlo, -.5, -.5+nlo);
93 bins.Apply(fHistLo);
94}
95
96Bool_t MHFadcPix::Fill(const MRawEvtData &evt)
97{
98 MRawEvtPixelIter pixel((MRawEvtData*)&evt);
99
100 if (!pixel.Jump(fPixId))
101 return kTRUE;
102
103 const Int_t nhisamples = evt.GetNumHiGainSamples();
104
105 if (fType==kValue)
106 for (Int_t i=0; i<nhisamples; i++)
107 fHistHi.Fill(pixel.GetHiGainSamples()[i]);
108 else
109 for (Int_t i=0; i<nhisamples; i++)
110 fHistHi.Fill(i, pixel.GetHiGainSamples()[i]);
111
112 if (!pixel.HasLoGain())
113 return kTRUE;
114
115 const Int_t nlosamples = evt.GetNumLoGainSamples();
116
117 if (fType==kValue)
118 for (Int_t i=0; i<nlosamples; i++)
119 fHistLo.Fill(pixel.GetLoGainSamples()[i]);
120 else
121 for (Int_t i=0; i<nlosamples; i++)
122 fHistLo.Fill(i, pixel.GetLoGainSamples()[i]);
123
124 return kTRUE;
125}
126
127// --------------------------------------------------------------------------
128void MHFadcPix::DrawHi()
129{
130 fHistHi.Draw();
131}
132
133// --------------------------------------------------------------------------
134void MHFadcPix::DrawLo()
135{
136 fHistLo.Draw();
137}
138
139// --------------------------------------------------------------------------
140//
141// We need our own clone function to get rid of the histogram in any
142// directory
143//
144TObject *MHFadcPix::Clone(const char *) const
145{
146 MHFadcPix &pix = *(MHFadcPix*)TObject::Clone();
147
148 pix.fHistHi.SetDirectory(NULL);
149 pix.fHistLo.SetDirectory(NULL);
150
151 return &pix;
152}
153
154// --------------------------------------------------------------------------
155void MHFadcPix::Draw(Option_t *)
156{
157 if (!gPad)
158 {
159 const char *name = StrDup(fPixId ? Form("Pixel #%d", fPixId) : "Pixel");
160 MH::MakeDefCanvas(name, fPixId ? Form("%s FADC Samples", name) : "FADC Samples");
161 delete [] name;
162 }
163
164 gPad->Divide(1, 2);
165
166 gPad->cd(1);
167 fHistHi.Draw();
168
169 gPad->cd(2);
170 fHistLo.Draw();
171}
Note: See TracBrowser for help on using the repository browser.