CMS 3D CMS Logo

DTTTrigSyncTOFCorr Class Reference

Concrete implementation of a DTTTrigBaseSync. More...

#include <CalibMuon/DTDigiSync/src/DTTTrigSyncTOFCorr.h>

Inheritance diagram for DTTTrigSyncTOFCorr:

DTTTrigBaseSync

List of all members.

Public Member Functions

 DTTTrigSyncTOFCorr (const edm::ParameterSet &config)
 Constructor.
virtual double offset (const DTWireId &wireId)
 Time (ns) to be subtracted to the digi time.
virtual double offset (const DTLayer *layer, const DTWireId &wireId, const GlobalPoint &globPos, double &tTrig, double &wirePropCorr, double &tofCorr)
 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.

virtual void setES (const edm::EventSetup &setup)
 Pass the Event Setup to the algo at each event.
virtual ~DTTTrigSyncTOFCorr ()
 Destructor.

Static Private Attributes

static bool debug
static int theTOFCorrType
static double theTTrig
static 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:

Date
2007/02/19 11:45:21
Revision
1.1
Author:
G. Cerminara - INFN Torino

Definition at line 42 of file DTTTrigSyncTOFCorr.h.


Constructor & Destructor Documentation

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

Constructor.

Definition at line 21 of file DTTTrigSyncTOFCorr.cc.

References debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), theTOFCorrType, theTTrig, and theVPropWire.

00021                                                                    {
00022   theTTrig = config.getParameter<double>("tTrig"); // FIXME: Default was 500 ns
00023   theVPropWire = config.getParameter<double>("vPropWire"); // FIXME: Default was 24.4 cm/ns
00024   theTOFCorrType = config.getParameter<int>("tofCorrType"); // FIXME: Default was 1
00025   debug = config.getUntrackedParameter<bool>("debug");
00026 
00027 }

DTTTrigSyncTOFCorr::~DTTTrigSyncTOFCorr (  )  [virtual]

Destructor.

Definition at line 31 of file DTTTrigSyncTOFCorr.cc.

00031 {}


Member Function Documentation

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 99 of file DTTTrigSyncTOFCorr.cc.

References theTTrig.

00099                                                         {
00100   // Correction for the float to int conversion
00101   // (half a bin on average) in DTDigi constructor
00102   static const float f2i_convCorr = (25./64.); // ns
00103   //FIXME: This should be considered only if the Digi is constructed from a float....
00104 
00105 
00106   // The tTrig is taken from a parameter
00107   return theTTrig - f2i_convCorr;
00108 }

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:

Implements DTTTrigBaseSync.

Definition at line 35 of file DTTTrigSyncTOFCorr.cc.

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

00040                                                    {
00041   tTrig = offset(wireId);
00042 
00043   //Compute the time spent in signal propagation along wire.
00044   // NOTE: the FE is always at y>0
00045   float halfL     = layer->specificTopology().cellLenght()/2;
00046   float wireCoord = layer->toLocal(globPos).y();
00047   float propgL    = halfL - wireCoord;
00048   wirePropCorr = propgL/theVPropWire;
00049 
00050 
00051   // Compute TOF correction treating it accordingly to
00052   // the tofCorrType card
00053   float flightToHit = globPos.mag();
00054   static const float cSpeed = 29.9792458; // cm/ns
00055   tofCorr = 0.;
00056   switch(theTOFCorrType) {
00057   case 0: {
00058     // In this mode the subtraction of the TOF from IP to
00059     // estimate 3D hit digi position is done here
00060     // (No correction is needed anymore)
00061     tofCorr = -flightToHit/cSpeed;
00062     break;
00063   }
00064   case 1: {
00065     // Correction for TOF from the center of the chamber to hit position
00066     const DTChamber * chamber = layer->chamber();
00067     double flightToChamber = chamber->surface().position().mag();
00068     tofCorr = (flightToChamber-flightToHit)/cSpeed;
00069     break;
00070   }
00071   case 2: {
00072     // TOF from 3D center of the wire to hit position
00073     float flightToWire =
00074       layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
00075     tofCorr = (flightToWire-flightToHit)/cSpeed;
00076     break;
00077   }
00078   default: {
00079     throw
00080       cms::Exception("[DTTTrigSyncTOFCorr]") << " Invalid parameter: tofCorrType = "
00081                                              << theTOFCorrType 
00082                                              << std::endl;
00083     break;
00084   }
00085   }
00086 
00087   if(debug) {
00088     cout << "[DTTTrigSyncTOFCorr] Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
00089          << "      various contributions are: " << endl
00090          << "      tTrig (ns):   " << tTrig << endl
00091          << "      Propagation along wire delay (ns): " <<  wirePropCorr << endl
00092          << "      TOF correction (ns): " << tofCorr << endl
00093          << endl;
00094   }
00095   //The global offset is the sum of various contributions
00096   return tTrig + wirePropCorr - tofCorr;
00097 }

virtual void DTTTrigSyncTOFCorr::setES ( const edm::EventSetup setup  )  [inline, virtual]

Pass the Event Setup to the algo at each event.

Implements DTTTrigBaseSync.

Definition at line 53 of file DTTTrigSyncTOFCorr.h.

00053 {}


Member Data Documentation

bool DTTTrigSyncTOFCorr::debug [static, private]

Definition at line 91 of file DTTTrigSyncTOFCorr.h.

Referenced by DTTTrigSyncTOFCorr(), and offset().

int DTTTrigSyncTOFCorr::theTOFCorrType [static, private]

Definition at line 88 of file DTTTrigSyncTOFCorr.h.

Referenced by DTTTrigSyncTOFCorr(), and offset().

double DTTTrigSyncTOFCorr::theTTrig [static, private]

Definition at line 75 of file DTTTrigSyncTOFCorr.h.

Referenced by DTTTrigSyncTOFCorr(), and offset().

double DTTTrigSyncTOFCorr::theVPropWire [static, private]

Definition at line 80 of file DTTTrigSyncTOFCorr.h.

Referenced by DTTTrigSyncTOFCorr(), and offset().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:19:10 2009 for CMSSW by  doxygen 1.5.4