CMS 3D CMS Logo

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

#include <EventFilter/RctRawToDigi/src/RctRawToDigi.cc>

Inheritance diagram for RctRawToDigi:
edm::stream::EDProducer<>

Public Member Functions

 RctRawToDigi (const edm::ParameterSet &)
 
 ~RctRawToDigi () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Private Member Functions

void checkHeaders ()
 
bool decodeLinkID (const uint32_t inputValue, uint32_t &crateNumber, uint32_t &linkNumber, bool &even)
 
virtual void endJob ()
 method called at job end - use to print summary report More...
 
bool printAll (const unsigned char *data, const unsigned size)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
void unpack (const FEDRawData &d, edm::Event &e, RctUnpackCollections *const colls)
 Unpacks the raw data. More...
 
void unpackCTP7 (const uint32_t *data, const unsigned block_id, const unsigned size, RctUnpackCollections *const colls)
 

Private Attributes

std::vector< unsigned > errorCounters_
 Counts number of errors for each code (index) More...
 
L1TriggerErrorCollectionerrors_
 pointer to error collection More...
 
int fedId_
 RCT FED ID. More...
 
edm::InputTag inputLabel_
 FED collection label. More...
 
unsigned unpackFailures_
 To count the total number of RCT unpack failures. More...
 
const bool verbose_
 If true, then debug print out for each event. More...
 

Static Private Attributes

static const unsigned amc13HeaderSize_ = 16
 
static const unsigned amc13TrailerSize_ = 8
 
static const unsigned ctp7HeaderSize_ = 8
 
static const unsigned ctp7TrailerSize_ = 8
 
static const unsigned MAX_DATA = 4680
 The maximum number of blocks we will try to unpack before thinking something is wrong. More...
 
static const unsigned MAX_ERR_CODE = 6
 
static const unsigned MIN_DATA = 900
 The minimum number of blocks we will try to unpack before thinking something is wrong (really this should be 920, to be tested) More...
 
static const unsigned sLinkHeaderSize_ = 8
 
static const unsigned sLinkTrailerSize_ = 8
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Description: Produce RCT digis from raw data

Implementation: <Notes on="" implementation>="">

Definition at line 40 of file RctRawToDigi.h.

Constructor & Destructor Documentation

RctRawToDigi::RctRawToDigi ( const edm::ParameterSet iConfig)
explicit

Register Products

Definition at line 36 of file RctRawToDigi.cc.

References fedId_, inputLabel_, and LogDebug.

36  :
37  inputLabel_(iConfig.getParameter<edm::InputTag>("inputLabel")),
38  fedId_(iConfig.getUntrackedParameter<int>("rctFedId", FEDNumbering::MINRCTFEDID)),
39  verbose_(iConfig.getUntrackedParameter<bool>("verbose",false))
40 {
41  LogDebug("RCT") << "RctRawToDigi will unpack FED Id " << fedId_;
43  // RCT collections
44  produces<L1CaloEmCollection>();
45  produces<L1CaloRegionCollection>();
46 
47  // Error collection
48  consumes<FEDRawDataCollection>(inputLabel_);
49 
50 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag inputLabel_
FED collection label.
Definition: RctRawToDigi.h:98
int fedId_
RCT FED ID.
Definition: RctRawToDigi.h:99
const bool verbose_
If true, then debug print out for each event.
Definition: RctRawToDigi.h:101
RctRawToDigi::~RctRawToDigi ( )
override

Definition at line 53 of file RctRawToDigi.cc.

54 {
55  // do anything here that needs to be done at destruction time
56  // (e.g. close files, deallocate resources etc.)
57  //delete formatTranslator_;
58 }

Member Function Documentation

void RctRawToDigi::checkHeaders ( )
private

Looks at the firmware version header in the S-Link packet and instantiates relevant format translator. check block headers for consistency

bool RctRawToDigi::decodeLinkID ( const uint32_t  inputValue,
uint32_t &  crateNumber,
uint32_t &  linkNumber,
bool &  even 
)
private
void RctRawToDigi::endJob ( void  )
privatevirtual

method called at job end - use to print summary report

Definition at line 348 of file RctRawToDigi.cc.

References DEFINE_FWK_MODULE, errorCounters_, mps_fire::i, MAX_ERR_CODE, pileupDistInMC::total, and verbose_.

349 {
350  unsigned total=0;
351  std::ostringstream os;
352 
353  for (unsigned i=0 ; i <= MAX_ERR_CODE ; ++i) {
354  total+=errorCounters_.at(i);
355  os << "Error " << i << " (" << errorCounters_.at(i) << ")";
356  if(i < MAX_ERR_CODE) { os << ", "; }
357  }
358 
359  if (total>0 && verbose_) {
360  edm::LogError("RCT") << "Encountered " << total << " unpacking errors: " << os.str();
361  }
362 }
std::vector< unsigned > errorCounters_
Counts number of errors for each code (index)
Definition: RctRawToDigi.h:109
static const unsigned MAX_ERR_CODE
Definition: RctRawToDigi.h:107
const bool verbose_
If true, then debug print out for each event.
Definition: RctRawToDigi.h:101
bool RctRawToDigi::printAll ( const unsigned char *  data,
const unsigned  size 
)
private

Definition at line 336 of file RctRawToDigi.cc.

References gather_cfg::cout, mps_fire::i, and findQualityFiles::size.

337 {
338  for(unsigned i = 0; i < size; i ++){
339  std::cout << data[i] << " ";
340  if(i%6==5)
341  std::cout<<std::endl;
342  }
343  return true;
344 }
size
Write out results.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
void RctRawToDigi::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 61 of file RctRawToDigi.cc.

References amc13HeaderSize_, amc13TrailerSize_, TauDecayModes::dec, FEDRawDataCollection::FEDData(), fedId_, edm::Event::get(), edm::Event::getByLabel(), inputLabel_, edm::HandleBase::isValid(), LogDebug, MIN_DATA, FEDRawData::size(), sLinkHeaderSize_, sLinkTrailerSize_, and unpack().

Referenced by JSONExport.JsonExport::export(), HTMLExport.HTMLExport::export(), and HTMLExport.HTMLExportStatic::export().

62 {
63 
64  using namespace edm;
65 
66  // Instantiate all the collections the unpacker needs; puts them in event when this object goes out of scope.
67  std::unique_ptr<RctUnpackCollections> colls(new RctUnpackCollections(iEvent));
68 
69  // get raw data collection
71  iEvent.getByLabel(inputLabel_, feds);
72 
73  // if raw data collection is present, check the headers and do the unpacking
74  if (feds.isValid()) {
75 
76  const FEDRawData& rctRcd = feds->FEDData(fedId_);
77 
78  //Check FED size
79  LogDebug("RCT") << "Upacking FEDRawData of size " << std::dec << rctRcd.size();
80 
81  //check header size
83  if (rctRcd.size() > 0) {
84  LogError("L1T") << "Cannot unpack: empty/invalid L1T raw data (size = "
85  << rctRcd.size() << ") for ID " << fedId_ << ". Returning empty collections!";
86  } else {
87  // there's no MC packer for this payload, so totally expected that sometimes it will be absent... no warning issued.
88  }
89  //continue;
90  return;
91  }
92 
93  // do the unpacking
94  unpack(rctRcd, iEvent, colls.get());
95 
96  }
97  else{
98 
99  LogError("L1T") << "Cannot unpack: no collection found";
100 
101  return;
102  }
103 
104 }
#define LogDebug(id)
static const unsigned sLinkTrailerSize_
Definition: RctRawToDigi.h:77
static const unsigned sLinkHeaderSize_
Definition: RctRawToDigi.h:74
static const unsigned amc13TrailerSize_
Definition: RctRawToDigi.h:83
static const unsigned MIN_DATA
The minimum number of blocks we will try to unpack before thinking something is wrong (really this sh...
Definition: RctRawToDigi.h:95
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:326
bool isValid() const
Definition: HandleBase.h:74
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:480
edm::InputTag inputLabel_
FED collection label.
Definition: RctRawToDigi.h:98
HLT enums.
void unpack(const FEDRawData &d, edm::Event &e, RctUnpackCollections *const colls)
Unpacks the raw data.
int fedId_
RCT FED ID.
Definition: RctRawToDigi.h:99
static const unsigned amc13HeaderSize_
Definition: RctRawToDigi.h:80
void RctRawToDigi::unpack ( const FEDRawData d,
edm::Event e,
RctUnpackCollections *const  colls 
)
private

Unpacks the raw data.

Parameters
invalidDataFlag- if true, then won't attempt unpack but just output empty collecions.

Definition at line 106 of file RctRawToDigi.cc.

References FEDHeader::bxID(), FEDHeader::check(), FEDTrailer::check(), FEDTrailer::crc(), FEDRawData::data(), data, FEDTrailer::evtStatus(), FEDTrailer::fragmentLength(), RecoTauValidation_cfi::header, LogDebug, FEDHeader::lvl1ID(), FEDRawData::size(), sLinkTrailerSize_, FEDHeader::sourceID(), FEDHeader::triggerType(), FEDTrailer::ttsBits(), unpackCTP7(), and FEDHeader::version().

Referenced by produce().

107 {
108 
109  const unsigned char * data = d.data(); // The 8-bit wide raw-data array.
110  // Data offset - starts at 16 as there is a 64-bit S-Link header followed
111  // by a 64-bit software-controlled header (for pipeline format version
112  // info that is not yet used).
113 
114 
115  FEDHeader header(data);
116 
117  if (header.check()) {
118  LogDebug("L1T") << "Found SLink header:"
119  << " Trigger type " << header.triggerType()
120  << " L1 event ID " << header.lvl1ID()
121  << " BX Number " << header.bxID()
122  << " FED source " << header.sourceID()
123  << " FED version " << header.version();
124  } else {
125  LogWarning("L1T") << "Did not find a valid SLink header!";
126  }
127 
128  FEDTrailer trailer(data + (d.size() - sLinkTrailerSize_));
129 
130  if (trailer.check()) {
131  LogDebug("L1T") << "Found SLink trailer:"
132  << " Length " << trailer.fragmentLength()
133  << " CRC " << trailer.crc()
134  << " Status " << trailer.evtStatus()
135  << " Throttling bits " << trailer.ttsBits();
136  } else {
137  LogWarning("L1T") << "Did not find a SLink trailer!";
138  }
139 
140  unpackCTP7(reinterpret_cast<const uint32_t*>(data), 0, sizeof(data), colls);
141 
142 }
#define LogDebug(id)
static const unsigned sLinkTrailerSize_
Definition: RctRawToDigi.h:77
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
void unpackCTP7(const uint32_t *data, const unsigned block_id, const unsigned size, RctUnpackCollections *const colls)
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28
void RctRawToDigi::unpackCTP7 ( const uint32_t *  data,
const unsigned  block_id,
const unsigned  size,
RctUnpackCollections *const  colls 
)
private

Definition at line 146 of file RctRawToDigi.cc.

References gen::k, LogDebug, funct::m, connectstrParser::o, haddnano::of, PFRecoTauDiscriminationByIsolation_cfi::offset, lumiQueryAPI::q, RctUnpackCollections::rctCalo(), RctUnpackCollections::rctEm(), L1CaloEmCand::setBx(), L1CaloRegion::setBx(), findQualityFiles::size, protons_cff::t, tmp, parallelization::uint(), and globals_cff::x1.

Referenced by unpack().

147 {
148  //offset from 6 link header words
149  uint32_t of = 6;
150  LogDebug("L1T") << "Block ID = " << block_id << " size = " << size;
151 
152  CTP7Format ctp7Format;
153  RctDataDecoder rctDataDecoder;
154  uint32_t nBXTemp = 0;
155  uint32_t ctp7FWVersion;
156  uint32_t L1ID, L1aBCID;
157  std::vector<RCTInfo> allCrateRCTInfo[5];
158 
159  L1ID = data[1+of]; // extract the L1 ID number
160  L1aBCID = data[5+of] & 0x00000FFF; // extract the BCID number of L1A
161  nBXTemp = (data[5+of] & 0x00FF0000) >> 16; // extract number of BXs readout per L1A
162  ctp7FWVersion = data[4+of];
163 
164  if(nBXTemp != 1 && nBXTemp != 3 && nBXTemp != 5)
165  nBXTemp = 1;
166 
167  const uint32_t nBX = nBXTemp;
168 
169  LogDebug("L1T") << "CTP7 L1ID = " << L1ID << " L1A BCID = " << L1aBCID << " BXs in capture = " << nBX << " CTP7 DAQ FW = " << ctp7FWVersion;
170 
171  struct link_data{
172  bool even;
173  unsigned int crateID;
174  unsigned int ctp7LinkNumber;
175  std::vector <unsigned int> uint;
176  };
177 
178  //nBX max 5, nLinks max 36 [nBX][nLinks]
179  link_data allLinks[5][36];
180  const uint32_t NLinks = ctp7Format.NLINKS;
181  assert(NLinks <= 36);
182 
183  //change this implementation
184  uint32_t iDAQBuffer = 0;
185 
186  //Step 1: Grab all data from ctp7 buffer and put into link format
187  for(unsigned int iLink = 0; iLink < NLinks; iLink++ ){
188  iDAQBuffer = of +
189  ctp7Format.EVENT_HEADER_WORDS + iLink * (ctp7Format.CHANNEL_HEADER_WORDS +
190  nBX * ctp7Format.CHANNEL_DATA_WORDS_PER_BX);
191 
192  //first decode linkID
193  uint32_t linkID = data[iDAQBuffer++];
194  uint32_t tmp = data[iDAQBuffer++];
195  uint32_t CRCErrCnt = tmp & 0x0000FFFF;
196  //uint32_t linkStatus = (tmp & 0xFFFF0000) >> 16;
197 
198  uint32_t crateID = 0; uint32_t expectedCrateID = 0;
199  bool even = false; bool expectedEven = false;
200 
201  //getExpected Link ID
202  rctDataDecoder.getExpectedLinkID(iLink, expectedCrateID, expectedEven);
203  //getDecodedLink ID
204  rctDataDecoder.decodeLinkID(linkID, crateID, even);
205 
206  //Perform a check to see if the link ID is as expected, if not then report an error but continue unpacking
207  if(expectedCrateID!=crateID || even!=expectedEven ){
208  LogError("L1T") <<"Expected Crate ID "<< expectedCrateID <<" expectedEven "<< expectedEven
209  <<"does not match actual Crate ID "<<crateID<<" even "<<even;
210  }
211 
212  if(CRCErrCnt!=0)
213  LogError("L1T")<<"WARNING CRC ErrorFound linkID "<< linkID<<" expected crateID "<< expectedCrateID;
214 
215  // Loop over multiple BX
216  for (uint32_t iBX=0; iBX<nBX; iBX++){
217  allLinks[iBX][iLink].uint.reserve(6);
218  allLinks[iBX][iLink].ctp7LinkNumber = iLink;
219  allLinks[iBX][iLink].crateID = expectedCrateID;
220  allLinks[iBX][iLink].even = expectedEven;
221 
222  //Notice 6 words per BX
223  for(unsigned int iWord = 0; iWord < 6 ; iWord++ ){
224  allLinks[iBX][iLink].uint.push_back(data[iDAQBuffer+iWord+iBX*6]);
225  }
226  }
227  }
228 
229  //Step 2: Dynamically match links and create RCTInfo Objects
230  uint32_t nCratesFound = 0;
231  for(unsigned int iCrate = 0; iCrate < 18 ; iCrate++){
232 
233  bool foundEven = false, foundOdd = false;
234  link_data even[5];
235  link_data odd[5];
236 
237  for(unsigned int iLink = 0; iLink < NLinks; iLink++){
238 
239  if( (allLinks[0][iLink].crateID==iCrate) && (allLinks[0][iLink].even == true) ){
240  foundEven = true;
241  for (unsigned int iBX=0; iBX<nBX; iBX++)
242  even[iBX] = allLinks[iBX][iLink];
243  }
244  else if( (allLinks[0][iLink].crateID==iCrate) && (allLinks[0][iLink].even == false) ){
245  foundOdd = true;
246  for (unsigned int iBX=0; iBX<nBX; iBX++)
247  odd[iBX] = allLinks[iBX][iLink];
248  }
249 
250  //if success then create RCTInfoVector and fill output object
251  if(foundEven && foundOdd){
252  nCratesFound++;
253 
254  //fill rctInfoVector for all BX read out
255  for (unsigned int iBX=0; iBX<nBX; iBX++){
256  //RCTInfoFactory rctInfoFactory;
257  std::vector <RCTInfo> rctInfoData;
258  rctDataDecoder.decodeLinks(even[iBX].uint, odd[iBX].uint, rctInfoData);
259  rctDataDecoder.setRCTInfoCrateID(rctInfoData, iCrate);
260  allCrateRCTInfo[iBX].push_back(rctInfoData.at(0));
261  }
262  break;
263  }
264  }
265  }
266 
267  if(nCratesFound != 18)
268  LogError("L1T") << "Warning -- only found "<< nCratesFound << " valid crates";
269 
270  //start assuming 1 BX readout
271  int32_t startBX = 0;
272  if(nBX == 1)
273  startBX = 0;
274  else if(nBX == 3)
275  startBX = -1;
276  else if(nBX == 5)
277  startBX = -2;
278 
279  //Step 3: Create Collections from RCTInfo Objects
280  //Notice, iBX used for grabbing data from array, startBX used for storing in Collection
281  for (uint32_t iBX=0; iBX<nBX; iBX++, startBX++){
282 
283  for(unsigned int iCrate = 0; iCrate < nCratesFound; iCrate++ ){
284 
285  RCTInfo rctInfo = allCrateRCTInfo[iBX].at(iCrate);
286  //Use Crate ID to identify eta/phi of candidate
287  for(int j = 0; j < 4; j++) {
288  L1CaloEmCand em = L1CaloEmCand(rctInfo.neRank[j],
289  rctInfo.neRegn[j],
290  rctInfo.neCard[j],
291  rctInfo.crateID,
292  false);
293  em.setBx(startBX);
294  colls->rctEm()->push_back(em);
295  }
296 
297  for(int j = 0; j < 4; j++) {
298  L1CaloEmCand em = L1CaloEmCand(rctInfo.ieRank[j],
299  rctInfo.ieRegn[j],
300  rctInfo.ieCard[j],
301  rctInfo.crateID,
302  true);
303  em.setBx(startBX);
304  colls->rctEm()->push_back(em);
305  }
306 
307  for(int j = 0; j < 7; j++) {
308  for(int k = 0; k < 2; k++) {
309  bool o = (((rctInfo.oBits >> (j * 2 + k)) & 0x1) == 0x1);
310  bool t = (((rctInfo.tBits >> (j * 2 + k)) & 0x1) == 0x1);
311  bool m = (((rctInfo.mBits >> (j * 2 + k)) & 0x1) == 0x1);
312  bool q = (((rctInfo.qBits >> (j * 2 + k)) & 0x1) == 0x1);
313  L1CaloRegion rgn = L1CaloRegion(rctInfo.rgnEt[j][k], o, t, m, q, rctInfo.crateID , j, k);
314  rgn.setBx(startBX);
315  colls->rctCalo()->push_back(rgn);
316  }
317  }
318 
319  for(int k = 0; k < 4; k++) {
320  for(int j = 0; j < 2; j++) {
321  // 0 1 4 5 2 3 6 7
322  uint32_t offset = j*2 + k%2 + (k/2)*4;
323  bool fg=(((rctInfo.hfQBits >> offset) & 0x1) == 0x1);
324  L1CaloRegion rgn = L1CaloRegion(rctInfo.hfEt[j][k], fg, rctInfo.crateID , (j * 4 + k));
325  rgn.setBx(startBX);
326  colls->rctCalo()->push_back(rgn);
327  }
328  }
329  }
330  }
331 }
#define LogDebug(id)
size
Write out results.
L1CaloRegionCollection *const rctCalo() const
Input calo regions from the RCT to the RCT.
void setBx(int16_t bx)
set bx
void setBx(int16_t bx)
set BX
Definition: L1CaloEmCand.cc:69
Level-1 Region Calorimeter Trigger EM candidate.
Definition: L1CaloEmCand.h:18
L1CaloEmCollection *const rctEm() const
Input electrons from the RCT to the RCT.
int k[5][pyjets_maxn]
def uint(string)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
A calorimeter trigger region (sum of 4x4 trigger towers)
Definition: L1CaloRegion.h:22

Member Data Documentation

const unsigned RctRawToDigi::amc13HeaderSize_ = 16
staticprivate

Definition at line 80 of file RctRawToDigi.h.

Referenced by produce().

const unsigned RctRawToDigi::amc13TrailerSize_ = 8
staticprivate

Definition at line 83 of file RctRawToDigi.h.

Referenced by produce().

const unsigned RctRawToDigi::ctp7HeaderSize_ = 8
staticprivate

Definition at line 86 of file RctRawToDigi.h.

const unsigned RctRawToDigi::ctp7TrailerSize_ = 8
staticprivate

Definition at line 89 of file RctRawToDigi.h.

std::vector<unsigned> RctRawToDigi::errorCounters_
private

Counts number of errors for each code (index)

Definition at line 109 of file RctRawToDigi.h.

Referenced by endJob().

L1TriggerErrorCollection* RctRawToDigi::errors_
private

pointer to error collection

Definition at line 108 of file RctRawToDigi.h.

int RctRawToDigi::fedId_
private

RCT FED ID.

Definition at line 99 of file RctRawToDigi.h.

Referenced by produce(), and RctRawToDigi().

edm::InputTag RctRawToDigi::inputLabel_
private

FED collection label.

Definition at line 98 of file RctRawToDigi.h.

Referenced by produce(), and RctRawToDigi().

const unsigned RctRawToDigi::MAX_DATA = 4680
staticprivate

The maximum number of blocks we will try to unpack before thinking something is wrong.

Definition at line 92 of file RctRawToDigi.h.

const unsigned RctRawToDigi::MAX_ERR_CODE = 6
staticprivate

Definition at line 107 of file RctRawToDigi.h.

Referenced by endJob().

const unsigned RctRawToDigi::MIN_DATA = 900
staticprivate

The minimum number of blocks we will try to unpack before thinking something is wrong (really this should be 920, to be tested)

Definition at line 95 of file RctRawToDigi.h.

Referenced by produce().

const unsigned RctRawToDigi::sLinkHeaderSize_ = 8
staticprivate

Definition at line 74 of file RctRawToDigi.h.

Referenced by produce().

const unsigned RctRawToDigi::sLinkTrailerSize_ = 8
staticprivate

Definition at line 77 of file RctRawToDigi.h.

Referenced by produce(), and unpack().

unsigned RctRawToDigi::unpackFailures_
private

To count the total number of RCT unpack failures.

Definition at line 110 of file RctRawToDigi.h.

const bool RctRawToDigi::verbose_
private

If true, then debug print out for each event.

Definition at line 101 of file RctRawToDigi.h.

Referenced by endJob().