ATLAS Offline Software
Loading...
Searching...
No Matches
MdtCalibrationTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
23
24#include "GeoModelKernel/throwExcept.h"
25
28using ToolSettings = MdtCalibrationTool::ToolSettings;
29
30
32 ToolSettings settings{};
33 using Property = ToolSettings::Property;
34 settings.setBit(Property::TofCorrection, m_doTof);
35 settings.setBit(Property::PropCorrection, m_doProp);
36 settings.setBit(Property::TempCorrection, m_doTemp);
37 settings.setBit(Property::MagFieldCorrection, m_doField);
38 settings.setBit(Property::SlewCorrection, m_doSlew);
39 settings.setBit(Property::BackgroundCorrection, m_doBkg);
40 settings.window = static_cast<timeWindowMode>(m_windowSetting.value());
41 return settings;
42}
44 ATH_MSG_DEBUG( "Initializing" );
45
46 switch(m_windowSetting.value()) {
47 case timeWindowMode::UserDefined:
48 ATH_MSG_DEBUG("Use predefined user values of "<<m_timeWindowLowerBound<<" & "<<m_timeWindowUpperBound);
49 break;
50 case timeWindowMode::Default:
51 ATH_MSG_DEBUG("Use 1000. & 2000. as the lower and upper time window values ");
54 break;
55 case timeWindowMode::CollisionG4:
56 ATH_MSG_DEBUG("Use Geant4 collision time window of 20-30");
59 break;
60 case timeWindowMode::CollisionData:
61 ATH_MSG_DEBUG("Use collision data time window of 10 to 30");
64 break;
65 case timeWindowMode::CollisionFitT0:
66 ATH_MSG_DEBUG("Use collision data time window of 50 to 100 to fit T0 in the end");
69 break;
70 default:
71 ATH_MSG_FATAL("Unknown time window setting "<<m_windowSetting<<" provided.");
72 return StatusCode::FAILURE;
73 };
74
75 ATH_CHECK(m_idHelperSvc.retrieve());
78 ATH_CHECK(m_calibDbKey.initialize());
80 ATH_CHECK(m_t0ShiftTool.retrieve(EnableTool{m_doT0Shift}));
81 ATH_CHECK(m_tMaxShiftTool.retrieve(EnableTool{m_doTMaxShift}));
82 ATH_MSG_DEBUG("Initialization finalized "<<std::endl
83 <<" TimeWindow: ["<< m_timeWindowLowerBound.value()<<";"<<m_timeWindowUpperBound.value()<<"]"<<std::endl
84 <<" Correct time of flight "<<(m_doTof ? "yay" : "nay")<<std::endl
85 <<" Correct propagation time "<<(m_doProp ? "si" : "no")<<std::endl
86 <<" Correct temperature "<<(m_doTemp ? "si" : "no")<<std::endl
87 <<" Correct magnetic field "<<(m_doField ? "si" : "no")<<std::endl
88 <<" Correct time slew "<<(m_doSlew ? "si" : "no")<<std::endl
89 <<" Correct background "<<(m_doBkg ? "si" : "no"));
90 return StatusCode::SUCCESS;
91} //end MdtCalibrationTool::initialize
92
93
95 const Identifier& channelId) const {
96 const MuonCalib::MdtCalibDataContainer* calibData{nullptr};
97 return SG::get(calibData, m_calibDbKey, ctx).isSuccess() ?
98 calibData->getCalibData(channelId, msgStream()) : nullptr;
99}
101 const MdtCalibInput& calibIn,
102 bool resolFromRtrack) const {
103
104 const Identifier& id{calibIn.identify()};
105
106 const MuonCalib::MdtCalibDataContainer* calibData{nullptr};
107 if (ATH_UNLIKELY(!SG::get(calibData, m_calibDbKey, ctx).isSuccess())) {
108 THROW_EXCEPTION("Failed to retrieve the Mdt calibration constants "<<m_calibDbKey.fullKey());
109 }
110
111 const MuonCalib::MdtFullCalibData* calibConstants = calibData->getCalibData(id, msgStream());
112 if (!calibConstants) {
113 ATH_MSG_WARNING("Could not find calibration data for channel "<<m_idHelperSvc->toString(id));
114 return MdtCalibOutput{};
115 }
116
117 // require at least the MdtRtRelation to be available
118 const RtRelationPtr& rtRelation{calibConstants->rtRelation};
119 // Hardcoded MDT tube radius 14.6mm here - not correct for sMDT
120 // on the other hand it should be rare that a tube does not have an RT
121 if(!rtRelation) {
122 ATH_MSG_WARNING("No rtRelation found, cannot calibrate tube "<<m_idHelperSvc->toString(id));
123 return MdtCalibOutput{};
124 }
125 if (!calibConstants->tubeCalib) {
126 ATH_MSG_WARNING("Cannot extract the single tube calibrations for tube "<<m_idHelperSvc->toString(id));
127 return MdtCalibOutput{};
128 }
130 const SingleTubeCalib* singleTubeData = calibConstants->tubeCalib->getCalib(id);
131 if (!singleTubeData) {
132 ATH_MSG_WARNING("Failed to access tubedata for " << m_idHelperSvc->toString(id));
133 return MdtCalibOutput{};
134 }
135 const float invPropSpeed = calibData->inversePropSpeed();
136
137 MdtCalibOutput calibResult{};
138 // correct for global t0 of rt-region
139
140 calibResult.setTubeT0(singleTubeData->t0);
141 calibResult.setMeanAdc(singleTubeData->adcCal);
142
143 // set propagation delay
144 if (m_doProp) {
145 const double propagationDistance = calibIn.signalPropagationDistance();
146 ATH_MSG_VERBOSE("Calibration of "<<m_idHelperSvc->toString(id)<<", propagation distance: "<<propagationDistance<<" -> "
147 <<(invPropSpeed * propagationDistance));
148 calibResult.setPropagationTime(invPropSpeed * propagationDistance);
149 }
150
152 const double driftTime = calibIn.tdc() * tdcBinSize
153 - (m_doTof ? calibIn.timeOfFlight() : 0.)
154 - calibIn.triggerTime()
155 - calibResult.tubeT0()
156 - calibResult.signalPropagationTime();
157
158 calibResult.setDriftTime(driftTime);
159 // apply corrections
160 double corrTime{0.};
161 const bool doCorrections = m_doField || m_doTemp || m_doBkg;
162 if (doCorrections) {
163 const CorrectionPtr& corrections{calibConstants->corrections};
164 const RtRelationPtr& rtRelation{calibConstants->rtRelation};
165 const MuonCalib::IRtRelation* rt = rtRelation->rt();
166 ATH_MSG_VERBOSE("There are correction functions.");
168 if (m_doSlew && corrections->slewing()) {
169 double slewTime=corrections->slewing()->correction(calibResult.driftTime(), calibIn.adc());
170 corrTime -= slewTime;
171 calibResult.setSlewingTime(slewTime);
172 }
173
174
175 if (m_doField && corrections->bField()) {
176 MagField::AtlasFieldCache fieldCache{};
177 const AtlasFieldCacheCondObj* bFieldCondCache{nullptr};
178 if (ATH_UNLIKELY(!SG::get(bFieldCondCache, m_fieldCacheCondObjInputKey, ctx).isSuccess())) {
179 THROW_EXCEPTION("calibrate: Failed to retrieve AtlasFieldCacheCondObj with key " << m_fieldCacheCondObjInputKey.key());
180 }
181 bFieldCondCache->getInitializedCache(fieldCache);
182
183 Amg::Vector3D globalB{Amg::Vector3D::Zero()};
184 fieldCache.getField(calibIn.closestApproach().data(), globalB.data());
185 const Amg::Vector2D locBField = calibIn.projectMagneticField(globalB);
186 using BFieldComp = MdtCalibInput::BFieldComp;
187 calibResult.setLorentzTime(corrections->bField()->correction(calibResult.driftTime(),
188 locBField[static_cast<int>(BFieldComp::alongWire)],
189 locBField[static_cast<int>(BFieldComp::alongTrack)]));
190 corrTime -= calibResult.lorentzTime();
191 }
192 if(m_doTemp && rt && rt->hasTmaxDiff()) {
193 const MdtIdHelper& id_helper{m_idHelperSvc->mdtIdHelper()};
194 const int mL = id_helper.multilayer(id);
195 const double tempTime = MuonCalib::RtScaleFunction(calibResult.driftTime(), mL == 2, *rt);
196 calibResult.setTemperatureTime(tempTime);
197 corrTime-=calibResult.temperatureTime();
198 }
199 // background corrections (I guess this is never active)
200 if (m_doBkg && corrections->background()) {
201 double bgLevel{0.};
202 calibResult.setBackgroundTime(corrections->background()->correction(calibResult.driftTime(), bgLevel ));
203 corrTime += calibResult.backgroundTime();
204 }
205 }
206
207 calibResult.setDriftTime(calibResult.driftTime() + corrTime);
208
209 // calculate drift radius + error
210 double r{0.}, reso{0.};
211 double t = calibResult.driftTime();
212 double t_inrange = t;
213 Muon::MdtDriftCircleStatus timeStatus = driftTimeStatus(t, *rtRelation);
214
215 assert(rtRelation->rt() != nullptr);
216 r = rtRelation->rt()->radius(t);
217 // apply tUpper gshift
218 if (m_doTMaxShift) {
219 float tShift = m_tMaxShiftTool->getValue(id);
220 r = rtRelation->rt()->radius( t * (1 + tShift) );
221 }
222 // check whether drift times are within range, if not fix them to the min/max range
223 if ( t < rtRelation->rt()->tLower()) {
224 t_inrange = rtRelation->rt()->tLower();
225 double rmin = rtRelation->rt()->radius( t_inrange );
226 double drdt = rtRelation->rt()->driftVelocity( t_inrange);
228 if (timeStatus == Muon::MdtStatusBeforeSpectrum) {
229 t = rtRelation->rt()->tLower() - m_timeWindowLowerBound;
230 }
232 r = rmin + drdt*(t-t_inrange);
233 } else if( t > rtRelation->rt()->tUpper() ) {
234 t_inrange = rtRelation->rt()->tUpper();
235 double rmax = rtRelation->rt()->radius( t_inrange );
236 double drdt = rtRelation->rt()->driftVelocity(t_inrange);
237 // now check whether we are outside the time window
238 if ( timeStatus == Muon::MdtStatusAfterSpectrum ) {
239 t = rtRelation->rt()->tUpper() + m_timeWindowUpperBound;
240 }
242 r = rmax + drdt*(t-t_inrange);
243 }
244
245 assert(rtRelation->rtRes() != nullptr);
246
249 ATH_MSG_VERBOSE("Calibrated radius below physical limit "<<r<<" vs. "<<calibIn.innerTubeR());
252 reso = std::pow(calibIn.innerTubeR(), 2);
253 } else if (r > calibIn.innerTubeR()) {
254 ATH_MSG_VERBOSE("Calibrated radius outside tube ."<<r<<" vs. "<<calibIn.innerTubeR());
255 timeStatus = Muon::MdtStatusAfterSpectrum;
256 r = calibIn.innerTubeR();
257 reso = std::pow(calibIn.innerTubeR(), 2);
258 } else if (!resolFromRtrack) {
259 reso = rtRelation->rtRes()->resolution( t_inrange );
260 } else {
261 const std::optional<double> tFromR = rtRelation->tr()->driftTime(std::abs(calibIn.distanceToTrack()));
262 reso = rtRelation->rtRes()->resolution(tFromR.value_or(0.));
263 }
264
265 if (m_doPropUncert && !calibIn.trackDirHasPhi()) {
266 assert(rtRelation->rt() != nullptr);
267 const double driftTimeUp = std::min(rtRelation->rt()->tUpper(),
268 calibIn.tdc() * tdcBinSize
269 - (m_doTof ? calibIn.timeOfFlight() : 0.)
270 - calibIn.triggerTime()
271 - calibResult.tubeT0());
272
273 const double driftTimeDn = std::max(rtRelation->rt()->tLower(),
274 calibIn.tdc() * tdcBinSize
275 - (m_doTof ? calibIn.timeOfFlight() : 0.)
276 - calibIn.triggerTime()
277 - calibResult.tubeT0()
278 - calibIn.tubeLength() * invPropSpeed);
279
280 const double radiusUp = rtRelation->rt()->radius(driftTimeUp);
281 const double radiusDn = rtRelation->rt()->radius(driftTimeDn);
282 ATH_MSG_VERBOSE("Measurement "<<m_idHelperSvc->toString(calibIn.identify())
283 <<" nominal drift time "<<calibResult.driftTime()<<", down: "<<driftTimeDn<<", up: "<<driftTimeUp
284 <<" --> driftRadius: "<<r<<" pm "<<reso<<", prop-up: "<<radiusUp<<", prop-dn: "<<radiusDn
285 <<" delta: "<<(radiusUp-radiusDn));
286 calibResult.setDriftUncertSigProp(0.5*std::abs(radiusUp - radiusDn));
287 reso = std::hypot(reso, calibResult.driftUncertSigProp());
288 }
289
290 calibResult.setDriftRadius(r, reso);
291 calibResult.setStatus(timeStatus);
292 // summary
293 ATH_MSG_VERBOSE( "Calibration for tube " << m_idHelperSvc->toString(id)
294 <<" passed. "<<std::endl<<"Input: "<<calibIn<<std::endl<<"Extracted calib constants: "<<calibResult<<std::endl);
295 return calibResult;
296} //end MdtCalibrationTool::calibrate
297
299 MdtCalibInput&& primHit,
300 MdtCalibInput&& twinHit) const {
301
302 MdtCalibOutput primResult = calibrate(ctx, primHit);
303 MdtCalibOutput twinResult = calibrate(ctx, twinHit);
304
305 // get Identifier and MdtReadOutElement for twin tubes
306 // get 'raw' drifttimes of twin pair; we don't use timeofFlight or propagationTime cause they are irrelevant for twin coordinate
307 double primdriftTime = primHit.tdc()*tdcBinSize - primResult.tubeT0();
308 double twinDriftTime = twinHit.tdc()*tdcBinSize - twinResult.tubeT0();
310 if (primdriftTime >= twinDriftTime) {
311 ATH_MSG_VERBOSE("Swap "<<m_idHelperSvc->toString(primHit.identify())<<" & "<<m_idHelperSvc->toString(twinHit.identify())
312 <<" primDriftTime: "<<primdriftTime<<", secondTime: "<<twinDriftTime);
313 std::swap(primdriftTime, twinDriftTime);
314 std::swap(primResult, twinResult);
315 std::swap(primHit, twinHit);
316 }
317 const Identifier& primId = primHit.identify();
318 const Identifier& twinId = twinHit.identify();
319
320 // get calibration constants from DbTool
321 const MuonCalib::MdtCalibDataContainer* constantHandle{nullptr};
322 if (ATH_UNLIKELY(!SG::get(constantHandle, m_calibDbKey, ctx).isSuccess())) {
323 THROW_EXCEPTION("Failed to retrieve the Mdt calibration constants "<<m_calibDbKey.fullKey());
324 }
325
326
327 const MuonCalib::MdtFullCalibData* data1st = constantHandle->getCalibData(primId, msgStream());
328 const MuonCalib::MdtFullCalibData* data2nd = constantHandle->getCalibData(twinId, msgStream());
329 if (!data1st || !data2nd) {
330 ATH_MSG_WARNING(__FILE__<<":"<<__LINE__<<" Failed to access calibration constants for tubes "<<
331 m_idHelperSvc->toString(primId)<<" & "<<m_idHelperSvc->toString(twinId));
332 return MdtCalibTwinOutput{};
333 }
334 const SingleTubeCalib* calibSingleTube1st = data1st->tubeCalib->getCalib(primId);
335 const SingleTubeCalib* calibSingleTube2nd = data2nd->tubeCalib->getCalib(twinId);
336 const double invPropSpeed = constantHandle->inversePropSpeed();
337 if (!calibSingleTube1st || !calibSingleTube2nd) {
338 ATH_MSG_WARNING(__FILE__<<":"<<__LINE__<<" Failed to access calibration constants for tubes "<<
339 m_idHelperSvc->toString(primId)<<" & "<<m_idHelperSvc->toString(twinId));
340 return MdtCalibTwinOutput{};
341 }
342
343 // get tubelength and set HV-delay (~6ns)
344 constexpr double HVdelay = 6.;
345
347 double twin_timedif = twinDriftTime - primdriftTime - invPropSpeed * twinHit.tubeLength() - HVdelay;
350
351 const double tubeHalfLength = 0.5*primHit.tubeLength();
352 const double zTwin = std::clamp(0.5* primHit.readOutSide()* twin_timedif / invPropSpeed, -tubeHalfLength, tubeHalfLength);
353 const double errZTwin = m_resTwin / invPropSpeed;
354
355 ATH_MSG_VERBOSE( "Twin calibration - tube: " << m_idHelperSvc->toString(primId)<< " twintube: " << m_idHelperSvc->toString(twinId)<<endmsg
356 << " prompthit tdc = " << primHit.tdc() << " twinhit tdc = " << twinHit.tdc()
357 << " tube driftTime = " << primResult<< " second tube driftTime = " << twinResult<<endmsg
358 << " Time difference =" << twin_timedif << " zTwin=" << zTwin<<", errorZ="<<errZTwin);
359 MdtCalibTwinOutput calibResult{primHit,twinHit,primResult, twinResult};
360
361
362 calibResult.setLocZ(zTwin, errZTwin);
363 return calibResult;
364}
365
367 const MuonCalib::MdtRtRelation& rtRelation ) const {
368 if (rtRelation.rt()) {
369 if(driftTime < rtRelation.rt()->tLower() - m_timeWindowLowerBound) {
370 ATH_MSG_VERBOSE( " drift time outside time window "
371 << driftTime << ". Mininum time = "
372 << rtRelation.rt()->tLower() - m_timeWindowLowerBound );
374 } else if (driftTime > rtRelation.rt()->tUpper() + m_timeWindowUpperBound) {
375 ATH_MSG_VERBOSE( " drift time outside time window "
376 << driftTime << ". Maximum time = "
377 << rtRelation.rt()->tUpper() + m_timeWindowUpperBound);
379 }
380 } else {
381 ATH_MSG_WARNING( "No valid rt relation supplied for driftTimeStatus method" );
383 }
385}
386double MdtCalibrationTool::getResolutionFromRt(const EventContext& ctx, const Identifier& moduleID, const double time) const {
387
388 const MuonCalib::MdtFullCalibData* moduleConstants = getCalibConstants(ctx, moduleID);
389 if (!moduleConstants){
390 THROW_EXCEPTION("Failed to retrieve set of calibration constants for "<<m_idHelperSvc->toString(moduleID));
391 }
392 const RtRelationPtr& rtRel{moduleConstants->rtRelation};
393 if (!rtRel) {
394 THROW_EXCEPTION("No rt-relation found for "<<m_idHelperSvc->toString(moduleID));
395 }
396 const double t = std::min(std::max(time, rtRel->rt()->tLower()), rtRel->rt()->tUpper());
397 return rtRel->rtRes()->resolution(t);
398}
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define ATH_UNLIKELY(x)
MdtCalibOutput::MdtDriftCircleStatus MdtDriftCircleStatus
MuonCalib::MdtTubeCalibContainer::SingleTubeCalib SingleTubeCalib
MdtCalibrationTool::ToolSettings ToolSettings
void getInitializedCache(MagField::AtlasFieldCache &cache) const
get B field cache for evaluation as a function of 2-d or 3-d position.
Local cache for magnetic field (based on MagFieldServices/AtlasFieldSvcTLS.h)
void getField(const double *ATH_RESTRICT xyz, double *ATH_RESTRICT bxyz, double *ATH_RESTRICT deriv=nullptr)
get B field value at given position xyz[3] is in mm, bxyz[3] is in kT if deriv[9] is given,...
double timeOfFlight() const
Returns the time of flight.
int16_t tdc() const
Returns the tdc counts of the hit.
bool trackDirHasPhi() const
Returns whether the track has a phi constaint or not.
double distanceToTrack() const
Returns the distance to track (signed)
double triggerTime() const
Returns the trigger offset time.
double signalPropagationDistance() const
Calculates the distance that the signal has to travel along the wire.
int16_t adc() const
Returns the amount of accumulated charge.
const Amg::Vector3D & closestApproach() const
Returns the point of closest approach to the wire.
const Identifier & identify() const
Returns the Identifier of the hit.
double tubeLength() const
Returns the tube length.
double innerTubeR() const
Returns the inner tube radius.
Amg::Vector2D projectMagneticField(const Amg::Vector3D &fieldInGlob) const
Splits the B-field into the components that point along the transverse track direction & along the tu...
double backgroundTime() const
Return the time correction arising from background processes.
void setSlewingTime(const double slewTime)
Sets the slewing time.
void setDriftTime(const double driftTime)
Sets the drift time.
void setBackgroundTime(const double bkgTime)
Sets the background time correction.
void setStatus(const MdtDriftCircleStatus stat)
double driftTime() const
Returns the drift time inside the tube.
void setTemperatureTime(const double tempTime)
Sets the temperature time correction.
void setDriftUncertSigProp(const double uncert)
Sets the uncertainty on the drift radius arising from the unknown position along the wires.
void setTubeT0(const double T0)
Sets the tube T0.
double driftUncertSigProp() const
Returns the uncertainty on the drift radius arising from the unknown position along the wire.
void setLorentzTime(const double time)
Sets the Lorentz time.
void setDriftRadius(const double radius, const double uncert)
Sets the charge drift radius and its associated uncertainty.
Muon::MdtDriftCircleStatus MdtDriftCircleStatus
void setMeanAdc(const double adc)
Sets the mean tube adc.
double temperatureTime() const
Returns the time corrections stemming from temperature & pressure corrections.
double signalPropagationTime() const
Returns the signal propagation time.
double lorentzTime() const
Returns the time corrections from the signal propgation inside a magnetic field.
double tubeT0() const
Returns the point in time where the muon typically enters the chamber.
void setPropagationTime(const double T0)
Sets the signal propagation time in the tube wire.
void setLocZ(const double locZ, const double locZuncert)
ToolHandle< MuonCalib::IShiftMapTools > m_tMaxShiftTool
Gaudi::Property< bool > m_doPropUncert
Gaudi::Property< bool > m_doSlew
Gaudi::Property< double > m_timeWindowLowerBound
Gaudi::Property< bool > m_doTMaxShift
Gaudi::Property< int > m_windowSetting
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
Gaudi::Property< bool > m_doProp
SG::ReadCondHandleKey< MuonCalib::MdtCalibDataContainer > m_calibDbKey
Gaudi::Property< double > m_unphysicalHitRadiusLowerBound
Muon::MdtDriftCircleStatus driftTimeStatus(double driftTime, const MuonCalib::MdtRtRelation &rtRelation) const
virtual double getResolutionFromRt(const EventContext &ctx, const Identifier &module, const double time) const override final
Gaudi::Property< double > m_timeWindowUpperBound
ToolHandle< MuonCalib::IShiftMapTools > m_t0ShiftTool
Gaudi::Property< bool > m_doBkg
virtual const MuonCalib::MdtFullCalibData * getCalibConstants(const EventContext &ctx, const Identifier &channelId) const override final
virtual MdtCalibOutput calibrate(const EventContext &ctx, const MdtCalibInput &hit, bool resolFromRtrack=false) const override final
Convert the raw MDT time (+charge) into a drift radius + error.
virtual ToolSettings getSettings() const override final
Gaudi::Property< double > m_resTwin
Gaudi::Property< bool > m_doTemp
MuonCalib::MdtFullCalibData::RtRelationPtr RtRelationPtr
virtual MdtCalibTwinOutput calibrateTwinTubes(const EventContext &ctx, MdtCalibInput &&hit, MdtCalibInput &&twinHit) const override final
Convert the raw MDT times of two twin hits into a Twin position (coordinate along tube) It returns wh...
virtual StatusCode initialize() override final
initialization
Gaudi::Property< bool > m_doField
MuonCalib::MdtFullCalibData::CorrectionPtr CorrectionPtr
Gaudi::Property< bool > m_doTof
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCacheCondObjInputKey
int multilayer(const Identifier &id) const
Access to components of the ID.
generic interface for a rt-relation
Definition IRtRelation.h:19
virtual double tLower() const =0
Returns the lower time covered by the r-t.
bool hasTmaxDiff() const
Definition IRtRelation.h:42
virtual double tUpper() const =0
Returns the upper time covered by the r-t.
const MdtFullCalibData * getCalibData(const Identifier &measId, MsgStream &msg) const
Returns the calibration data associated with this station.
class which holds calibration constants per rt-region
const IRtRelation * rt() const
rt relation
Support class for PropertyMgr.
Definition Property.h:23
int r
Definition globals.cxx:22
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
float RtScaleFunction(const float t, const bool ml2, const IRtRelation &rtrelation)
MdtDriftCircleStatus
Enum to represent the 'status' of Mdt measurements e.g.
@ MdtStatusAfterSpectrum
The tube produced a hit that is inconsistent with the drift time spectrum, the drift time is larger t...
@ MdtStatusBeforeSpectrum
The tube produced a hit that is inconsistent with the drift time spectrum, the drift time is smaller ...
@ MdtStatusUnDefined
Undefined.
@ MdtStatusDriftTime
The tube produced a vaild measurement.
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
void swap(ElementLinkVector< DOBJ > &lhs, ElementLinkVector< DOBJ > &rhs)
class which holds the full set of calibration constants for a given tube
float adcCal
quality flag for the SingleTubeCalib constants: 0 all ok, 1 no hits found, 2 too few hits,...
#define THROW_EXCEPTION(MESSAGE)
Definition throwExcept.h:10