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
DTDigiSyncTOFCorr Class Reference

#include <DTDigiSyncTOFCorr.h>

Inheritance diagram for DTDigiSyncTOFCorr:
DTDigiSyncBase

Public Member Functions

virtual double digitizerOffset (const DTWireId *id, const DTLayer *layer=0) const
 Delays to be added to digi times during digitization, in ns. More...
 
 DTDigiSyncTOFCorr (const edm::ParameterSet &)
 Constructor. More...
 
virtual double emulatorOffset (const DTWireId *id) const
 Offset to obtain "raw" TDCs for the L1 emulator from digis. More...
 
virtual ~DTDigiSyncTOFCorr ()
 Destructor. More...
 
- Public Member Functions inherited from DTDigiSyncBase
 DTDigiSyncBase ()
 Constructor. More...
 
virtual ~DTDigiSyncBase ()
 Destructor. More...
 

Private Attributes

int corrType
 
double theOffset
 

Detailed Description

Digi offset computed as:
t0 = Tcommon - aTOF

where Tcommon is a fixed offset defined in
DTDigiSyncTOFCorr:offset (in ORCA the default was = 500 ns)

and aTOF is set according to MuBarDigiSyncTOFCorr:TOFCorrection:
0: no TOF correction (aTOF=0)
1: aTOF = the TOF of an infinite-momentum particle travelling from the nominal IP to the 3D center of the chamber
2: ditto, but for a particle travelling to the 3D center of the wire. (This mode is avaliable for comparison with older data which were produced in this way) 3: aTOF = the TOF of an infinite-momentum particle travelling from the nominal IP to the 3D center of the SL. Use this, unless you really know what you are doing.

Author
N. Amapane, R. Bellan - INFN Torino

Definition at line 31 of file DTDigiSyncTOFCorr.h.

Constructor & Destructor Documentation

DTDigiSyncTOFCorr::DTDigiSyncTOFCorr ( const edm::ParameterSet pSet)

Constructor.

Definition at line 18 of file DTDigiSyncTOFCorr.cc.

References edm::ParameterSet::getParameter().

18  {
19 
20  theOffset = pSet.getParameter<double>("offset"); //500ns
21  corrType = pSet.getParameter<int>("TOFCorrection"); //1
22 }
T getParameter(std::string const &) const
DTDigiSyncTOFCorr::~DTDigiSyncTOFCorr ( )
virtual

Destructor.

Definition at line 25 of file DTDigiSyncTOFCorr.cc.

25 {}

Member Function Documentation

double DTDigiSyncTOFCorr::digitizerOffset ( const DTWireId id,
const DTLayer layer = 0 
) const
virtual

Delays to be added to digi times during digitization, in ns.

Implements DTDigiSyncBase.

Definition at line 29 of file DTDigiSyncTOFCorr.cc.

References DTLayer::chamber(), gather_cfg::cout, PV3DBase< T, PVType, FrameType >::mag(), hltrates_dqm_sourceclient-live_cfg::offset, GloballyPositioned< T >::position(), DTLayer::specificTopology(), DTLayer::superLayer(), GeomDet::surface(), Surface::toGlobal(), DTWireId::wire(), and DTTopology::wirePosition().

29  {
30 
31  double offset = theOffset;
32  const double cSpeed = 29.9792458; // cm/ns
33 
34  if (corrType==1) {
35  // Subtraction of assumed TOF, per CHAMBER
36  double flightL = layer->chamber()->surface().position().mag();
37  offset -= flightL/cSpeed;
38  } else if (corrType==2) {
39  // Subtraction of assumed TOF, per WIRE
40 
41  // Position of the wire in the Layer's reference frame
42  float localXPos = layer->specificTopology().wirePosition(id->wire());
43  LocalPoint localPos(localXPos,0,0);
44 
45  // Distance of the wire to the CMS's I.P.
46  double flightL = layer->surface().toGlobal(localPos).mag();
47 
48  offset -= flightL/cSpeed;
49 
50  } else if (corrType==3) {
51  // Subtraction of assumed TOF, per SUPERLAYER
52  double flightL = layer->superLayer()->surface().position().mag();
53  offset -= flightL/cSpeed;
54  } else if (corrType!=0){
55  cout << "ERROR: SimMuon:DTDigitizer:DTDigiSyncTOFCorr:TOFCorrection = " << corrType
56  << "is not defined " << endl;
57  }
58  return offset;
59 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:86
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
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
const DTSuperLayer * superLayer() const
Definition: DTLayer.cc:54
tuple cout
Definition: gather_cfg.py:121
const PositionType & position() const
double DTDigiSyncTOFCorr::emulatorOffset ( const DTWireId id) const
virtual

Offset to obtain "raw" TDCs for the L1 emulator from digis.

Implements DTDigiSyncBase.

Definition at line 63 of file DTDigiSyncTOFCorr.cc.

63  {
64  return theOffset;
65 }

Member Data Documentation

int DTDigiSyncTOFCorr::corrType
private

Definition at line 47 of file DTDigiSyncTOFCorr.h.

double DTDigiSyncTOFCorr::theOffset
private

Definition at line 46 of file DTDigiSyncTOFCorr.h.