CMS 3D CMS Logo

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

#include <DTDigitizer.h>

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

Classes

struct  hitLessT
 

Public Member Functions

 DTDigitizer (const edm::ParameterSet &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Types

typedef std::map< DTWireId, std::vector< const PSimHit * > > DTWireIdMap
 
typedef DTWireIdMap::const_iterator DTWireIdMapConstIter
 
typedef DTWireIdMap::iterator DTWireIdMapIter
 
typedef std::pair< const PSimHit *, float > hitAndT
 
typedef std::vector< hitAndTTDContainer
 

Private Member Functions

float asymGausSmear (double mean, double sigmaLeft, double sigmaRight, CLHEP::HepRandomEngine *) const
 
std::pair< float, bool > computeTime (const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit, const LocalVector &BLoc, CLHEP::HepRandomEngine *)
 
std::pair< float, bool > driftTimeFromParametrization (float x, float alpha, float By, float Bz, CLHEP::HepRandomEngine *) const
 
std::pair< float, bool > driftTimeFromTimeMap () const
 
void dumpHit (const PSimHit *hit, float xEntry, float xExit, const DTTopology &topo)
 
float externalDelays (const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit) const
 
void storeDigis (DTWireId &wireId, TDContainer &hits, DTDigiCollection &output, DTDigiSimLinkCollection &outputLinks)
 

Private Attributes

int base
 
edm::EDGetTokenT< CrossingFrame< PSimHit > > cf_token
 
std::string collection_for_XF
 
float deadTime
 
bool debug
 
std::string geometryType
 
bool IdealModel
 
bool interpolate
 
float LinksTimeWindow
 
edm::ESGetToken< MagneticField, IdealMagneticFieldRecordmagnField_token
 
std::string mix_
 
bool MultipleLinks
 
edm::ESGetToken< DTGeometry, MuonGeometryRecordmuonGeom_token
 
bool onlyMuHits
 
float smearing
 
std::string syncName
 
float theConstVDrift
 
std::unique_ptr< DTDigiSyncBasetheSync
 
double vPropWire
 

Friends

class DTDigitizerAnalysis
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Digitize the muon drift tubes. The parametrisation function in DTDriftTimeParametrization from P.G.Abia, J.Puerta is used in all cases where it is applicable.

Authors
: G. Bevilacqua, N. Amapane, G. Cerminara, R. Bellan

Definition at line 52 of file DTDigitizer.h.

Member Typedef Documentation

◆ DTWireIdMap

typedef std::map<DTWireId, std::vector<const PSimHit *> > DTDigitizer::DTWireIdMap
private

Definition at line 62 of file DTDigitizer.h.

◆ DTWireIdMapConstIter

typedef DTWireIdMap::const_iterator DTDigitizer::DTWireIdMapConstIter
private

Definition at line 64 of file DTDigitizer.h.

◆ DTWireIdMapIter

typedef DTWireIdMap::iterator DTDigitizer::DTWireIdMapIter
private

Definition at line 63 of file DTDigitizer.h.

◆ hitAndT

typedef std::pair<const PSimHit *, float> DTDigitizer::hitAndT
private

Definition at line 59 of file DTDigitizer.h.

◆ TDContainer

typedef std::vector<hitAndT> DTDigitizer::TDContainer
private

Definition at line 60 of file DTDigitizer.h.

Constructor & Destructor Documentation

◆ DTDigitizer()

DTDigitizer::DTDigitizer ( const edm::ParameterSet conf_)
explicit

Definition at line 51 of file DTDigitizer.cc.

References get, and edm::ParameterSet::getParameter().

52  : // Sync Algo
53  theSync{DTDigiSyncFactory::get()->create(conf_.getParameter<string>("SyncName"),
54  conf_.getParameter<ParameterSet>("pset"))} {
55  // Set verbose output
56  debug = conf_.getUntrackedParameter<bool>("debug");
57 
58  if (debug)
59  LogPrint("DTDigitizer") << "Creating a DTDigitizer" << endl;
60 
61  // register the Producer with a label
62  // produces<DTDigiCollection>("MuonDTDigis"); // FIXME: Do I pass it by
63  // ParameterSet?
64  produces<DTDigiCollection>(); // FIXME: Do I pass it by ParameterSet?
65  // produces<DTDigiSimLinkCollection>("MuonDTDigiSimLinks");
66  produces<DTDigiSimLinkCollection>();
67 
68  // Parameters:
69 
70  // build digis only for mu hits (for debug purposes)
71  onlyMuHits = conf_.getParameter<bool>("onlyMuHits");
72 
73  // interpolate parametrization function
74  interpolate = conf_.getParameter<bool>("interpolate");
75 
76  // Velocity of signal propagation along the wire (cm/ns)
77  // For the default value
78  // cfr. CMS-IN 2000-021: (2.56+-0.17)x1e8 m/s
79  // CMS NOTE 2003-17: (0.244) m/ns
80  vPropWire = conf_.getParameter<double>("vPropWire"); // 24.4
81 
82  // Dead time for signals on the same wire (number from M. Pegoraro)
83  deadTime = conf_.getParameter<double>("deadTime"); // 150
84 
85  // further configurable smearing
86  smearing = conf_.getParameter<double>("Smearing"); // 3.
87 
88  // Debug flag to switch to the Ideal model
89  // it uses a constant drift velocity and doesn't set any external delay
90  IdealModel = conf_.getParameter<bool>("IdealModel");
91 
92  // Flag to specify that we want digis in phase-2 units (25./30. ns)
93  // instead that the old units (25./32.)
94  if (conf_.getParameter<bool>("phase2Digis"))
95  base = 30;
96  else
97  base = 32;
98 
99  // Constant drift velocity needed by the above flag
100  if (IdealModel)
101  theConstVDrift = conf_.getParameter<double>("IdealModelConstantDriftVelocity"); // 55 um/ns
102  else
103  theConstVDrift = 55.;
104 
105  // get random engine
107  if (!rng.isAvailable()) {
108  throw cms::Exception("Configuration") << "RandomNumberGeneratorService for DTDigitizer missing in cfg file";
109  }
110 
111  // MultipleLinks=false ==> one-to-one correspondence between digis and SimHits
112  MultipleLinks = conf_.getParameter<bool>("MultipleLinks");
113  // MultipleLinks=true ==> association of SimHits within a time window
114  // LinksTimeWindow (of the order of the resolution)
115  LinksTimeWindow = conf_.getParameter<double>("LinksTimeWindow"); // (10 ns)
116 
117  // Name of Collection used for create the XF
118  mix_ = conf_.getParameter<std::string>("mixLabel");
119  collection_for_XF = conf_.getParameter<std::string>("InputCollection");
120  cf_token = consumes<CrossingFrame<PSimHit>>(edm::InputTag(mix_, collection_for_XF));
121  magnField_token = esConsumes<MagneticField, IdealMagneticFieldRecord>();
122 
123  // String to choose between ideal (the default) and (mis)aligned geometry
124  // for the digitization step
125  geometryType = conf_.getParameter<std::string>("GeometryType");
126  muonGeom_token = esConsumes<DTGeometry, MuonGeometryRecord>(edm::ESInputTag("", geometryType));
127 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
bool onlyMuHits
Definition: DTDigitizer.h:116
float theConstVDrift
Definition: DTDigitizer.h:126
float LinksTimeWindow
Definition: DTDigitizer.h:130
T getUntrackedParameter(std::string const &, T const &) const
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magnField_token
Definition: DTDigitizer.h:138
edm::EDGetTokenT< CrossingFrame< PSimHit > > cf_token
Definition: DTDigitizer.h:136
float smearing
Definition: DTDigitizer.h:113
bool interpolate
Definition: DTDigitizer.h:115
bool IdealModel
Definition: DTDigitizer.h:125
Log< level::Warning, true > LogPrint
std::unique_ptr< DTDigiSyncBase > theSync
Definition: DTDigitizer.h:120
std::string mix_
Definition: DTDigitizer.h:133
std::string geometryType
Definition: DTDigitizer.h:122
double vPropWire
Definition: DTDigitizer.h:111
bool isAvailable() const
Definition: Service.h:40
#define get
float deadTime
Definition: DTDigitizer.h:112
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeom_token
Definition: DTDigitizer.h:137
bool MultipleLinks
Definition: DTDigitizer.h:129
std::string collection_for_XF
Definition: DTDigitizer.h:134

Member Function Documentation

◆ asymGausSmear()

float DTDigitizer::asymGausSmear ( double  mean,
double  sigmaLeft,
double  sigmaRight,
CLHEP::HepRandomEngine *  engine 
) const
private

Definition at line 443 of file DTDigitizer.cc.

References f, SiStripPI::mean, and submitPVValidationJobs::t.

Referenced by driftTimeFromParametrization().

446  {
447  double f = sigmaLeft / (sigmaLeft + sigmaRight);
448  double t;
449 
450  if (CLHEP::RandFlat::shoot(engine) <= f) {
451  t = CLHEP::RandGaussQ::shoot(engine, mean, sigmaLeft);
452  t = mean - fabs(t - mean);
453  } else {
454  t = CLHEP::RandGaussQ::shoot(engine, mean, sigmaRight);
455  t = mean + fabs(t - mean);
456  }
457  return static_cast<float>(t);
458 }
double f[11][100]

◆ computeTime()

pair< float, bool > DTDigitizer::computeTime ( const DTLayer layer,
const DTWireId wireId,
const PSimHit hit,
const LocalVector BLoc,
CLHEP::HepRandomEngine *  engine 
)
private

Definition at line 225 of file DTDigitizer.cc.

References Vector3DBase< T, FrameTag >::cross(), ztail::d, debug, dumpMFGeometry_cfg::delta, DeadROC_duringRun::dir, Vector3DBase< T, FrameTag >::dot(), driftTimeFromParametrization(), driftTimeFromTimeMap(), dumpHit(), externalDelays(), IdealModel, nano_mu_digi_cff::layer, M_PI, PV3DBase< T, PVType, FrameType >::mag(), DTTopology::none, DTTopology::onWhichBorder(), DiDispStaMuonMonitor_cfi::pt, pv::pT, mathSSE::sqrt(), theConstVDrift, theta(), Vector3DBase< T, FrameTag >::unit(), DTWireId::wire(), DTTopology::wirePosition(), x, PV3DBase< T, PVType, FrameType >::x(), DTTopology::xMax, DTTopology::xMin, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by produce().

229  {
230  LocalPoint entryP = hit->entryPoint();
231  LocalPoint exitP = hit->exitPoint();
232  int partType = hit->particleType();
233 
234  const DTTopology &topo = layer->specificTopology();
235 
236  // Pay attention: in CMSSW the rf of the SimHit is in the layer's rf
237 
238  if (debug)
239  LogPrint("DTDigitizer") << "Hit local entry point: " << entryP << endl << "Hit local exit point: " << exitP << endl;
240 
241  float xwire = topo.wirePosition(wireId.wire());
242  float xEntry = entryP.x() - xwire;
243  float xExit = exitP.x() - xwire;
244 
245  if (debug)
246  LogPrint("DTDigitizer") << "wire position: " << xwire << " x entry in cell rf: " << xEntry
247  << " x exit in cell rf: " << xExit << endl;
248 
249  DTTopology::Side entrySide = topo.onWhichBorder(xEntry, entryP.y(), entryP.z());
250  DTTopology::Side exitSide = topo.onWhichBorder(xExit, exitP.y(), exitP.z());
251 
252  if (debug)
253  dumpHit(hit, xEntry, xExit, topo);
254 
255  // The bolean is used to flag the drift time computation
256  pair<float, bool> driftTime(0., false);
257 
258  // if delta in gas->ignore, since it is included in the parametrisation.
259  // FIXME: should check that it is actually a delta ray produced by a nearby
260  // muon hit.
261 
262  if (partType == 11 && entrySide == DTTopology::none) {
263  if (debug)
264  LogPrint("DTDigitizer") << " e- hit in gas; discarding " << endl;
265  return driftTime;
266  }
267 
268  float By = BLoc.y();
269  float Bz = BLoc.z();
270 
271  // Radius and sagitta according to direction of momentum
272  // (just for printing)
273  // NOTE: in cmsim, d is always taken // pHat!
274  LocalVector d = (exitP - entryP);
275  LocalVector pHat = hit->localDirection().unit();
276  LocalVector hHat = (d.cross(pHat.cross(d))).unit();
277  float cosAlpha = hHat.dot(pHat);
278  float sinAlpha = sqrt(1. - cosAlpha * cosAlpha);
279  float radius_P = (d.mag()) / (2. * cosAlpha);
280  float sagitta_P = radius_P * (1. - sinAlpha);
281 
282  // Radius, sagitta according to field bending
283  // (just for printing)
284  float halfd = d.mag() / 2.;
285  float BMag = BLoc.mag();
286  LocalVector pT = (pHat - (BLoc.unit() * pHat.dot(BLoc.unit()))) * (hit->pabs());
287  float radius_B = (pT.mag() / (0.3 * BMag)) * 100.;
288  float sagitta_B;
289  if (radius_B > halfd) {
290  sagitta_B = radius_B - sqrt(radius_B * radius_B - halfd * halfd);
291  } else {
292  sagitta_B = radius_B;
293  }
294 
295  // cos(delta), delta= angle between direction at entry and hit segment
296  // (just for printing)
297  float delta = pHat.dot(d.unit());
298  if (debug)
299  LogPrint("DTDigitizer") << " delta = " << delta << endl
300  << " cosAlpha = " << cosAlpha << endl
301  << " sinAlpha = " << sinAlpha << endl
302  << " pMag = " << pT.mag() << endl
303  << " bMag = " << BMag << endl
304  << " pT = " << pT << endl
305  << " halfd = " << halfd << endl
306  << " radius_P (cm) = " << radius_P << endl
307  << " sagitta_P (um) = " << sagitta_P * 10000. << endl
308  << " radius_B (cm) = " << radius_B << endl
309  << " sagitta_B (um) = " << sagitta_B * 10000. << endl;
310 
311  // Select cases where parametrization can not be used.
312  bool noParametrisation = ((entrySide == DTTopology::none || exitSide == DTTopology::none) // case # 2,3,8,9 or 11
313  || (entrySide == exitSide) // case # 4 or 10
314  || ((entrySide == DTTopology::xMin && exitSide == DTTopology::xMax) ||
315  (entrySide == DTTopology::xMax && exitSide == DTTopology::xMin)) // Hit is case # 7
316  );
317 
318  // FIXME: now, debug warning only; consider treating those
319  // with TM algo.
320  if (delta < 0.99996 // Track is not straight. FIXME: use sagitta?
321  && (noParametrisation == false)) {
322  if (debug)
323  LogPrint("DTDigitizer") << "*** WARNING: hit is not straight, type = " << partType << endl;
324  }
325 
326  //************ 5A ***************
327 
328  if (!noParametrisation) {
329  LocalVector dir = hit->momentumAtEntry(); // ex Measurement3DVector dir =
330  // hit->measurementDirection(); //FIXME?
331  float theta = atan(dir.x() / -dir.z()) * 180 / M_PI;
332 
333  // FIXME: use dir if M.S. is included as GARFIELD option...
334  // otherwise use hit segment dirction.
335  // LocalVector dir0 = (exitP-entryP).unit();
336  // float theta = atan(dir0.x()/-dir0.z())*180/M_PI;
337  float x;
338 
339  Local3DPoint pt = hit->localPosition(); // ex Measurement3DPoint pt =
340  // hit->measurementPosition(); // FIXME?
341 
342  if (fabs(pt.z()) < 0.002) {
343  // hit center within 20 um from z=0, no need to extrapolate.
344  x = pt.x() - xwire;
345  } else {
346  x = xEntry - (entryP.z() * (xExit - xEntry)) / (exitP.z() - entryP.z());
347  }
348 
349  if (IdealModel)
350  return make_pair(fabs(x) / theConstVDrift, true);
351  else
352  driftTime = driftTimeFromParametrization(x, theta, By, Bz, engine);
353  }
354 
355  if ((driftTime.second) == false) {
356  // Parametrisation not applicable, or failed. Use time map.
357  driftTime = driftTimeFromTimeMap();
358  }
359 
360  //************ 5B ***************
361 
362  // Signal propagation, TOF etc.
363  if (driftTime.second) {
364  driftTime.first += externalDelays(layer, wireId, hit);
365  }
366  return driftTime;
367 }
void dumpHit(const PSimHit *hit, float xEntry, float xExit, const DTTopology &topo)
Definition: DTDigitizer.cc:561
std::pair< float, bool > driftTimeFromParametrization(float x, float alpha, float By, float Bz, CLHEP::HepRandomEngine *) const
Definition: DTDigitizer.cc:371
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:110
int wire() const
Return the wire number.
Definition: DTWireId.h:45
T z() const
Definition: PV3DBase.h:61
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
float theConstVDrift
Definition: DTDigitizer.h:126
Side onWhichBorder(float x, float y, float z) const
Definition: DTTopology.cc:82
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T sqrt(T t)
Definition: SSEVec.h:23
T mag() const
Definition: PV3DBase.h:64
Side
Sides of the cell.
Definition: DTTopology.h:89
bool IdealModel
Definition: DTDigitizer.h:125
Log< level::Warning, true > LogPrint
d
Definition: ztail.py:151
#define M_PI
std::pair< float, bool > driftTimeFromTimeMap() const
Definition: DTDigitizer.cc:460
float externalDelays(const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit) const
Definition: DTDigitizer.cc:469
Vector3DBase unit() const
Definition: Vector3DBase.h:54
Geom::Theta< T > theta() const
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:63

◆ driftTimeFromParametrization()

pair< float, bool > DTDigitizer::driftTimeFromParametrization ( float  x,
float  alpha,
float  By,
float  Bz,
CLHEP::HepRandomEngine *  engine 
) const
private

Definition at line 371 of file DTDigitizer.cc.

References asymGausSmear(), debug, GeomDetEnumerators::DT, f, RemoveAddSevLevel::flag, interpolate, SiStripPI::max, DTDriftTimeParametrization::MB_DT_drift_time(), smearing, theta(), hcalRecHitTable_cff::time, and x.

Referenced by computeTime().

372  {
373  // Convert from CMSSW frame/units r.f. to parametrization ones.
374  x *= 10.; // cm -> mm
375 
376  // FIXME: Current parametrisation can extrapolate above 21 mm,
377  // however a detailed study is needed before using this.
378  if (fabs(x) > 21.) {
379  if (debug)
380  LogPrint("DTDigitizer") << "*** WARNING: parametrisation: x out of range = " << x << ", skipping" << endl;
381  return pair<float, bool>(0.f, false);
382  }
383 
384  // Different r.f. of the parametrization:
385  // X_par = X_ORCA; Y_par=Z_ORCA; Z_par = -Y_ORCA
386 
387  float By_par = Bz; // Bnorm
388  float Bz_par = -By; // Bwire
389  float theta_par = theta;
390 
391  // Parametrisation uses interpolation up to |theta|=45 deg,
392  // |Bwire|=0.4, |Bnorm|=0.75; extrapolation above.
393  if (fabs(theta_par) > 45.) {
394  if (debug)
395  LogPrint("DTDigitizer") << "*** WARNING: extrapolating theta > 45: " << theta << endl;
396  // theta_par = min(fabs(theta_par),45.f)*((theta_par<0.)?-1.:1.);
397  }
398  if (fabs(By_par) > 0.75) {
399  if (debug)
400  LogPrint("DTDigitizer") << "*** WARNING: extrapolating Bnorm > 0.75: " << By_par << endl;
401  // By_par = min(fabs(By_par),0.75f)*((By_par<0.)?-1.:1.);
402  }
403  if (fabs(Bz_par) > 0.4) {
404  if (debug)
405  LogPrint("DTDigitizer") << "*** WARNING: extrapolating Bwire >0.4: " << Bz_par << endl;
406  // Bz_par = min(fabs(Bz_par),0.4)*((Bz_par<0.)?-1.:1.);
407  }
408 
410  static const DTDriftTimeParametrization par;
411  unsigned short flag = par.MB_DT_drift_time(x, theta_par, By_par, Bz_par, 0, &DT, interpolate);
412 
413  if (debug) {
414  LogPrint("DTDigitizer") << " Parametrisation: x, theta, Bnorm, Bwire = " << x << " " << theta_par << " "
415  << By_par << " " << Bz_par << endl
416  << " time=" << DT.t_drift << " sigma_m=" << DT.t_width_m << " sigma_p=" << DT.t_width_p
417  << endl;
418  if (flag != 1) {
419  LogPrint("DTDigitizer") << "*** WARNING: call to parametrisation failed" << endl;
420  return pair<float, bool>(0.f, false);
421  }
422  }
423 
424  // Double half-gaussian smearing
425  float time = asymGausSmear(DT.t_drift, DT.t_width_m, DT.t_width_p, engine);
426 
427  // Do not allow the smearing to lead to negative values
428  time = max(time, 0.f);
429 
430  // Apply a Gaussian smearing to account for electronic effects (cf. 2004 TB
431  // analysis) The width of the Gaussian can be configured with the "Smearing"
432  // parameter
433 
434  double u = CLHEP::RandGaussQ::shoot(engine, 0., smearing);
435  time += u;
436 
437  if (debug)
438  LogPrint("DTDigitizer") << " drift time = " << time << endl;
439 
440  return pair<float, bool>(time, true);
441 }
unsigned short MB_DT_drift_time(double x, double alpha, double by, double bz, short ifl, drift_time *DT, short interpolate) const
Calculate drift time and spread.
float asymGausSmear(double mean, double sigmaLeft, double sigmaRight, CLHEP::HepRandomEngine *) const
Definition: DTDigitizer.cc:443
double f[11][100]
float smearing
Definition: DTDigitizer.h:113
bool interpolate
Definition: DTDigitizer.h:115
Log< level::Warning, true > LogPrint
Structure used to return output values.
Geom::Theta< T > theta() const

◆ driftTimeFromTimeMap()

pair< float, bool > DTDigitizer::driftTimeFromTimeMap ( ) const
private

Definition at line 460 of file DTDigitizer.cc.

References debug.

Referenced by computeTime().

460  {
461  // FIXME: not yet implemented.
462  if (debug)
463  LogPrint("DTDigitizer") << " TimeMap " << endl;
464  return pair<float, bool>(0., false);
465 }
Log< level::Warning, true > LogPrint

◆ dumpHit()

void DTDigitizer::dumpHit ( const PSimHit hit,
float  xEntry,
float  xExit,
const DTTopology topo 
)
private

Definition at line 561 of file DTDigitizer.cc.

References DTTopology::cellHeight(), DTTopology::cellLenght(), DTTopology::cellWidth(), createfilelist::int, mag(), DTTopology::onWhichBorder(), x, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by computeTime().

561  {
562  LocalPoint entryP = hit->entryPoint();
563  LocalPoint exitP = hit->exitPoint();
564 
565  DTTopology::Side entrySide = topo.onWhichBorder(xEntry, entryP.y(), entryP.z());
566  DTTopology::Side exitSide = topo.onWhichBorder(xExit, exitP.y(), exitP.z());
567  // ProcessTypeEnumerator pTypes;
568 
569  LogPrint("DTDigitizer") << endl
570  << "------- SimHit: " << endl
571  << " Particle type = " << hit->particleType() << endl
572  << " process type = " << hit->processType() << endl
573  << " process type = " << hit->processType()
574  << endl
575  // << " packedTrackId = " << hit->packedTrackId() << endl
576  << " trackId = " << hit->trackId()
577  << endl // new,is the same as the
578  // previous?? FIXME-Check
579  << " |p| = " << hit->pabs() << endl
580  << " Energy loss = " << hit->energyLoss()
581  << endl
582  // << " timeOffset = " << hit->timeOffset() << endl
583  // << " measurementPosition = " << hit->measurementPosition() << endl
584  // << " measurementDirection = " << hit->measurementDirection() << endl
585  // //FIXME
586  << " localDirection = " << hit->momentumAtEntry().unit()
587  << endl // FIXME is it a versor?
588  << " Entry point = " << entryP << " cell x = " << xEntry << endl
589  << " Exit point = " << exitP << " cell x = " << xExit << endl
590  << " DR = = " << (exitP - entryP).mag() << endl
591  << " Dx = = " << (exitP - entryP).x() << endl
592  << " Cell w,h,l = (" << topo.cellWidth() << " , " << topo.cellHeight() << " , "
593  << topo.cellLenght() << ") cm" << endl
594  << " DY entry from edge = " << topo.cellLenght() / 2. - fabs(entryP.y())
595  << " DY exit from edge = " << topo.cellLenght() / 2. - fabs(exitP.y())
596  << " entrySide = " << (int)entrySide << " ; exitSide = " << (int)exitSide << endl;
597 }
float cellLenght() const
Definition: DTTopology.h:74
T z() const
Definition: PV3DBase.h:61
Side onWhichBorder(float x, float y, float z) const
Definition: DTTopology.cc:82
T y() const
Definition: PV3DBase.h:60
float cellWidth() const
Returns the cell width.
Definition: DTTopology.h:69
Side
Sides of the cell.
Definition: DTTopology.h:89
Log< level::Warning, true > LogPrint
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float cellHeight() const
Returns the cell height.
Definition: DTTopology.h:71

◆ externalDelays()

float DTDigitizer::externalDelays ( const DTLayer layer,
const DTWireId wireId,
const PSimHit hit 
) const
private

Definition at line 469 of file DTDigitizer.cc.

References debug, nano_mu_digi_cff::layer, theSync, vPropWire, and hit::y.

Referenced by computeTime().

469  {
470  // Time of signal propagation along wire.
471 
472  float wireCoord = hit->localPosition().y();
473  float halfL = (layer->specificTopology().cellLenght()) / 2.;
474  float propgL = halfL - wireCoord; // the FE is always located at the pos coord.
475 
476  float propDelay = propgL / vPropWire;
477 
478  // Real TOF.
479  float tof = hit->tof();
480 
481  // Delays and t0 according to theSync
482 
483  double sync = theSync->digitizerOffset(&wireId, layer);
484 
485  if (debug) {
486  LogPrint("DTDigitizer") << " propDelay =" << propDelay << "; TOF=" << tof << "; sync= " << sync << endl;
487  }
488 
489  return propDelay + tof + sync;
490 }
Log< level::Warning, true > LogPrint
std::unique_ptr< DTDigiSyncBase > theSync
Definition: DTDigitizer.h:120
double vPropWire
Definition: DTDigitizer.h:111

◆ produce()

void DTDigitizer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 130 of file DTDigitizer.cc.

References cf_token, computeTime(), debug, edm::RandomNumberGenerator::getEngine(), edm::EventSetup::getHandle(), iEvent, MagneticField::inTesla(), nano_mu_digi_cff::layer, DTGeometry::layer(), DTWireId::layerId(), magnField_token, eostools::move(), muonGeom_token, convertSQLitetoXML_cfg::output, edm::Handle< T >::product(), rpcPointValidation_cfi::simHit, FastTrackerRecHitCombiner_cfi::simHits, storeDigis(), hcalRecHitTable_cff::time, and nano_mu_digi_cff::wire.

130  {
132  CLHEP::HepRandomEngine *engine = &rng->getEngine(iEvent.streamID());
133 
134  if (debug)
135  LogPrint("DTDigitizer") << "--- Run: " << iEvent.id().run() << " Event: " << iEvent.id().event() << endl;
136 
137  //************ 1 ***************
138  // create the container for the SimHits
139  // Handle<PSimHitContainer> simHits;
140  // iEvent.getByLabel("g4SimHits","MuonDTHits",simHits);
141 
142  // use MixCollection instead of the previous
144  iEvent.getByToken(cf_token, xFrame);
145 
146  unique_ptr<MixCollection<PSimHit>> simHits(new MixCollection<PSimHit>(xFrame.product()));
147 
148  // create the pointer to the Digi container
149  unique_ptr<DTDigiCollection> output(new DTDigiCollection());
150  // pointer to the DigiSimLink container
151  unique_ptr<DTDigiSimLinkCollection> outputLinks(new DTDigiSimLinkCollection());
152 
153  // Muon Geometry
154  ESHandle<DTGeometry> muonGeom = iSetup.getHandle(muonGeom_token);
155 
156  // Magnetic Field
158 
159  //************ 2 ***************
160 
161  // These are sorted by DetId, i.e. by layer and then by wire #
162  // map<DTDetId, vector<const PSimHit*> > wireMap;
163  DTWireIdMap wireMap;
164 
165  for (MixCollection<PSimHit>::MixItr simHit = simHits->begin(); simHit != simHits->end(); simHit++) {
166  // Create the id of the wire, the simHits in the DT known also the wireId
167 
168  DTWireId wireId((*simHit).detUnitId());
169  // Fill the map
170  wireMap[wireId].push_back(&(*simHit));
171  }
172 
173  pair<float, bool> time(0., false);
174 
175  //************ 3 ***************
176  // Loop over the wires
177  for (DTWireIdMapConstIter wire = wireMap.begin(); wire != wireMap.end(); wire++) {
178  // SimHit Container associated to the wire
179  const vector<const PSimHit *> &vhit = (*wire).second;
180  if (!vhit.empty()) {
181  TDContainer tdCont; // It is a vector<pair<const PSimHit*,float> >;
182 
183  //************ 4 ***************
184  DTWireId wireId = (*wire).first;
185 
186  // const DTLayer* layer = dynamic_cast< const DTLayer* >
187  // (muonGeom->idToDet(wireId.layerId()));
188  const DTLayer *layer = muonGeom->layer(wireId.layerId());
189 
190  // Loop on the hits of this wire
191  for (vector<const PSimHit *>::const_iterator hit = vhit.begin(); hit != vhit.end(); hit++) {
192  //************ 5 ***************
193  LocalPoint locPos = (*hit)->localPosition();
194 
195  const LocalVector BLoc = layer->surface().toLocal(magnField->inTesla(layer->surface().toGlobal(locPos)));
196 
197  time = computeTime(layer, wireId, *hit, BLoc, engine);
198 
199  //************ 6 ***************
200  if (time.second) {
201  tdCont.push_back(make_pair((*hit), time.first));
202  } else {
203  if (debug)
204  LogPrint("DTDigitizer") << "hit discarded" << endl;
205  }
206  }
207 
208  //************ 7 ***************
209 
210  // the loading must be done by layer but
211  // the digitization must be done by wire (in order to take into account
212  // the dead time)
213 
214  storeDigis(wireId, tdCont, *output, *outputLinks);
215  }
216  }
217 
218  //************ 8 ***************
219  // Load the Digi Container in the Event
220  // iEvent.put(std::move(output),"MuonDTDigis");
221  iEvent.put(std::move(output));
222  iEvent.put(std::move(outputLinks));
223 }
std::pair< float, bool > computeTime(const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit, const LocalVector &BLoc, CLHEP::HepRandomEngine *)
Definition: DTDigitizer.cc:225
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
T const * product() const
Definition: Handle.h:70
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
int iEvent
Definition: GenABIO.cc:224
void storeDigis(DTWireId &wireId, TDContainer &hits, DTDigiCollection &output, DTDigiSimLinkCollection &outputLinks)
Definition: DTDigitizer.cc:494
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magnField_token
Definition: DTDigitizer.h:138
edm::EDGetTokenT< CrossingFrame< PSimHit > > cf_token
Definition: DTDigitizer.h:136
DTWireIdMap::const_iterator DTWireIdMapConstIter
Definition: DTDigitizer.h:64
Log< level::Warning, true > LogPrint
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
std::vector< hitAndT > TDContainer
Definition: DTDigitizer.h:60
MuonDigiCollection< DTLayerId, DTDigi > DTDigiCollection
Definition: output.py:1
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:48
MuonDigiCollection< DTLayerId, DTDigiSimLink > DTDigiSimLinkCollection
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeom_token
Definition: DTDigitizer.h:137
std::map< DTWireId, std::vector< const PSimHit * > > DTWireIdMap
Definition: DTDigitizer.h:62
def move(src, dest)
Definition: eostools.py:511
const DTLayer * layer(const DTLayerId &id) const
Return a layer given its id.
Definition: DTGeometry.cc:96

◆ storeDigis()

void DTDigitizer::storeDigis ( DTWireId wireId,
TDContainer hits,
DTDigiCollection output,
DTDigiSimLinkCollection outputLinks 
)
private

Definition at line 494 of file DTDigitizer.cc.

References funct::abs(), base, deadTime, debug, hfClusterShapes_cfi::hits, MuonDigiCollection< IndexType, DigiType >::insertDigi(), DTWireId::layerId(), LinksTimeWindow, MultipleLinks, DTDigi::number(), onlyMuHits, jetUpdater_cfi::sort, hcalRecHitTable_cff::time, DTDigi::time(), and DTWireId::wire().

Referenced by produce().

497  {
498  //************ 7A ***************
499 
500  // sort signal times
501  sort(hits.begin(), hits.end(), hitLessT());
502 
503  //************ 7B ***************
504 
505  float wakeTime = -999999.0;
506  float resolTime = -999999.0;
507  int digiN = -1; // Digi identifier within the cell (for multiple digis)
508  DTDigi digi;
509 
510  // loop over signal times and drop signals inside dead time
511  for (TDContainer::const_iterator hit = hits.begin(); hit != hits.end(); hit++) {
512  if (onlyMuHits && abs((*hit).first->particleType()) != 13)
513  continue;
514 
515  //************ 7C ***************
516 
517  float time = (*hit).second;
518  if (time > wakeTime) {
519  // Note that digi is constructed with a float value (in ns)
520  int wireN = wireId.wire();
521  digiN++;
522  digi = DTDigi(wireN, time, digiN, base);
523 
524  // Add association between THIS digi and the corresponding SimTrack
525  unsigned int SimTrackId = (*hit).first->trackId();
526  EncodedEventId evId = (*hit).first->eventId();
527  DTDigiSimLink digisimLink(wireN, digiN, time, SimTrackId, evId, base);
528 
529  if (debug) {
530  LogPrint("DTDigitizer") << endl << "---- DTDigitizer ----" << endl;
531  LogPrint("DTDigitizer") << "wireId: " << wireId << endl;
532  LogPrint("DTDigitizer") << "sim. time = " << time << endl;
533  LogPrint("DTDigitizer") << "digi number = " << digi.number() << ", digi time = " << digi.time()
534  << ", linked to SimTrack Id = " << SimTrackId << endl;
535  }
536 
537  //************ 7D ***************
538  DTLayerId layerID = wireId.layerId(); // taking the layer of the wire
539  output.insertDigi(layerID, digi); // ordering Digis by layer
540  outputLinks.insertDigi(layerID, digisimLink);
541  wakeTime = time + deadTime;
542  resolTime = time + LinksTimeWindow;
543  } else if (MultipleLinks && time < resolTime) {
544  int wireN = wireId.wire();
545  unsigned int SimTrackId = (*hit).first->trackId();
546  EncodedEventId evId = (*hit).first->eventId();
547  DTDigiSimLink digisimLink(wireN, digiN, time, SimTrackId, evId, base);
548  DTLayerId layerID = wireId.layerId();
549  outputLinks.insertDigi(layerID, digisimLink);
550 
551  if (debug) {
552  LogPrint("DTDigitizer") << "\nAdded multiple link: \n"
553  << "digi number = " << digi.number() << ", digi time = " << digi.time()
554  << " (sim. time = " << time << ")"
555  << ", linked to SimTrack Id = " << SimTrackId << endl;
556  }
557  }
558  }
559 }
bool onlyMuHits
Definition: DTDigitizer.h:116
int wire() const
Return the wire number.
Definition: DTWireId.h:45
float LinksTimeWindow
Definition: DTDigitizer.h:130
void insertDigi(const IndexType &index, const DigiType &digi)
insert a digi for a given DetUnit
double time() const
Get time in ns.
Definition: DTDigi.cc:37
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: DTDigi.h:17
Log< level::Warning, true > LogPrint
int number() const
Identifies different digis within the same cell.
Definition: DTDigi.cc:43
Definition: output.py:1
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:48
float deadTime
Definition: DTDigitizer.h:112
bool MultipleLinks
Definition: DTDigitizer.h:129

Friends And Related Function Documentation

◆ DTDigitizerAnalysis

friend class DTDigitizerAnalysis
friend

Definition at line 108 of file DTDigitizer.h.

Member Data Documentation

◆ base

int DTDigitizer::base
private

Definition at line 117 of file DTDigitizer.h.

Referenced by storeDigis().

◆ cf_token

edm::EDGetTokenT<CrossingFrame<PSimHit> > DTDigitizer::cf_token
private

Definition at line 136 of file DTDigitizer.h.

Referenced by produce().

◆ collection_for_XF

std::string DTDigitizer::collection_for_XF
private

Definition at line 134 of file DTDigitizer.h.

◆ deadTime

float DTDigitizer::deadTime
private

Definition at line 112 of file DTDigitizer.h.

Referenced by storeDigis().

◆ debug

bool DTDigitizer::debug
private

◆ geometryType

std::string DTDigitizer::geometryType
private

Definition at line 122 of file DTDigitizer.h.

◆ IdealModel

bool DTDigitizer::IdealModel
private

Definition at line 125 of file DTDigitizer.h.

Referenced by computeTime().

◆ interpolate

bool DTDigitizer::interpolate
private

Definition at line 115 of file DTDigitizer.h.

Referenced by driftTimeFromParametrization().

◆ LinksTimeWindow

float DTDigitizer::LinksTimeWindow
private

Definition at line 130 of file DTDigitizer.h.

Referenced by storeDigis().

◆ magnField_token

edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> DTDigitizer::magnField_token
private

Definition at line 138 of file DTDigitizer.h.

Referenced by produce().

◆ mix_

std::string DTDigitizer::mix_
private

Definition at line 133 of file DTDigitizer.h.

◆ MultipleLinks

bool DTDigitizer::MultipleLinks
private

Definition at line 129 of file DTDigitizer.h.

Referenced by storeDigis().

◆ muonGeom_token

edm::ESGetToken<DTGeometry, MuonGeometryRecord> DTDigitizer::muonGeom_token
private

Definition at line 137 of file DTDigitizer.h.

Referenced by produce().

◆ onlyMuHits

bool DTDigitizer::onlyMuHits
private

Definition at line 116 of file DTDigitizer.h.

Referenced by storeDigis().

◆ smearing

float DTDigitizer::smearing
private

Definition at line 113 of file DTDigitizer.h.

Referenced by driftTimeFromParametrization().

◆ syncName

std::string DTDigitizer::syncName
private

Definition at line 119 of file DTDigitizer.h.

◆ theConstVDrift

float DTDigitizer::theConstVDrift
private

Definition at line 126 of file DTDigitizer.h.

Referenced by computeTime().

◆ theSync

std::unique_ptr<DTDigiSyncBase> DTDigitizer::theSync
private

Definition at line 120 of file DTDigitizer.h.

Referenced by externalDelays().

◆ vPropWire

double DTDigitizer::vPropWire
private

Definition at line 111 of file DTDigitizer.h.

Referenced by externalDelays().