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.

Date:
2010/07/16 14:44:04
Revision:
1.4
Author
N. Amapane, R. Bellan - INFN Torino

Definition at line 33 of file DTDigiSyncTOFCorr.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 20 of file DTDigiSyncTOFCorr.cc.

References edm::ParameterSet::getParameter().

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

Destructor.

Definition at line 27 of file DTDigiSyncTOFCorr.cc.

27 {}

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 31 of file DTDigiSyncTOFCorr.cc.

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

31  {
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 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:78
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:88
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
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
double DTDigiSyncTOFCorr::emulatorOffset ( const DTWireId id) const
virtual

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

Implements DTDigiSyncBase.

Definition at line 65 of file DTDigiSyncTOFCorr.cc.

65  {
66  return theOffset;
67 }

Member Data Documentation

int DTDigiSyncTOFCorr::corrType
private

Definition at line 49 of file DTDigiSyncTOFCorr.h.

double DTDigiSyncTOFCorr::theOffset
private

Definition at line 48 of file DTDigiSyncTOFCorr.h.