CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
ZdcTBAnalysis Class Reference

#include <ZdcTBAnalysis.h>

Public Member Functions

void analyze (const ZDCRecHitCollection &hf)
 
void analyze (const HcalTBTriggerData &trg)
 
void analyze (const HcalTBBeamCounters &bc)
 
void analyze (const HcalTBTiming &times)
 
void analyze (const HcalTBEventPosition &chpos)
 
void done ()
 
void fillTree ()
 
void setup (const std::string &histoFileName)
 
 ZdcTBAnalysis ()
 

Private Attributes

ADC adc
 
double beam_coincidence [5]
 
double BH1adc
 
double bh1hits [5]
 
double BH2adc
 
double bh2hits [5]
 
double BH3adc
 
double bh3hits [5]
 
double BH4adc
 
double bh4hits [5]
 
CHAMB chamb
 
double CK1adc
 
double CK2adc
 
double CK3adc
 
HcalZDCDetId detID
 
double Ecal7x7adc
 
double energy
 
int eventNumber
 
int ichannel
 
int idepth
 
bool isBeamTrigger
 
bool isCalibTrigger
 
int isection
 
bool isFakeTrigger
 
int iside
 
bool isInSpillPedestalTrigger
 
bool isLaserTrigger
 
bool isLedTrigger
 
bool isOutSpillPedestalTrigger
 
bool isSpillTrigger
 
double laser_flash
 
double m1hits [5]
 
double m2hits [5]
 
double m3hits [5]
 
TFile * outFile
 
double qie_phase
 
int runNumber
 
double S1adc
 
double s1hits [5]
 
double S2adc
 
double s2hits [5]
 
double S3adc
 
double s3hits [5]
 
double S4adc
 
double s4hits [5]
 
double Sci521adc
 
double Sci528adc
 
double SciVLEadc
 
TDC tdc
 
double TOF1_time
 
double TOF1adc
 
double TOF2_time
 
double TOF2adc
 
TRIGGER trigger
 
double trigger_time
 
double ttc_L1a_time
 
double V3adc
 
double V6adc
 
double VH1adc
 
double VH2adc
 
double VH3adc
 
double VH4adc
 
double VM1adc
 
double VM2adc
 
double VM3adc
 
double VM4adc
 
double VM5adc
 
double VM6adc
 
double VM7adc
 
double VM8adc
 
double VMadc
 
double VMBadc
 
double VMFadc
 
std::vector< double > wcax
 
std::vector< double > wcay
 
std::vector< double > wcbx
 
std::vector< double > wcby
 
std::vector< double > wccx
 
std::vector< double > wccy
 
std::vector< double > wcdx
 
std::vector< double > wcdy
 
std::vector< double > wcex
 
std::vector< double > wcey
 
std::vector< double > wcfx
 
std::vector< double > wcfy
 
std::vector< double > wcgx
 
std::vector< double > wcgy
 
std::vector< double > wchx
 
std::vector< double > wchy
 
TTree * ZdcAnalize
 
ZDCN zdcn
 
ZDCP zdcp
 

Detailed Description

Definition at line 132 of file ZdcTBAnalysis.h.

Constructor & Destructor Documentation

◆ ZdcTBAnalysis()

ZdcTBAnalysis::ZdcTBAnalysis ( )

Definition at line 6 of file ZdcTBAnalysis.cc.

6 { ; }

Member Function Documentation

◆ analyze() [1/5]

void ZdcTBAnalysis::analyze ( const ZDCRecHitCollection hf)

Definition at line 241 of file ZdcTBAnalysis.cc.

References edm::SortedCollection< T, SORT >::begin(), HcalZDCDetId::channel(), gather_cfg::cout, HcalZDCDetId::depth(), detID, edm::SortedCollection< T, SORT >::end(), energy, mps_fire::i, ichannel, idepth, isection, iside, HcalZDCDetId::section(), ZDCN::zdcEMMod1, ZDCP::zdcEMMod1, ZDCN::zdcEMMod2, ZDCP::zdcEMMod2, ZDCN::zdcEMMod3, ZDCP::zdcEMMod3, ZDCN::zdcEMMod4, ZDCP::zdcEMMod4, ZDCN::zdcEMMod5, ZDCP::zdcEMMod5, ZDCN::zdcHADMod1, ZDCP::zdcHADMod1, ZDCN::zdcHADMod2, ZDCP::zdcHADMod2, ZDCN::zdcHADMod3, ZDCP::zdcHADMod3, ZDCN::zdcHADMod4, ZDCP::zdcHADMod4, zdcn, zdcp, ZDCN::zdcScint1, ZDCP::zdcScint1, and HcalZDCDetId::zside().

Referenced by ZdcTBAnalyzer::analyze().

241  {
242  // zdc hits
243  std::cout << "****************************************************" << std::endl;
245  for (i = zdcHits.begin(); i != zdcHits.end(); i++) {
246  energy = i->energy();
247  detID = i->id();
248  iside = detID.zside();
249  isection = detID.section();
250  ichannel = detID.channel();
251  idepth = detID.depth();
252  std::cout << "energy: " << energy << " detID: " << detID << " side: " << iside << " section: " << isection
253  << " channel: " << ichannel << " depth: " << idepth << std::endl;
254 
255  if (iside > 0) {
256  if (ichannel == 1 && isection == 1)
258  if (ichannel == 2 && isection == 1)
260  if (ichannel == 3 && isection == 1)
262  if (ichannel == 4 && isection == 1)
264  if (ichannel == 5 && isection == 1)
266  if (ichannel == 1 && isection == 2)
268  if (ichannel == 2 && isection == 2)
270  if (ichannel == 3 && isection == 2)
272  if (ichannel == 4 && isection == 2)
274  if (ichannel == 1 && isection == 3)
276  }
277  if (iside < 0) {
278  if (ichannel == 1 && isection == 1)
280  if (ichannel == 2 && isection == 1)
282  if (ichannel == 3 && isection == 1)
284  if (ichannel == 4 && isection == 1)
286  if (ichannel == 5 && isection == 1)
288  if (ichannel == 1 && isection == 2)
290  if (ichannel == 2 && isection == 2)
292  if (ichannel == 3 && isection == 2)
294  if (ichannel == 4 && isection == 2)
296  if (ichannel == 1 && isection == 3)
298  }
299  }
300 }
double zdcEMMod2
Section section() const
get the section
Definition: HcalZDCDetId.cc:44
double zdcEMMod3
double zdcHADMod3
double zdcEMMod3
std::vector< T >::const_iterator const_iterator
double zdcEMMod5
double zdcEMMod1
double zdcScint1
double zdcEMMod1
double zdcEMMod2
double zdcHADMod3
double zdcHADMod2
double zdcScint1
int depth() const
get the depth (1 for EM, channel + 1 for HAD, 2 for RPD, not sure yet for LUM, leave as default) ...
Definition: HcalZDCDetId.cc:51
double zdcHADMod4
double zdcEMMod4
HcalZDCDetId detID
double zdcHADMod1
double zdcHADMod2
double zdcHADMod4
double zdcHADMod1
double zdcEMMod5
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalZDCDetId.h:39
int channel() const
get the channel
Definition: HcalZDCDetId.cc:63
double zdcEMMod4

◆ analyze() [2/5]

void ZdcTBAnalysis::analyze ( const HcalTBTriggerData trg)

Definition at line 55 of file ZdcTBAnalysis.cc.

References HcalTBTriggerData::eventNumber(), eventNumber, isBeamTrigger, isCalibTrigger, isFakeTrigger, isInSpillPedestalTrigger, isLaserTrigger, isLedTrigger, isOutSpillPedestalTrigger, isSpillTrigger, HcalTBTriggerData::runNumber(), runNumber, HcalTBTriggerData::wasBeamTrigger(), HcalTBTriggerData::wasFakeTrigger(), HcalTBTriggerData::wasInSpill(), HcalTBTriggerData::wasInSpillPedestalTrigger(), HcalTBTriggerData::wasLaserTrigger(), HcalTBTriggerData::wasLEDTrigger(), HcalTBTriggerData::wasOutSpillPedestalTrigger(), and HcalTBTriggerData::wasSpillIgnorantPedestalTrigger().

55  {
56  // trigger
57  trigger.runNum = runNumber = trg.runNumber();
58  trigger.eventNum = eventNumber = trg.eventNumber();
66  isSpillTrigger = trg.wasInSpill();
67 
68  trigger.beamTrigger = trigger.fakeTrigger = trigger.calibTrigger = trigger.outSpillPedestalTrigger =
69  trigger.inSpillPedestalTrigger = trigger.laserTrigger = trigger.ledTrigger = trigger.spillTrigger = 0;
70 
71  if (isBeamTrigger)
72  trigger.beamTrigger = 1;
73  if (isFakeTrigger)
74  trigger.fakeTrigger = 1;
75  if (isCalibTrigger)
76  trigger.calibTrigger = 1;
78  trigger.outSpillPedestalTrigger = 1;
80  trigger.inSpillPedestalTrigger = 1;
81  if (isLaserTrigger)
82  trigger.laserTrigger = 1;
83  if (isLedTrigger)
84  trigger.ledTrigger = 1;
85  if (isSpillTrigger)
86  trigger.spillTrigger = 1;
87 }
uint16_t eventNumber() const
Returns the event number of this trigger.
uint32_t runNumber() const
Returns the current run number.
bool wasLEDTrigger() const
returns true if this was a LED trigger
bool wasSpillIgnorantPedestalTrigger() const
returns true if this trigger was a calibration trigger
bool wasInSpillPedestalTrigger() const
returns true if this was an in-spill pedestal trigger
bool wasOutSpillPedestalTrigger() const
returns true if this was an out-of-spill pedestal trigger
bool wasFakeTrigger() const
returns true if this trigger was fake (from a non-H2 manager)
bool wasInSpill() const
returns true if the "spill" bit was set
bool wasLaserTrigger() const
returns true if this was a laser trigger
bool isOutSpillPedestalTrigger
bool wasBeamTrigger() const
returns true if this trigger came from beam data
bool isInSpillPedestalTrigger

◆ analyze() [3/5]

void ZdcTBAnalysis::analyze ( const HcalTBBeamCounters bc)

Definition at line 139 of file ZdcTBAnalysis.cc.

References adc, ADC::BH1, HcalTBBeamCounters::BH1adc(), BH1adc, ADC::BH2, HcalTBBeamCounters::BH2adc(), BH2adc, ADC::BH3, HcalTBBeamCounters::BH3adc(), BH3adc, ADC::BH4, HcalTBBeamCounters::BH4adc(), BH4adc, ADC::CK1, HcalTBBeamCounters::CK1adc(), CK1adc, ADC::CK2, HcalTBBeamCounters::CK2adc(), CK2adc, ADC::CK3, HcalTBBeamCounters::CK3adc(), CK3adc, HcalTBBeamCounters::Ecal7x7(), ADC::Ecal7x7, Ecal7x7adc, ADC::S1, HcalTBBeamCounters::S1adc(), S1adc, ADC::S2, HcalTBBeamCounters::S2adc(), S2adc, ADC::S3, HcalTBBeamCounters::S3adc(), S3adc, ADC::S4, HcalTBBeamCounters::S4adc(), S4adc, ADC::Sci521, HcalTBBeamCounters::Sci521adc(), Sci521adc, ADC::Sci528, HcalTBBeamCounters::Sci528adc(), Sci528adc, ADC::SciVLE, HcalTBBeamCounters::SciVLEadc(), SciVLEadc, ADC::TOF1, TOF1adc, HcalTBBeamCounters::TOF1Sadc(), ADC::TOF2, TOF2adc, HcalTBBeamCounters::TOF2Sadc(), ADC::V3, HcalTBBeamCounters::V3adc(), V3adc, ADC::V6, HcalTBBeamCounters::V6adc(), V6adc, ADC::VH1, HcalTBBeamCounters::VH1adc(), VH1adc, ADC::VH2, HcalTBBeamCounters::VH2adc(), VH2adc, ADC::VH3, HcalTBBeamCounters::VH3adc(), VH3adc, ADC::VH4, HcalTBBeamCounters::VH4adc(), VH4adc, ADC::VM, ADC::VM1, HcalTBBeamCounters::VM1adc(), VM1adc, ADC::VM2, HcalTBBeamCounters::VM2adc(), VM2adc, ADC::VM3, HcalTBBeamCounters::VM3adc(), VM3adc, ADC::VM4, HcalTBBeamCounters::VM4adc(), VM4adc, ADC::VM5, HcalTBBeamCounters::VM5adc(), VM5adc, ADC::VM6, HcalTBBeamCounters::VM6adc(), VM6adc, ADC::VM7, HcalTBBeamCounters::VM7adc(), VM7adc, ADC::VM8, HcalTBBeamCounters::VM8adc(), VM8adc, HcalTBBeamCounters::VMadc(), VMadc, ADC::VMB, HcalTBBeamCounters::VMBadc(), VMBadc, ADC::VMF, HcalTBBeamCounters::VMFadc(), and VMFadc.

139  {
140  //beam counters
141  adc.VM = VMadc = bc.VMadc();
142  adc.V3 = V3adc = bc.V3adc();
143  adc.V6 = V6adc = bc.V6adc();
144  adc.VH1 = VH1adc = bc.VH1adc();
145  adc.VH2 = VH2adc = bc.VH2adc();
146  adc.VH3 = VH3adc = bc.VH3adc();
147  adc.VH4 = VH4adc = bc.VH4adc();
148  adc.Ecal7x7 = Ecal7x7adc = bc.Ecal7x7();
149  adc.Sci521 = Sci521adc = bc.Sci521adc();
150  adc.Sci528 = Sci528adc = bc.Sci528adc();
151  adc.CK1 = CK1adc = bc.CK1adc();
152  adc.CK2 = CK2adc = bc.CK2adc();
153  adc.CK3 = CK3adc = bc.CK3adc();
154  adc.SciVLE = SciVLEadc = bc.SciVLEadc();
155  adc.S1 = S1adc = bc.S1adc();
156  adc.S2 = S2adc = bc.S2adc();
157  adc.S3 = S3adc = bc.S3adc();
158  adc.S4 = S4adc = bc.S4adc();
159  adc.VMF = VMFadc = bc.VMFadc();
160  adc.VMB = VMBadc = bc.VMBadc();
161  adc.VM1 = VM1adc = bc.VM1adc();
162  adc.VM2 = VM2adc = bc.VM2adc();
163  adc.VM3 = VM3adc = bc.VM3adc();
164  adc.VM4 = VM4adc = bc.VM4adc();
165  adc.VM5 = VM5adc = bc.VM5adc();
166  adc.VM6 = VM6adc = bc.VM6adc();
167  adc.VM7 = VM7adc = bc.VM7adc();
168  adc.VM8 = VM8adc = bc.VM8adc();
169  adc.TOF1 = TOF1adc = bc.TOF1Sadc();
170  adc.TOF2 = TOF2adc = bc.TOF2Sadc();
171  adc.BH1 = BH1adc = bc.BH1adc();
172  adc.BH2 = BH2adc = bc.BH2adc();
173  adc.BH3 = BH3adc = bc.BH3adc();
174  adc.BH4 = BH4adc = bc.BH4adc();
175 }
double VM4
Definition: ZdcTBAnalysis.h:70
double CK1adc() const
double VH2adc() const
double VH3adc() const
double S2
Definition: ZdcTBAnalysis.h:62
double VMBadc() const
double S2adc() const
double CK3adc() const
double CK1
Definition: ZdcTBAnalysis.h:57
double Sci528
Definition: ZdcTBAnalysis.h:56
double VM
Definition: ZdcTBAnalysis.h:47
double VH2
Definition: ZdcTBAnalysis.h:51
double Sci521adc() const
double VH3
Definition: ZdcTBAnalysis.h:52
double VM5
Definition: ZdcTBAnalysis.h:71
double VMF
Definition: ZdcTBAnalysis.h:65
double VM3adc() const
double VM7adc() const
double Sci528adc() const
double S1
Definition: ZdcTBAnalysis.h:61
double V3adc() const
double TOF2
Definition: ZdcTBAnalysis.h:76
double VM2
Definition: ZdcTBAnalysis.h:68
double VH1
Definition: ZdcTBAnalysis.h:50
double BH2adc() const
double VM8adc() const
double VM1
Definition: ZdcTBAnalysis.h:67
double VM8
Definition: ZdcTBAnalysis.h:74
double V6
Definition: ZdcTBAnalysis.h:49
double SciVLE
Definition: ZdcTBAnalysis.h:60
double VM3
Definition: ZdcTBAnalysis.h:69
double BH4adc() const
double S3adc() const
double Sci521
Definition: ZdcTBAnalysis.h:55
double VH4
Definition: ZdcTBAnalysis.h:53
double S4adc() const
double VMFadc() const
double TOF2Sadc() const
double V6adc() const
double S3
Definition: ZdcTBAnalysis.h:63
double CK3
Definition: ZdcTBAnalysis.h:59
double VM6
Definition: ZdcTBAnalysis.h:72
double VH1adc() const
double VMadc() const
Muon Veto adc.
double SciVLEadc() const
double VM5adc() const
double BH3adc() const
double Ecal7x7
Definition: ZdcTBAnalysis.h:54
double S1adc() const
double VMB
Definition: ZdcTBAnalysis.h:66
double CK2adc() const
double BH4
Definition: ZdcTBAnalysis.h:80
double BH1adc() const
double TOF1
Definition: ZdcTBAnalysis.h:75
double S4
Definition: ZdcTBAnalysis.h:64
double CK2
Definition: ZdcTBAnalysis.h:58
double VM1adc() const
double BH3
Definition: ZdcTBAnalysis.h:79
double Ecal7x7() const
double V3
Definition: ZdcTBAnalysis.h:48
double VM7
Definition: ZdcTBAnalysis.h:73
double VM4adc() const
double VH4adc() const
double VM2adc() const
double TOF1Sadc() const
double BH2
Definition: ZdcTBAnalysis.h:78
double VM6adc() const
double BH1
Definition: ZdcTBAnalysis.h:77

◆ analyze() [4/5]

void ZdcTBAnalysis::analyze ( const HcalTBTiming times)

Definition at line 89 of file ZdcTBAnalysis.cc.

References beam_coincidence, TDC::beamCoincidence, HcalTBTiming::BeamCoincidenceCount(), HcalTBTiming::BeamCoincidenceHits(), TDC::bh1, HcalTBTiming::BH1Count(), HcalTBTiming::BH1Hits(), bh1hits, TDC::bh2, HcalTBTiming::BH2Count(), HcalTBTiming::BH2Hits(), bh2hits, TDC::bh3, HcalTBTiming::BH3Count(), HcalTBTiming::BH3Hits(), bh3hits, TDC::bh4, HcalTBTiming::BH4Count(), HcalTBTiming::BH4Hits(), bh4hits, laser_flash, TDC::laserFlash, HcalTBTiming::laserFlash(), TDC::m1, HcalTBTiming::M1Count(), HcalTBTiming::M1Hits(), m1hits, TDC::m2, HcalTBTiming::M2Count(), HcalTBTiming::M2Hits(), m2hits, TDC::m3, HcalTBTiming::M3Count(), HcalTBTiming::M3Hits(), m3hits, qie_phase, TDC::qiePhase, HcalTBTiming::qiePhase(), TDC::s1, HcalTBTiming::S1Count(), HcalTBTiming::S1Hits(), s1hits, TDC::s2, HcalTBTiming::S2Count(), HcalTBTiming::S2Hits(), s2hits, TDC::s3, HcalTBTiming::S3Count(), HcalTBTiming::S3Hits(), s3hits, TDC::s4, HcalTBTiming::S4Count(), HcalTBTiming::S4Hits(), s4hits, tdc, TDC::TOF1, TOF1_time, HcalTBTiming::TOF1Stime(), TDC::TOF2, TOF2_time, HcalTBTiming::TOF2Stime(), TDC::trigger, trigger_time, HcalTBTiming::triggerTime(), ttc_L1a_time, TDC::ttcL1, and HcalTBTiming::ttcL1Atime().

89  {
90  //times
91  tdc.trigger = trigger_time = times.triggerTime();
92  tdc.ttcL1 = ttc_L1a_time = times.ttcL1Atime();
94  tdc.qiePhase = qie_phase = times.qiePhase();
95  tdc.TOF1 = TOF1_time = times.TOF1Stime();
96  tdc.TOF2 = TOF2_time = times.TOF2Stime();
97 
98  // just take 5 first hits of multihit tdc (5 tick cycles)
99  int indx = 0;
100  int indTop = 5;
101  for (indx = 0; indx < times.BeamCoincidenceCount(); indx++)
102  if (indx < indTop)
103  tdc.beamCoincidence[indx] = beam_coincidence[indx] = times.BeamCoincidenceHits(indx);
104  for (indx = 0; indx < times.M1Count(); indx++)
105  if (indx < indTop)
106  tdc.m1[indx] = m1hits[indx] = times.M1Hits(indx);
107  for (indx = 0; indx < times.M2Count(); indx++)
108  if (indx < indTop)
109  tdc.m2[indx] = m2hits[indx] = times.M2Hits(indx);
110  for (indx = 0; indx < times.M3Count(); indx++)
111  if (indx < indTop)
112  tdc.m3[indx] = m3hits[indx] = times.M3Hits(indx);
113  for (indx = 0; indx < times.S1Count(); indx++)
114  if (indx < indTop)
115  tdc.s1[indx] = s1hits[indx] = times.S1Hits(indx);
116  for (indx = 0; indx < times.S2Count(); indx++)
117  if (indx < indTop)
118  tdc.s2[indx] = s2hits[indx] = times.S2Hits(indx);
119  for (indx = 0; indx < times.S3Count(); indx++)
120  if (indx < indTop)
121  tdc.s3[indx] = s3hits[indx] = times.S3Hits(indx);
122  for (indx = 0; indx < times.S4Count(); indx++)
123  if (indx < indTop)
124  tdc.s4[indx] = s4hits[indx] = times.S4Hits(indx);
125  for (indx = 0; indx < times.BH1Count(); indx++)
126  if (indx < indTop)
127  tdc.bh1[indx] = bh1hits[indx] = times.BH1Hits(indx);
128  for (indx = 0; indx < times.BH2Count(); indx++)
129  if (indx < indTop)
130  tdc.bh2[indx] = bh2hits[indx] = times.BH2Hits(indx);
131  for (indx = 0; indx < times.BH3Count(); indx++)
132  if (indx < indTop)
133  tdc.bh3[indx] = bh3hits[indx] = times.BH3Hits(indx);
134  for (indx = 0; indx < times.BH4Count(); indx++)
135  if (indx < indTop)
136  tdc.bh4[indx] = bh4hits[indx] = times.BH4Hits(indx);
137 }
double M3Hits(int index) const
Returns the indexed hit time from muon veto scintillator M3.
Definition: HcalTBTiming.h:77
double s2hits[5]
double BH3Hits(int index) const
Returns the indexed hit time from beam halo counter BEAM RIGHT FROM PARTICLE&#39;S VIEW.
Definition: HcalTBTiming.h:93
double laser_flash
int BH1Count() const
Returns the number of hits from beam halo counter up horizontal.
Definition: HcalTBTiming.h:62
double S3Hits(int index) const
Returns the indexed hit time from scintillator S3, which is 2cm x 2cm.
Definition: HcalTBTiming.h:84
double M1Hits(int index) const
Returns the indexed hit time from muon veto scintillator M1.
Definition: HcalTBTiming.h:73
double trigger
Definition: ZdcTBAnalysis.h:26
double bh1hits[5]
int S3Count() const
Returns the number of hits from scintillator S3, which is 2cm x 2cm.
Definition: HcalTBTiming.h:57
double s4hits[5]
double BeamCoincidenceHits(int index) const
Returns the indexed hit time from Beam Coincidence.
Definition: HcalTBTiming.h:71
double bh4hits[5]
double triggerTime() const
Returns the trigger time in ns.
Definition: HcalTBTiming.h:24
double s1[5]
Definition: ZdcTBAnalysis.h:36
double beamCoincidence[5]
Definition: ZdcTBAnalysis.h:28
double trigger_time
double beam_coincidence[5]
double laserFlash() const
Returns the laser activation time in ns.
Definition: HcalTBTiming.h:30
double bh2hits[5]
double ttcL1Atime() const
Returns the Level 1 Accept time in ns.
Definition: HcalTBTiming.h:27
double BH4Hits(int index) const
Returns the indexed hit time from beam halo counter DOWN HORZINTAL.
Definition: HcalTBTiming.h:95
double m2[5]
Definition: ZdcTBAnalysis.h:34
double m1[5]
Definition: ZdcTBAnalysis.h:33
double s2[5]
Definition: ZdcTBAnalysis.h:37
int BeamCoincidenceCount() const
Returns the number of hits from Beam Coincidence.
Definition: HcalTBTiming.h:44
double m3hits[5]
double s4[5]
Definition: ZdcTBAnalysis.h:39
double s1hits[5]
int M2Count() const
Returns the number of hits from muon veto scintillator M2.
Definition: HcalTBTiming.h:48
double TOF2Stime() const
Returns the TOF2S time (zero otherwise)
Definition: HcalTBTiming.h:39
double bh3hits[5]
double ttc_L1a_time
double TOF1Stime() const
Returns the TOF1S time (zero otherwise)
Definition: HcalTBTiming.h:35
double BH2Hits(int index) const
Returns the indexed hit time from from beam halo counter BEAM LEFT FROM PARTICLE&#39;S VIEW...
Definition: HcalTBTiming.h:91
double bh4[5]
Definition: ZdcTBAnalysis.h:43
double bh1[5]
Definition: ZdcTBAnalysis.h:40
double BH1Hits(int index) const
Returns the indexed hit time from beam halo counter UP HORIZONTAL.
Definition: HcalTBTiming.h:89
double s3[5]
Definition: ZdcTBAnalysis.h:38
double S2Hits(int index) const
Returns the indexed hit time from scintillator S2, which is 4cm x 4cm.
Definition: HcalTBTiming.h:82
int S2Count() const
Returns the number of hits from scintillator S2, which is 4cm x 4cm.
Definition: HcalTBTiming.h:55
double m3[5]
Definition: ZdcTBAnalysis.h:35
double bh3[5]
Definition: ZdcTBAnalysis.h:42
double m1hits[5]
double bh2[5]
Definition: ZdcTBAnalysis.h:41
int M1Count() const
Returns the number of hits from muon veto scintillator M1.
Definition: HcalTBTiming.h:46
double TOF1
Definition: ZdcTBAnalysis.h:31
int M3Count() const
Returns the number of hits from muon veto scintillator M3.
Definition: HcalTBTiming.h:50
double s3hits[5]
int BH4Count() const
Returns the number of hits from beam halo counter down horizontal.
Definition: HcalTBTiming.h:68
int S1Count() const
Returns the number of hits from scintillator S1, which is 12cm x 12cm.
Definition: HcalTBTiming.h:53
double qiePhase() const
Returns the QIE phase for 2003 testbeam data (zero otherwise)
Definition: HcalTBTiming.h:32
double TOF2
Definition: ZdcTBAnalysis.h:32
int BH2Count() const
Returns the number of hits from beam halo counter left from particle view.
Definition: HcalTBTiming.h:64
int S4Count() const
Returns the number of hits from scintillator S4, which is 12cm x 12cm.
Definition: HcalTBTiming.h:59
double laserFlash
Definition: ZdcTBAnalysis.h:29
double M2Hits(int index) const
Returns the indexed hit time from muon veto scintillator M2.
Definition: HcalTBTiming.h:75
double ttcL1
Definition: ZdcTBAnalysis.h:27
int BH3Count() const
Returns the number of hits from beam halo counter right from particle view.
Definition: HcalTBTiming.h:66
double S1Hits(int index) const
Returns the indexed hit time from scintillator S1, which is 12cm x 12cm.
Definition: HcalTBTiming.h:80
double m2hits[5]
double S4Hits(int index) const
Returns the indexed hit time from scintillator S4, which is 12cm x 12cm.
Definition: HcalTBTiming.h:86
double qiePhase
Definition: ZdcTBAnalysis.h:30

◆ analyze() [5/5]

void ZdcTBAnalysis::analyze ( const HcalTBEventPosition chpos)

Definition at line 177 of file ZdcTBAnalysis.cc.

References chamb, HcalTBEventPosition::getChamberHits(), CHAMB::WCAx, wcax, CHAMB::WCAy, wcay, CHAMB::WCBx, wcbx, CHAMB::WCBy, wcby, CHAMB::WCCx, wccx, CHAMB::WCCy, wccy, CHAMB::WCDx, wcdx, CHAMB::WCDy, wcdy, CHAMB::WCEx, wcex, CHAMB::WCEy, wcey, CHAMB::WCFx, wcfx, CHAMB::WCFy, wcfy, CHAMB::WCGx, wcgx, CHAMB::WCGy, wcgy, CHAMB::WCHx, wchx, CHAMB::WCHy, and wchy.

177  {
178  //chambers position
179  chpos.getChamberHits('A', wcax, wcay);
180  chpos.getChamberHits('B', wcbx, wcby);
181  chpos.getChamberHits('C', wccx, wccy);
182  chpos.getChamberHits('D', wcdx, wcdy);
183  chpos.getChamberHits('E', wcex, wcey);
184  chpos.getChamberHits('F', wcfx, wcfy);
185  chpos.getChamberHits('G', wcgx, wcgy);
186  chpos.getChamberHits('H', wchx, wchy);
187 
188  // just take 5 first hits of chambers (5 tick cycles)
189  unsigned int indTop = 5;
190  unsigned int indx = 0;
191  for (indx = 0; indx < wcax.size(); indx++)
192  if (indx < indTop)
193  chamb.WCAx[indx] = wcax[indx];
194  for (indx = 0; indx < wcay.size(); indx++)
195  if (indx < indTop)
196  chamb.WCAy[indx] = wcay[indx];
197  for (indx = 0; indx < wcbx.size(); indx++)
198  if (indx < indTop)
199  chamb.WCBx[indx] = wcbx[indx];
200  for (indx = 0; indx < wcby.size(); indx++)
201  if (indx < indTop)
202  chamb.WCBy[indx] = wcby[indx];
203  for (indx = 0; indx < wccx.size(); indx++)
204  if (indx < indTop)
205  chamb.WCCx[indx] = wccx[indx];
206  for (indx = 0; indx < wccy.size(); indx++)
207  if (indx < indTop)
208  chamb.WCCy[indx] = wccy[indx];
209  for (indx = 0; indx < wcdx.size(); indx++)
210  if (indx < indTop)
211  chamb.WCDx[indx] = wcdx[indx];
212  for (indx = 0; indx < wcdy.size(); indx++)
213  if (indx < indTop)
214  chamb.WCDy[indx] = wcdy[indx];
215  for (indx = 0; indx < wcdx.size(); indx++)
216  if (indx < indTop)
217  chamb.WCEx[indx] = wcex[indx];
218  for (indx = 0; indx < wcey.size(); indx++)
219  if (indx < indTop)
220  chamb.WCEy[indx] = wcey[indx];
221  for (indx = 0; indx < wcfx.size(); indx++)
222  if (indx < indTop)
223  chamb.WCFx[indx] = wcfx[indx];
224  for (indx = 0; indx < wcfy.size(); indx++)
225  if (indx < indTop)
226  chamb.WCFy[indx] = wcfy[indx];
227  for (indx = 0; indx < wcgx.size(); indx++)
228  if (indx < indTop)
229  chamb.WCGx[indx] = wcgx[indx];
230  for (indx = 0; indx < wcgy.size(); indx++)
231  if (indx < indTop)
232  chamb.WCGy[indx] = wcgy[indx];
233  for (indx = 0; indx < wchx.size(); indx++)
234  if (indx < indTop)
235  chamb.WCHx[indx] = wchx[indx];
236  for (indx = 0; indx < wchy.size(); indx++)
237  if (indx < indTop)
238  chamb.WCHy[indx] = wchy[indx];
239 }
double WCFx[5]
Definition: ZdcTBAnalysis.h:94
double WCGy[5]
Definition: ZdcTBAnalysis.h:97
std::vector< double > wcfx
std::vector< double > wcgx
std::vector< double > wcex
std::vector< double > wchx
double WCBy[5]
Definition: ZdcTBAnalysis.h:87
std::vector< double > wcax
double WCEy[5]
Definition: ZdcTBAnalysis.h:93
std::vector< double > wccx
std::vector< double > wchy
double WCCx[5]
Definition: ZdcTBAnalysis.h:88
std::vector< double > wcbx
double WCHx[5]
Definition: ZdcTBAnalysis.h:98
std::vector< double > wccy
double WCBx[5]
Definition: ZdcTBAnalysis.h:86
std::vector< double > wcgy
std::vector< double > wcey
double WCAx[5]
Definition: ZdcTBAnalysis.h:84
double WCEx[5]
Definition: ZdcTBAnalysis.h:92
std::vector< double > wcdy
double WCDx[5]
Definition: ZdcTBAnalysis.h:90
std::vector< double > wcfy
double WCFy[5]
Definition: ZdcTBAnalysis.h:95
double WCAy[5]
Definition: ZdcTBAnalysis.h:85
std::vector< double > wcdx
double WCGx[5]
Definition: ZdcTBAnalysis.h:96
std::vector< double > wcby
double WCHy[5]
Definition: ZdcTBAnalysis.h:99
void getChamberHits(char chamberch, std::vector< double > &xvec, std::vector< double > &yvec) const
Get the wire chamber hits for the specified chamber For HB/HE/HO running, chambers A...
double WCDy[5]
Definition: ZdcTBAnalysis.h:91
double WCCy[5]
Definition: ZdcTBAnalysis.h:89
std::vector< double > wcay

◆ done()

void ZdcTBAnalysis::done ( )

Definition at line 304 of file ZdcTBAnalysis.cc.

References outFile, and ZdcAnalize.

Referenced by ZdcTBAnalyzer::endJob().

304  {
305  ZdcAnalize->Print();
306  outFile->cd();
307  ZdcAnalize->Write();
308  outFile->Close();
309 }
TTree * ZdcAnalize

◆ fillTree()

void ZdcTBAnalysis::fillTree ( )

Definition at line 302 of file ZdcTBAnalysis.cc.

References ZdcAnalize.

Referenced by ZdcTBAnalyzer::analyze().

302 { ZdcAnalize->Fill(); }
TTree * ZdcAnalize

◆ setup()

void ZdcTBAnalysis::setup ( const std::string &  histoFileName)

Definition at line 8 of file ZdcTBAnalysis.cc.

References adc, chamb, outFile, DeadROC_duringRun::outFileName, dataset::outName, tdc, ZdcAnalize, zdcn, and zdcp.

Referenced by ZdcTBAnalyzer::ZdcTBAnalyzer().

8  {
9  TString outName = outFileName;
10  outFile = new TFile(outName, "RECREATE");
11  ZdcAnalize = new TTree("ZdcAnaTree", "ZdcAnaTree");
12  ZdcAnalize->Branch("Trigger",
13  0,
14  "run/I:event/I:beamTrigger/I:fakeTrigger/I:"
15  "pedestalTrigger/I:outSpillPedestalTrigger/I:inSpillPedestalTrigger/I:"
16  "laserTrigger/I:laserTrigger/I:ledTrigger/I:spillTrigger/I");
17  ZdcAnalize->Branch("TDC",
18  0,
19  "trigger/D:ttcL1/D:beamCoincidence/D:laserFlash/D:qiePhase/D:"
20  "TTOF1/D:TTOF2/D:m1[5]/D:m2[5]/D:m3[5]/D:"
21  "s1[5]/D:s2[5]/D:s3[5]/D:s4[5]/D:"
22  "bh1[5]/D:bh2[5]/D:bh3[5]/D:bh4[5]/D");
23  ZdcAnalize->Branch("ADC",
24  0,
25  "VM/D:V3/D:V6/D:VH1/D:VH2/D:VH3/D:VH4/D:Ecal7x7/D:"
26  "Sci521/D:Sci528/D:CK1/D:CK2/D:CK3/D:SciVLE/D:S1/D:S2/D:S3/D:S4/D:"
27  "VMF/D:VMB/D:VM1/D:VM2/D:VM3/D:VM4/D:VM5/D:VM6/D:VM7/D:VM8/D:"
28  "TOF1/D:TOF2/D:BH1/D:BH2/D:BH3/BH4/D");
29  ZdcAnalize->Branch("Chamb",
30  0,
31  "WCAx[5]/D:WCAy[5]/D:WCBx[5]/D:WCBy[5]/D:"
32  "WCCx[5]/D:WCCy[5]/D:WCDx[5]/D:WCDy[5]/D:WCEx[5]/D:WCEy[5]/D:"
33  "WCFx[5]/D:WCFy[5]/D:WCGx[5]/D:WCGy[5]/D:WCHx[5]/D:WCHy[5]/D");
34  ZdcAnalize->Branch("ZDCP",
35  0,
36  "zdcHAD1/D:zdcHAD2/D:zdcHAD3/D:zdcHAD4/D:"
37  "zdcEM1/D:zdcEM2/D:zdcEM3/D:zdcEM4/D:zdcEM5/D:"
38  "zdcScint1/D:zdcScint2/D:"
39  "zdcExtras[7]/D");
40  ZdcAnalize->Branch("ZDCN",
41  0,
42  "zdcHAD1/D:zdcHAD2/D:zdcHAD3/D:zdcHAD4/D:"
43  "zdcEM1/D:zdcEM2/D:zdcEM3/D:zdcEM4/D:zdcEM5/D:"
44  "zdcScint1/D:zdcScint2/D:"
45  "zdcExtras[7]/D");
46  ZdcAnalize->GetBranch("Trigger")->SetAddress(&trigger);
47  ZdcAnalize->GetBranch("TDC")->SetAddress(&tdc);
48  ZdcAnalize->GetBranch("ADC")->SetAddress(&adc);
49  ZdcAnalize->GetBranch("Chamb")->SetAddress(&chamb);
50  ZdcAnalize->GetBranch("ZDCP")->SetAddress(&zdcp);
51  ZdcAnalize->GetBranch("ZDCN")->SetAddress(&zdcn);
52  ZdcAnalize->SetAutoSave();
53 }
TTree * ZdcAnalize

Member Data Documentation

◆ adc

ADC ZdcTBAnalysis::adc
private

Definition at line 239 of file ZdcTBAnalysis.h.

Referenced by analyze(), and setup().

◆ beam_coincidence

double ZdcTBAnalysis::beam_coincidence[5]
private

Definition at line 165 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ BH1adc

double ZdcTBAnalysis::BH1adc
private

Definition at line 215 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ bh1hits

double ZdcTBAnalysis::bh1hits[5]
private

Definition at line 178 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ BH2adc

double ZdcTBAnalysis::BH2adc
private

Definition at line 216 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ bh2hits

double ZdcTBAnalysis::bh2hits[5]
private

Definition at line 179 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ BH3adc

double ZdcTBAnalysis::BH3adc
private

Definition at line 217 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ bh3hits

double ZdcTBAnalysis::bh3hits[5]
private

Definition at line 180 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ BH4adc

double ZdcTBAnalysis::BH4adc
private

Definition at line 218 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ bh4hits

double ZdcTBAnalysis::bh4hits[5]
private

Definition at line 181 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ chamb

CHAMB ZdcTBAnalysis::chamb
private

Definition at line 240 of file ZdcTBAnalysis.h.

Referenced by analyze(), and setup().

◆ CK1adc

double ZdcTBAnalysis::CK1adc
private

Definition at line 194 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ CK2adc

double ZdcTBAnalysis::CK2adc
private

Definition at line 195 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ CK3adc

double ZdcTBAnalysis::CK3adc
private

Definition at line 196 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ detID

HcalZDCDetId ZdcTBAnalysis::detID
private

Definition at line 150 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ Ecal7x7adc

double ZdcTBAnalysis::Ecal7x7adc
private

Definition at line 190 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ energy

double ZdcTBAnalysis::energy
private

Definition at line 149 of file ZdcTBAnalysis.h.

Referenced by analyze(), and Jet.Jet::rawEnergy().

◆ eventNumber

int ZdcTBAnalysis::eventNumber
private

Definition at line 153 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ ichannel

int ZdcTBAnalysis::ichannel
private

Definition at line 147 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ idepth

int ZdcTBAnalysis::idepth
private

Definition at line 148 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isBeamTrigger

bool ZdcTBAnalysis::isBeamTrigger
private

Definition at line 154 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isCalibTrigger

bool ZdcTBAnalysis::isCalibTrigger
private

Definition at line 156 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isection

int ZdcTBAnalysis::isection
private

Definition at line 146 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isFakeTrigger

bool ZdcTBAnalysis::isFakeTrigger
private

Definition at line 155 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ iside

int ZdcTBAnalysis::iside
private

Definition at line 145 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isInSpillPedestalTrigger

bool ZdcTBAnalysis::isInSpillPedestalTrigger
private

Definition at line 158 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isLaserTrigger

bool ZdcTBAnalysis::isLaserTrigger
private

Definition at line 159 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isLedTrigger

bool ZdcTBAnalysis::isLedTrigger
private

Definition at line 160 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isOutSpillPedestalTrigger

bool ZdcTBAnalysis::isOutSpillPedestalTrigger
private

Definition at line 157 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ isSpillTrigger

bool ZdcTBAnalysis::isSpillTrigger
private

Definition at line 161 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ laser_flash

double ZdcTBAnalysis::laser_flash
private

Definition at line 166 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ m1hits

double ZdcTBAnalysis::m1hits[5]
private

Definition at line 171 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ m2hits

double ZdcTBAnalysis::m2hits[5]
private

Definition at line 172 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ m3hits

double ZdcTBAnalysis::m3hits[5]
private

Definition at line 173 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ outFile

TFile* ZdcTBAnalysis::outFile
private

◆ qie_phase

double ZdcTBAnalysis::qie_phase
private

Definition at line 167 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ runNumber

int ZdcTBAnalysis::runNumber
private

Definition at line 152 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ S1adc

double ZdcTBAnalysis::S1adc
private

Definition at line 198 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ s1hits

double ZdcTBAnalysis::s1hits[5]
private

Definition at line 174 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ S2adc

double ZdcTBAnalysis::S2adc
private

Definition at line 199 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ s2hits

double ZdcTBAnalysis::s2hits[5]
private

Definition at line 175 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ S3adc

double ZdcTBAnalysis::S3adc
private

Definition at line 200 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ s3hits

double ZdcTBAnalysis::s3hits[5]
private

Definition at line 176 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ S4adc

double ZdcTBAnalysis::S4adc
private

Definition at line 201 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ s4hits

double ZdcTBAnalysis::s4hits[5]
private

Definition at line 177 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ Sci521adc

double ZdcTBAnalysis::Sci521adc
private

Definition at line 191 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ Sci528adc

double ZdcTBAnalysis::Sci528adc
private

Definition at line 192 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ SciVLEadc

double ZdcTBAnalysis::SciVLEadc
private

Definition at line 197 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ tdc

TDC ZdcTBAnalysis::tdc
private

Definition at line 238 of file ZdcTBAnalysis.h.

Referenced by analyze(), and setup().

◆ TOF1_time

double ZdcTBAnalysis::TOF1_time
private

Definition at line 168 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ TOF1adc

double ZdcTBAnalysis::TOF1adc
private

Definition at line 213 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ TOF2_time

double ZdcTBAnalysis::TOF2_time
private

Definition at line 169 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ TOF2adc

double ZdcTBAnalysis::TOF2adc
private

Definition at line 214 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ trigger

TRIGGER ZdcTBAnalysis::trigger
private

Definition at line 237 of file ZdcTBAnalysis.h.

◆ trigger_time

double ZdcTBAnalysis::trigger_time
private

Definition at line 163 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ ttc_L1a_time

double ZdcTBAnalysis::ttc_L1a_time
private

Definition at line 164 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ V3adc

double ZdcTBAnalysis::V3adc
private

Definition at line 184 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ V6adc

double ZdcTBAnalysis::V6adc
private

Definition at line 185 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VH1adc

double ZdcTBAnalysis::VH1adc
private

Definition at line 186 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VH2adc

double ZdcTBAnalysis::VH2adc
private

Definition at line 187 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VH3adc

double ZdcTBAnalysis::VH3adc
private

Definition at line 188 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VH4adc

double ZdcTBAnalysis::VH4adc
private

Definition at line 189 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VM1adc

double ZdcTBAnalysis::VM1adc
private

Definition at line 205 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VM2adc

double ZdcTBAnalysis::VM2adc
private

Definition at line 206 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VM3adc

double ZdcTBAnalysis::VM3adc
private

Definition at line 207 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VM4adc

double ZdcTBAnalysis::VM4adc
private

Definition at line 208 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VM5adc

double ZdcTBAnalysis::VM5adc
private

Definition at line 209 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VM6adc

double ZdcTBAnalysis::VM6adc
private

Definition at line 210 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VM7adc

double ZdcTBAnalysis::VM7adc
private

Definition at line 211 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VM8adc

double ZdcTBAnalysis::VM8adc
private

Definition at line 212 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VMadc

double ZdcTBAnalysis::VMadc
private

Definition at line 183 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VMBadc

double ZdcTBAnalysis::VMBadc
private

Definition at line 204 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ VMFadc

double ZdcTBAnalysis::VMFadc
private

Definition at line 203 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcax

std::vector<double> ZdcTBAnalysis::wcax
private

Definition at line 220 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcay

std::vector<double> ZdcTBAnalysis::wcay
private

Definition at line 221 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcbx

std::vector<double> ZdcTBAnalysis::wcbx
private

Definition at line 222 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcby

std::vector<double> ZdcTBAnalysis::wcby
private

Definition at line 223 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wccx

std::vector<double> ZdcTBAnalysis::wccx
private

Definition at line 224 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wccy

std::vector<double> ZdcTBAnalysis::wccy
private

Definition at line 225 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcdx

std::vector<double> ZdcTBAnalysis::wcdx
private

Definition at line 226 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcdy

std::vector<double> ZdcTBAnalysis::wcdy
private

Definition at line 227 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcex

std::vector<double> ZdcTBAnalysis::wcex
private

Definition at line 228 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcey

std::vector<double> ZdcTBAnalysis::wcey
private

Definition at line 229 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcfx

std::vector<double> ZdcTBAnalysis::wcfx
private

Definition at line 230 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcfy

std::vector<double> ZdcTBAnalysis::wcfy
private

Definition at line 231 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcgx

std::vector<double> ZdcTBAnalysis::wcgx
private

Definition at line 232 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wcgy

std::vector<double> ZdcTBAnalysis::wcgy
private

Definition at line 233 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wchx

std::vector<double> ZdcTBAnalysis::wchx
private

Definition at line 234 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ wchy

std::vector<double> ZdcTBAnalysis::wchy
private

Definition at line 235 of file ZdcTBAnalysis.h.

Referenced by analyze().

◆ ZdcAnalize

TTree* ZdcTBAnalysis::ZdcAnalize
private

Definition at line 245 of file ZdcTBAnalysis.h.

Referenced by done(), fillTree(), and setup().

◆ zdcn

ZDCN ZdcTBAnalysis::zdcn
private

Definition at line 242 of file ZdcTBAnalysis.h.

Referenced by analyze(), and setup().

◆ zdcp

ZDCP ZdcTBAnalysis::zdcp
private

Definition at line 241 of file ZdcTBAnalysis.h.

Referenced by analyze(), and setup().