CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes
ME0ReDigiProducer Class Reference

#include <ME0ReDigiProducer.h>

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

Classes

class  TemporaryGeometry
 

Public Member Functions

void beginRun (const edm::Run &, const edm::EventSetup &) override
 
void buildDigis (const ME0DigiPreRecoCollection &, ME0DigiPreRecoCollection &, ME0DigiPreRecoMap &, CLHEP::HepRandomEngine *engine)
 
 ME0ReDigiProducer (const edm::ParameterSet &ps)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~ME0ReDigiProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Private Types

typedef std::map< DigiIndicies, unsigned int > ChamberDigiMap
 
typedef std::tuple< unsigned int, unsigned int, unsigned int > DigiIndicies
 

Private Member Functions

void fillCentralTOFs ()
 
unsigned int fillDigiMap (ChamberDigiMap &chDigiMap, unsigned int bx, unsigned int part, unsigned int strip, unsigned int currentIDX) const
 
int getCustomStripProperties (const ME0DetId &detId, const ME0DigiPreReco *inDigi, float &tof, int &strip, LocalPoint &digiLocalPoint, LocalError &digiLocalError) const
 
void getStripProperties (const ME0EtaPartition *etaPart, const ME0DigiPreReco *inDigi, float &tof, int &strip, LocalPoint &digiLocalPoint, LocalError &digiLocalError) const
 

Private Attributes

const float bxWidth
 
const ME0Geometrygeometry
 
std::vector< int > layerReadout
 
int maxBXReadout
 
bool mergeDigis
 
int minBXReadout
 
double neutronAcceptance
 
unsigned int numberOfPartitions
 
unsigned int numberOfStrips
 
TemporaryGeometrytempGeo
 
double timeResolution
 
std::vector< std::vector< double > > tofs
 
edm::EDGetTokenT< ME0DigiPreRecoCollectiontoken
 
bool useBuiltinGeo
 
bool useCusGeoFor1PartGeo
 
bool usePads
 

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

Definition at line 29 of file ME0ReDigiProducer.h.

Member Typedef Documentation

typedef std::map<DigiIndicies,unsigned int> ME0ReDigiProducer::ChamberDigiMap
private

Definition at line 77 of file ME0ReDigiProducer.h.

typedef std::tuple<unsigned int, unsigned int, unsigned int> ME0ReDigiProducer::DigiIndicies
private

Definition at line 76 of file ME0ReDigiProducer.h.

Constructor & Destructor Documentation

ME0ReDigiProducer::ME0ReDigiProducer ( const edm::ParameterSet ps)
explicit

Definition at line 136 of file ME0ReDigiProducer.cc.

References Exception, neutronAcceptance, numberOfPartitions, numberOfStrips, tempGeo, useBuiltinGeo, useCusGeoFor1PartGeo, and usePads.

136  :
137  bxWidth (25.0),
138  useCusGeoFor1PartGeo(ps.getParameter<bool>("useCusGeoFor1PartGeo")),
139  usePads(ps.getParameter<bool>("usePads")),
140  numberOfStrips (ps.getParameter<unsigned int>("numberOfStrips")),
141  numberOfPartitions (ps.getParameter<unsigned int>("numberOfPartitions")),
142  neutronAcceptance (ps.getParameter<double>("neutronAcceptance")),
143  timeResolution (ps.getParameter<double>("timeResolution")),
144  minBXReadout (ps.getParameter<int>("minBXReadout")),
145  maxBXReadout (ps.getParameter<int>("maxBXReadout")),
146  layerReadout (ps.getParameter<std::vector<int>>("layerReadout")),
147  mergeDigis (ps.getParameter<bool>("mergeDigis")),
148  token(consumes<ME0DigiPreRecoCollection>(edm::InputTag(ps.getParameter<std::string>("inputCollection"))))
149 {
150  produces<ME0DigiPreRecoCollection>();
151  produces<ME0DigiPreRecoMap>();
152 
154  if (!rng.isAvailable()){
155  throw cms::Exception("Configuration")
156  << "ME0ReDigiProducer::ME0PreRecoDigiProducer() - RandomNumberGeneratorService is not present in configuration file.\n"
157  << "Add the service in the configuration file or remove the modules that require it.";
158  }
159  geometry = nullptr;
160  tempGeo = nullptr;
161  useBuiltinGeo = true;
162 
164  if (usePads)
166  if(numberOfStrips == 0)
167  throw cms::Exception("Setup") << "ME0ReDigiProducer::ME0PreRecoDigiProducer() - Must have at least one strip if using custom geometry.";
168  if(numberOfPartitions == 0)
169  throw cms::Exception("Setup") << "ME0ReDigiProducer::ME0PreRecoDigiProducer() - Must have at least one partition if using custom geometry.";
170  }
171 
172  if(neutronAcceptance < 0 )
173  throw cms::Exception("Setup") << "ME0ReDigiProducer::ME0PreRecoDigiProducer() - neutronAcceptance must be >= 0.";
174 }
T getParameter(std::string const &) const
edm::EDGetTokenT< ME0DigiPreRecoCollection > token
std::vector< int > layerReadout
TemporaryGeometry * tempGeo
unsigned int numberOfPartitions
unsigned int numberOfStrips
ME0ReDigiProducer::~ME0ReDigiProducer ( )
override

Definition at line 177 of file ME0ReDigiProducer.cc.

References tempGeo.

178 {
179  if(tempGeo) delete tempGeo;
180 }
TemporaryGeometry * tempGeo

Member Function Documentation

void ME0ReDigiProducer::beginRun ( const edm::Run ,
const edm::EventSetup eventSetup 
)
override

Definition at line 183 of file ME0ReDigiProducer.cc.

References chambers, Exception, fillCentralTOFs(), edm::EventSetup::get(), layerReadout, LogDebug, MuonTCMETValueMapProducer_cff::nLayers, numberOfPartitions, numberOfStrips, ME0ReDigiProducer::TemporaryGeometry::numLayers(), tempGeo, useBuiltinGeo, and useCusGeoFor1PartGeo.

184 {
185  // set geometry
187  eventSetup.get<MuonGeometryRecord>().get(hGeom);
188  geometry= &*hGeom;
189 
190  const auto& chambers = geometry->chambers();
191  if(chambers.empty())
192  throw cms::Exception("Setup") << "ME0ReDigiProducer::beginRun() - No ME0Chambers in geometry.";
193 
194  const unsigned int nLayers = chambers.front()->nLayers();
195  if(!nLayers) throw cms::Exception("Setup") << "ME0ReDigiProducer::beginRun() - No layers in ME0 geometry.";
196 
197  const unsigned int nPartitions = chambers.front()->layers()[0]->nEtaPartitions();
198 
199  if(useCusGeoFor1PartGeo && nPartitions == 1){
200  useBuiltinGeo = false;
201  }
202 
203  if(useBuiltinGeo){
204  if(nLayers != layerReadout.size() )
205  throw cms::Exception("Configuration") << "ME0ReDigiProducer::beginRun() - The geometry has "<<nLayers
206  << " layers, but the readout of "<<layerReadout.size() << " were specified with the layerReadout parameter." ;
207  fillCentralTOFs();
208  } else {
209  LogDebug("ME0ReDigiProducer")
210  << "Building temporary geometry:" << std::endl;
211  tempGeo = new TemporaryGeometry(geometry,numberOfStrips,numberOfPartitions);
212  LogDebug("ME0ReDigiProducer")
213  << "Done building temporary geometry!" << std::endl;
214 
215  if(tempGeo->numLayers() != layerReadout.size() )
216  throw cms::Exception("Configuration") << "ME0ReDigiProducer::beginRun() - The geometry has "<<tempGeo->numLayers()
217  << " layers, but the readout of "<<layerReadout.size() << " were specified with the layerReadout parameter." ;
218  }
219 }
#define LogDebug(id)
std::vector< int > layerReadout
TemporaryGeometry * tempGeo
unsigned int numberOfPartitions
unsigned int numberOfStrips
T get() const
Definition: EventSetup.h:71
static char chambers[264][20]
Definition: ReadPGInfo.cc:243
void ME0ReDigiProducer::buildDigis ( const ME0DigiPreRecoCollection input_digis,
ME0DigiPreRecoCollection output_digis,
ME0DigiPreRecoMap output_digimap,
CLHEP::HepRandomEngine *  engine 
)

Definition at line 245 of file ME0ReDigiProducer.cc.

References bxWidth, fillDigiMap(), getCustomStripProperties(), getStripProperties(), MuonDigiCollection< IndexType, DigiType >::insertDigi(), layerReadout, LogDebug, LogTrace, maxBXReadout, mergeDigis, minBXReadout, neutronAcceptance, fftjetvertexadder_cfi::sigmaX, fftjetvertexadder_cfi::sigmaY, mathSSE::sqrt(), digitizers_cfi::strip, timeResolution, useBuiltinGeo, PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), LocalError::xy(), PV3DBase< T, PVType, FrameType >::y(), and LocalError::yy().

Referenced by produce().

249 {
250 
251  /*
252  Starting form the incoming pseudo-digi, which has perfect time and position resolution:
253  1A. Smear time using sigma_t by some value
254  1B. Correct the smeared time with the central arrival time for partition
255  1C. Apply discretization: if the smeared time is outside the BX window (-12.5ns;+12.5ns),
256  the hit should be assigned to the next (or previous) BX
257 
258  2A. Find strip that the digi belongs to
259  2B. Get the center of this strip and the error on the position assuming the geometry
260 
261  3A. Filter event if a digi at this partition/strip/BX already exists
262  3B. Add to collection
263  */
264 
265  LogDebug("ME0ReDigiProducer::buildDigis") << "Begin building digis."<<std::endl;
267  for (me0dgIt = input_digis.begin(); me0dgIt != input_digis.end();
268  ++me0dgIt){
269 
270  const auto& me0Id = (*me0dgIt).first;
271  LogTrace("ME0ReDigiProducer::buildDigis") << "Starting with roll: "<< me0Id<<std::endl;
272 
273  //setup map for this chamber/eta partition
274  ChamberDigiMap chDigiMap;
275 
276  int newDigiIdx = 0;
277  const ME0DigiPreRecoCollection::Range& range = (*me0dgIt).second;
278  for (ME0DigiPreRecoCollection::const_iterator digi = range.first;
279  digi != range.second;digi++) {
280  LogTrace("ME0ReDigiProducer::buildDigis") << std::endl<< "(" <<digi->x() <<","<< digi->y()<<","<<digi->tof()<<","<<digi->pdgid()<<","<<digi->prompt()<<")-> ";
281 
282  //If we don't readout this layer skip
283  if(!layerReadout[me0Id.layer() -1 ]) {
284  output_digimap.insertDigi(me0Id, -1);
285  continue;
286  }
287 
288  //if neutron and we are filtering skip
289  if(!digi->prompt() && neutronAcceptance < 1.0 )
290  if (CLHEP::RandFlat::shoot(engine) > neutronAcceptance){
291  output_digimap.insertDigi(me0Id, -1);
292  continue;
293  }
294 
295  //smear TOF if necessary
296  float tof = digi->tof() + (timeResolution < 0 ? 0.0 : CLHEP::RandGaussQ::shoot(engine, 0, timeResolution));
297 
298  //Values used to fill objet
299  int mapPartIDX = me0Id.roll() -1;
300  int strip = 0;
301  LocalPoint digiLocalPoint;
302  LocalError digiLocalError;
303  if(useBuiltinGeo){
304  getStripProperties(geometry->etaPartition(me0Id),&*digi,tof,strip,digiLocalPoint,digiLocalError);
305  } else {
306  mapPartIDX = getCustomStripProperties(me0Id,&*digi,tof,strip,digiLocalPoint,digiLocalError);
307 
308  }
309 
310  //filter if outside of readout window
311  const int bxIdx = std::round(tof/bxWidth);
312  LogTrace("ME0ReDigiProducer::buildDigis") << tof <<"("<<bxIdx<<") ";
313  if(bxIdx < minBXReadout) {output_digimap.insertDigi(me0Id, -1); continue; }
314  if(bxIdx > maxBXReadout) {output_digimap.insertDigi(me0Id, -1); continue; }
315  tof = bxIdx*bxWidth;
316 
317 
318  //If we are merging check to see if it already exists
319  LogTrace("ME0ReDigiProducer::buildDigis") << "("<<bxIdx<<","<<mapPartIDX<<","<<strip<<") ";
320  if(mergeDigis){
321  int matchIDX = fillDigiMap(chDigiMap, bxIdx,mapPartIDX,strip,newDigiIdx);
322  if(matchIDX >= 0){
323  output_digimap.insertDigi(me0Id, matchIDX);
324  continue;
325  }
326  }
327 
328  //Digis store sigmaX,sigmaY, correlationCoef
329  const float sigmaX = std::sqrt(digiLocalError.xx());
330  const float sigmaY = std::sqrt(digiLocalError.yy());
331  const float corrCoef = digiLocalError.xy() /(sigmaX*sigmaY);
332 
333  //Fill in the new collection
334  ME0DigiPreReco out_digi(digiLocalPoint.x(), digiLocalPoint.y(),
335  sigmaX, sigmaY, corrCoef, tof, digi->pdgid(), digi->prompt());
336  output_digis.insertDigi(me0Id, out_digi);
337 
338  // store index of previous detid and digi
339  output_digimap.insertDigi(me0Id, newDigiIdx);
340  newDigiIdx++;
341 
342  LogTrace("ME0ReDigiProducer::buildDigis") << "("<<digiLocalPoint.x()<<","<<digiLocalPoint.y()<<","<<sigmaX<<","<<sigmaY<<","<< tof<<") ";
343  }
344 
345  chDigiMap.clear();
346 
347 
348  }
349 
350 }
#define LogDebug(id)
float xx() const
Definition: LocalError.h:24
int getCustomStripProperties(const ME0DetId &detId, const ME0DigiPreReco *inDigi, float &tof, int &strip, LocalPoint &digiLocalPoint, LocalError &digiLocalError) const
std::vector< int > layerReadout
T y() const
Definition: PV3DBase.h:63
void insertDigi(const IndexType &index, const DigiType &digi)
insert a digi for a given DetUnit
unsigned int fillDigiMap(ChamberDigiMap &chDigiMap, unsigned int bx, unsigned int part, unsigned int strip, unsigned int currentIDX) const
float xy() const
Definition: LocalError.h:25
float yy() const
Definition: LocalError.h:26
T sqrt(T t)
Definition: SSEVec.h:18
std::map< DigiIndicies, unsigned int > ChamberDigiMap
#define LogTrace(id)
std::vector< DigiType >::const_iterator const_iterator
std::pair< const_iterator, const_iterator > Range
T x() const
Definition: PV3DBase.h:62
void getStripProperties(const ME0EtaPartition *etaPart, const ME0DigiPreReco *inDigi, float &tof, int &strip, LocalPoint &digiLocalPoint, LocalError &digiLocalError) const
void ME0ReDigiProducer::fillCentralTOFs ( )
private

Definition at line 352 of file ME0ReDigiProducer.cc.

References Exception, LogDebug, MuonTCMETValueMapProducer_cff::nLayers, and tofs.

Referenced by beginRun().

352  {
353  const auto* mainChamber = geometry->chambers().front();
354  const unsigned int nLayers = mainChamber->nLayers();
355  //Get TOF at center of each partition
356  tofs.clear();
357  tofs.resize(nLayers);
358  LogDebug("ME0ReDigiProducer::fillCentralTOFs()") << "TOF numbers [layer][partition]: " ;
359  for(unsigned int iL = 0; iL < nLayers; ++iL){
360  const auto* layer = mainChamber->layers()[iL];
361  const unsigned int mapLayIDX = layer->id().layer() -1;
362  const unsigned int nPartitions = layer->nEtaPartitions();
363  if(!nPartitions)
364  throw cms::Exception("Setup") << "ME0ReDigiProducer::fillCentralTOFs() - ME0Layer has no partitions.";
365  tofs[mapLayIDX].resize(nPartitions);
366  for(unsigned int iP = 0; iP < nPartitions; ++iP){
367  const unsigned int mapPartIDX = layer->etaPartitions()[iP]->id().roll() -1;
368  const GlobalPoint centralGP(layer->etaPartitions()[iP]->position());
369  tofs[mapLayIDX][mapPartIDX] = (centralGP.mag() / (CLHEP::c_light/CLHEP::cm)); //speed of light [cm/ns]
370  LogDebug("ME0ReDigiProducer::fillCentralTOFs()") << "["<<mapLayIDX<<"]["<<mapPartIDX<<"]="<< tofs[mapLayIDX][mapPartIDX] <<" "<<std::endl;
371  }
372  }
373 }
#define LogDebug(id)
std::vector< std::vector< double > > tofs
unsigned int ME0ReDigiProducer::fillDigiMap ( ChamberDigiMap chDigiMap,
unsigned int  bx,
unsigned int  part,
unsigned int  strip,
unsigned int  currentIDX 
) const
private

Definition at line 417 of file ME0ReDigiProducer.cc.

Referenced by buildDigis().

417  {
418  DigiIndicies newIDX(bx,part,strip);
419  auto it1 = chDigiMap.find(newIDX);
420  if (it1 == chDigiMap.end()){
421  chDigiMap[newIDX] = currentIDX;
422  return -1;
423  }
424  return it1->second;
425 }
std::tuple< unsigned int, unsigned int, unsigned int > DigiIndicies
part
Definition: HCALResponse.h:20
int ME0ReDigiProducer::getCustomStripProperties ( const ME0DetId detId,
const ME0DigiPreReco inDigi,
float &  tof,
int &  strip,
LocalPoint digiLocalPoint,
LocalError digiLocalError 
) const
private

Definition at line 374 of file ME0ReDigiProducer.cc.

References ME0ReDigiProducer::TemporaryGeometry::findEtaPartition(), objects.autophobj::float, ME0ReDigiProducer::TemporaryGeometry::getCentralTOF(), ME0ReDigiProducer::TemporaryGeometry::getPartCenter(), ME0ReDigiProducer::TemporaryGeometry::getTopo(), LogTrace, partIdx(), mathSSE::sqrt(), tempGeo, ME0DigiPreReco::x(), PV3DBase< T, PVType, FrameType >::x(), ME0DigiPreReco::y(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by buildDigis().

374  {
375  const unsigned int partIdx = tempGeo->findEtaPartition(inDigi->y());
376  LogTrace("ME0ReDigiProducer::buildDigis") << partIdx <<" ";
377  const float partMeanTof = tempGeo->getCentralTOF(detId,partIdx);
378 
379  //convert to relative to partition
380  tof -= partMeanTof;
381 
382  //get coordinates and errors
383  const float partCenter = tempGeo->getPartCenter(partIdx);
384  const auto* topo = tempGeo->getTopo(partIdx);
385 
386  //find channel
387  const LocalPoint partLocalPoint(inDigi->x(), inDigi->y() - partCenter ,0.);
388  strip = topo->channel(partLocalPoint);
389  const float stripF = float(strip)+0.5;
390 
391  //get digitized location
392  LocalPoint digiPartLocalPoint = topo->localPosition(stripF);
393  digiLocalError = topo->localError(stripF, 1./sqrt(12.)); //std dev. flat distribution with length L is L/sqrt(12). The strip topology expects the error in units of strips.
394  digiLocalPoint = LocalPoint(digiPartLocalPoint.x(),digiPartLocalPoint.y() + partCenter,0.0);
395  return partIdx;
396 
397 
398 }
static unsigned int partIdx(const InputGenJetsParticleSelector::ParticleVector &p, const reco::Candidate *particle)
float y() const
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
float getCentralTOF(const ME0DetId &me0Id, unsigned int partIdx) const
T y() const
Definition: PV3DBase.h:63
T sqrt(T t)
Definition: SSEVec.h:18
#define LogTrace(id)
unsigned int findEtaPartition(float locY) const
float getPartCenter(const unsigned int partIdx) const
TemporaryGeometry * tempGeo
float x() const
T x() const
Definition: PV3DBase.h:62
const TrapezoidalStripTopology * getTopo(const unsigned int partIdx) const
void ME0ReDigiProducer::getStripProperties ( const ME0EtaPartition etaPart,
const ME0DigiPreReco inDigi,
float &  tof,
int &  strip,
LocalPoint digiLocalPoint,
LocalError digiLocalError 
) const
private

Definition at line 399 of file ME0ReDigiProducer.cc.

References objects.autophobj::float, ME0EtaPartition::id(), ME0DetId::layer(), TrapezoidalStripTopology::nstrips(), TrapezoidalStripTopology::pitch(), TrapezoidalStripTopology::radius(), ME0DetId::roll(), ME0EtaPartition::specificTopology(), mathSSE::sqrt(), TrapezoidalStripTopology::stripLength(), tofs, usePads, ME0DigiPreReco::x(), and ME0DigiPreReco::y().

Referenced by buildDigis().

399  {
400  //convert to relative to partition
401  tof -= tofs[etaPart->id().layer()-1][etaPart->id().roll() -1];
402 
403  const TrapezoidalStripTopology * origTopo = (const TrapezoidalStripTopology*)(&etaPart->specificTopology());
404  TrapezoidalStripTopology padTopo(origTopo->nstrips()/2,origTopo->pitch()*2,origTopo->stripLength(),origTopo->radius());
405  const auto & topo = usePads ? padTopo : etaPart->specificTopology();
406 
407  //find channel
408  const LocalPoint partLocalPoint(inDigi->x(), inDigi->y(),0.);
409  strip = topo.channel(partLocalPoint);
410  const float stripF = float(strip)+0.5;
411 
412  //get digitized location
413  digiLocalPoint = topo.localPosition(stripF);
414  digiLocalError = topo.localError(stripF, 1./sqrt(12.));
415 }
const StripTopology & specificTopology() const
float y() const
T sqrt(T t)
Definition: SSEVec.h:18
ME0DetId id() const
float stripLength() const override
det heigth (strip length in the middle)
float x() const
int roll() const
Definition: ME0DetId.h:62
std::vector< std::vector< double > > tofs
int layer() const
Layer id: each chamber has six layers of chambers: layer 1 is the inner layer and layer 6 is the oute...
Definition: ME0DetId.h:56
void ME0ReDigiProducer::produce ( edm::Event e,
const edm::EventSetup eventSetup 
)
override

Definition at line 222 of file ME0ReDigiProducer.cc.

References buildDigis(), edm::Event::getByToken(), edm::RandomNumberGenerator::getEngine(), eostools::move(), edm::Handle< T >::product(), edm::Event::put(), edm::Event::streamID(), and token.

223 {
225  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
226 
228  e.getByToken(token, input_digis);
229 
230  std::unique_ptr<ME0DigiPreRecoCollection> output_digis(new ME0DigiPreRecoCollection());
231  std::unique_ptr<ME0DigiPreRecoMap> output_digimap(new ME0DigiPreRecoMap());
232 
233  // build the digis
234  buildDigis(*(input_digis.product()),
235  *output_digis,
236  *output_digimap,
237  engine);
238 
239  // store them in the event
240  e.put(std::move(output_digis));
241  e.put(std::move(output_digimap));
242 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< ME0DigiPreRecoCollection > token
void buildDigis(const ME0DigiPreRecoCollection &, ME0DigiPreRecoCollection &, ME0DigiPreRecoMap &, CLHEP::HepRandomEngine *engine)
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
MuonDigiCollection< ME0DetId, int > ME0DigiPreRecoMap
MuonDigiCollection< ME0DetId, ME0DigiPreReco > ME0DigiPreRecoCollection
T const * product() const
Definition: Handle.h:74
StreamID streamID() const
Definition: Event.h:95
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

const float ME0ReDigiProducer::bxWidth
private

Definition at line 82 of file ME0ReDigiProducer.h.

Referenced by buildDigis().

const ME0Geometry* ME0ReDigiProducer::geometry
private
std::vector<int> ME0ReDigiProducer::layerReadout
private

Definition at line 91 of file ME0ReDigiProducer.h.

Referenced by beginRun(), and buildDigis().

int ME0ReDigiProducer::maxBXReadout
private

Definition at line 90 of file ME0ReDigiProducer.h.

Referenced by buildDigis().

bool ME0ReDigiProducer::mergeDigis
private

Definition at line 92 of file ME0ReDigiProducer.h.

Referenced by buildDigis().

int ME0ReDigiProducer::minBXReadout
private

Definition at line 89 of file ME0ReDigiProducer.h.

Referenced by buildDigis().

double ME0ReDigiProducer::neutronAcceptance
private

Definition at line 87 of file ME0ReDigiProducer.h.

Referenced by buildDigis(), and ME0ReDigiProducer().

unsigned int ME0ReDigiProducer::numberOfPartitions
private
unsigned int ME0ReDigiProducer::numberOfStrips
private
TemporaryGeometry* ME0ReDigiProducer::tempGeo
private
double ME0ReDigiProducer::timeResolution
private

Definition at line 88 of file ME0ReDigiProducer.h.

Referenced by buildDigis().

std::vector<std::vector<double> > ME0ReDigiProducer::tofs
private

Definition at line 98 of file ME0ReDigiProducer.h.

Referenced by fillCentralTOFs(), and getStripProperties().

edm::EDGetTokenT<ME0DigiPreRecoCollection> ME0ReDigiProducer::token
private

Definition at line 93 of file ME0ReDigiProducer.h.

Referenced by produce().

bool ME0ReDigiProducer::useBuiltinGeo
private

Definition at line 95 of file ME0ReDigiProducer.h.

Referenced by beginRun(), buildDigis(), and ME0ReDigiProducer().

bool ME0ReDigiProducer::useCusGeoFor1PartGeo
private

Definition at line 83 of file ME0ReDigiProducer.h.

Referenced by beginRun(), and ME0ReDigiProducer().

bool ME0ReDigiProducer::usePads
private

Definition at line 84 of file ME0ReDigiProducer.h.

Referenced by getStripProperties(), and ME0ReDigiProducer().