#include <CalibMuon/DTDigiSync/src/DTTTrigSyncFromDB.h>
Public Member Functions | |
DTTTrigSyncFromDB (const edm::ParameterSet &config) | |
Constructor. | |
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). | |
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 | kFactor |
int | theTOFCorrType |
double | theVPropWire |
int | theWirePropCorrType |
const DTTtrig * | tTrigMap |
const DTT0 * | tZeroMap |
Static Private Attributes | |
static bool | debug = false |
This class define the offset for RecHit building of data and simulation. The offset is computes as:
offset = t0 + tTrig + wirePropCorr - tofCorr
where:
Definition at line 47 of file DTTTrigSyncFromDB.h.
DTTTrigSyncFromDB::DTTTrigSyncFromDB | ( | const edm::ParameterSet & | config | ) |
Constructor.
Definition at line 29 of file DTTTrigSyncFromDB.cc.
References debug, doT0Correction, doTOFCorrection, doWirePropCorrection, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), kFactor, theTOFCorrType, theVPropWire, and theWirePropCorrType.
00029 { 00030 debug = config.getUntrackedParameter<bool>("debug"); 00031 // The ttrig is defined as mean + kFactor * sigma 00032 kFactor = config.getParameter<double>("kFactor"); 00033 // The velocity of signal propagation along the wire (cm/ns) 00034 theVPropWire = config.getParameter<double>("vPropWire"); 00035 // Switch on/off the T0 correction from pulses 00036 doT0Correction = config.getParameter<bool>("doT0Correction"); 00037 // Switch on/off the TOF correction for particles from IP 00038 doTOFCorrection = config.getParameter<bool>("doTOFCorrection"); 00039 theTOFCorrType = config.getParameter<int>("tofCorrType"); 00040 // Switch on/off the correction for the signal propagation along the wire 00041 doWirePropCorrection = config.getParameter<bool>("doWirePropCorrection"); 00042 theWirePropCorrType = config.getParameter<int>("wirePropCorrType"); 00043 }
DTTTrigSyncFromDB::~DTTTrigSyncFromDB | ( | ) | [virtual] |
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 170 of file DTTTrigSyncFromDB.cc.
References DTT0::cellT0(), doT0Correction, kFactor, DTTimeUnits::ns, DTTtrig::slTtrig(), DTLayerId::superlayerId(), tTrigMap, and tZeroMap.
00170 { 00171 float t0 = 0; 00172 float t0rms = 0; 00173 if(doT0Correction) 00174 { 00175 // Read the t0 from pulses for this wire (ns) 00176 tZeroMap->cellT0(wireId, 00177 t0, 00178 t0rms, 00179 DTTimeUnits::ns); 00180 00181 } 00182 00183 // Read the ttrig for this wire 00184 float ttrigMean = 0; 00185 float ttrigSigma = 0;//FIXME: should use this! 00186 tTrigMap->slTtrig(wireId.superlayerId(), 00187 ttrigMean, 00188 ttrigSigma, 00189 DTTimeUnits::ns); 00190 00191 return t0 + ttrigMean + kFactor * ttrigSigma; 00192 }
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 73 of file DTTTrigSyncFromDB.cc.
References DTTopology::cellLenght(), GenMuonPlsPt100GeV_cfg::cout, debug, doTOFCorrection, doWirePropCorrection, lat::endl(), Exception, PV3DBase< T, PVType, FrameType >::mag(), muonGeometry::mag(), GloballyPositioned< T >::position(), sl, DTLayer::specificTopology(), DTLayer::superLayer(), GeomDet::surface(), theTOFCorrType, theVPropWire, theWirePropCorrType, GeomDet::toGlobal(), GeomDet::toLocal(), DTWireId::wire(), DTTopology::wirePosition(), and PV3DBase< T, PVType, FrameType >::y().
00078 { 00079 // Correction for the float to int conversion while writeing the ttrig in ns into an int variable 00080 // (half a bin on average) 00081 // FIXME: this should disappear as soon as the ttrig object will become a float 00082 // static const float f2i_convCorr = (25./64.); // ns //FIXME: check how the conversion is performed 00083 00084 tTrig = offset(wireId); 00085 00086 // Compute the time spent in signal propagation along wire. 00087 // NOTE: the FE is always at y>0 00088 wirePropCorr = 0; 00089 if(doWirePropCorrection) { 00090 switch(theWirePropCorrType){ 00091 // The ttrig computed from the timebox accounts on average for the signal propagation time 00092 // from the center of the wire to the frontend. Here we just have to correct for 00093 // the distance of the hit from the wire center. 00094 case 0: { 00095 float wireCoord = layer->toLocal(globPos).y(); 00096 wirePropCorr = -wireCoord/theVPropWire; 00097 break; 00098 // FIXME: What if hits used for the time box are not distributed uniformly along the wire? 00099 } 00100 //On simulated data you need to subtract the total propagation time 00101 case 1: { 00102 float halfL = layer->specificTopology().cellLenght()/2; 00103 float wireCoord = layer->toLocal(globPos).y(); 00104 float propgL = halfL - wireCoord; 00105 wirePropCorr = propgL/theVPropWire; 00106 break; 00107 } 00108 default: { 00109 throw 00110 cms::Exception("[DTTTrigSyncFromDB]") << " Invalid parameter: wirePropCorrType = " 00111 << theWirePropCorrType 00112 << std::endl; 00113 break; 00114 } 00115 } 00116 } 00117 00118 // Compute TOF correction: 00119 tofCorr = 0.; 00120 // TOF Correction can be switched off with appropriate parameter 00121 if(doTOFCorrection) { 00122 float flightToHit = globPos.mag(); 00123 static const float cSpeed = 29.9792458; // cm/ns 00124 switch(theTOFCorrType) { 00125 case 0: { 00126 // The ttrig computed from the real data accounts on average for the TOF correction 00127 // Depending on the granularity used for the ttrig computation we just have to correct for the 00128 // TOF from the center of the chamber, SL, layer or wire to the hit position. 00129 // At the moment only SL granularity is considered 00130 // Correction for TOF from the center of the SL to hit position 00131 const DTSuperLayer *sl = layer->superLayer(); 00132 double flightToSL = sl->surface().position().mag(); 00133 tofCorr = (flightToSL-flightToHit)/cSpeed; 00134 break; 00135 } 00136 case 1: { 00137 // On simulated data you need to consider only the TOF from 3D center of the wire to hit position 00138 // (because the TOF from the IP to the wire has been already subtracted in the digitization: 00139 // SimMuon/DTDigitizer/DTDigiSyncTOFCorr.cc corrType=2) 00140 float flightToWire = 00141 layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag(); 00142 tofCorr = (flightToWire-flightToHit)/cSpeed; 00143 break; 00144 } 00145 default: { 00146 throw 00147 cms::Exception("[DTTTrigSyncFromDB]") << " Invalid parameter: tofCorrType = " 00148 << theTOFCorrType 00149 << std::endl; 00150 break; 00151 } 00152 } 00153 } 00154 00155 00156 if(debug) { 00157 cout << "[DTTTrigSyncFromDB] Channel: " << wireId << endl 00158 << " Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl 00159 << " various contributions are: " << endl 00160 << " tTrig + t0 (ns): " << tTrig << endl 00161 //<< " tZero (ns): " << t0 << endl 00162 << " Propagation along wire delay (ns): " << wirePropCorr << endl 00163 << " TOF correction (ns): " << tofCorr << endl 00164 << endl; 00165 } 00166 //The global offset is the sum of various contributions 00167 return tTrig + wirePropCorr - tofCorr; 00168 }
void DTTTrigSyncFromDB::setES | ( | const edm::EventSetup & | setup | ) | [virtual] |
Pass the Event Setup to the algo at each event.
Implements DTTTrigBaseSync.
Definition at line 51 of file DTTTrigSyncFromDB.cc.
References GenMuonPlsPt100GeV_cfg::cout, debug, doT0Correction, lat::endl(), edm::EventSetup::get(), tTrigMap, tZeroMap, DTTtrig::version(), and DTT0::version().
00051 { 00052 if(doT0Correction) 00053 { 00054 // Get the map of t0 from pulses from the Setup 00055 ESHandle<DTT0> t0Handle; 00056 setup.get<DTT0Rcd>().get(t0Handle); 00057 tZeroMap = &*t0Handle; 00058 if(debug) { 00059 cout << "[DTTTrigSyncFromDB] t0 version: " << tZeroMap->version()<<endl; 00060 } 00061 } 00062 00063 // Get the map of ttrig from the Setup 00064 ESHandle<DTTtrig> ttrigHandle; 00065 setup.get<DTTtrigRcd>().get(ttrigHandle); 00066 tTrigMap = &*ttrigHandle; 00067 if(debug) { 00068 cout << "[DTTTrigSyncFromDB] ttrig version: " << tTrigMap->version() << endl; 00069 } 00070 }
bool DTTTrigSyncFromDB::debug = false [static, private] |
Definition at line 81 of file DTTTrigSyncFromDB.h.
Referenced by DTTTrigSyncFromDB(), offset(), and setES().
bool DTTTrigSyncFromDB::doT0Correction [private] |
Definition at line 87 of file DTTTrigSyncFromDB.h.
Referenced by DTTTrigSyncFromDB(), offset(), and setES().
bool DTTTrigSyncFromDB::doTOFCorrection [private] |
bool DTTTrigSyncFromDB::doWirePropCorrection [private] |
double DTTTrigSyncFromDB::kFactor [private] |
int DTTTrigSyncFromDB::theTOFCorrType [private] |
double DTTTrigSyncFromDB::theVPropWire [private] |
int DTTTrigSyncFromDB::theWirePropCorrType [private] |
const DTTtrig* DTTTrigSyncFromDB::tTrigMap [private] |
const DTT0* DTTTrigSyncFromDB::tZeroMap [private] |