Public Member Functions | Private Attributes | List of all members
AnalysisExample::RawEVD35t Class Reference
Inheritance diagram for AnalysisExample::RawEVD35t:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Public Member Functions

 RawEVD35t (fhicl::ParameterSet const &pset)
 
virtual ~RawEVD35t ()
 
void beginJob ()
 
void beginRun (const art::Run &run)
 
void reconfigure (fhicl::ParameterSet const &pset)
 
void analyze (const art::Event &evt)
 
unsigned int getAPAindex (unsigned int apa)
 
unsigned int getTPCindex (unsigned int tpc)
 
- Public Member Functions inherited from art::EDAnalyzer
 EDAnalyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDAnalyzer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Analyzer
virtual ~Analyzer () noexcept
 
 Analyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 Analyzer (Table< Config > const &config)
 
void doBeginJob ()
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Observer
 ~Observer () noexcept
 
 Observer (Observer const &)=delete
 
 Observer (Observer &&)=delete
 
Observeroperator= (Observer const &)=delete
 
Observeroperator= (Observer &&)=delete
 
void registerProducts (ProductDescriptions &, ModuleDescription const &)
 
void fillDescriptions (ModuleDescription const &)
 
std::string const & processName () const
 
bool wantAllEvents () const
 
bool wantEvent (Event const &e)
 
fhicl::ParameterSetID selectorConfig () const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Attributes

std::string fRawDigitLabel
 
unsigned int fNticks
 
bool fUncompressWithPed
 
unsigned int fUChanMin
 
unsigned int fUChanMax
 
unsigned int fVChanMin
 
unsigned int fVChanMax
 
unsigned int fZ0ChanMin
 
unsigned int fZ0ChanMax
 
unsigned int fZ1ChanMin
 
unsigned int fZ1ChanMax
 
unsigned int fNofAPA
 
unsigned int fChansPerAPA
 
art::ServiceHandle< geo::GeometryfGeom
 
std::vector< TH2I * > fTimeChanU
 
std::vector< TH2I * > fTimeChanV
 
std::vector< TH2I * > fTimeChanZ0
 
std::vector< TH2I * > fTimeChanZ1
 
std::vector< TH2I * > fTimeChanThumbU
 
std::vector< TH2I * > fTimeChanThumbV
 
std::vector< TH2I * > fTimeChanThumbZ0
 
std::vector< TH2I * > fTimeChanThumbZ1
 
TH2I * fChargeSumU
 
TH2I * fChargeSumV
 
TH2I * fChargeSumZ
 
unsigned int UVPlane [4] ={3,2,1,0}
 
unsigned int ZPlane [8] ={7,6,5,4,3,2,1,0}
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &paths, fhicl::ParameterSet const &config)
 
detail::ProcessAndEventSelectorsprocessAndEventSelectors ()
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 56 of file RawEVD35t_module.cc.

Constructor & Destructor Documentation

AnalysisExample::RawEVD35t::RawEVD35t ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 134 of file RawEVD35t_module.cc.

135  : EDAnalyzer(parameterSet)
136  {
137  this->reconfigure(parameterSet);
138  }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:27
void reconfigure(fhicl::ParameterSet const &pset)
AnalysisExample::RawEVD35t::~RawEVD35t ( )
virtual

Definition at line 155 of file RawEVD35t_module.cc.

155  {
156  }

Member Function Documentation

void AnalysisExample::RawEVD35t::analyze ( const art::Event evt)

Definition at line 320 of file RawEVD35t_module.cc.

320  {
321 
322  unsigned int tpcid, cryoid;
323  std::stringstream thumbnameZ0, thumbnameZ1;
324 
325  // get the objects holding all of the raw data information
327  std::cout << "raw digit label check: " << fRawDigitLabel << std::endl;
328  event.getByLabel(fRawDigitLabel, Raw);
329 
330  // put it in a more easily usable form
331  std::vector< art::Ptr<raw::RawDigit> > Digits;
332  art::fill_ptr_vector(Digits, Raw);
333 
334  //loop through all RawDigits (over entire channels)
335  for(size_t d = 0; d < Digits.size(); d++){
337  digit=Digits.at(d);
338 
339  // get the channel number for this digit
340  uint32_t chan = digit->Channel();
341  unsigned int apa = std::floor( chan/fChansPerAPA );
342  tpcid=fGeom->ChannelToWire(chan)[0].TPC;
343  cryoid=fGeom->ChannelToWire(chan)[0].Cryostat;
344 
345  int nSamples = digit->Samples();
346  std::vector<short> uncompressed(nSamples);
347  // raw::Uncompress(digit->ADCs(), uncompressed, digit->Compression());
348  int pedestal = (int)digit->GetPedestal();
349  // std::cout << "channel " << chan << " pedestal " << pedestal << std::endl;
350  // uncompress the data
351  if (fUncompressWithPed){
352  raw::Uncompress(digit->ADCs(), uncompressed, pedestal, digit->Compression());
353  }
354  else{
355  raw::Uncompress(digit->ADCs(), uncompressed, digit->Compression());
356  }
357 
358  // subtract pedestals
359  std::vector<short> ladc(nSamples);
360  for (int i=0; i<nSamples; i++) ladc[i]=uncompressed[i]-pedestal;
361  // for (int i=0; i<nSamples; i++) {ladc[i]=uncompressed[i]-pedestal;
362  // if (i<10) std::cout << uncompressed[i] << " " << ladc[i] << std::endl;
363  // }
364  if( fGeom->View(chan) == geo::kU ){
365  for(unsigned int l=0;l<ladc.size();l++) {
366  if(ladc.at(l)!=0){
367  fTimeChanU[apa]->Fill(chan,l, ladc.at(l));
368  if(ladc.at(l)>0) fTimeChanThumbU[apa]->Fill(chan,l, ladc.at(l));
369  fChargeSumU->Fill(1+getAPAindex(apa)%2, 2-(getAPAindex(apa)/2),std::abs(ladc.at(l))/2);
370  }
371  }
372  }
373  if( fGeom->View(chan) == geo::kV ){
374  for(unsigned int l=0;l<ladc.size();l++) {
375  if(ladc.at(l)!=0){
376  fTimeChanV[apa]->Fill(chan,l, ladc.at(l));
377  if(ladc.at(l)>0) fTimeChanThumbV[apa]->Fill(chan,l, ladc.at(l));
378  fChargeSumV->Fill(1+getAPAindex(apa)%2, 2-(getAPAindex(apa)/2),std::abs(ladc.at(l))/2);
379  }
380  }
381  }
382  if ( fGeom->View(chan) == geo::kZ && fGeom->ChannelToWire(chan)[0].TPC % 2 == 0 ){
383  for(unsigned int l=0;l<ladc.size();l++) {
384  if(ladc.at(l)!=0){
385  fTimeChanZ0[apa]->Fill(chan,l, ladc.at(l));
386  if(ladc.at(l)>0) fTimeChanThumbZ0[apa]->Fill(chan,l, ladc.at(l));
387  fChargeSumZ->Fill(1+getTPCindex(cryoid*4+tpcid)%4, 4-(getTPCindex(cryoid*4+tpcid)/4),ladc.at(l));
388  }
389  }
390  }
391  if ( fGeom->View(chan) == geo::kZ && fGeom->ChannelToWire(chan)[0].TPC % 2 == 1 ){
392  for(unsigned int l=0;l<ladc.size();l++) {
393  if(ladc.at(l)!=0){
394  fTimeChanZ1[apa]->Fill(chan,l, ladc.at(l));
395  if(ladc.at(l)>0) fTimeChanThumbZ1[apa]->Fill(chan,l, ladc.at(l));
396  fChargeSumZ->Fill(1+getTPCindex(cryoid*4+tpcid)%4, 4-(getTPCindex(cryoid*4+tpcid)/4),ladc.at(l));
397  }
398  }
399  }
400 
401  } // end RawDigit loop
402 
403  return;
404  }
float GetPedestal() const
Definition: RawDigit.h:214
const ADCvector_t & ADCs() const
Reference to the compressed ADC count vector.
Definition: RawDigit.h:210
ULong64_t Samples() const
Number of samples in the uncompressed ADC data.
Definition: RawDigit.h:213
art::ServiceHandle< geo::Geometry > fGeom
Planes which measure V.
Definition: geo_types.h:126
ChannelID_t Channel() const
DAQ channel this raw data was read from.
Definition: RawDigit.h:212
std::vector< geo::WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
Planes which measure Z direction.
Definition: geo_types.h:128
static QStrList * l
Definition: config.cpp:1044
std::vector< TH2I * > fTimeChanU
unsigned int getTPCindex(unsigned int tpc)
T abs(T value)
Planes which measure U.
Definition: geo_types.h:125
std::vector< TH2I * > fTimeChanThumbU
std::vector< TH2I * > fTimeChanV
unsigned int uint32_t
Definition: stdint.h:126
unsigned int getAPAindex(unsigned int apa)
std::vector< TH2I * > fTimeChanThumbV
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
raw::Compress_t Compression() const
Compression algorithm used to store the ADC counts.
Definition: RawDigit.h:216
std::vector< TH2I * > fTimeChanThumbZ1
std::vector< TH2I * > fTimeChanZ1
type * at(uint i)
Definition: qinternallist.h:81
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:291
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:755
std::vector< TH2I * > fTimeChanZ0
QTextStream & endl(QTextStream &s)
std::vector< TH2I * > fTimeChanThumbZ0
void AnalysisExample::RawEVD35t::beginJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 160 of file RawEVD35t_module.cc.

160  {
161  // Access ART's TFileService, which will handle creating and writing
162  // histograms and n-tuples for us.
164 
165  //Histogram names and titles
166  std::stringstream name, title;
167 
168  unsigned int UChMin;
169  unsigned int UChMax;
170  unsigned int VChMin;
171  unsigned int VChMax;
172  unsigned int Z0ChMin;
173  unsigned int Z0ChMax;
174  unsigned int Z1ChMin;
175  unsigned int Z1ChMax;
176  TH2I* TempHisto;
177 
178  std::ofstream outfile;
179  outfile.open("msglog.txt");
180 
181  // loop through channels in the first APA to find the channel boundaries for each view
182  // will adjust for desired APA after
183  fUChanMin = 0;
186  fZ1ChanMax = fChansPerAPA - 1;
187  for ( unsigned int c = fUChanMin + 1; c < fZ1ChanMax; c++ ){
188  if ( fGeom->View(c) == geo::kV && fGeom->View(c-1) == geo::kU ){
189  fVChanMin = c;
190  fUChanMax = c - 1;
191  }
192  if ( fGeom->View(c) == geo::kZ && fGeom->View(c-1) == geo::kV ){
193  fZ0ChanMin = c;
194  fVChanMax = c-1;
195  }
196  if ( fGeom->View(c) == geo::kZ && fGeom->ChannelToWire(c)[0].TPC == fGeom->ChannelToWire(c-1)[0].TPC + 1 ){
197  fZ1ChanMin = c;
198  fZ0ChanMax = c-1;
199  }
200  }
201 
202 outfile<<fChansPerAPA<<" "<<fGeom->Ncryostats()<<" "<<fNofAPA<<std::endl;
203 
204  unsigned int minT = 0;
205  unsigned int maxT = 0;
206  minT = 0;
207  maxT = fNticks;
208  unsigned int binT = (maxT-minT);
209 
210  for(unsigned int i=0;i<fNofAPA;i++){
211  UChMin=fUChanMin + i*fChansPerAPA;
212  UChMax=fUChanMax + i*fChansPerAPA;
213  VChMin=fVChanMin + i*fChansPerAPA;
214  VChMax=fVChanMax + i*fChansPerAPA;
215  Z0ChMin=fZ0ChanMin + i*fChansPerAPA;
216  Z0ChMax=fZ0ChanMax + i*fChansPerAPA;
217  Z1ChMin=fZ1ChanMin + i*fChansPerAPA;
218  Z1ChMax=fZ1ChanMax + i*fChansPerAPA;
219 
220  // construct the histograms; TH2 constructors: ("Name", "Title", NxBin, xMin, xMax, NyBin, yMax, yMin)
221  name.str("");
222  name << "fTimeChanU";
223  name << i;
224  title.str("");
225  title << "Time vs Channel(Plane U, APA";
226  title << i<<")";
227  TempHisto = tfs->make<TH2I>(name.str().c_str(),title.str().c_str(), UChMax - UChMin + 1, UChMin, UChMax, binT, minT, maxT);
228  fTimeChanU.push_back(TempHisto);
229 
230  name.str("");
231  name << "fTimeChanThumbU";
232  name << i;
233  title.str("");
234  title << "Time vs Channel(Plane U, APA";
235  title << i<<")";
236  TempHisto = tfs->make<TH2I>(name.str().c_str(),title.str().c_str(), 138, UChMin, UChMax, 138, minT, maxT);
237  fTimeChanThumbU.push_back(TempHisto);
238 
239  name.str("");
240  name << "fTimeChanV";
241  name << i;
242  title.str("");
243  title << "Time vs Channel(Plane V, APA";
244  title << i<<")";
245  TempHisto = tfs->make<TH2I>(name.str().c_str(),title.str().c_str(), VChMax - VChMin + 1, VChMin, VChMax, binT, minT, maxT);
246  fTimeChanV.push_back(TempHisto);
247 
248  name.str("");
249  name << "fTimeChanThumbV";
250  name << i;
251  title.str("");
252  title << "Time vs Channel(Plane V, APA";
253  title << i<<")";
254  TempHisto = tfs->make<TH2I>(name.str().c_str(),title.str().c_str(), 138, VChMin, VChMax, 138, minT, maxT);
255  fTimeChanThumbV.push_back(TempHisto);
256 
257  name.str("");
258  name << "fTimeChanZ0";
259  name << i;
260  title.str("");
261  title << "Time vs Channel(Plane Z0, APA";
262  title <<i<<")";
263  TempHisto = tfs->make<TH2I>(name.str().c_str(),title.str().c_str(), Z0ChMax - Z0ChMin + 1, Z0ChMin, Z0ChMax, binT, minT, maxT);
264  fTimeChanZ0.push_back(TempHisto);
265 
266  name.str("");
267  name << "fTimeChanThumbZ0";
268  name << i;
269  title.str("");
270  title << "Time vs Channel(Plane Z0, APA";
271  title <<i<<")";
272  TempHisto = tfs->make<TH2I>(name.str().c_str(),title.str().c_str(), 111, Z0ChMin, Z0ChMax, 111, minT, maxT);
273  fTimeChanThumbZ0.push_back(TempHisto);
274 
275  name.str("");
276  name << "fTimeChanZ1";
277  name << i;
278  title.str("");
279  title << "Time vs Channel(Plane Z1, APA";
280  title << i<<")";
281  TempHisto = tfs->make<TH2I>(name.str().c_str(),title.str().c_str(), Z1ChMax - Z1ChMin + 1, Z1ChMin, Z1ChMax, binT, minT, maxT);
282  fTimeChanZ1.push_back(TempHisto);
283 
284  name.str("");
285  name << "fTimeChanThumbZ1";
286  name << i;
287  title.str("");
288  title << "Time vs Channel(Plane Z1, APA";
289  title << i<<")";
290  TempHisto = tfs->make<TH2I>(name.str().c_str(),title.str().c_str(), 111, Z1ChMin, Z1ChMax, 111, minT, maxT);
291  fTimeChanThumbZ1.push_back(TempHisto);
292 
293  fTimeChanU[i]->SetStats(0); fTimeChanV[i]->SetStats(0);
294  fTimeChanZ0[i]->SetStats(0); fTimeChanZ1[i]->SetStats(0);
295  fTimeChanThumbU[i]->SetStats(0); fTimeChanThumbV[i]->SetStats(0);
296  fTimeChanThumbZ0[i]->SetStats(0); fTimeChanThumbZ1[i]->SetStats(0);
297 
298  fTimeChanU[i]->GetXaxis()->SetTitle("Channel"); fTimeChanU[i]->GetYaxis()->SetTitle("TDC");
299  fTimeChanV[i]->GetXaxis()->SetTitle("Channel"); fTimeChanV[i]->GetYaxis()->SetTitle("TDC");
300  fTimeChanZ0[i]->GetXaxis()->SetTitle("Channel"); fTimeChanZ0[i]->GetYaxis()->SetTitle("TDC");
301  fTimeChanZ1[i]->GetXaxis()->SetTitle("Channel"); fTimeChanZ1[i]->GetYaxis()->SetTitle("TDC");
302  }
303 
304  fChargeSumU= tfs->make<TH2I>("hChargeSumU","Charge Sum on U-planes", 2,0.5,2.5,2,0.5,2.5 );
305  fChargeSumV= tfs->make<TH2I>("hChargeSumV","Charge Sum on V-planes", 2,0.5,2.5,2,0.5,2.5 );
306  fChargeSumZ= tfs->make<TH2I>("hChargeSumZ","Charge Sum on Z-planes", 4,0.5,4.5,4,0.5,4.5);
307  fChargeSumU->SetStats(0); fChargeSumV->SetStats(0); fChargeSumZ->SetStats(0);
308 
309  }
static QCString name
Definition: declinfo.cpp:673
art::ServiceHandle< geo::Geometry > fGeom
Planes which measure V.
Definition: geo_types.h:126
std::vector< geo::WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
Planes which measure Z direction.
Definition: geo_types.h:128
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
std::vector< TH2I * > fTimeChanU
unsigned int Nchannels() const
Returns the number of TPC readout channels in the detector.
Planes which measure U.
Definition: geo_types.h:125
std::vector< TH2I * > fTimeChanThumbU
std::vector< TH2I * > fTimeChanV
std::vector< TH2I * > fTimeChanThumbV
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
std::vector< TH2I * > fTimeChanThumbZ1
std::vector< TH2I * > fTimeChanZ1
std::vector< TH2I * > fTimeChanZ0
QTextStream & endl(QTextStream &s)
std::vector< TH2I * > fTimeChanThumbZ0
void AnalysisExample::RawEVD35t::beginRun ( const art::Run run)

Definition at line 313 of file RawEVD35t_module.cc.

313  {
314 
315  }
unsigned int AnalysisExample::RawEVD35t::getAPAindex ( unsigned int  apa)
inline

Definition at line 70 of file RawEVD35t_module.cc.

70  {
71  for (int k=0; k<4; k++){if(apa==UVPlane[k]) {return k; }}
72  return -1;
73  }
unsigned int AnalysisExample::RawEVD35t::getTPCindex ( unsigned int  tpc)
inline

Definition at line 75 of file RawEVD35t_module.cc.

75  {
76  for (int k=0; k<8; k++){if(tpc==ZPlane[k]) {return k; }}
77  return -1;
78  }
void AnalysisExample::RawEVD35t::reconfigure ( fhicl::ParameterSet const &  pset)

Definition at line 143 of file RawEVD35t_module.cc.

143  {
144  fRawDigitLabel = p.get< std::string >("RawDigitLabel");
145  // auto const* fDetProp = lar::providerFrom<detinfo::DetectorPropertiesService>();
146  // fNticks = fDetProp->NumberTimeSamples();
147  fNticks = (unsigned int) p.get< int >("TicksToDraw");
148  fUncompressWithPed = p.get< bool >("UncompressWithPed", true);
149  return;
150  }
std::string string
Definition: nybbler.cc:12
p
Definition: test.py:223

Member Data Documentation

unsigned int AnalysisExample::RawEVD35t::fChansPerAPA
private

Definition at line 101 of file RawEVD35t_module.cc.

TH2I* AnalysisExample::RawEVD35t::fChargeSumU
private

Definition at line 122 of file RawEVD35t_module.cc.

TH2I* AnalysisExample::RawEVD35t::fChargeSumV
private

Definition at line 123 of file RawEVD35t_module.cc.

TH2I* AnalysisExample::RawEVD35t::fChargeSumZ
private

Definition at line 124 of file RawEVD35t_module.cc.

art::ServiceHandle<geo::Geometry> AnalysisExample::RawEVD35t::fGeom
private

Definition at line 110 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fNofAPA
private

Definition at line 100 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fNticks
private

Definition at line 87 of file RawEVD35t_module.cc.

std::string AnalysisExample::RawEVD35t::fRawDigitLabel
private

Definition at line 83 of file RawEVD35t_module.cc.

std::vector<TH2I*> AnalysisExample::RawEVD35t::fTimeChanThumbU
private

Definition at line 117 of file RawEVD35t_module.cc.

std::vector<TH2I*> AnalysisExample::RawEVD35t::fTimeChanThumbV
private

Definition at line 118 of file RawEVD35t_module.cc.

std::vector<TH2I*> AnalysisExample::RawEVD35t::fTimeChanThumbZ0
private

Definition at line 119 of file RawEVD35t_module.cc.

std::vector<TH2I*> AnalysisExample::RawEVD35t::fTimeChanThumbZ1
private

Definition at line 120 of file RawEVD35t_module.cc.

std::vector<TH2I*> AnalysisExample::RawEVD35t::fTimeChanU
private

Definition at line 112 of file RawEVD35t_module.cc.

std::vector<TH2I*> AnalysisExample::RawEVD35t::fTimeChanV
private

Definition at line 113 of file RawEVD35t_module.cc.

std::vector<TH2I*> AnalysisExample::RawEVD35t::fTimeChanZ0
private

Definition at line 114 of file RawEVD35t_module.cc.

std::vector<TH2I*> AnalysisExample::RawEVD35t::fTimeChanZ1
private

Definition at line 115 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fUChanMax
private

Definition at line 92 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fUChanMin
private

Definition at line 91 of file RawEVD35t_module.cc.

bool AnalysisExample::RawEVD35t::fUncompressWithPed
private

Definition at line 88 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fVChanMax
private

Definition at line 94 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fVChanMin
private

Definition at line 93 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fZ0ChanMax
private

Definition at line 96 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fZ0ChanMin
private

Definition at line 95 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fZ1ChanMax
private

Definition at line 98 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::fZ1ChanMin
private

Definition at line 97 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::UVPlane[4] ={3,2,1,0}
private

Definition at line 126 of file RawEVD35t_module.cc.

unsigned int AnalysisExample::RawEVD35t::ZPlane[8] ={7,6,5,4,3,2,1,0}
private

Definition at line 127 of file RawEVD35t_module.cc.


The documentation for this class was generated from the following file: