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::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Classes

struct  hitLessT
 

Public Member Functions

 DTDigitizer (const edm::ParameterSet &)
 
virtual void produce (edm::Event &, const edm::EventSetup &)
 
 ~DTDigitizer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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
 
DTDigiSyncBasetheSync
 
double vPropWire
 

Friends

class DTDigitizerAnalysis
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

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 44 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 56 of file DTDigitizer.cc.

References gather_cfg::cout, debug, edm::hlt::Exception, reco::get(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), edm::Service< T >::isAvailable(), and AlCaHLTBitMon_QueryRunRegistry::string.

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

Definition at line 126 of file DTDigitizer.cc.

126  {
127 }

Member Function Documentation

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

Definition at line 445 of file DTDigitizer.cc.

References f, and edmStreamStallGrapher::t.

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

References gather_cfg::cout, Vector3DBase< T, FrameTag >::cross(), debug, delta, dir, Vector3DBase< T, FrameTag >::dot(), PSimHit::entryPoint(), PSimHit::exitPoint(), PSimHit::localDirection(), PSimHit::localPosition(), M_PI, PV3DBase< T, PVType, FrameType >::mag(), PSimHit::momentumAtEntry(), DTTopology::none, DTTopology::onWhichBorder(), PSimHit::pabs(), PSimHit::particleType(), RecoTauCleanerPlugins::pt, DTLayer::specificTopology(), mathSSE::sqrt(), 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().

234  {
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) cout<<"Hit local entry point: "<<entryP<<endl
245  <<"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) cout<<"wire position: "<<xwire
252  <<" 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) dumpHit(hit, xEntry, xExit,topo);
259 
260  // The bolean is used to flag the drift time computation
261  pair<float,bool> driftTime(0.,false);
262 
263  // if delta in gas->ignore, since it is included in the parametrisation.
264  // FIXME: should check that it is actually a delta ray produced by a nearby
265  // muon hit.
266 
267  if (partType == 11 && entrySide == DTTopology::none) {
268  if (debug) cout << " e- hit in gas; discarding " << endl;
269  return driftTime;
270  }
271 
272  float By = BLoc.y();
273  float Bz = BLoc.z();
274 
275  // Radius and sagitta according to direction of momentum
276  // (just for printing)
277  // NOTE: in cmsim, d is always taken // pHat!
278  LocalVector d = (exitP-entryP);
279  LocalVector pHat = hit->localDirection().unit();
280  LocalVector hHat = (d.cross(pHat.cross(d))).unit();
281  float cosAlpha = hHat.dot(pHat);
282  float sinAlpha = sqrt(1.-cosAlpha*cosAlpha);
283  float radius_P = (d.mag())/(2.*cosAlpha);
284  float sagitta_P = radius_P*(1.-sinAlpha);
285 
286  // Radius, sagitta according to field bending
287  // (just for printing)
288  float halfd = d.mag()/2.;
289  float BMag = BLoc.mag();
290  LocalVector pT = (pHat - (BLoc.unit()*pHat.dot(BLoc.unit())))*(hit->pabs());
291  float radius_B = (pT.mag()/(0.3*BMag))*100.;
292  float sagitta_B;
293  if (radius_B > halfd) {
294  sagitta_B = radius_B - sqrt(radius_B*radius_B - halfd*halfd);
295  } else {
296  sagitta_B = radius_B;
297  }
298 
299  // cos(delta), delta= angle between direction at entry and hit segment
300  // (just for printing)
301  float delta = pHat.dot(d.unit());
302  if (debug) cout << " delta = " << delta << endl
303  << " cosAlpha = " << cosAlpha << endl
304  << " sinAlpha = " << sinAlpha << endl
305  << " pMag = " << pT.mag() << endl
306  << " bMag = " << BMag << endl
307  << " pT = " << pT << endl
308  << " halfd = " << halfd << endl
309  << " radius_P (cm) = " << radius_P << endl
310  << " sagitta_P (um) = " << sagitta_P*10000. << endl
311  << " radius_B (cm) = " << radius_B << endl
312  << " sagitta_B (um) = " << sagitta_B*10000. << endl;
313 
314  // Select cases where parametrization can not be used.
315  bool noParametrisation =
316  ( ( entrySide == DTTopology::none || exitSide == DTTopology::none ) // case # 2,3,8,9 or 11
317  || (entrySide == exitSide) // case # 4 or 10
318  || ((entrySide == DTTopology::xMin && exitSide == DTTopology::xMax) ||
319  (entrySide == DTTopology::xMax && exitSide == DTTopology::xMin)) // Hit is case # 7
320  );
321 
322  // FIXME: now, debug warning only; consider treating those
323  // with TM algo.
324  if ( delta < 0.99996 // Track is not straight. FIXME: use sagitta?
325  && (noParametrisation == false)) {
326  if (debug) cout << "*** WARNING: hit is not straight, type = " << partType << endl;
327  }
328 
329  //************ 5A ***************
330 
331  if (!noParametrisation) {
332 
333  LocalVector dir = hit->momentumAtEntry(); // ex Measurement3DVector dir = hit->measurementDirection(); //FIXME?
334  float theta = atan(dir.x()/-dir.z())*180/M_PI;
335 
336  // FIXME: use dir if M.S. is included as GARFIELD option...
337  // otherwise use hit segment dirction.
338  // LocalVector dir0 = (exitP-entryP).unit();
339  // float theta = atan(dir0.x()/-dir0.z())*180/M_PI;
340  float x;
341 
342  Local3DPoint pt = hit->localPosition(); //ex Measurement3DPoint pt = hit->measurementPosition(); // FIXME?
343 
344  if(fabs(pt.z()) < 0.002) {
345  // hit center within 20 um from z=0, no need to extrapolate.
346  x = pt.x() - xwire;
347  } else {
348  x = xEntry - (entryP.z()*(xExit-xEntry))/(exitP.z()-entryP.z());
349  }
350 
351  if(IdealModel) return make_pair(fabs(x)/theConstVDrift,true);
352  else driftTime = driftTimeFromParametrization(x, theta, By, Bz, engine);
353 
354  }
355 
356 
357  if ((driftTime.second)==false) {
358  // Parametrisation not applicable, or failed. Use time map.
359  driftTime = driftTimeFromTimeMap();
360  }
361 
362  //************ 5B ***************
363 
364  // Signal propagation, TOF etc.
365  if (driftTime.second) {
366  driftTime.first += externalDelays(layer,wireId,hit);
367  }
368  return driftTime;
369 }
dbl * delta
Definition: mlp_gen.cc:36
void dumpHit(const PSimHit *hit, float xEntry, float xExit, const DTTopology &topo)
Definition: DTDigitizer.cc:576
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:47
float theConstVDrift
Definition: DTDigitizer.h:121
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:38
T mag() const
Definition: PV3DBase.h:67
Local3DPoint localPosition() const
Definition: PSimHit.h:44
T sqrt(T t)
Definition: SSEVec.h:48
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:63
Side
Sides of the cell.
Definition: DTTopology.h:88
bool IdealModel
Definition: DTDigitizer.h:120
#define M_PI
LocalVector localDirection() const
Obsolete. Same as momentumAtEntry().unit(), for backward compatibility.
Definition: PSimHit.h:52
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:373
int particleType() const
Definition: PSimHit.h:85
float externalDelays(const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit) const
Definition: DTDigitizer.cc:469
tuple cout
Definition: gather_cfg.py:121
dbl *** dir
Definition: mlp_gen.cc:35
Definition: DDAxes.h:10
T x() const
Definition: PV3DBase.h:62
std::pair< float, bool > driftTimeFromTimeMap() const
Definition: DTDigitizer.cc:461
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
pair< float, bool > DTDigitizer::driftTimeFromParametrization ( float  x,
float  alpha,
float  By,
float  Bz,
CLHEP::HepRandomEngine *  engine 
) const
private

Definition at line 373 of file DTDigitizer.cc.

References gather_cfg::cout, debug, GeomDetEnumerators::DT, f, archive::flag, edm::max(), DTDriftTimeParametrization::MB_DT_drift_time(), DTDriftTimeParametrization::drift_time::t_drift, DTDriftTimeParametrization::drift_time::t_width_m, DTDriftTimeParametrization::drift_time::t_width_p, theta(), and cond::rpcobgas::time.

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

Definition at line 461 of file DTDigitizer.cc.

References gather_cfg::cout, and debug.

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

Definition at line 576 of file DTDigitizer.cc.

References DTTopology::cellHeight(), DTTopology::cellLenght(), DTTopology::cellWidth(), gather_cfg::cout, PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), 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().

578  {
579 
580  LocalPoint entryP = hit->entryPoint();
581  LocalPoint exitP = hit->exitPoint();
582 
583  DTTopology::Side entrySide = topo.onWhichBorder(xEntry,entryP.y(),entryP.z());
584  DTTopology::Side exitSide = topo.onWhichBorder(xExit,exitP.y(),exitP.z());
585  // ProcessTypeEnumerator pTypes;
586 
587  cout << endl
588  << "------- SimHit: " << endl
589  << " Particle type = " << hit->particleType() << endl
590  << " process type = " << hit->processType() << endl
591  << " process type = " << hit->processType() << endl
592  // << " packedTrackId = " << hit->packedTrackId() << endl
593  << " trackId = " << hit->trackId() << endl // new,is the same as the
594  // previous?? FIXME-Check
595  << " |p| = " << hit->pabs() << endl
596  << " Energy loss = " << hit->energyLoss() << endl
597  // << " timeOffset = " << hit->timeOffset() << endl
598  // << " measurementPosition = " << hit->measurementPosition() << endl
599  // << " measurementDirection = " << hit->measurementDirection() << endl //FIXME
600  << " localDirection = " << hit->momentumAtEntry().unit() << endl //FIXME is it a versor?
601  << " Entry point = " << entryP << " cell x = " << xEntry << endl
602  << " Exit point = " << exitP << " cell x = " << xExit << endl
603  << " DR = = " << (exitP-entryP).mag() << endl
604  << " Dx = = " << (exitP-entryP).x() << endl
605  << " Cell w,h,l = (" << topo.cellWidth()
606  << " , " << topo.cellHeight()
607  << " , " << topo.cellLenght() << ") cm" << endl
608  << " DY entry from edge = " << topo.cellLenght()/2.-fabs(entryP.y())
609  << " DY exit from edge = " << topo.cellLenght()/2.-fabs(exitP.y())
610  << " entrySide = " << (int)entrySide
611  << " ; exitSide = " << (int)exitSide << endl;
612 
613 }
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:47
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:38
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:63
Side
Sides of the cell.
Definition: DTTopology.h:88
Vector3DBase unit() const
Definition: Vector3DBase.h:57
unsigned short processType() const
Definition: PSimHit.h:118
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
int particleType() const
Definition: PSimHit.h:85
unsigned int trackId() const
Definition: PSimHit.h:102
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10
float cellLenght() const
Definition: DTTopology.h:73
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
float DTDigitizer::externalDelays ( const DTLayer layer,
const DTWireId wireId,
const PSimHit hit 
) const
private

Definition at line 469 of file DTDigitizer.cc.

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

471  {
472 
473  // Time of signal propagation along wire.
474 
475  float wireCoord = hit->localPosition().y();
476  float halfL = (layer->specificTopology().cellLenght())/2.;
477  float propgL = halfL - wireCoord; // the FE is always located at the pos coord.
478 
479  float propDelay = propgL/vPropWire;
480 
481  // Real TOF.
482  float tof = hit->tof();
483 
484  // Delays and t0 according to theSync
485 
486  double sync= theSync->digitizerOffset(&wireId,layer);
487 
488 
489  if (debug) {
490  cout << " propDelay =" << propDelay
491  << "; TOF=" << tof
492  << "; sync= " << sync
493  << endl;
494  }
495 
496  return propDelay + tof + sync;
497 }
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:72
T y() const
Definition: PV3DBase.h:63
DTDigiSyncBase * theSync
Definition: DTDigitizer.h:115
const DTTopology & specificTopology() const
Definition: DTLayer.cc:42
Local3DPoint localPosition() const
Definition: PSimHit.h:44
virtual double digitizerOffset(const DTWireId *id, const DTLayer *layer) const =0
Delays to be added to digi times during digitization, in ns.
double vPropWire
Definition: DTDigitizer.h:107
tuple cout
Definition: gather_cfg.py:121
float cellLenght() const
Definition: DTTopology.h:73
void DTDigitizer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDProducer.

Definition at line 130 of file DTDigitizer.cc.

References gather_cfg::cout, debug, edm::EventID::event(), edm::EventSetup::get(), edm::Event::getByToken(), edm::RandomNumberGenerator::getEngine(), edm::EventBase::id(), DTWireId::layerId(), convertSQLitetoXML_cfg::output, edm::Handle< T >::product(), edm::Event::put(), edm::EventID::run(), trackerHits::simHits, edm::Event::streamID(), GeomDet::surface(), cond::rpcobgas::time, Surface::toGlobal(), and GloballyPositioned< T >::toLocal().

130  {
131 
133  CLHEP::HepRandomEngine* engine = &rng->getEngine(iEvent.streamID());
134 
135  if(debug)
136  cout << "--- Run: " << iEvent.id().run()
137  << " Event: " << iEvent.id().event() << endl;
138 
139  //************ 1 ***************
140  // create the container for the SimHits
141  // Handle<PSimHitContainer> simHits;
142  // iEvent.getByLabel("g4SimHits","MuonDTHits",simHits);
143 
144  // use MixCollection instead of the previous
146  iEvent.getByToken(cf_token, xFrame);
147 
148  auto_ptr<MixCollection<PSimHit> >
149  simHits( new MixCollection<PSimHit>(xFrame.product()) );
150 
151  // create the pointer to the Digi container
152  auto_ptr<DTDigiCollection> output(new DTDigiCollection());
153  // pointer to the DigiSimLink container
154  auto_ptr<DTDigiSimLinkCollection> outputLinks(new DTDigiSimLinkCollection());
155 
156  // Muon Geometry
157  ESHandle<DTGeometry> muonGeom;
158  iSetup.get<MuonGeometryRecord>().get(geometryType,muonGeom);
159 
160  // Magnetic Field
161  ESHandle<MagneticField> magnField;
162  iSetup.get<IdealMagneticFieldRecord>().get(magnField);
163 
164  //************ 2 ***************
165 
166  // These are sorted by DetId, i.e. by layer and then by wire #
167  // map<DTDetId, vector<const PSimHit*> > wireMap;
168  DTWireIdMap wireMap;
169 
170  for(MixCollection<PSimHit>::MixItr simHit = simHits->begin();
171  simHit != simHits->end(); simHit++){
172 
173  // Create the id of the wire, the simHits in the DT known also the wireId
174 
175  DTWireId wireId( (*simHit).detUnitId() );
176  // Fill the map
177  wireMap[wireId].push_back(&(*simHit));
178  }
179 
180  pair<float,bool> time(0.,false);
181 
182  //************ 3 ***************
183  // Loop over the wires
184  for(DTWireIdMapConstIter wire = wireMap.begin(); wire!=wireMap.end(); wire++){
185  // SimHit Container associated to the wire
186  const vector<const PSimHit*> & vhit = (*wire).second;
187  if(vhit.size()!=0) {
188  TDContainer tdCont; // It is a vector<pair<const PSimHit*,float> >;
189 
190  //************ 4 ***************
191  DTWireId wireId = (*wire).first;
192 
193  //const DTLayer* layer = dynamic_cast< const DTLayer* > (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();
198  hit != vhit.end(); hit++){
199  //************ 5 ***************
200  LocalPoint locPos = (*hit)->localPosition();
201 
202  const LocalVector BLoc=layer->surface().toLocal(magnField->inTesla(layer->surface().toGlobal(locPos)));
203 
204  time = computeTime(layer, wireId, *hit, BLoc, engine);
205 
206  //************ 6 ***************
207  if (time.second) {
208  tdCont.push_back(make_pair((*hit),time.first));
209  } else {
210  if (debug) cout << "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 the dead time)
218 
219  storeDigis(wireId,tdCont,*output,*outputLinks);
220  }
221 
222  }
223 
224  //************ 8 ***************
225  // Load the Digi Container in the Event
226  //iEvent.put(output,"MuonDTDigis");
227  iEvent.put(output);
228  iEvent.put(outputLinks);
229 
230 }
RunNumber_t run() const
Definition: EventID.h:42
std::pair< float, bool > computeTime(const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit, const LocalVector &BLoc, CLHEP::HepRandomEngine *)
Definition: DTDigitizer.cc:232
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
EventNumber_t event() const
Definition: EventID.h:44
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
void storeDigis(DTWireId &wireId, TDContainer &hits, DTDigiCollection &output, DTDigiSimLinkCollection &outputLinks)
Definition: DTDigitizer.cc:502
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
LocalPoint toLocal(const GlobalPoint &gp) const
edm::EDGetTokenT< CrossingFrame< PSimHit > > cf_token
Definition: DTDigitizer.h:131
DTWireIdMap::const_iterator DTWireIdMapConstIter
Definition: DTDigitizer.h:58
std::string geometryType
Definition: DTDigitizer.h:117
const T & get() const
Definition: EventSetup.h:55
tuple simHits
Definition: trackerHits.py:16
T const * product() const
Definition: Handle.h:81
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &) const =0
Use this engine in event methods.
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:56
MuonDigiCollection< DTLayerId, DTDigi > DTDigiCollection
StreamID streamID() const
Definition: Event.h:75
tuple cout
Definition: gather_cfg.py:121
MuonDigiCollection< DTLayerId, DTDigiSimLink > DTDigiSimLinkCollection
std::map< DTWireId, std::vector< const PSimHit * > > DTWireIdMap
Definition: DTDigitizer.h:56
void DTDigitizer::storeDigis ( DTWireId wireId,
TDContainer hits,
DTDigiCollection output,
DTDigiSimLinkCollection outputLinks 
)
private

Definition at line 502 of file DTDigitizer.cc.

References funct::abs(), DTDigi::countsTDC(), gather_cfg::cout, debug, MuonDigiCollection< IndexType, DigiType >::insertDigi(), DTWireId::layerId(), DTDigi::number(), funct::pow(), python.multivaluedict::sort(), DTDigi::time(), cond::rpcobgas::time, and DTWireId::wire().

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

Friends And Related Function Documentation

friend class DTDigitizerAnalysis
friend

Definition at line 104 of file DTDigitizer.h.

Member Data Documentation

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

Definition at line 131 of file DTDigitizer.h.

std::string DTDigitizer::collection_for_XF
private

Definition at line 129 of file DTDigitizer.h.

float DTDigitizer::deadTime
private

Definition at line 108 of file DTDigitizer.h.

bool DTDigitizer::debug
private

Definition at line 110 of file DTDigitizer.h.

std::string DTDigitizer::geometryType
private

Definition at line 117 of file DTDigitizer.h.

bool DTDigitizer::IdealModel
private

Definition at line 120 of file DTDigitizer.h.

bool DTDigitizer::interpolate
private

Definition at line 111 of file DTDigitizer.h.

float DTDigitizer::LinksTimeWindow
private

Definition at line 125 of file DTDigitizer.h.

std::string DTDigitizer::mix_
private

Definition at line 128 of file DTDigitizer.h.

bool DTDigitizer::MultipleLinks
private

Definition at line 124 of file DTDigitizer.h.

bool DTDigitizer::onlyMuHits
private

Definition at line 112 of file DTDigitizer.h.

float DTDigitizer::smearing
private

Definition at line 109 of file DTDigitizer.h.

std::string DTDigitizer::syncName
private

Definition at line 114 of file DTDigitizer.h.

float DTDigitizer::theConstVDrift
private

Definition at line 121 of file DTDigitizer.h.

DTDigiSyncBase* DTDigitizer::theSync
private

Definition at line 115 of file DTDigitizer.h.

double DTDigitizer::vPropWire
private

Definition at line 107 of file DTDigitizer.h.