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 ()
 
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 ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) 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) const
 
std::pair< float, bool > computeTime (const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit, const LocalVector &BLoc)
 
std::pair< float, bool > driftTimeFromParametrization (float x, float alpha, float By, float Bz) 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

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
 
CLHEP::RandFlat * theFlatDistribution
 
CLHEP::RandGaussQ * theGaussianDistribution
 
DTDigiSyncBasetheSync
 
double vPropWire
 

Friends

class DTDigitizerAnalysis
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- 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::EDProducer
CurrentProcessingContext const * currentContext () const
 
- 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.

Date:
2013/03/02 12:37:47
Revision:
1.13
Authors
: G. Bevilacqua, N. Amapane, G. Cerminara, R. Bellan

Definition at line 39 of file DTDigitizer.h.

Member Typedef Documentation

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

Definition at line 51 of file DTDigitizer.h.

typedef DTWireIdMap::const_iterator DTDigitizer::DTWireIdMapConstIter
private

Definition at line 53 of file DTDigitizer.h.

typedef DTWireIdMap::iterator DTDigitizer::DTWireIdMapIter
private

Definition at line 52 of file DTDigitizer.h.

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

Definition at line 48 of file DTDigitizer.h.

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

Definition at line 49 of file DTDigitizer.h.

Constructor & Destructor Documentation

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

Definition at line 63 of file DTDigitizer.cc.

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

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

Definition at line 133 of file DTDigitizer.cc.

133  {
135  delete theFlatDistribution;
136 }
CLHEP::RandGaussQ * theGaussianDistribution
Definition: DTDigitizer.h:118
CLHEP::RandFlat * theFlatDistribution
Definition: DTDigitizer.h:119

Member Function Documentation

float DTDigitizer::asymGausSmear ( double  mean,
double  sigmaLeft,
double  sigmaRight 
) const
private

Definition at line 448 of file DTDigitizer.cc.

References f, and lumiQTWidget::t.

448  {
449 
450  double f = sigmaLeft/(sigmaLeft+sigmaRight);
451  double t;
452 
453  if (theFlatDistribution->fire() <= f) {
454  t = theGaussianDistribution->fire(mean,sigmaLeft);
455  t = mean - fabs(t - mean);
456  } else {
457  t = theGaussianDistribution->fire(mean,sigmaRight);
458  t = mean + fabs(t - mean);
459  }
460  return static_cast<float>(t);
461 }
CLHEP::RandGaussQ * theGaussianDistribution
Definition: DTDigitizer.h:118
double f[11][100]
CLHEP::RandFlat * theFlatDistribution
Definition: DTDigitizer.h:119
pair< float, bool > DTDigitizer::computeTime ( const DTLayer layer,
const DTWireId wireId,
const PSimHit hit,
const LocalVector BLoc 
)
private

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

238  {
239 
240  LocalPoint entryP = hit->entryPoint();
241  LocalPoint exitP = hit->exitPoint();
242  int partType = hit->particleType();
243 
244  const DTTopology &topo = layer->specificTopology();
245 
246  // Pay attention: in CMSSW the rf of the SimHit is in the layer's rf
247 
248  if(debug) cout<<"Hit local entry point: "<<entryP<<endl
249  <<"Hit local exit point: "<<exitP<<endl;
250 
251  float xwire = topo.wirePosition(wireId.wire());
252  float xEntry = entryP.x() - xwire;
253  float xExit = exitP.x() - xwire;
254 
255  if(debug) cout<<"wire position: "<<xwire
256  <<" x entry in cell rf: "<<xEntry
257  <<" x exit in cell rf: "<<xExit<<endl;
258 
259  DTTopology::Side entrySide = topo.onWhichBorder(xEntry,entryP.y(),entryP.z());
260  DTTopology::Side exitSide = topo.onWhichBorder(xExit,exitP.y(),exitP.z());
261 
262  if (debug) dumpHit(hit, xEntry, xExit,topo);
263 
264  // The bolean is used to flag the drift time computation
265  pair<float,bool> driftTime(0.,false);
266 
267  // if delta in gas->ignore, since it is included in the parametrisation.
268  // FIXME: should check that it is actually a delta ray produced by a nearby
269  // muon hit.
270 
271  if (partType == 11 && entrySide == DTTopology::none) {
272  if (debug) cout << " e- hit in gas; discarding " << endl;
273  return driftTime;
274  }
275 
276  float By = BLoc.y();
277  float Bz = BLoc.z();
278 
279  // Radius and sagitta according to direction of momentum
280  // (just for printing)
281  // NOTE: in cmsim, d is always taken // pHat!
282  LocalVector d = (exitP-entryP);
283  LocalVector pHat = hit->localDirection().unit();
284  LocalVector hHat = (d.cross(pHat.cross(d))).unit();
285  float cosAlpha = hHat.dot(pHat);
286  float sinAlpha = sqrt(1.-cosAlpha*cosAlpha);
287  float radius_P = (d.mag())/(2.*cosAlpha);
288  float sagitta_P = radius_P*(1.-sinAlpha);
289 
290  // Radius, sagitta according to field bending
291  // (just for printing)
292  float halfd = d.mag()/2.;
293  float BMag = BLoc.mag();
294  LocalVector pT = (pHat - (BLoc.unit()*pHat.dot(BLoc.unit())))*(hit->pabs());
295  float radius_B = (pT.mag()/(0.3*BMag))*100.;
296  float sagitta_B;
297  if (radius_B > halfd) {
298  sagitta_B = radius_B - sqrt(radius_B*radius_B - halfd*halfd);
299  } else {
300  sagitta_B = radius_B;
301  }
302 
303  // cos(delta), delta= angle between direction at entry and hit segment
304  // (just for printing)
305  float delta = pHat.dot(d.unit());
306  if (debug) cout << " delta = " << delta << endl
307  << " cosAlpha = " << cosAlpha << endl
308  << " sinAlpha = " << sinAlpha << endl
309  << " pMag = " << pT.mag() << endl
310  << " bMag = " << BMag << endl
311  << " pT = " << pT << endl
312  << " halfd = " << halfd << endl
313  << " radius_P (cm) = " << radius_P << endl
314  << " sagitta_P (um) = " << sagitta_P*10000. << endl
315  << " radius_B (cm) = " << radius_B << endl
316  << " sagitta_B (um) = " << sagitta_B*10000. << endl;
317 
318  // Select cases where parametrization can not be used.
319  bool noParametrisation =
320  ( ( entrySide == DTTopology::none || exitSide == DTTopology::none ) // case # 2,3,8,9 or 11
321  || (entrySide == exitSide) // case # 4 or 10
322  || ((entrySide == DTTopology::xMin && exitSide == DTTopology::xMax) ||
323  (entrySide == DTTopology::xMax && exitSide == DTTopology::xMin)) // Hit is case # 7
324  );
325 
326  // FIXME: now, debug warning only; consider treating those
327  // with TM algo.
328  if ( delta < 0.99996 // Track is not straight. FIXME: use sagitta?
329  && (noParametrisation == false)) {
330  if (debug) cout << "*** WARNING: hit is not straight, type = " << partType << endl;
331  }
332 
333  //************ 5A ***************
334 
335  if (!noParametrisation) {
336 
337  LocalVector dir = hit->momentumAtEntry(); // ex Measurement3DVector dir = hit->measurementDirection(); //FIXME?
338  float theta = atan(dir.x()/-dir.z())*180/M_PI;
339 
340  // FIXME: use dir if M.S. is included as GARFIELD option...
341  // otherwise use hit segment dirction.
342  // LocalVector dir0 = (exitP-entryP).unit();
343  // float theta = atan(dir0.x()/-dir0.z())*180/M_PI;
344  float x;
345 
346  Local3DPoint pt = hit->localPosition(); //ex Measurement3DPoint pt = 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) return make_pair(fabs(x)/theConstVDrift,true);
356  else driftTime = driftTimeFromParametrization(x, theta, By, Bz);
357 
358  }
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 }
dbl * delta
Definition: mlp_gen.cc:36
void dumpHit(const PSimHit *hit, float xEntry, float xExit, const DTTopology &topo)
Definition: DTDigitizer.cc:578
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:88
Side onWhichBorder(float x, float y, float z) const
Definition: DTTopology.cc:111
LocalVector momentumAtEntry() const
The momentum of the track that produced the hit, at entry point.
Definition: PSimHit.h:47
float theConstVDrift
Definition: DTDigitizer.h:115
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
std::pair< float, bool > driftTimeFromParametrization(float x, float alpha, float By, float Bz) const
Definition: DTDigitizer.cc:377
const DTTopology & specificTopology() const
Definition: DTLayer.cc:44
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:90
bool IdealModel
Definition: DTDigitizer.h:114
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:58
Vector3DBase unit() const
Definition: Vector3DBase.h:57
#define M_PI
Definition: BFit3D.cc:3
int particleType() const
Definition: PSimHit.h:85
float externalDelays(const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit) const
Definition: DTDigitizer.cc:471
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:463
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 
) const
private

Definition at line 377 of file DTDigitizer.cc.

References gather_cfg::cout, debug, GeomDetEnumerators::DT, f, 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.

377  {
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) cout << "*** WARNING: parametrisation: x out of range = "
386  << 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) cout << "*** WARNING: extrapolating theta > 45: "
401  << 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) cout << "*** WARNING: extrapolating Bnorm > 0.75: "
406  << 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) cout << "*** WARNING: extrapolating Bwire >0.4: "
411  << Bz_par << endl;
412  // Bz_par = min(fabs(Bz_par),0.4)*((Bz_par<0.)?-1.:1.);
413  }
414 
416  static 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  cout << " Parametrisation: x, theta, Bnorm, Bwire = "
421  << x << " " << theta_par << " " << By_par << " " << Bz_par << endl
422  << " time=" << DT.t_drift
423  << " sigma_m=" << DT.t_width_m
424  << " sigma_p=" << DT.t_width_p << endl;
425  if (flag!=1) {
426  cout << "*** WARNING: call to parametrisation failed" << endl;
427  return pair<float,bool>(0.f,false);
428  }
429  }
430 
431  // Double half-gaussian smearing
432  float time = asymGausSmear(DT.t_drift, DT.t_width_m, DT.t_width_p);
433 
434  // Do not allow the smearing to lead to negative values
435  time = max(time,0.f);
436 
437  // Apply a Gaussian smearing to account for electronic effects (cf. 2004 TB analysis)
438  // The width of the Gaussian can be configured with the "Smearing" parameter
439 
440  double u = theGaussianDistribution->fire(0.,smearing);
441  time += u;
442 
443  if (debug) cout << " drift time = " << time << endl;
444 
445  return pair<float,bool>(time,true);
446 }
long int flag
Definition: mlp_lapack.h:47
float asymGausSmear(double mean, double sigmaLeft, double sigmaRight) const
Definition: DTDigitizer.cc:448
CLHEP::RandGaussQ * theGaussianDistribution
Definition: DTDigitizer.h:118
Geom::Theta< T > theta() const
const T & max(const T &a, const T &b)
double f[11][100]
float smearing
Definition: DTDigitizer.h:103
bool interpolate
Definition: DTDigitizer.h:105
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
pair< float, bool > DTDigitizer::driftTimeFromTimeMap ( ) const
private

Definition at line 463 of file DTDigitizer.cc.

References gather_cfg::cout, and debug.

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

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

580  {
581 
582  LocalPoint entryP = hit->entryPoint();
583  LocalPoint exitP = hit->exitPoint();
584 
585  DTTopology::Side entrySide = topo.onWhichBorder(xEntry,entryP.y(),entryP.z());
586  DTTopology::Side exitSide = topo.onWhichBorder(xExit,exitP.y(),exitP.z());
587  // ProcessTypeEnumerator pTypes;
588 
589  cout << endl
590  << "------- SimHit: " << endl
591  << " Particle type = " << hit->particleType() << endl
592  << " process type = " << hit->processType() << endl
593  << " process type = " << hit->processType() << endl
594  // << " packedTrackId = " << hit->packedTrackId() << endl
595  << " trackId = " << hit->trackId() << endl // new,is the same as the
596  // previous?? FIXME-Check
597  << " |p| = " << hit->pabs() << endl
598  << " Energy loss = " << hit->energyLoss() << endl
599  // << " timeOffset = " << hit->timeOffset() << endl
600  // << " measurementPosition = " << hit->measurementPosition() << endl
601  // << " measurementDirection = " << hit->measurementDirection() << endl //FIXME
602  << " localDirection = " << hit->momentumAtEntry().unit() << endl //FIXME is it a versor?
603  << " Entry point = " << entryP << " cell x = " << xEntry << endl
604  << " Exit point = " << exitP << " cell x = " << xExit << endl
605  << " DR = = " << (exitP-entryP).mag() << endl
606  << " Dx = = " << (exitP-entryP).x() << endl
607  << " Cell w,h,l = (" << topo.cellWidth()
608  << " , " << topo.cellHeight()
609  << " , " << topo.cellLenght() << ") cm" << endl
610  << " DY entry from edge = " << topo.cellLenght()/2.-fabs(entryP.y())
611  << " DY exit from edge = " << topo.cellLenght()/2.-fabs(exitP.y())
612  << " entrySide = " << (int)entrySide
613  << " ; exitSide = " << (int)exitSide << endl;
614 
615 }
Side onWhichBorder(float x, float y, float z) const
Definition: DTTopology.cc:111
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:70
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:72
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
Side
Sides of the cell.
Definition: DTTopology.h:90
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:75
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 471 of file DTDigitizer.cc.

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

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

Implements edm::EDProducer.

Definition at line 139 of file DTDigitizer.cc.

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

139  {
140  if(debug)
141  cout << "--- Run: " << iEvent.id().run()
142  << " Event: " << iEvent.id().event() << endl;
143 
144  //************ 1 ***************
145  // create the container for the SimHits
146  // Handle<PSimHitContainer> simHits;
147  // iEvent.getByLabel("g4SimHits","MuonDTHits",simHits);
148 
149  // use MixCollection instead of the previous
151  iEvent.getByLabel(mix_,collection_for_XF,xFrame);
152 
153  auto_ptr<MixCollection<PSimHit> >
154  simHits( new MixCollection<PSimHit>(xFrame.product()) );
155 
156  // create the pointer to the Digi container
157  auto_ptr<DTDigiCollection> output(new DTDigiCollection());
158  // pointer to the DigiSimLink container
159  auto_ptr<DTDigiSimLinkCollection> outputLinks(new DTDigiSimLinkCollection());
160 
161  // Muon Geometry
162  ESHandle<DTGeometry> muonGeom;
163  iSetup.get<MuonGeometryRecord>().get(geometryType,muonGeom);
164 
165  // Magnetic Field
166  ESHandle<MagneticField> magnField;
167  iSetup.get<IdealMagneticFieldRecord>().get(magnField);
168 
169  //************ 2 ***************
170 
171  // These are sorted by DetId, i.e. by layer and then by wire #
172  // map<DTDetId, vector<const PSimHit*> > wireMap;
173  DTWireIdMap wireMap;
174 
175  for(MixCollection<PSimHit>::MixItr simHit = simHits->begin();
176  simHit != simHits->end(); simHit++){
177 
178  // Create the id of the wire, the simHits in the DT known also the wireId
179 
180  DTWireId wireId( (*simHit).detUnitId() );
181  // Fill the map
182  wireMap[wireId].push_back(&(*simHit));
183  }
184 
185  pair<float,bool> time(0.,false);
186 
187  //************ 3 ***************
188  // Loop over the wires
189  for(DTWireIdMapConstIter wire = wireMap.begin(); wire!=wireMap.end(); wire++){
190  // SimHit Container associated to the wire
191  const vector<const PSimHit*> & vhit = (*wire).second;
192  if(vhit.size()!=0) {
193  TDContainer tdCont; // It is a vector<pair<const PSimHit*,float> >;
194 
195  //************ 4 ***************
196  DTWireId wireId = (*wire).first;
197 
198  //const DTLayer* layer = dynamic_cast< const DTLayer* > (muonGeom->idToDet(wireId.layerId()));
199  const DTLayer* layer = muonGeom->layer(wireId.layerId());
200 
201  // Loop on the hits of this wire
202  for (vector<const PSimHit*>::const_iterator hit=vhit.begin();
203  hit != vhit.end(); hit++){
204  //************ 5 ***************
205  LocalPoint locPos = (*hit)->localPosition();
206 
207  const LocalVector BLoc=layer->surface().toLocal(magnField->inTesla(layer->surface().toGlobal(locPos)));
208 
209  time = computeTime(layer, wireId, *hit, BLoc);
210 
211  //************ 6 ***************
212  if (time.second) {
213  tdCont.push_back(make_pair((*hit),time.first));
214  } else {
215  if (debug) cout << "hit discarded" << endl;
216  }
217  }
218 
219  //************ 7 ***************
220 
221  // the loading must be done by layer but
222  // the digitization must be done by wire (in order to take into account the dead time)
223 
224  storeDigis(wireId,tdCont,*output,*outputLinks);
225  }
226 
227  }
228 
229  //************ 8 ***************
230  // Load the Digi Container in the Event
231  //iEvent.put(output,"MuonDTDigis");
232  iEvent.put(output);
233  iEvent.put(outputLinks);
234 
235 }
RunNumber_t run() const
Definition: EventID.h:42
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
EventNumber_t event() const
Definition: EventID.h:44
std::pair< float, bool > computeTime(const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit, const LocalVector &BLoc)
Definition: DTDigitizer.cc:237
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:504
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
LocalPoint toLocal(const GlobalPoint &gp) const
DTWireIdMap::const_iterator DTWireIdMapConstIter
Definition: DTDigitizer.h:53
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
std::string mix_
Definition: DTDigitizer.h:126
std::string geometryType
Definition: DTDigitizer.h:111
const T & get() const
Definition: EventSetup.h:55
tuple simHits
Definition: trackerHits.py:16
T const * product() const
Definition: Handle.h:74
std::vector< hitAndT > TDContainer
Definition: DTDigitizer.h:49
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:64
edm::EventID id() const
Definition: EventBase.h:56
MuonDigiCollection< DTLayerId, DTDigi > DTDigiCollection
tuple cout
Definition: gather_cfg.py:121
MuonDigiCollection< DTLayerId, DTDigiSimLink > DTDigiSimLinkCollection
std::map< DTWireId, std::vector< const PSimHit * > > DTWireIdMap
Definition: DTDigitizer.h:51
std::string collection_for_XF
Definition: DTDigitizer.h:127
void DTDigitizer::storeDigis ( DTWireId wireId,
TDContainer hits,
DTDigiCollection output,
DTDigiSimLinkCollection outputLinks 
)
private

Definition at line 504 of file DTDigitizer.cc.

References 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().

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

Member Data Documentation

std::string DTDigitizer::collection_for_XF
private

Definition at line 127 of file DTDigitizer.h.

float DTDigitizer::deadTime
private

Definition at line 102 of file DTDigitizer.h.

bool DTDigitizer::debug
private

Definition at line 104 of file DTDigitizer.h.

std::string DTDigitizer::geometryType
private

Definition at line 111 of file DTDigitizer.h.

bool DTDigitizer::IdealModel
private

Definition at line 114 of file DTDigitizer.h.

bool DTDigitizer::interpolate
private

Definition at line 105 of file DTDigitizer.h.

float DTDigitizer::LinksTimeWindow
private

Definition at line 123 of file DTDigitizer.h.

std::string DTDigitizer::mix_
private

Definition at line 126 of file DTDigitizer.h.

bool DTDigitizer::MultipleLinks
private

Definition at line 122 of file DTDigitizer.h.

bool DTDigitizer::onlyMuHits
private

Definition at line 106 of file DTDigitizer.h.

float DTDigitizer::smearing
private

Definition at line 103 of file DTDigitizer.h.

std::string DTDigitizer::syncName
private

Definition at line 108 of file DTDigitizer.h.

float DTDigitizer::theConstVDrift
private

Definition at line 115 of file DTDigitizer.h.

CLHEP::RandFlat* DTDigitizer::theFlatDistribution
private

Definition at line 119 of file DTDigitizer.h.

CLHEP::RandGaussQ* DTDigitizer::theGaussianDistribution
private

Definition at line 118 of file DTDigitizer.h.

DTDigiSyncBase* DTDigitizer::theSync
private

Definition at line 109 of file DTDigitizer.h.

double DTDigitizer::vPropWire
private

Definition at line 101 of file DTDigitizer.h.