#include <DTTTrigSyncFromDB.h>
Public Member Functions | |
DTTTrigSyncFromDB (const edm::ParameterSet &config) | |
Constructor. | |
virtual double | emulatorOffset (const DTWireId &wireId, double &tTrig, double &t0cell) |
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 | ~DTTTrigSyncFromDB () |
Destructor. | |
Private Attributes | |
bool | doT0Correction |
bool | doTOFCorrection |
bool | doWirePropCorrection |
double | theBXspace |
int | theTOFCorrType |
std::string | thetTrigLabel |
double | theVPropWire |
int | theWirePropCorrType |
const DTTtrig * | tTrigMap |
const DTT0 * | tZeroMap |
Static Private Attributes | |
static bool | debug = false |
Concrete implementation of a DTTTrigBaseSync. This class define the offset for RecHit building of data and simulation. The offset is computes as:
offset = t0 + tTrig + wirePropCorr - tofCorr
where:
The emulatorOffset is computed as:
offset = int(ttrig/BXspace)*BXspace + t0
where:
NOTE: this should approximate what is seen online by the BTI
Definition at line 56 of file DTTTrigSyncFromDB.h.
DTTTrigSyncFromDB::DTTTrigSyncFromDB | ( | const edm::ParameterSet & | config | ) |
Constructor.
Definition at line 29 of file DTTTrigSyncFromDB.cc.
References debug, edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().
{ debug = config.getUntrackedParameter<bool>("debug"); // The velocity of signal propagation along the wire (cm/ns) theVPropWire = config.getParameter<double>("vPropWire"); // Switch on/off the T0 correction from pulses doT0Correction = config.getParameter<bool>("doT0Correction"); // Switch on/off the TOF correction for particles from IP doTOFCorrection = config.getParameter<bool>("doTOFCorrection"); theTOFCorrType = config.getParameter<int>("tofCorrType"); // Switch on/off the correction for the signal propagation along the wire doWirePropCorrection = config.getParameter<bool>("doWirePropCorrection"); theWirePropCorrType = config.getParameter<int>("wirePropCorrType"); // spacing of BX in ns theBXspace = config.getUntrackedParameter<double>("bxSpace", 25.); thetTrigLabel = config.getParameter<string>("tTrigLabel"); }
DTTTrigSyncFromDB::~DTTTrigSyncFromDB | ( | ) | [virtual] |
double DTTTrigSyncFromDB::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:
Implements DTTTrigBaseSync.
Definition at line 208 of file DTTTrigSyncFromDB.cc.
References gather_cfg::cout, DTTimeUnits::ns, and DTLayerId::superlayerId().
{ float t0 = 0; float t0rms = 0; if(doT0Correction) { // Read the t0 from pulses for this wire (ns) tZeroMap->get(wireId, t0, t0rms, DTTimeUnits::ns); } // Read the ttrig for this wire float ttrigMean = 0; float ttrigSigma = 0; float kFactor = 0; // FIXME: should check the return value of the DTTtrigRcd::get(..) method if(tTrigMap->get(wireId.superlayerId(), ttrigMean, ttrigSigma, kFactor, DTTimeUnits::ns) != 0) { cout << "[DTTTrigSyncFromDB]*Error: ttrig not found for SL: " << wireId.superlayerId() << endl; // FIXME: LogError..... } tTrig = ttrigMean + kFactor * ttrigSigma; t0cell = t0; return int(tTrig/theBXspace)*theBXspace + t0cell; }
double DTTTrigSyncFromDB::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)
Implements DTTTrigBaseSync.
Definition at line 75 of file DTTTrigSyncFromDB.cc.
References DTTopology::cellLenght(), gather_cfg::cout, debug, Exception, PV3DBase< T, PVType, FrameType >::mag(), mag(), evf::evtn::offset(), GloballyPositioned< T >::position(), DTLayer::specificTopology(), DTLayer::superLayer(), GeomDet::surface(), GeomDet::toGlobal(), GeomDet::toLocal(), DTWireId::wire(), DTTopology::wirePosition(), and PV3DBase< T, PVType, FrameType >::y().
{ // Correction for the float to int conversion while writeing the ttrig in ns into an int variable // (half a bin on average) // FIXME: this should disappear as soon as the ttrig object will become a float // static const float f2i_convCorr = (25./64.); // ns //FIXME: check how the conversion is performed tTrig = offset(wireId); // Compute the time spent in signal propagation along wire. // NOTE: the FE is always at y>0 wirePropCorr = 0; if(doWirePropCorrection) { switch(theWirePropCorrType){ // The ttrig computed from the timebox accounts on average for the signal propagation time // from the center of the wire to the frontend. Here we just have to correct for // the distance of the hit from the wire center. case 0: { float wireCoord = layer->toLocal(globPos).y(); wirePropCorr = -wireCoord/theVPropWire; break; // FIXME: What if hits used for the time box are not distributed uniformly along the wire? } //On simulated data you need to subtract the total propagation time case 1: { float halfL = layer->specificTopology().cellLenght()/2; float wireCoord = layer->toLocal(globPos).y(); float propgL = halfL - wireCoord; wirePropCorr = propgL/theVPropWire; break; } default: { throw cms::Exception("[DTTTrigSyncFromDB]") << " Invalid parameter: wirePropCorrType = " << theWirePropCorrType << std::endl; break; } } } // Compute TOF correction: tofCorr = 0.; // TOF Correction can be switched off with appropriate parameter if(doTOFCorrection) { float flightToHit = globPos.mag(); static const float cSpeed = 29.9792458; // cm/ns switch(theTOFCorrType) { case 0: { // The ttrig computed from the real data accounts on average for the TOF correction // Depending on the granularity used for the ttrig computation we just have to correct for the // TOF from the center of the chamber, SL, layer or wire to the hit position. // At the moment only SL granularity is considered // Correction for TOF from the center of the SL to hit position const DTSuperLayer *sl = layer->superLayer(); double flightToSL = sl->surface().position().mag(); tofCorr = (flightToSL-flightToHit)/cSpeed; break; } case 1: { // On simulated data you need to consider only the TOF from 3D center of the wire to hit position // (because the TOF from the IP to the wire has been already subtracted in the digitization: // SimMuon/DTDigitizer/DTDigiSyncTOFCorr.cc corrType=2) float flightToWire = layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag(); tofCorr = (flightToWire-flightToHit)/cSpeed; break; } default: { throw cms::Exception("[DTTTrigSyncFromDB]") << " Invalid parameter: tofCorrType = " << theTOFCorrType << std::endl; break; } } } if(debug) { cout << "[DTTTrigSyncFromDB] Channel: " << wireId << endl << " Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl << " various contributions are: " << endl << " tTrig + t0 (ns): " << tTrig << endl //<< " tZero (ns): " << t0 << 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 DTTTrigSyncFromDB::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 172 of file DTTTrigSyncFromDB.cc.
References gather_cfg::cout, DTTimeUnits::ns, and DTLayerId::superlayerId().
{ float t0 = 0; float t0rms = 0; if(doT0Correction) { // Read the t0 from pulses for this wire (ns) tZeroMap->get(wireId, t0, t0rms, DTTimeUnits::ns); } // Read the ttrig for this wire float ttrigMean = 0; float ttrigSigma = 0; float kFactor = 0; // FIXME: should check the return value of the DTTtrigRcd::get(..) method if(tTrigMap->get(wireId.superlayerId(), ttrigMean, ttrigSigma, kFactor, DTTimeUnits::ns) != 0) { cout << "[DTTTrigSyncFromDB]*Error: ttrig not found for SL: " << wireId.superlayerId() << endl; // FIXME: LogError..... } return t0 + ttrigMean + kFactor * ttrigSigma; }
void DTTTrigSyncFromDB::setES | ( | const edm::EventSetup & | setup | ) | [virtual] |
Pass the Event Setup to the algo at each event.
Implements DTTTrigBaseSync.
Definition at line 53 of file DTTTrigSyncFromDB.cc.
References gather_cfg::cout, debug, and edm::EventSetup::get().
{ if(doT0Correction) { // Get the map of t0 from pulses from the Setup ESHandle<DTT0> t0Handle; setup.get<DTT0Rcd>().get(t0Handle); tZeroMap = &*t0Handle; if(debug) { cout << "[DTTTrigSyncFromDB] t0 version: " << tZeroMap->version()<<endl; } } // Get the map of ttrig from the Setup ESHandle<DTTtrig> ttrigHandle; setup.get<DTTtrigRcd>().get(thetTrigLabel,ttrigHandle); tTrigMap = &*ttrigHandle; if(debug) { cout << "[DTTTrigSyncFromDB] ttrig version: " << tTrigMap->version() << endl; } }
bool DTTTrigSyncFromDB::debug = false [static, private] |
Definition at line 101 of file DTTTrigSyncFromDB.h.
bool DTTTrigSyncFromDB::doT0Correction [private] |
Definition at line 105 of file DTTTrigSyncFromDB.h.
bool DTTTrigSyncFromDB::doTOFCorrection [private] |
Definition at line 107 of file DTTTrigSyncFromDB.h.
bool DTTTrigSyncFromDB::doWirePropCorrection [private] |
Definition at line 110 of file DTTTrigSyncFromDB.h.
double DTTTrigSyncFromDB::theBXspace [private] |
Definition at line 113 of file DTTTrigSyncFromDB.h.
int DTTTrigSyncFromDB::theTOFCorrType [private] |
Definition at line 108 of file DTTTrigSyncFromDB.h.
std::string DTTTrigSyncFromDB::thetTrigLabel [private] |
Definition at line 115 of file DTTTrigSyncFromDB.h.
double DTTTrigSyncFromDB::theVPropWire [private] |
Definition at line 103 of file DTTTrigSyncFromDB.h.
int DTTTrigSyncFromDB::theWirePropCorrType [private] |
Definition at line 111 of file DTTTrigSyncFromDB.h.
const DTTtrig* DTTTrigSyncFromDB::tTrigMap [private] |
Definition at line 99 of file DTTTrigSyncFromDB.h.
const DTT0* DTTTrigSyncFromDB::tZeroMap [private] |
Definition at line 98 of file DTTTrigSyncFromDB.h.