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
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () 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

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 cf_token, collection_for_XF, deadTime, debug, Exception, geometryType, reco::get(), edm::ParameterSet::getParameter(), IdealModel, 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  // Constant drift velocity needed by the above flag
99  if (IdealModel)
100  theConstVDrift = conf_.getParameter<double>("IdealModelConstantDriftVelocity"); // 55 um/ns
101  else
102  theConstVDrift = 55.;
103 
104  // get random engine
106  if (!rng.isAvailable()) {
107  throw cms::Exception("Configuration") << "RandomNumberGeneratorService for DTDigitizer missing in cfg file";
108  }
109 
110  // MultipleLinks=false ==> one-to-one correspondence between digis and SimHits
111  MultipleLinks = conf_.getParameter<bool>("MultipleLinks");
112  // MultipleLinks=true ==> association of SimHits within a time window
113  // LinksTimeWindow (of the order of the resolution)
114  LinksTimeWindow = conf_.getParameter<double>("LinksTimeWindow"); // (10 ns)
115 
116  // Name of Collection used for create the XF
117  mix_ = conf_.getParameter<std::string>("mixLabel");
118  collection_for_XF = conf_.getParameter<std::string>("InputCollection");
119  cf_token = consumes<CrossingFrame<PSimHit>>(edm::InputTag(mix_, collection_for_XF));
120 
121  // String to choice between ideal (the deafult) and (mis)aligned geometry for
122  // the digitization step
123  geometryType = conf_.getParameter<std::string>("GeometryType");
124 }
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:119
float LinksTimeWindow
Definition: DTDigitizer.h:123
bool isAvailable() const
Definition: Service.h:40
edm::EDGetTokenT< CrossingFrame< PSimHit > > cf_token
Definition: DTDigitizer.h:129
float smearing
Definition: DTDigitizer.h:107
bool interpolate
Definition: DTDigitizer.h:109
bool IdealModel
Definition: DTDigitizer.h:118
std::unique_ptr< DTDigiSyncBase > theSync
Definition: DTDigitizer.h:113
std::string mix_
Definition: DTDigitizer.h:126
std::string geometryType
Definition: DTDigitizer.h:115
double vPropWire
Definition: DTDigitizer.h:105
float deadTime
Definition: DTDigitizer.h:106
T get(const Candidate &c)
Definition: component.h:55
bool MultipleLinks
Definition: DTDigitizer.h:122
std::string collection_for_XF
Definition: DTDigitizer.h:127

Member Function Documentation

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

Definition at line 442 of file DTDigitizer.cc.

References f, and lumiQTWidget::t.

Referenced by driftTimeFromParametrization().

445  {
446  double f = sigmaLeft / (sigmaLeft + sigmaRight);
447  double t;
448 
449  if (CLHEP::RandFlat::shoot(engine) <= f) {
450  t = CLHEP::RandGaussQ::shoot(engine, mean, sigmaLeft);
451  t = mean - fabs(t - mean);
452  } else {
453  t = CLHEP::RandGaussQ::shoot(engine, mean, sigmaRight);
454  t = mean + fabs(t - mean);
455  }
456  return static_cast<float>(t);
457 }
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 224 of file DTDigitizer.cc.

References Vector3DBase< T, FrameTag >::cross(), edmIntegrityCheck::d, debug, delta, 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(), EnergyCorrector::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().

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

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

Definition at line 459 of file DTDigitizer.cc.

References debug.

Referenced by computeTime().

459  {
460  // FIXME: not yet implemented.
461  if (debug)
462  LogPrint("DTDigitizer") << " TimeMap " << endl;
463  return pair<float, bool>(0., false);
464 }
void DTDigitizer::dumpHit ( const PSimHit hit,
float  xEntry,
float  xExit,
const DTTopology topo 
)
private

Definition at line 564 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().

564  {
565  LocalPoint entryP = hit->entryPoint();
566  LocalPoint exitP = hit->exitPoint();
567 
568  DTTopology::Side entrySide = topo.onWhichBorder(xEntry, entryP.y(), entryP.z());
569  DTTopology::Side exitSide = topo.onWhichBorder(xExit, exitP.y(), exitP.z());
570  // ProcessTypeEnumerator pTypes;
571 
572  LogPrint("DTDigitizer") << endl
573  << "------- SimHit: " << endl
574  << " Particle type = " << hit->particleType() << endl
575  << " process type = " << hit->processType() << endl
576  << " process type = " << hit->processType()
577  << endl
578  // << " packedTrackId = " << hit->packedTrackId() << endl
579  << " trackId = " << hit->trackId()
580  << endl // new,is the same as the
581  // previous?? FIXME-Check
582  << " |p| = " << hit->pabs() << endl
583  << " Energy loss = " << hit->energyLoss()
584  << endl
585  // << " timeOffset = " << hit->timeOffset() << endl
586  // << " measurementPosition = " << hit->measurementPosition() << endl
587  // << " measurementDirection = " << hit->measurementDirection() << endl
588  // //FIXME
589  << " localDirection = " << hit->momentumAtEntry().unit()
590  << endl // FIXME is it a versor?
591  << " Entry point = " << entryP << " cell x = " << xEntry << endl
592  << " Exit point = " << exitP << " cell x = " << xExit << endl
593  << " DR = = " << (exitP - entryP).mag() << endl
594  << " Dx = = " << (exitP - entryP).x() << endl
595  << " Cell w,h,l = (" << topo.cellWidth() << " , " << topo.cellHeight() << " , "
596  << topo.cellLenght() << ") cm" << endl
597  << " DY entry from edge = " << topo.cellLenght() / 2. - fabs(entryP.y())
598  << " DY exit from edge = " << topo.cellLenght() / 2. - fabs(exitP.y())
599  << " entrySide = " << (int)entrySide << " ; exitSide = " << (int)exitSide << endl;
600 }
Side onWhichBorder(float x, float y, float z) const
Definition: DTTopology.cc:109
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:63
float cellWidth() const
Returns the cell width.
Definition: DTTopology.h:68
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
T z() const
Definition: PV3DBase.h:64
float cellHeight() const
Returns the cell height.
Definition: DTTopology.h:70
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
Side
Sides of the cell.
Definition: DTTopology.h:88
Vector3DBase unit() const
Definition: Vector3DBase.h:57
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:73
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 468 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().

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

Definition at line 127 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, trackerHits::simHits, storeDigis(), edm::Event::streamID(), GeomDet::surface(), ntuplemaker::time, Surface::toGlobal(), and GloballyPositioned< T >::toLocal().

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

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

Referenced by produce().

496  {
497  //************ 7A ***************
498 
499  // sort signal times
500  sort(hits.begin(), hits.end(), hitLessT());
501 
502  //************ 7B ***************
503 
504  float wakeTime = -999999.0;
505  float resolTime = -999999.0;
506  int digiN = -1; // Digi identifier within the cell (for multiple digis)
507  DTDigi digi;
508 
509  // loop over signal times and drop signals inside dead time
510  for (TDContainer::const_iterator hit = hits.begin(); hit != hits.end(); hit++) {
511  if (onlyMuHits && abs((*hit).first->particleType()) != 13)
512  continue;
513 
514  //************ 7C ***************
515 
516  float time = (*hit).second;
517  if (time > wakeTime) {
518  // Note that digi is constructed with a float value (in ns)
519  int wireN = wireId.wire();
520  digiN++;
521  digi = DTDigi(wireN, time, digiN);
522 
523  // Add association between THIS digi and the corresponding SimTrack
524  unsigned int SimTrackId = (*hit).first->trackId();
525  EncodedEventId evId = (*hit).first->eventId();
526  DTDigiSimLink digisimLink(wireN, digiN, time, SimTrackId, evId);
527 
528  if (debug) {
529  LogPrint("DTDigitizer") << endl << "---- DTDigitizer ----" << endl;
530  LogPrint("DTDigitizer") << "wireId: " << wireId << endl;
531  LogPrint("DTDigitizer") << "sim. time = " << time << endl;
532  LogPrint("DTDigitizer") << "digi number = " << digi.number() << ", digi time = " << digi.time()
533  << ", linked to SimTrack Id = " << SimTrackId << endl;
534  }
535 
536  //************ 7D ***************
537  if (digi.countsTDC() < pow(2., 16)) {
538  DTLayerId layerID = wireId.layerId(); // taking the layer in which reside 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 {
544  digiN--;
545  }
546  } else if (MultipleLinks && time < resolTime) {
547  int wireN = wireId.wire();
548  unsigned int SimTrackId = (*hit).first->trackId();
549  EncodedEventId evId = (*hit).first->eventId();
550  DTDigiSimLink digisimLink(wireN, digiN, time, SimTrackId, evId);
551  DTLayerId layerID = wireId.layerId();
552  outputLinks.insertDigi(layerID, digisimLink);
553 
554  if (debug) {
555  LogPrint("DTDigitizer") << "\nAdded multiple link: \n"
556  << "digi number = " << digi.number() << ", digi time = " << digi.time()
557  << " (sim. time = " << time << ")"
558  << ", linked to SimTrack Id = " << SimTrackId << endl;
559  }
560  }
561  }
562 }
bool onlyMuHits
Definition: DTDigitizer.h:110
float LinksTimeWindow
Definition: DTDigitizer.h:123
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:63
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: DTDigi.h:17
uint32_t countsTDC() const
Get raw TDC count.
Definition: DTDigi.cc:65
int wire() const
Return the wire number.
Definition: DTWireId.h:56
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:62
int number() const
Identifies different digis within the same.
Definition: DTDigi.cc:69
float deadTime
Definition: DTDigitizer.h:106
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
bool MultipleLinks
Definition: DTDigitizer.h:122

Friends And Related Function Documentation

friend class DTDigitizerAnalysis
friend

Definition at line 102 of file DTDigitizer.h.

Member Data Documentation

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

Definition at line 129 of file DTDigitizer.h.

Referenced by DTDigitizer(), and produce().

std::string DTDigitizer::collection_for_XF
private

Definition at line 127 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 115 of file DTDigitizer.h.

Referenced by DTDigitizer(), and produce().

bool DTDigitizer::IdealModel
private

Definition at line 118 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 123 of file DTDigitizer.h.

Referenced by DTDigitizer(), and storeDigis().

std::string DTDigitizer::mix_
private

Definition at line 126 of file DTDigitizer.h.

Referenced by DTDigitizer().

bool DTDigitizer::MultipleLinks
private

Definition at line 122 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 112 of file DTDigitizer.h.

float DTDigitizer::theConstVDrift
private

Definition at line 119 of file DTDigitizer.h.

Referenced by computeTime(), and DTDigitizer().

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

Definition at line 113 of file DTDigitizer.h.

Referenced by externalDelays().

double DTDigitizer::vPropWire
private

Definition at line 105 of file DTDigitizer.h.

Referenced by DTDigitizer(), and externalDelays().