CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTDigitizer.h
Go to the documentation of this file.
1 #ifndef SimMuon_DTDigitizer_h
2 #define SimMuon_DTDigitizer_h
3 
15 
19 
21 
22 #include <vector>
23 
24 namespace CLHEP {
25  class RandGaussQ;
26  class RandFlat;
27 }
28 
29 class DTLayer;
30 class PSimHit;
31 class DTWireType;
32 class DTBaseDigiSync;
33 class DTTopology;
34 class DTDigiSyncBase;
35 
36 
37 namespace edm {class ParameterSet; class Event; class EventSetup;}
38 
39 class DTDigitizer : public edm::EDProducer {
40 
41  public:
42 
43  explicit DTDigitizer(const edm::ParameterSet&);
44  ~DTDigitizer();
45  virtual void produce(edm::Event&, const edm::EventSetup&);
46 
47  private:
48  typedef std::pair<const PSimHit*,float> hitAndT; // hit & corresponding time
49  typedef std::vector<hitAndT> TDContainer; // hits & times for one wire
50 
51  typedef std::map<DTWireId, std::vector<const PSimHit*> > DTWireIdMap;
52  typedef DTWireIdMap::iterator DTWireIdMapIter;
53  typedef DTWireIdMap::const_iterator DTWireIdMapConstIter;
54 
55  // Sort hits container by time.
56  struct hitLessT {
57  bool operator()(const hitAndT & h1, const hitAndT & h2) {
58  if (h1.second < h2.second) return true;
59  return false;
60  }
61  };
62 
63  // Calculate the drift time for one hit.
64  // if status flag == false, hit has to be discarded.
65  std::pair<float,bool> computeTime(const DTLayer* layer,const DTWireId &wireId,
66  const PSimHit *hit,
67  const LocalVector &BLoc); //FIXME??
68 
69  // Calculate the drift time using the GARFIELD cell parametrization,
70  // taking care of all conversions from CMSSW local coordinates
71  // to the conventions used for the parametrization.
72  std::pair<float,bool> driftTimeFromParametrization(float x, float alpha, float By,
73  float Bz) const;
74 
75  // Calculate the drift time for the cases where it is not possible
76  // to use the GARFIELD cell parametrization.
77  std::pair<float,bool> driftTimeFromTimeMap() const;
78 
79  // Add all delays other than drift times (signal propagation along the wire,
80  // TOF etc.; subtract calibration time.
81  float externalDelays(const DTLayer* layer,
82  const DTWireId &wireId,
83  const PSimHit *hit) const;
84 
85  // Store digis for one wire, taking into account the dead time.
86  //FiXME put alias for the map.
87  void storeDigis(DTWireId &wireId,
88  TDContainer &hits,
90 
91  // Debug output
92  void dumpHit(const PSimHit * hit, float xEntry, float xExit, const DTTopology &topo);
93 
94  // Double half-gaussian smearing.
95  float asymGausSmear(double mean, double sigmaLeft, double sigmaRight) const;
96 
97  // Allow debugging and testing.
98  friend class DTDigitizerAnalysis;
99 
100  //Its Atributes:
101  double vPropWire;
102  float deadTime;
103  float smearing;
104  bool debug;
107 
110 
112 
113  // Ideal model. Used for debug
116 
117  // the random generator
118  CLHEP::RandGaussQ* theGaussianDistribution;
119  CLHEP::RandFlat* theFlatDistribution;
120 
121  // to configure the creation of Digi-Sim links
124 
125  //Name of Collection use for create the XF
128 
129 };
130 #endif
void dumpHit(const PSimHit *hit, float xEntry, float xExit, const DTTopology &topo)
Definition: DTDigitizer.cc:578
float alpha
Definition: AMPTWrapper.h:95
std::pair< float, bool > computeTime(const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit, const LocalVector &BLoc)
Definition: DTDigitizer.cc:237
bool onlyMuHits
Definition: DTDigitizer.h:106
float asymGausSmear(double mean, double sigmaLeft, double sigmaRight) const
Definition: DTDigitizer.cc:448
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
std::pair< float, bool > driftTimeFromParametrization(float x, float alpha, float By, float Bz) const
Definition: DTDigitizer.cc:377
std::pair< const PSimHit *, float > hitAndT
Definition: DTDigitizer.h:48
void storeDigis(DTWireId &wireId, TDContainer &hits, DTDigiCollection &output, DTDigiSimLinkCollection &outputLinks)
Definition: DTDigitizer.cc:504
DTDigitizer(const edm::ParameterSet &)
Definition: DTDigitizer.cc:63
float smearing
Definition: DTDigitizer.h:103
bool interpolate
Definition: DTDigitizer.h:105
DTWireIdMap::const_iterator DTWireIdMapConstIter
Definition: DTDigitizer.h:53
bool IdealModel
Definition: DTDigitizer.h:114
bool operator()(const hitAndT &h1, const hitAndT &h2)
Definition: DTDigitizer.h:57
std::string mix_
Definition: DTDigitizer.h:126
std::string geometryType
Definition: DTDigitizer.h:111
std::vector< hitAndT > TDContainer
Definition: DTDigitizer.h:49
std::string syncName
Definition: DTDigitizer.h:108
double vPropWire
Definition: DTDigitizer.h:101
float externalDelays(const DTLayer *layer, const DTWireId &wireId, const PSimHit *hit) const
Definition: DTDigitizer.cc:471
friend class DTDigitizerAnalysis
Definition: DTDigitizer.h:98
virtual void produce(edm::Event &, const edm::EventSetup &)
Definition: DTDigitizer.cc:139
Definition: DDAxes.h:10
float deadTime
Definition: DTDigitizer.h:102
CLHEP::RandFlat * theFlatDistribution
Definition: DTDigitizer.h:119
DTWireIdMap::iterator DTWireIdMapIter
Definition: DTDigitizer.h:52
std::pair< float, bool > driftTimeFromTimeMap() const
Definition: DTDigitizer.cc:463
std::map< DTWireId, std::vector< const PSimHit * > > DTWireIdMap
Definition: DTDigitizer.h:51
bool MultipleLinks
Definition: DTDigitizer.h:122
std::string collection_for_XF
Definition: DTDigitizer.h:127
A container for a generic type of digis indexed by some index, implemented with a map&lt;IndexType...