CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Attributes
cms::DigitizerFP420 Class Reference

#include <DigitizerFP420.h>

Inheritance diagram for cms::DigitizerFP420:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 DigitizerFP420 (const edm::ParameterSet &conf)
 
virtual void produce (edm::Event &e, const edm::EventSetup &c)
 
virtual ~DigitizerFP420 ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Types

typedef std::map< unsigned int,
std::vector< PSimHit >
, std::less< unsigned int > > 
simhit_map
 
typedef simhit_map::iterator simhit_map_iterator
 
typedef std::vector< std::string > vstring
 

Private Attributes

std::vector< HDigiFP420collector
 
edm::ParameterSet conf_
 
int dn0
 
int numStrips
 
int pn0
 
int rn0
 
simhit_map SimHitMap
 
int sn0
 
FP420DigiMainstripDigitizer_
 
FP420NumberingSchemetheFP420NumberingScheme
 
vstring trackerContainers
 
int verbosity
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Definition at line 42 of file DigitizerFP420.h.

Member Typedef Documentation

typedef std::map<unsigned int, std::vector<PSimHit>,std::less<unsigned int> > cms::DigitizerFP420::simhit_map
private

Definition at line 63 of file DigitizerFP420.h.

typedef simhit_map::iterator cms::DigitizerFP420::simhit_map_iterator
private

Definition at line 64 of file DigitizerFP420.h.

typedef std::vector<std::string> cms::DigitizerFP420::vstring
private

Definition at line 62 of file DigitizerFP420.h.

Constructor & Destructor Documentation

cms::DigitizerFP420::DigitizerFP420 ( const edm::ParameterSet conf)
explicit

Definition at line 96 of file DigitizerFP420.cc.

References conf_, gather_cfg::cout, dn0, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), pn0, rn0, sn0, theFP420NumberingScheme, trackerContainers, and verbosity.

96  :conf_(conf),stripDigitizer_(new FP420DigiMain(conf)) {
97 
98 
99  std::string alias ( conf.getParameter<std::string>("@module_label") );
100 
101  // produces<edm::DetSetVector<HDigiFP420> >().setBranchAlias( alias );
102  // produces<edm::DetSetVector<HDigiFP420SimLink> >().setBranchAlias ( alias + "hDigiFP420SimLink");
103  produces<DigiCollectionFP420>().setBranchAlias( alias );
104 
105  trackerContainers.clear();
106  trackerContainers = conf.getParameter<std::vector<std::string> >("ROUList");
107 
108  verbosity = conf_.getUntrackedParameter<int>("VerbosityLevel");
109  dn0 = conf_.getParameter<int>("NumberFP420Detectors");
110  sn0 = conf_.getParameter<int>("NumberFP420Stations");
111  pn0 = conf_.getParameter<int>("NumberFP420SPlanes");
112  rn0 = 7;
113  //rn0 = 3;
115 
116  // produces<DigiCollectionFP420>();
117 
118  // produces<StripDigiCollection>();
119  // produces<HDigiFP420>();
120  // produces<edm::DetSetVector<HDigiFP420> >().setBranchAlias( alias );
121 
122  // produces<DigiCollectionFP420>();
123  // produces<DigiCollectionFP420>("HDigiFP420");
124 
125  // produces<edm::DigiCollectionFP420>();
126 
127  // produces<edm::DetSetVector<DigiCollectionFP420> >();
128 
129  if(verbosity>0) {
130  std::cout << "Creating a DigitizerFP420" << std::endl;
131  std::cout << "DigitizerFP420: dn0=" << dn0 << " sn0=" << sn0 << " pn0=" << pn0 << " rn0=" << rn0 << std::endl;
132  std::cout << "DigitizerFP420:trackerContainers.size()=" << trackerContainers.size() << std::endl;
133 
134  }
135  }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
FP420DigiMain * stripDigitizer_
edm::ParameterSet conf_
FP420NumberingScheme * theFP420NumberingScheme
tuple cout
Definition: gather_cfg.py:121
cms::DigitizerFP420::~DigitizerFP420 ( )
virtual

Definition at line 138 of file DigitizerFP420.cc.

References gather_cfg::cout, stripDigitizer_, and verbosity.

138  {
139  if(verbosity>0) {
140  std::cout << "Destroying a DigitizerFP420" << std::endl;
141  }
142  delete stripDigitizer_;
143 
144  }
FP420DigiMain * stripDigitizer_
tuple cout
Definition: gather_cfg.py:121

Member Function Documentation

void cms::DigitizerFP420::produce ( edm::Event e,
const edm::EventSetup c 
)
virtual

Implements edm::EDProducer.

Definition at line 149 of file DigitizerFP420.cc.

References HDigiFP420::adc(), HDigiFP420::channel(), collector, gather_cfg::cout, PSimHit::detUnitId(), dn0, first, edm::Event::getByLabel(), i, getHLTprescales::index, estimatePileup::inputRange, convertSQLitetoXML_cfg::output, pn0, edm::Event::put(), rn0, FP420DigiMain::run(), SimHitMap, sn0, HDigiFP420::strip(), stripDigitizer_, HDigiFP420::stripV(), HDigiFP420::stripVW(), theFP420NumberingScheme, trackerContainers, FP420NumberingScheme::unpackFP420Index(), and verbosity.

149  {
150  // be lazy and include the appropriate namespaces
151  using namespace edm;
152  using namespace std;
153 
154  if(verbosity>0) {
155  std::cout <<" ===" << std::endl;
156  std::cout <<" ============== DigitizerFP420: start produce= " << std::endl;
157  std::cout <<" ===" << std::endl;
158  }
159  // Get input
160  // std::cout << "DigitizerFP420 start produce" << std::endl;
161  // edm::ESHandle < ParticleDataTable > pdt;
162  // iSetup.getData( pdt );
163 
164  // Step A: Get Inputs for allTrackerHits
165 
166 
168  std::vector<const CrossingFrame<PSimHit> *> cf_simhitvec;
169  for(uint32_t i = 0; i< trackerContainers.size();i++){
170  iEvent.getByLabel("mix",trackerContainers[i],cf_simhit);
171  cf_simhitvec.push_back(cf_simhit.product()); }
172  std::auto_ptr<MixCollection<PSimHit> > allTrackerHits(new MixCollection<PSimHit>(cf_simhitvec));
173 
174  // use instead of the previous
175  /*
176  std::cout <<" ============== DigitizerFP420: start loop 1 " << std::endl;
177  edm::Handle<CrossingFrame<PSimHit> > xFrame;
178  std::cout <<" ============== DigitizerFP420: start loop 2 " << std::endl;
179  iEvent.getByLabel("mix","FP420SI",xFrame);
180  std::cout <<" ============== DigitizerFP420: start loop 3 " << std::endl;
181  std::auto_ptr<MixCollection<PSimHit> > allTrackerHits( new MixCollection<PSimHit>(xFrame.product()) );
182  std::cout <<" ============== DigitizerFP420: start loop 4 " << std::endl;
183  */
184 
185  // use instead of the previous
186  /*
187  edm::Handle<CrossingFrame<PSimHit> > crossingFrame;
188  const std::string FP420HitsName("FP420SI");
189  bool isHit = true;
190  iEvent.getByLabel("mix",FP420HitsName,crossingFrame);
191  MixCollection<PSimHit> * FP420Hits = 0 ;
192  std::cout <<" ============== DigitizerFP420: start loop 1 " << std::endl;
193  // std::auto_ptr<MixCollection<PSimHit> > allTrackerHits(new MixCollection<PSimHit>(crossingFrame.product()));
194  FP420Hits = new MixCollection<PSimHit>(crossingFrame.product());
195  std::cout <<" ============== DigitizerFP420: start loop 2 " << std::endl;
196  // if ( ! FP420Hits->inRegistry() ) isHit = false;
197  // if ( isHit ) {
198  std::auto_ptr<MixCollection<PSimHit> > allTrackerHits( FP420Hits );
199  std::cout <<" ============== DigitizerFP420: start loop 3 " << std::endl;
200  // }
201  */
202 
203  // std::cout << "DigitizerFP420 Step A done" << std::endl;
204 
205  //Loop on PSimHit
206 
207 
209  // Step C: create empty output collection
210  std::auto_ptr<DigiCollectionFP420> output(new DigiCollectionFP420);
211  // std::auto_ptr<edm::DetSetVector<HDigiFP420> > outputfinal(new edm::DetSetVector<HDigiFP420>(output) );
212  // std::auto_ptr<edm::DetSetVector<HDigiFP420> > outputfinal(new edm::DetSetVector<HDigiFP420>(output) );
213  // std::auto_ptr<edm::DetSetVector<HDigiFP420SimLink> > outputlink(new edm::DetSetVector<HDigiFP420SimLink>(output) );
214 
215  SimHitMap.clear();
216 
217  // ==================================
218  if(verbosity>0) {
219  std::cout <<" ===" << std::endl;
220  std::cout <<" ============== DigitizerFP420: MixCollection treatment= " << std::endl;
221  std::cout <<" ===" << std::endl;
222  }
223 
225  for (isim=allTrackerHits->begin(); isim!= allTrackerHits->end();isim++) {
226  unsigned int unitID = (*isim).detUnitId();
227  int det, zside, sector, zmodule;
228  FP420NumberingScheme::unpackFP420Index(unitID, det, zside, sector, zmodule);
229  // below, the continues plane index should be (for even different sensor index zside)
230  // unsigned int intindex = packMYIndex(rn0, pn0, sn0, det, zside, sector, zmodule);
231  unsigned int intindex = theFP420NumberingScheme->FP420NumberingScheme::packMYIndex(rn0, pn0, sn0, det, zside, sector, zmodule);
232  // int zScale=(rn0-1), sScale = (rn0-1)*(pn0-1), dScale = (rn0-1)*(pn0-1)*(sn0-1);
233  // unsigned int intindex = dScale*(det - 1)+sScale*(sector - 1)+zScale*(zmodule - 1)+zside;
234 
235  if(verbosity>0) {
236  double losenergy = (*isim).energyLoss();
237  std::cout <<" ===" << std::endl;
238  std::cout <<" ============== DigitizerFP420: losenergy= " << losenergy << std::endl;
239  std::cout <<" === for intindex = " << intindex << std::endl;
240  }
241  // does not matter which index is used: intindex or unitID - mainly to collect hits under every index
242  SimHitMap[intindex].push_back((*isim));
243  // for development later one( cal be used another index):
244  // SimHitMap[unitID].push_back((*isim));
245  }
246  //============================================================================================================================
247 
248  if(verbosity>0) {
249  std::cout <<" ===" << std::endl;
250  std::cout <<" ============== DigitizerFP420: put zero to container " << std::endl;
251  std::cout <<" ===" << std::endl;
252  }
253  // put zero to container info from the beginning (important! because not any detID is updated with coming of new event !!!!!!
254  // clean info of container from previous event
255  for (int det=1; det<dn0; det++) {
256  for (int sector=1; sector<sn0; sector++) {
257  for (int zmodule=1; zmodule<pn0; zmodule++) {
258  for (int zside=1; zside<rn0; zside++) {
259  // intindex is a continues numbering of FP420
260  unsigned int detID = theFP420NumberingScheme->FP420NumberingScheme::packMYIndex(rn0, pn0, sn0, det, zside, sector, zmodule);
261  // int zScale=(rn0-1), sScale = (rn0-1)*(pn0-1), dScale = (rn0-1)*(pn0-1)*(sn0-1);
262  // unsigned int detID = dScale*(det - 1)+sScale*(sector - 1)+zScale*(zmodule - 1)+zside;
263  std::vector<HDigiFP420> collector;
264  collector.clear();
266  inputRange.first = collector.begin();
267  inputRange.second = collector.end();
268  output->putclear(inputRange,detID);
269  }//for
270  }//for
271  }//for
272  }//for
273  // !!!!!!
274  // if we want to keep Digi container/Collection for one event uncomment the line below and vice versa
275  output->clear(); //container_.clear() --> start from the beginning of the container
276 
277  //============================================================================================================================================
278 
279 
280  if(verbosity>0) {
281  std::cout <<" ===" << std::endl;
282  std::cout <<" ============== DigitizerFP420: start loop over det iu " << std::endl;
283  std::cout <<" ============== DigitizerFP420: SimHitMap.size()= " << SimHitMap.size() << std::endl;
284  std::cout <<" ===" << std::endl;
285  }
286  bool first = true;
287 
289  /*
290  if(verbosity>0) std::cout <<"======= DigitizerFP420: SimHitMap size = " << SimHitMap.size() << std::endl;
291  for(unsigned int i = 0; i < SimHitMap.size(); i++ ) {
292  // std::cout <<" ====== DigitizerFP420: i= " << i << std::endl;
293  vector<PSimHit>::const_iterator simHitIter = SimHitMap[i].begin();
294  vector<PSimHit>::const_iterator simHitIterEnd = SimHitMap[i].end();
295  for (;simHitIter != simHitIterEnd; ++simHitIter) {
296  const PSimHit ihit = *simHitIter;
297  unsigned int unitID = ihit.detUnitId();
298  if(verbosity>0) std::cout <<" ====== DigitizerFP420: unitID= " << unitID << " i= " << i << std::endl;
299  int det, zside, sector, zmodule;
300  FP420NumberingScheme::unpackFP420Index(unitID, det, zside, sector, zmodule);
301  int zScale=(rn0-1), sScale = (rn0-1)*(pn0-1), dScale = (rn0-1)*(pn0-1)*(sn0-1);
302  unsigned int iu = dScale*(det - 1)+sScale*(sector - 1)+zScale*(zmodule - 1)+zside;
303  }
304  }
305 */
307  //============================================================================================================================================
308  // new: <------
309  for(unsigned int i = 0; i < SimHitMap.size(); i++ ) {
310  vector<PSimHit>::const_iterator simHitIter = SimHitMap[i].begin();
311  vector<PSimHit>::const_iterator simHitIterEnd = SimHitMap[i].end();
312  for (;simHitIter != simHitIterEnd; ++simHitIter) {
313  const PSimHit ihit = *simHitIter;
314  unsigned int unitID = ihit.detUnitId();
315  if(verbosity>0 || verbosity==-50) std::cout <<" ====== DigitizerFP420: unitID= " << unitID << "Hit number i= " << i << std::endl;
316  int det, zside, sector, zmodule;
317  FP420NumberingScheme::unpackFP420Index(unitID, det, zside, sector, zmodule);
318  // <------
319  // old: <------
320  // for (int det=1; det<dn0; det++) {
321  // for (int sector=1; sector<sn0; sector++) {
322  // for (int zmodule=1; zmodule<pn0; zmodule++) {
323  // for (int zside=1; zside<rn0; zside++) {
324  // <------
325 
326 
327 
328 
329 
330  unsigned int iu = theFP420NumberingScheme->FP420NumberingScheme::packMYIndex(rn0, pn0, sn0, det, zside, sector, zmodule);
331  if(verbosity>0 || verbosity==-50) std::cout <<"for Hits iu = " << iu <<" sector = " << sector <<" zmodule = " << zmodule <<" zside = " << zside << " det=" << det << std::endl;
332  // int zScale=(rn0-1), sScale = (rn0-1)*(pn0-1), dScale = (rn0-1)*(pn0-1)*(sn0-1);
333  // unsigned int iu = dScale*(det - 1)+sScale*(sector - 1)+zScale*(zmodule - 1)+zside;
334 
335  if(verbosity>0) {
336  unsigned int index = theFP420NumberingScheme->FP420NumberingScheme::packFP420Index(det, zside, sector, zmodule);
337  std::cout << " DigitizerFP420: index = " << index << " iu = " << iu << std::endl;
338  }
339 
340  // GlobalVector bfield=pSetup->inTesla((*iu)->surface().position());
341  // CLHEP::Hep3Vector Bfieldloc=bfield();
342  G4ThreeVector bfield( 0., 0., 0.0 );
343  // G4ThreeVector bfield( 0.5, 0.5, 1.0 );
344 
345 
346  if(verbosity>0) {
347  std::cout <<" ===" << std::endl;
348  std::cout <<" ============== DigitizerFP420: call run for iu= " << iu << std::endl;
349  std::cout <<" ===" << std::endl;
350  }
351  collector.clear();
352 
353  collector= stripDigitizer_->run( SimHitMap[iu],
354  bfield,
355  iu
356  ); // stripDigitizer_.run... return
357  // ,sScale
358 
359 
360 
361  if(verbosity>0) {
362  std::cout <<" ===" << std::endl;
363  std::cout <<" ===" << std::endl;
364  std::cout <<"======= DigitizerFP420: collector size = " << collector.size() << std::endl;
365  std::cout <<" ===" << std::endl;
366  std::cout <<" ===" << std::endl;
367  }
368  /*
369 
370  std::vector<HDigiFP420> collector;
371  collector.clear();
372  DigiCollectionFP420::Range inputRange;
373  inputRange.first = collector.begin();
374  inputRange.second = collector.end();
375  output->putclear(inputRange,detID);
376  */
377  if (collector.size()>0){
378  if(verbosity>0) {
379  std::cout <<" ============= DigitizerFP420:collector start!!!!!!!!!!!!!!" << std::endl;
380  }
381  DigiCollectionFP420::Range outputRange;
382  outputRange.first = collector.begin();
383  outputRange.second = collector.end();
384 
385  if ( first ) {
386  // use it only if ClusterCollectionFP420 is the ClusterCollection of one event, otherwise, do not use (loose 1st cl. of 1st event only)
387  first = false;
388  unsigned int detID0= 0;
389  output->put(outputRange,detID0); // !!! put into adress 0 for detID which will not be used never
390  } //if ( first )
391 
392  // put !!!
393  output->put(outputRange,iu);
394 
395  } // if(collector.size()>0
396 
397  // } // for
398  // } // for
399  // } // for
400  // } // for
401 
402  }//for
403  }//for
404 
405  // END
406 
407  /*
408  if(verbosity>0) {
409  std::vector<HDigiFP420> theAllDigis;
410  theAllDigis.clear();
411  DigiCollectionFP420::Range outputRange;
412  DigiCollectionFP420::ContainerIterator sort_begin = outputRange.first;
413  DigiCollectionFP420::ContainerIterator sort_end = outputRange.second;
414  theAllDigis.insert(theAllDigis.end(), sort_begin, sort_end);
415  std::cout <<"====== theAllDigis size = " << theAllDigis.size() << std::endl;
416  for (std::vector<HDigiFP420>::iterator isim = theAllDigis.begin();
417  isim != theAllDigis.end(); ++isim){
418  const HDigiFP420 istrip = *isim;
419  std::cout << "*******************************************DigitizerFP420:check1" << std::endl;
420  std::cout << " strip number=" << istrip.strip() << " adc=" << istrip.adc() << std::endl;
421  std::cout <<" channel =" << istrip.channel() <<" V " << istrip.stripV() <<" VW " << istrip.stripVW() << std::endl;
422  std::cout <<" ===" << std::endl;
423  std::cout <<" ===" << std::endl;
424  std::cout <<" =======================" << std::endl;
425  }// for
426  }
427 */
428  if(verbosity==-50) {
429  // check of access to the collector:
430  for (int det=1; det<dn0; det++) {
431  for (int sector=1; sector<sn0; sector++) {
432  for (int zmodule=1; zmodule<pn0; zmodule++) {
433  for (int zside=1; zside<rn0; zside++) {
434  unsigned int iu = theFP420NumberingScheme->FP420NumberingScheme::packMYIndex(rn0, pn0, sn0, det, zside, sector, zmodule);
435  int layer = theFP420NumberingScheme->FP420NumberingScheme::unpackLayerIndex(rn0,zside);
436  int orient = theFP420NumberingScheme->FP420NumberingScheme::unpackOrientation(rn0,zside);
437  std::cout << "****DigitizerFP420:check2" << std::endl;
438  // std::cout <<" iu = " << iu <<" sector = " << sector <<" zmodule = " << zmodule <<" zside = " << zside << " det=" << det << std::endl;
439  // std::cout <<" layer = " << layer <<" orient = " << orient << std::endl;
440  int newdet, newzside, newsector, newzmodule;
441  theFP420NumberingScheme->FP420NumberingScheme::unpackMYIndex(iu, rn0, pn0, sn0, newdet, newzside, newsector, newzmodule);
442  std::cout <<" newdet = " << newdet <<" newsector = " << newsector <<" newzmodule = " << newzmodule <<" newzside = " << newzside << std::endl;
443 
444  collector.clear();
445  DigiCollectionFP420::Range outputRange;
446  // outputRange = output->get(iu);
447  outputRange = output->get(iu);
448 
449  // fill output in collector vector (for may be sorting? or other checks)
450  std::vector<HDigiFP420> collector;
451  // collector.clear();
452  DigiCollectionFP420::ContainerIterator sort_begin = outputRange.first;
453  DigiCollectionFP420::ContainerIterator sort_end = outputRange.second;
454  for ( ;sort_begin != sort_end; ++sort_begin ) {
455  collector.push_back(*sort_begin);
456  } // for
457  // std::sort(collector.begin(),collector.end());
458  std::cout <<" ===" << std::endl;
459  std::cout <<"====== collector size = " << collector.size() << std::endl;
460  if(collector.size()>0) {
461  std::cout <<" iu = " << iu <<" sector = " << sector <<" zmodule = " << zmodule <<" zside = " << zside << " det=" << det <<" layer = " << layer <<" orient = " << orient << std::endl;
462  std::cout <<" ===" << std::endl;
463  }
464  vector<HDigiFP420>::const_iterator simHitIter = collector.begin();
465  vector<HDigiFP420>::const_iterator simHitIterEnd = collector.end();
466  for (;simHitIter != simHitIterEnd; ++simHitIter) {
467  const HDigiFP420 istrip = *simHitIter;
468  std::cout << " strip number=" << istrip.strip() << " adc=" << istrip.adc() << std::endl;
469  std::cout <<" channel =" << istrip.channel() <<" V " << istrip.stripV() <<" VW " << istrip.stripVW() << std::endl;
470  std::cout <<" ===" << std::endl;
471  std::cout <<" ===" << std::endl;
472  std::cout <<" ===================================================" << std::endl;
473  }
474 
475  //==================================
476 
477  } // for
478  } // for
479  } // for
480  } // for
481 
482  // end of check of access to the strip collection
483 
484  }// if(verbosity
485  //
486 
487 
488  // Step D: write output to file
489  // iEvent.put(output);
490 
491  if(verbosity>0) {
492  std::cout << "DigitizerFP420 recoutput" << std::endl;
493  }
494  // Step D: write output to file
495  iEvent.put(output);
496  // iEvent.put(outputlink);
497  // iEvent.put(pDigis);
498 
499  // Step D: write output to file
500  // iEvent.put(output);
501  // iEvent.put(outputlink);
502  //-------------------------------------------------------------------
503  // std::cout << "DigitizerFP420 recoutput" << std::endl;
504  // iEvent.put(pDigis);
505 
506 
507 
508 
509  }//produce
int channel() const
Definition: HDigiFP420.h:22
int i
Definition: DBlmapReader.cc:9
inputRange
Get input source.
std::vector< HDigiFP420 > run(const std::vector< PSimHit > &input, G4ThreeVector, unsigned int)
FP420DigiMain * stripDigitizer_
int stripVW() const
Definition: HDigiFP420.h:24
std::vector< HDigiFP420 >::const_iterator ContainerIterator
static void unpackFP420Index(const unsigned int &idx, int &det, int &zside, int &station, int &superplane)
std::vector< HDigiFP420 > collector
int strip() const
Definition: HDigiFP420.h:20
int iEvent
Definition: GenABIO.cc:243
int adc() const
Definition: HDigiFP420.h:21
bool first
Definition: L1TdeRCT.cc:94
std::pair< ContainerIterator, ContainerIterator > Range
FP420NumberingScheme * theFP420NumberingScheme
tuple cout
Definition: gather_cfg.py:121
unsigned int detUnitId() const
Definition: PSimHit.h:93
int stripV() const
Definition: HDigiFP420.h:26

Member Data Documentation

std::vector<HDigiFP420> cms::DigitizerFP420::collector
private

Definition at line 80 of file DigitizerFP420.h.

Referenced by produce().

edm::ParameterSet cms::DigitizerFP420::conf_
private

Definition at line 67 of file DigitizerFP420.h.

Referenced by DigitizerFP420().

int cms::DigitizerFP420::dn0
private

Definition at line 77 of file DigitizerFP420.h.

Referenced by DigitizerFP420(), and produce().

int cms::DigitizerFP420::numStrips
private

Definition at line 75 of file DigitizerFP420.h.

int cms::DigitizerFP420::pn0
private

Definition at line 77 of file DigitizerFP420.h.

Referenced by DigitizerFP420(), and produce().

int cms::DigitizerFP420::rn0
private

Definition at line 77 of file DigitizerFP420.h.

Referenced by DigitizerFP420(), and produce().

simhit_map cms::DigitizerFP420::SimHitMap
private

Definition at line 65 of file DigitizerFP420.h.

Referenced by produce().

int cms::DigitizerFP420::sn0
private

Definition at line 77 of file DigitizerFP420.h.

Referenced by DigitizerFP420(), and produce().

FP420DigiMain* cms::DigitizerFP420::stripDigitizer_
private

Definition at line 72 of file DigitizerFP420.h.

Referenced by produce(), and ~DigitizerFP420().

FP420NumberingScheme* cms::DigitizerFP420::theFP420NumberingScheme
private

Definition at line 73 of file DigitizerFP420.h.

Referenced by DigitizerFP420(), and produce().

vstring cms::DigitizerFP420::trackerContainers
private

Definition at line 68 of file DigitizerFP420.h.

Referenced by DigitizerFP420(), and produce().

int cms::DigitizerFP420::verbosity
private

Definition at line 77 of file DigitizerFP420.h.

Referenced by DigitizerFP420(), produce(), and ~DigitizerFP420().