CMS 3D CMS Logo

DTTTrigSyncTOFCorr.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author G. Cerminara - INFN Torino
5  */
6 
7 #include "DTTTrigSyncTOFCorr.h"
8 
14 #include <iostream>
15 
16 using namespace std;
17 
19  : // The fixed t0 (or t_trig) to be subtracted to digi time (ns)
20  theTTrig(config.getParameter<double>("tTrig")), // FIXME: Default was 500 ns
21  // Velocity of signal propagation along the wire (cm/ns)
22  theVPropWire(config.getParameter<double>("vPropWire")), // FIXME: Default was 24.4 cm/ns
23  // Select the mode for TOF correction:
24  // 0: tofCorr = TOF from IP to 3D Hit position (globPos)
25  // 1: tofCorr = TOF correction for distance difference
26  // between 3D center of the chamber and hit position
27  // (default)
28  // 2: tofCorr = TOF correction for distance difference
29  // between 3D center of the wire and hit position
30  // (This mode in available for backward compatibility)
31  theTOFCorrType(config.getParameter<int>("tofCorrType")), // FIXME: Default was 1
32  debug(config.getUntrackedParameter<bool>("debug")),
33  theBXspace(config.getUntrackedParameter<double>("bxSpace", 25.)) // spacing of BX in ns
34 {}
35 
37 
38 double DTTTrigSyncTOFCorr::offset(const DTLayer* layer,
39  const DTWireId& wireId,
40  const GlobalPoint& globPos,
41  double& tTrig,
42  double& wirePropCorr,
43  double& tofCorr) const {
44  tTrig = offset(wireId);
45 
46  //Compute the time spent in signal propagation along wire.
47  // NOTE: the FE is always at y>0
48  float halfL = layer->specificTopology().cellLenght() / 2;
49  float wireCoord = layer->toLocal(globPos).y();
50  float propgL = halfL - wireCoord;
51  wirePropCorr = propgL / theVPropWire;
52 
53  // Compute TOF correction treating it accordingly to
54  // the tofCorrType card
55  float flightToHit = globPos.mag();
56  static const float cSpeed = 29.9792458; // cm/ns
57  tofCorr = 0.;
58  switch (theTOFCorrType) {
59  case 0: {
60  // In this mode the subtraction of the TOF from IP to
61  // estimate 3D hit digi position is done here
62  // (No correction is needed anymore)
63  tofCorr = -flightToHit / cSpeed;
64  break;
65  }
66  case 1: {
67  // Correction for TOF from the center of the chamber to hit position
68  const DTChamber* chamber = layer->chamber();
69  double flightToChamber = chamber->surface().position().mag();
70  tofCorr = (flightToChamber - flightToHit) / cSpeed;
71  break;
72  }
73  case 2: {
74  // TOF from 3D center of the wire to hit position
75  float flightToWire =
76  layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
77  tofCorr = (flightToWire - flightToHit) / cSpeed;
78  break;
79  }
80  default: {
81  throw cms::Exception("[DTTTrigSyncTOFCorr]")
82  << " Invalid parameter: tofCorrType = " << theTOFCorrType << std::endl;
83  break;
84  }
85  }
86 
87  if (debug) {
88  cout << "[DTTTrigSyncTOFCorr] Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
89  << " various contributions are: " << endl
90  << " tTrig (ns): " << tTrig << endl
91  << " Propagation along wire delay (ns): " << wirePropCorr << endl
92  << " TOF correction (ns): " << tofCorr << endl
93  << endl;
94  }
95  //The global offset is the sum of various contributions
96  return tTrig + wirePropCorr - tofCorr;
97 }
98 
99 double DTTTrigSyncTOFCorr::offset(const DTWireId& wireId) const {
100  // Correction for the float to int conversion
101  // (half a bin on average) in DTDigi constructor
102  static const float f2i_convCorr = (25. / 64.); // ns
103  //FIXME: This should be considered only if the Digi is constructed from a float....
104 
105  // The tTrig is taken from a parameter
106  return theTTrig - f2i_convCorr;
107 }
108 
109 double DTTTrigSyncTOFCorr::emulatorOffset(const DTWireId& wireId, double& tTrig, double& t0cell) const {
110  tTrig = theTTrig;
111  t0cell = 0.;
112 
113  return int(tTrig / theBXspace) * theBXspace;
114 }
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:86
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
double emulatorOffset(const DTWireId &wireId, double &tTrig, double &t0cell) const override
const double theVPropWire
T y() const
Definition: PV3DBase.h:63
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:69
Definition: config.py:1
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
DTTTrigSyncTOFCorr(const edm::ParameterSet &config)
Constructor.
const DTTopology & specificTopology() const
Definition: DTLayer.cc:42
T mag() const
Definition: PV3DBase.h:67
int wire() const
Return the wire number.
Definition: DTWireId.h:56
#define debug
Definition: HDRShower.cc:19
~DTTTrigSyncTOFCorr() override
Destructor.
const DTChamber * chamber() const
Definition: DTLayer.cc:58
double offset(const DTLayer *layer, const DTWireId &wireId, const GlobalPoint &globPos, double &tTrig, double &wirePropCorr, double &tofCorr) const override
float cellLenght() const
Definition: DTTopology.h:73
const PositionType & position() const