CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTDigiSyncTOFCorr.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * $Date: 2010/07/16 14:44:04 $
5  * $Revision: 1.8 $
6  * \author N. Amapane, R. Bellan - INFN Torino
7  */
8 
12 
16 #include <iostream>
17 
18 using namespace std;
19 
21 
22  theOffset = pSet.getParameter<double>("offset"); //500ns
23  corrType = pSet.getParameter<int>("TOFCorrection"); //1
24 }
25 
26 
28 
29 
30 // Delays to be added to digi times during digitization, in ns.
31 double DTDigiSyncTOFCorr::digitizerOffset(const DTWireId * id, const DTLayer* layer) const {
32 
33  double offset = theOffset;
34  const double cSpeed = 29.9792458; // cm/ns
35 
36  if (corrType==1) {
37  // Subtraction of assumed TOF, per CHAMBER
38  double flightL = layer->chamber()->surface().position().mag();
39  offset -= flightL/cSpeed;
40  } else if (corrType==2) {
41  // Subtraction of assumed TOF, per WIRE
42 
43  // Position of the wire in the Layer's reference frame
44  float localXPos = layer->specificTopology().wirePosition(id->wire());
45  LocalPoint localPos(localXPos,0,0);
46 
47  // Distance of the wire to the CMS's I.P.
48  double flightL = layer->surface().toGlobal(localPos).mag();
49 
50  offset -= flightL/cSpeed;
51 
52  } else if (corrType==3) {
53  // Subtraction of assumed TOF, per SUPERLAYER
54  double flightL = layer->superLayer()->surface().position().mag();
55  offset -= flightL/cSpeed;
56  } else if (corrType!=0){
57  cout << "ERROR: SimMuon:DTDigitizer:DTDigiSyncTOFCorr:TOFCorrection = " << corrType
58  << "is not defined " << endl;
59  }
60  return offset;
61 }
62 
63 
64 // Offset to obtain "raw" TDCs for the L1 emulator from digis.
65 double DTDigiSyncTOFCorr::emulatorOffset(const DTWireId * id) const {
66  return theOffset;
67 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:78
T getParameter(std::string const &) const
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:88
DTDigiSyncTOFCorr(const edm::ParameterSet &)
Constructor.
virtual double digitizerOffset(const DTWireId *id, const DTLayer *layer=0) const
Delays to be added to digi times during digitization, in ns.
const DTTopology & specificTopology() const
Definition: DTLayer.cc:44
T mag() const
Definition: PV3DBase.h:66
unsigned int offset(bool)
int wire() const
Return the wire number.
Definition: DTWireId.h:58
virtual ~DTDigiSyncTOFCorr()
Destructor.
const DTChamber * chamber() const
Definition: DTLayer.cc:60
const DTSuperLayer * superLayer() const
Definition: DTLayer.cc:56
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
tuple cout
Definition: gather_cfg.py:121
const PositionType & position() const
virtual double emulatorOffset(const DTWireId *id) const
Offset to obtain &quot;raw&quot; TDCs for the L1 emulator from digis.