CMS 3D CMS Logo

Public Member Functions | Private Attributes | Static Private Attributes

DTTTrigSyncTOFCorr Class Reference

#include <DTTTrigSyncTOFCorr.h>

Inheritance diagram for DTTTrigSyncTOFCorr:
DTTTrigBaseSync

List of all members.

Public Member Functions

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

Private Attributes

double theBXspace

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:

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

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

Date:
2009/10/21 17:05:47
Revision:
1.2
Author:
G. Cerminara - INFN Torino

Definition at line 52 of file DTTTrigSyncTOFCorr.h.


Constructor & Destructor Documentation

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

Constructor.

Definition at line 22 of file DTTTrigSyncTOFCorr.cc.

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

                                                                   {
  theTTrig = config.getParameter<double>("tTrig"); // FIXME: Default was 500 ns
  theVPropWire = config.getParameter<double>("vPropWire"); // FIXME: Default was 24.4 cm/ns
  theTOFCorrType = config.getParameter<int>("tofCorrType"); // FIXME: Default was 1
  debug = config.getUntrackedParameter<bool>("debug");
  // spacing of BX in ns
  theBXspace  = config.getUntrackedParameter<double>("bxSpace", 25.);

}
DTTTrigSyncTOFCorr::~DTTTrigSyncTOFCorr ( ) [virtual]

Destructor.

Definition at line 34 of file DTTTrigSyncTOFCorr.cc.

{}

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

                                                         {
  tTrig = theTTrig;
  t0cell = 0.;

  return int(tTrig/theBXspace)*theBXspace;
}
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 38 of file DTTTrigSyncTOFCorr.cc.

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

                                                   {
  tTrig = offset(wireId);

  //Compute the time spent in signal propagation along wire.
  // NOTE: the FE is always at y>0
  float halfL     = layer->specificTopology().cellLenght()/2;
  float wireCoord = layer->toLocal(globPos).y();
  float propgL    = halfL - wireCoord;
  wirePropCorr = propgL/theVPropWire;


  // Compute TOF correction treating it accordingly to
  // the tofCorrType card
  float flightToHit = globPos.mag();
  static const float cSpeed = 29.9792458; // cm/ns
  tofCorr = 0.;
  switch(theTOFCorrType) {
  case 0: {
    // In this mode the subtraction of the TOF from IP to
    // estimate 3D hit digi position is done here
    // (No correction is needed anymore)
    tofCorr = -flightToHit/cSpeed;
    break;
  }
  case 1: {
    // Correction for TOF from the center of the chamber to hit position
    const DTChamber * chamber = layer->chamber();
    double flightToChamber = chamber->surface().position().mag();
    tofCorr = (flightToChamber-flightToHit)/cSpeed;
    break;
  }
  case 2: {
    // TOF from 3D center of the wire to hit position
    float flightToWire =
      layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
    tofCorr = (flightToWire-flightToHit)/cSpeed;
    break;
  }
  default: {
    throw
      cms::Exception("[DTTTrigSyncTOFCorr]") << " Invalid parameter: tofCorrType = "
                                             << theTOFCorrType 
                                             << std::endl;
    break;
  }
  }

  if(debug) {
    cout << "[DTTTrigSyncTOFCorr] Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
         << "      various contributions are: " << endl
         << "      tTrig (ns):   " << tTrig << endl
         << "      Propagation along wire delay (ns): " <<  wirePropCorr << endl
         << "      TOF correction (ns): " << tofCorr << endl
         << endl;
  }
  //The global offset is the sum of various contributions
  return tTrig + wirePropCorr - tofCorr;
}
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 102 of file DTTTrigSyncTOFCorr.cc.

                                                        {
  // Correction for the float to int conversion
  // (half a bin on average) in DTDigi constructor
  static const float f2i_convCorr = (25./64.); // ns
  //FIXME: This should be considered only if the Digi is constructed from a float....


  // The tTrig is taken from a parameter
  return theTTrig - f2i_convCorr;
}
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 63 of file DTTTrigSyncTOFCorr.h.

{}

Member Data Documentation

bool DTTTrigSyncTOFCorr::debug [static, private]

Definition at line 110 of file DTTTrigSyncTOFCorr.h.

Definition at line 112 of file DTTTrigSyncTOFCorr.h.

int DTTTrigSyncTOFCorr::theTOFCorrType [static, private]

Definition at line 107 of file DTTTrigSyncTOFCorr.h.

double DTTTrigSyncTOFCorr::theTTrig [static, private]

Definition at line 94 of file DTTTrigSyncTOFCorr.h.

double DTTTrigSyncTOFCorr::theVPropWire [static, private]

Definition at line 99 of file DTTTrigSyncTOFCorr.h.