CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
DTTTrigSyncTOFCorr Class Reference

#include <DTTTrigSyncTOFCorr.h>

Inheritance diagram for DTTTrigSyncTOFCorr:
DTTTrigBaseSync

Public Member Functions

 DTTTrigSyncTOFCorr (const edm::ParameterSet &config)
 Constructor. More...
 
virtual double emulatorOffset (const DTWireId &wireId, double &tTrig, double &t0cell)
 
virtual double offset (const DTLayer *layer, const DTWireId &wireId, const GlobalPoint &globPos, double &tTrig, double &wirePropCorr, double &tofCorr)
 
virtual double offset (const DTWireId &wireId)
 
virtual void setES (const edm::EventSetup &setup)
 Pass the Event Setup to the algo at each event. More...
 
virtual ~DTTTrigSyncTOFCorr ()
 Destructor. More...
 
- Public Member Functions inherited from DTTTrigBaseSync
 DTTTrigBaseSync ()
 Constructor. More...
 
virtual double emulatorOffset (const DTWireId &wireId)
 
double offset (const DTLayer *layer, const DTWireId &wireId, const GlobalPoint &globalPos)
 
virtual ~DTTTrigBaseSync ()
 Destructor. More...
 

Private Attributes

const bool debug
 
double theBXspace
 
const int theTOFCorrType
 
const double theTTrig
 
const double theVPropWire
 

Detailed Description

Concrete implementation of a DTTTrigBaseSync. This class define the offsets for RecHit building coherently to the digitization realized with the DTDigiSyncTOFCorr module. The offset is computes as:
offset = tTrig + wirePropCorr - tofCorr
where:

The emulatorOffset is computed as:
offset = int(ttrig/BXspace)*BXspace
where:

NOTE: this should approximate what is seen online by the BTI

Author
G. Cerminara - INFN Torino

Definition at line 50 of file DTTTrigSyncTOFCorr.h.

Constructor & Destructor Documentation

DTTTrigSyncTOFCorr::DTTTrigSyncTOFCorr ( const edm::ParameterSet config)

Constructor.

Definition at line 18 of file DTTTrigSyncTOFCorr.cc.

19 :
20  // The fixed t0 (or t_trig) to be subtracted to digi time (ns)
21  theTTrig(config.getParameter<double>("tTrig")), // FIXME: Default was 500 ns
22  // Velocity of signal propagation along the wire (cm/ns)
23  theVPropWire(config.getParameter<double>("vPropWire")), // FIXME: Default was 24.4 cm/ns
24  // Select the mode for TOF correction:
25  // 0: tofCorr = TOF from IP to 3D Hit position (globPos)
26  // 1: tofCorr = TOF correction for distance difference
27  // between 3D center of the chamber and hit position
28  // (default)
29  // 2: tofCorr = TOF correction for distance difference
30  // between 3D center of the wire and hit position
31  // (This mode in available for backward compatibility)
32  theTOFCorrType(config.getParameter<int>("tofCorrType")), // FIXME: Default was 1
33  debug(config.getUntrackedParameter<bool>("debug")),
34  theBXspace(config.getUntrackedParameter<double>("bxSpace", 25.)) // spacing of BX in ns
35 {
36 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const double theVPropWire
DTTTrigSyncTOFCorr::~DTTTrigSyncTOFCorr ( )
virtual

Destructor.

Definition at line 40 of file DTTTrigSyncTOFCorr.cc.

40 {}

Member Function Documentation

double DTTTrigSyncTOFCorr::emulatorOffset ( const DTWireId wireId,
double &  tTrig,
double &  t0cell 
)
virtual

Time (ns) to be subtracted to the digi time for emulation purposes It does not take into account TOF and signal propagation along the wire It also returns the different contributions separately:

  • tTrig is the offset (t_trig)
  • t0cell is the t0 from pulses (always 0 in this case)

Implements DTTTrigBaseSync.

Definition at line 120 of file DTTTrigSyncTOFCorr.cc.

References theBXspace, and theTTrig.

122  {
123  tTrig = theTTrig;
124  t0cell = 0.;
125 
126  return int(tTrig/theBXspace)*theBXspace;
127 }
double DTTTrigSyncTOFCorr::offset ( const DTLayer layer,
const DTWireId wireId,
const GlobalPoint globPos,
double &  tTrig,
double &  wirePropCorr,
double &  tofCorr 
)
virtual

Time (ns) to be subtracted to the digi time, Parameters are the layer and the wireId to which the digi is referred and the estimation of the 3D hit position (globPos) It also returns the different contributions separately:

  • tTrig is the offset (t_trig)
  • wirePropCorr is the delay for signal propagation along the wire
  • tofCorr is the correction due to the particle TOF

Implements DTTTrigBaseSync.

Definition at line 44 of file DTTTrigSyncTOFCorr.cc.

References DTTopology::cellLenght(), DTLayer::chamber(), gather_cfg::cout, debug, Exception, PV3DBase< T, PVType, FrameType >::mag(), mag(), GloballyPositioned< T >::position(), DTLayer::specificTopology(), GeomDet::surface(), theTOFCorrType, theVPropWire, GeomDet::toGlobal(), GeomDet::toLocal(), DTWireId::wire(), DTTopology::wirePosition(), and PV3DBase< T, PVType, FrameType >::y().

49  {
50  tTrig = offset(wireId);
51 
52  //Compute the time spent in signal propagation along wire.
53  // NOTE: the FE is always at y>0
54  float halfL = layer->specificTopology().cellLenght()/2;
55  float wireCoord = layer->toLocal(globPos).y();
56  float propgL = halfL - wireCoord;
57  wirePropCorr = propgL/theVPropWire;
58 
59 
60  // Compute TOF correction treating it accordingly to
61  // the tofCorrType card
62  float flightToHit = globPos.mag();
63  static const float cSpeed = 29.9792458; // cm/ns
64  tofCorr = 0.;
65  switch(theTOFCorrType) {
66  case 0: {
67  // In this mode the subtraction of the TOF from IP to
68  // estimate 3D hit digi position is done here
69  // (No correction is needed anymore)
70  tofCorr = -flightToHit/cSpeed;
71  break;
72  }
73  case 1: {
74  // Correction for TOF from the center of the chamber to hit position
75  const DTChamber * chamber = layer->chamber();
76  double flightToChamber = chamber->surface().position().mag();
77  tofCorr = (flightToChamber-flightToHit)/cSpeed;
78  break;
79  }
80  case 2: {
81  // TOF from 3D center of the wire to hit position
82  float flightToWire =
83  layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
84  tofCorr = (flightToWire-flightToHit)/cSpeed;
85  break;
86  }
87  default: {
88  throw
89  cms::Exception("[DTTTrigSyncTOFCorr]") << " Invalid parameter: tofCorrType = "
90  << theTOFCorrType
91  << std::endl;
92  break;
93  }
94  }
95 
96  if(debug) {
97  cout << "[DTTTrigSyncTOFCorr] Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
98  << " various contributions are: " << endl
99  << " tTrig (ns): " << tTrig << endl
100  << " Propagation along wire delay (ns): " << wirePropCorr << endl
101  << " TOF correction (ns): " << tofCorr << endl
102  << endl;
103  }
104  //The global offset is the sum of various contributions
105  return tTrig + wirePropCorr - tofCorr;
106 }
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:86
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
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
virtual double offset(const DTLayer *layer, const DTWireId &wireId, const GlobalPoint &globPos, double &tTrig, double &wirePropCorr, double &tofCorr)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
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
const DTChamber * chamber() const
Definition: DTLayer.cc:58
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
tuple cout
Definition: gather_cfg.py:145
float cellLenght() const
Definition: DTTopology.h:73
const PositionType & position() const
double DTTTrigSyncTOFCorr::offset ( const DTWireId wireId)
virtual

Time (ns) to be subtracted to the digi time. It does not take into account TOF and signal propagation along the wire

Implements DTTTrigBaseSync.

Definition at line 108 of file DTTTrigSyncTOFCorr.cc.

References theTTrig.

108  {
109  // Correction for the float to int conversion
110  // (half a bin on average) in DTDigi constructor
111  static const float f2i_convCorr = (25./64.); // ns
112  //FIXME: This should be considered only if the Digi is constructed from a float....
113 
114 
115  // The tTrig is taken from a parameter
116  return theTTrig - f2i_convCorr;
117 }
virtual void DTTTrigSyncTOFCorr::setES ( const edm::EventSetup setup)
inlinevirtual

Pass the Event Setup to the algo at each event.

Implements DTTTrigBaseSync.

Definition at line 61 of file DTTTrigSyncTOFCorr.h.

61 {}

Member Data Documentation

const bool DTTTrigSyncTOFCorr::debug
private
double DTTTrigSyncTOFCorr::theBXspace
private

Definition at line 110 of file DTTTrigSyncTOFCorr.h.

Referenced by emulatorOffset().

const int DTTTrigSyncTOFCorr::theTOFCorrType
private

Definition at line 105 of file DTTTrigSyncTOFCorr.h.

Referenced by offset().

const double DTTTrigSyncTOFCorr::theTTrig
private

Definition at line 92 of file DTTTrigSyncTOFCorr.h.

Referenced by emulatorOffset(), and offset().

const double DTTTrigSyncTOFCorr::theVPropWire
private

Definition at line 97 of file DTTTrigSyncTOFCorr.h.

Referenced by offset().