CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

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
 
std::string mix_
 
bool MultipleLinks
 
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<>
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

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 46 of file DTDigitizer.h.

Member Typedef Documentation

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

Definition at line 56 of file DTDigitizer.h.

typedef DTWireIdMap::const_iterator DTDigitizer::DTWireIdMapConstIter
private

Definition at line 58 of file DTDigitizer.h.

typedef DTWireIdMap::iterator DTDigitizer::DTWireIdMapIter
private

Definition at line 57 of file DTDigitizer.h.

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

Definition at line 53 of file DTDigitizer.h.

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

Definition at line 54 of file DTDigitizer.h.

Constructor & Destructor Documentation

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

Definition at line 57 of file DTDigitizer.cc.

References base, cf_token, collection_for_XF, deadTime, debug, Exception, geometryType, timingPdfMaker::get, edm::ParameterSet::getParameter(), IdealModel, HLT_2018_cff::InputTag, interpolate, edm::Service< T >::isAvailable(), LinksTimeWindow, mix_, MultipleLinks, onlyMuHits, smearing, AlCaHLTBitMon_QueryRunRegistry::string, theConstVDrift, and vPropWire.

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

Member Function Documentation

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

Definition at line 449 of file DTDigitizer.cc.

References f, and OrderedSet::t.

Referenced by driftTimeFromParametrization().

452  {
453  double f = sigmaLeft / (sigmaLeft + sigmaRight);
454  double t;
455 
456  if (CLHEP::RandFlat::shoot(engine) <= f) {
457  t = CLHEP::RandGaussQ::shoot(engine, mean, sigmaLeft);
458  t = mean - fabs(t - mean);
459  } else {
460  t = CLHEP::RandGaussQ::shoot(engine, mean, sigmaRight);
461  t = mean + fabs(t - mean);
462  }
463  return static_cast<float>(t);
464 }
double f[11][100]
pair< float, bool > DTDigitizer::computeTime ( const DTLayer layer,
const DTWireId wireId,
const PSimHit hit,
const LocalVector BLoc,
CLHEP::HepRandomEngine *  engine 
)
private

Definition at line 231 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(), PSimHit::entryPoint(), PSimHit::exitPoint(), externalDelays(), IdealModel, PSimHit::localDirection(), PSimHit::localPosition(), M_PI, PV3DBase< T, PVType, FrameType >::mag(), PSimHit::momentumAtEntry(), DTTopology::none, DTTopology::onWhichBorder(), PSimHit::pabs(), PSimHit::particleType(), DiDispStaMuonMonitor_cfi::pt, PVValHelper::pT, DTLayer::specificTopology(), 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().

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

Definition at line 377 of file DTDigitizer.cc.

References asymGausSmear(), debug, GeomDetEnumerators::DT, f, RemoveAddSevLevel::flag, interpolate, SiStripPI::max, DTDriftTimeParametrization::MB_DT_drift_time(), smearing, DTDriftTimeParametrization::drift_time::t_drift, DTDriftTimeParametrization::drift_time::t_width_m, DTDriftTimeParametrization::drift_time::t_width_p, theta(), and ntuplemaker::time.

Referenced by computeTime().

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

Definition at line 466 of file DTDigitizer.cc.

References debug.

Referenced by computeTime().

466  {
467  // FIXME: not yet implemented.
468  if (debug)
469  LogPrint("DTDigitizer") << " TimeMap " << endl;
470  return pair<float, bool>(0., false);
471 }
void DTDigitizer::dumpHit ( const PSimHit hit,
float  xEntry,
float  xExit,
const DTTopology topo 
)
private

Definition at line 567 of file DTDigitizer.cc.

References DTTopology::cellHeight(), DTTopology::cellLenght(), DTTopology::cellWidth(), PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), createfilelist::int, mag(), PSimHit::momentumAtEntry(), DTTopology::onWhichBorder(), PSimHit::pabs(), PSimHit::particleType(), PSimHit::processType(), PSimHit::trackId(), Vector3DBase< T, FrameTag >::unit(), x, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by computeTime().

567  {
568  LocalPoint entryP = hit->entryPoint();
569  LocalPoint exitP = hit->exitPoint();
570 
571  DTTopology::Side entrySide = topo.onWhichBorder(xEntry, entryP.y(), entryP.z());
572  DTTopology::Side exitSide = topo.onWhichBorder(xExit, exitP.y(), exitP.z());
573  // ProcessTypeEnumerator pTypes;
574 
575  LogPrint("DTDigitizer") << endl
576  << "------- SimHit: " << endl
577  << " Particle type = " << hit->particleType() << endl
578  << " process type = " << hit->processType() << endl
579  << " process type = " << hit->processType()
580  << endl
581  // << " packedTrackId = " << hit->packedTrackId() << endl
582  << " trackId = " << hit->trackId()
583  << endl // new,is the same as the
584  // previous?? FIXME-Check
585  << " |p| = " << hit->pabs() << endl
586  << " Energy loss = " << hit->energyLoss()
587  << endl
588  // << " timeOffset = " << hit->timeOffset() << endl
589  // << " measurementPosition = " << hit->measurementPosition() << endl
590  // << " measurementDirection = " << hit->measurementDirection() << endl
591  // //FIXME
592  << " localDirection = " << hit->momentumAtEntry().unit()
593  << endl // FIXME is it a versor?
594  << " Entry point = " << entryP << " cell x = " << xEntry << endl
595  << " Exit point = " << exitP << " cell x = " << xExit << endl
596  << " DR = = " << (exitP - entryP).mag() << endl
597  << " Dx = = " << (exitP - entryP).x() << endl
598  << " Cell w,h,l = (" << topo.cellWidth() << " , " << topo.cellHeight() << " , "
599  << topo.cellLenght() << ") cm" << endl
600  << " DY entry from edge = " << topo.cellLenght() / 2. - fabs(entryP.y())
601  << " DY exit from edge = " << topo.cellLenght() / 2. - fabs(exitP.y())
602  << " entrySide = " << (int)entrySide << " ; exitSide = " << (int)exitSide << endl;
603 }
Side onWhichBorder(float x, float y, float z) const
Definition: DTTopology.cc:78
LocalVector momentumAtEntry() const
The momentum of the track that produced the hit, at entry point.
Definition: PSimHit.h:55
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
T y() const
Definition: PV3DBase.h:60
float cellWidth() const
Returns the cell width.
Definition: DTTopology.h:69
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
T z() const
Definition: PV3DBase.h:61
float cellHeight() const
Returns the cell height.
Definition: DTTopology.h:71
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
Side
Sides of the cell.
Definition: DTTopology.h:89
Vector3DBase unit() const
Definition: Vector3DBase.h:54
unsigned short processType() const
Definition: PSimHit.h:120
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:79
int particleType() const
Definition: PSimHit.h:89
unsigned int trackId() const
Definition: PSimHit.h:106
float cellLenght() const
Definition: DTTopology.h:74
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
float DTDigitizer::externalDelays ( const DTLayer layer,
const DTWireId wireId,
const PSimHit hit 
) const
private

Definition at line 475 of file DTDigitizer.cc.

References DTTopology::cellLenght(), debug, PSimHit::localPosition(), DTLayer::specificTopology(), theSync, PSimHit::tof(), vPropWire, and PV3DBase< T, PVType, FrameType >::y().

Referenced by computeTime().

475  {
476  // Time of signal propagation along wire.
477 
478  float wireCoord = hit->localPosition().y();
479  float halfL = (layer->specificTopology().cellLenght()) / 2.;
480  float propgL = halfL - wireCoord; // the FE is always located at the pos coord.
481 
482  float propDelay = propgL / vPropWire;
483 
484  // Real TOF.
485  float tof = hit->tof();
486 
487  // Delays and t0 according to theSync
488 
489  double sync = theSync->digitizerOffset(&wireId, layer);
490 
491  if (debug) {
492  LogPrint("DTDigitizer") << " propDelay =" << propDelay << "; TOF=" << tof << "; sync= " << sync << endl;
493  }
494 
495  return propDelay + tof + sync;
496 }
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:76
T y() const
Definition: PV3DBase.h:60
const DTTopology & specificTopology() const
Definition: DTLayer.cc:37
Local3DPoint localPosition() const
Definition: PSimHit.h:52
std::unique_ptr< DTDigiSyncBase > theSync
Definition: DTDigitizer.h:114
double vPropWire
Definition: DTDigitizer.h:105
float cellLenght() const
Definition: DTTopology.h:74
void DTDigitizer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 134 of file DTDigitizer.cc.

References cf_token, computeTime(), debug, edm::EventID::event(), geometryType, edm::EventSetup::get(), edm::Event::getByToken(), edm::RandomNumberGenerator::getEngine(), edm::EventBase::id(), MagneticField::inTesla(), DTGeometry::layer(), DTWireId::layerId(), eostools::move(), convertSQLitetoXML_cfg::output, edm::Handle< T >::product(), edm::Event::put(), edm::EventID::run(), rpcPointValidation_cfi::simHit, FastTrackerRecHitCombiner_cfi::simHits, storeDigis(), edm::Event::streamID(), GeomDet::surface(), ntuplemaker::time, Surface::toGlobal(), and GloballyPositioned< T >::toLocal().

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

Definition at line 500 of file DTDigitizer.cc.

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

Referenced by produce().

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

Friends And Related Function Documentation

friend class DTDigitizerAnalysis
friend

Definition at line 102 of file DTDigitizer.h.

Member Data Documentation

int DTDigitizer::base
private

Definition at line 111 of file DTDigitizer.h.

Referenced by DTDigitizer(), and storeDigis().

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

Definition at line 130 of file DTDigitizer.h.

Referenced by DTDigitizer(), and produce().

std::string DTDigitizer::collection_for_XF
private

Definition at line 128 of file DTDigitizer.h.

Referenced by DTDigitizer().

float DTDigitizer::deadTime
private

Definition at line 106 of file DTDigitizer.h.

Referenced by DTDigitizer(), and storeDigis().

bool DTDigitizer::debug
private
std::string DTDigitizer::geometryType
private

Definition at line 116 of file DTDigitizer.h.

Referenced by DTDigitizer(), and produce().

bool DTDigitizer::IdealModel
private

Definition at line 119 of file DTDigitizer.h.

Referenced by computeTime(), and DTDigitizer().

bool DTDigitizer::interpolate
private

Definition at line 109 of file DTDigitizer.h.

Referenced by driftTimeFromParametrization(), and DTDigitizer().

float DTDigitizer::LinksTimeWindow
private

Definition at line 124 of file DTDigitizer.h.

Referenced by DTDigitizer(), and storeDigis().

std::string DTDigitizer::mix_
private

Definition at line 127 of file DTDigitizer.h.

Referenced by DTDigitizer().

bool DTDigitizer::MultipleLinks
private

Definition at line 123 of file DTDigitizer.h.

Referenced by DTDigitizer(), and storeDigis().

bool DTDigitizer::onlyMuHits
private

Definition at line 110 of file DTDigitizer.h.

Referenced by DTDigitizer(), and storeDigis().

float DTDigitizer::smearing
private

Definition at line 107 of file DTDigitizer.h.

Referenced by driftTimeFromParametrization(), and DTDigitizer().

std::string DTDigitizer::syncName
private

Definition at line 113 of file DTDigitizer.h.

float DTDigitizer::theConstVDrift
private

Definition at line 120 of file DTDigitizer.h.

Referenced by computeTime(), and DTDigitizer().

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

Definition at line 114 of file DTDigitizer.h.

Referenced by externalDelays().

double DTDigitizer::vPropWire
private

Definition at line 105 of file DTDigitizer.h.

Referenced by DTDigitizer(), and externalDelays().