ATLAS Offline Software
Loading...
Searching...
No Matches
MdtCalibrationTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 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 SG::ReadCondHandle constantHandle{m_calibDbKey, ctx};
97 if (!constantHandle.isValid()){
98 THROW_EXCEPTION("Failed to retrieve the Mdt calibration constants "<<m_calibDbKey.fullKey());
99 }
100 return constantHandle->getCalibData(channelId, msgStream());
101}
103 const MdtCalibInput& calibIn,
104 bool resolFromRtrack) const {
105
106 const Identifier& id{calibIn.identify()};
107
108 SG::ReadCondHandle constantHandle{m_calibDbKey, ctx};
109 if (!constantHandle.isValid()){
110 THROW_EXCEPTION("Failed to retrieve the Mdt calibration constants "<<m_calibDbKey.fullKey());
111 }
112
113 const MuonCalib::MdtFullCalibData* calibConstants = constantHandle->getCalibData(id, msgStream());
114 if (!calibConstants) {
115 ATH_MSG_WARNING("Could not find calibration data for channel "<<m_idHelperSvc->toString(id));
116 return MdtCalibOutput{};
117 }
118
119 // require at least the MdtRtRelation to be available
120 const RtRelationPtr& rtRelation{calibConstants->rtRelation};
121 // Hardcoded MDT tube radius 14.6mm here - not correct for sMDT
122 // on the other hand it should be rare that a tube does not have an RT
123 if(!rtRelation) {
124 ATH_MSG_WARNING("No rtRelation found, cannot calibrate tube "<<m_idHelperSvc->toString(id));
125 return MdtCalibOutput{};
126 }
127 if (!calibConstants->tubeCalib) {
128 ATH_MSG_WARNING("Cannot extract the single tube calibrations for tube "<<m_idHelperSvc->toString(id));
129 return MdtCalibOutput{};
130 }
132 const SingleTubeCalib* singleTubeData = calibConstants->tubeCalib->getCalib(id);
133 if (!singleTubeData) {
134 ATH_MSG_WARNING("Failed to access tubedata for " << m_idHelperSvc->toString(id));
135 return MdtCalibOutput{};
136 }
137 const float invPropSpeed = constantHandle->inversePropSpeed();
138
139 MdtCalibOutput calibResult{};
140 // correct for global t0 of rt-region
141
142 calibResult.setTubeT0(singleTubeData->t0);
143 calibResult.setMeanAdc(singleTubeData->adcCal);
144
145 // set propagation delay
146 if (m_doProp) {
147 const double propagationDistance = calibIn.signalPropagationDistance();
148 ATH_MSG_VERBOSE("Calibration of "<<m_idHelperSvc->toString(id)<<", propagation distance: "<<propagationDistance<<" -> "
149 <<(invPropSpeed * propagationDistance));
150 calibResult.setPropagationTime(invPropSpeed * propagationDistance);
151 }
152
154 const double driftTime = calibIn.tdc() * tdcBinSize
155 - (m_doTof ? calibIn.timeOfFlight() : 0.)
156 - calibIn.triggerTime()
157 - calibResult.tubeT0()
158 - calibResult.signalPropagationTime();
159
160 calibResult.setDriftTime(driftTime);
161 // apply corrections
162 double corrTime{0.};
163 const bool doCorrections = m_doField || m_doTemp || m_doBkg;
164 if (doCorrections) {
165 const CorrectionPtr& corrections{calibConstants->corrections};
166 const RtRelationPtr& rtRelation{calibConstants->rtRelation};
167 const MuonCalib::IRtRelation* rt = rtRelation->rt();
168 ATH_MSG_VERBOSE("There are correction functions.");
170 if (m_doSlew && corrections->slewing()) {
171 double slewTime=corrections->slewing()->correction(calibResult.driftTime(), calibIn.adc());
172 corrTime -= slewTime;
173 calibResult.setSlewingTime(slewTime);
174 }
175
176
177 if (m_doField && corrections->bField()) {
178 MagField::AtlasFieldCache fieldCache{};
179
181 if (!readHandle.isValid()) {
182 THROW_EXCEPTION("calibrate: Failed to retrieve AtlasFieldCacheCondObj with key " << m_fieldCacheCondObjInputKey.key());
183 }
184 readHandle->getInitializedCache(fieldCache);
185
186 Amg::Vector3D globalB{Amg::Vector3D::Zero()};
187 fieldCache.getField(calibIn.closestApproach().data(), globalB.data());
188 const Amg::Vector2D locBField = calibIn.projectMagneticField(globalB);
189 using BFieldComp = MdtCalibInput::BFieldComp;
190 calibResult.setLorentzTime(corrections->bField()->correction(calibResult.driftTime(),
191 locBField[static_cast<int>(BFieldComp::alongWire)],
192 locBField[static_cast<int>(BFieldComp::alongTrack)]));
193 corrTime -= calibResult.lorentzTime();
194 }
195 if(m_doTemp && rt && rt->hasTmaxDiff()) {
196 const MdtIdHelper& id_helper{m_idHelperSvc->mdtIdHelper()};
197 const int mL = id_helper.multilayer(id);
198 const double tempTime = MuonCalib::RtScaleFunction(calibResult.driftTime(), mL == 2, *rt);
199 calibResult.setTemperatureTime(tempTime);
200 corrTime-=calibResult.temperatureTime();
201 }
202 // background corrections (I guess this is never active)
203 if (m_doBkg && corrections->background()) {
204 double bgLevel{0.};
205 calibResult.setBackgroundTime(corrections->background()->correction(calibResult.driftTime(), bgLevel ));
206 corrTime += calibResult.backgroundTime();
207 }
208 }
209
210 calibResult.setDriftTime(calibResult.driftTime() + corrTime);
211
212 // calculate drift radius + error
213 double r{0.}, reso{0.};
214 double t = calibResult.driftTime();
215 double t_inrange = t;
216 Muon::MdtDriftCircleStatus timeStatus = driftTimeStatus(t, *rtRelation);
217
218 assert(rtRelation->rt() != nullptr);
219 r = rtRelation->rt()->radius(t);
220 // apply tUpper gshift
221 if (m_doTMaxShift) {
222 float tShift = m_tMaxShiftTool->getValue(id);
223 r = rtRelation->rt()->radius( t * (1 + tShift) );
224 }
225 // check whether drift times are within range, if not fix them to the min/max range
226 if ( t < rtRelation->rt()->tLower()) {
227 t_inrange = rtRelation->rt()->tLower();
228 double rmin = rtRelation->rt()->radius( t_inrange );
229 double drdt = rtRelation->rt()->driftVelocity( t_inrange);
231 if (timeStatus == Muon::MdtStatusBeforeSpectrum) {
232 t = rtRelation->rt()->tLower() - m_timeWindowLowerBound;
233 }
235 r = rmin + drdt*(t-t_inrange);
236 } else if( t > rtRelation->rt()->tUpper() ) {
237 t_inrange = rtRelation->rt()->tUpper();
238 double rmax = rtRelation->rt()->radius( t_inrange );
239 double drdt = rtRelation->rt()->driftVelocity(t_inrange);
240 // now check whether we are outside the time window
241 if ( timeStatus == Muon::MdtStatusAfterSpectrum ) {
242 t = rtRelation->rt()->tUpper() + m_timeWindowUpperBound;
243 }
245 r = rmax + drdt*(t-t_inrange);
246 }
247
248 assert(rtRelation->rtRes() != nullptr);
249
252 ATH_MSG_VERBOSE("Calibrated radius below physical limit "<<r<<" vs. "<<calibIn.innerTubeR());
255 reso = std::pow(calibIn.innerTubeR(), 2);
256 } else if (r > calibIn.innerTubeR()) {
257 ATH_MSG_VERBOSE("Calibrated radius outside tube ."<<r<<" vs. "<<calibIn.innerTubeR());
258 timeStatus = Muon::MdtStatusAfterSpectrum;
259 r = calibIn.innerTubeR();
260 reso = std::pow(calibIn.innerTubeR(), 2);
261 } else if (!resolFromRtrack) {
262 reso = rtRelation->rtRes()->resolution( t_inrange );
263 } else {
264 const std::optional<double> tFromR = rtRelation->tr()->driftTime(std::abs(calibIn.distanceToTrack()));
265 reso = rtRelation->rtRes()->resolution(tFromR.value_or(0.));
266 }
267
268 if (m_doPropUncert && !calibIn.trackDirHasPhi()) {
269 assert(rtRelation->rt() != nullptr);
270 const double driftTimeUp = std::min(rtRelation->rt()->tUpper(),
271 calibIn.tdc() * tdcBinSize
272 - (m_doTof ? calibIn.timeOfFlight() : 0.)
273 - calibIn.triggerTime()
274 - calibResult.tubeT0());
275
276 const double driftTimeDn = std::max(rtRelation->rt()->tLower(),
277 calibIn.tdc() * tdcBinSize
278 - (m_doTof ? calibIn.timeOfFlight() : 0.)
279 - calibIn.triggerTime()
280 - calibResult.tubeT0()
281 - calibIn.tubeLength() * invPropSpeed);
282
283 const double radiusUp = rtRelation->rt()->radius(driftTimeUp);
284 const double radiusDn = rtRelation->rt()->radius(driftTimeDn);
285 ATH_MSG_VERBOSE("Measurement "<<m_idHelperSvc->toString(calibIn.identify())
286 <<" nominal drift time "<<calibResult.driftTime()<<", down: "<<driftTimeDn<<", up: "<<driftTimeUp
287 <<" --> driftRadius: "<<r<<" pm "<<reso<<", prop-up: "<<radiusUp<<", prop-dn: "<<radiusDn
288 <<" delta: "<<(radiusUp-radiusDn));
289 calibResult.setDriftUncertSigProp(0.5*std::abs(radiusUp - radiusDn));
290 reso = std::hypot(reso, calibResult.driftUncertSigProp());
291 }
292
293 calibResult.setDriftRadius(r, reso);
294 calibResult.setStatus(timeStatus);
295 // summary
296 ATH_MSG_VERBOSE( "Calibration for tube " << m_idHelperSvc->toString(id)
297 <<" passed. "<<std::endl<<"Input: "<<calibIn<<std::endl<<"Extracted calib constants: "<<calibResult<<std::endl);
298 return calibResult;
299} //end MdtCalibrationTool::calibrate
300
302 MdtCalibInput&& primHit,
303 MdtCalibInput&& twinHit) const {
304
305 MdtCalibOutput primResult = calibrate(ctx, primHit);
306 MdtCalibOutput twinResult = calibrate(ctx, twinHit);
307
308 // get Identifier and MdtReadOutElement for twin tubes
309 // get 'raw' drifttimes of twin pair; we don't use timeofFlight or propagationTime cause they are irrelevant for twin coordinate
310 double primdriftTime = primHit.tdc()*tdcBinSize - primResult.tubeT0();
311 double twinDriftTime = twinHit.tdc()*tdcBinSize - twinResult.tubeT0();
313 if (primdriftTime >= twinDriftTime) {
314 ATH_MSG_VERBOSE("Swap "<<m_idHelperSvc->toString(primHit.identify())<<" & "<<m_idHelperSvc->toString(twinHit.identify())
315 <<" primDriftTime: "<<primdriftTime<<", secondTime: "<<twinDriftTime);
316 std::swap(primdriftTime, twinDriftTime);
317 std::swap(primResult, twinResult);
318 std::swap(primHit, twinHit);
319 }
320 const Identifier& primId = primHit.identify();
321 const Identifier& twinId = twinHit.identify();
322
323 // get calibration constants from DbTool
324 SG::ReadCondHandle constantHandle{m_calibDbKey, ctx};
325 if (!constantHandle.isValid()){
326 THROW_EXCEPTION("Failed to retrieve the Mdt calibration constants "<<m_calibDbKey.fullKey());
327 }
328
329
330 const MuonCalib::MdtFullCalibData* data1st = constantHandle->getCalibData(primId, msgStream());
331 const MuonCalib::MdtFullCalibData* data2nd = constantHandle->getCalibData(twinId, msgStream());
332 if (!data1st || !data2nd) {
333 ATH_MSG_WARNING(__FILE__<<":"<<__LINE__<<" Failed to access calibration constants for tubes "<<
334 m_idHelperSvc->toString(primId)<<" & "<<m_idHelperSvc->toString(twinId));
335 return MdtCalibTwinOutput{};
336 }
337 const SingleTubeCalib* calibSingleTube1st = data1st->tubeCalib->getCalib(primId);
338 const SingleTubeCalib* calibSingleTube2nd = data2nd->tubeCalib->getCalib(twinId);
339 const double invPropSpeed = constantHandle->inversePropSpeed();
340 if (!calibSingleTube1st || !calibSingleTube2nd) {
341 ATH_MSG_WARNING(__FILE__<<":"<<__LINE__<<" Failed to access calibration constants for tubes "<<
342 m_idHelperSvc->toString(primId)<<" & "<<m_idHelperSvc->toString(twinId));
343 return MdtCalibTwinOutput{};
344 }
345
346 // get tubelength and set HV-delay (~6ns)
347 constexpr double HVdelay = 6.;
348
350 double twin_timedif = twinDriftTime - primdriftTime - invPropSpeed * twinHit.tubeLength() - HVdelay;
353
354 const double tubeHalfLength = 0.5*primHit.tubeLength();
355 const double zTwin = std::clamp(0.5* primHit.readOutSide()* twin_timedif / invPropSpeed, -tubeHalfLength, tubeHalfLength);
356 const double errZTwin = m_resTwin / invPropSpeed;
357
358 ATH_MSG_VERBOSE( "Twin calibration - tube: " << m_idHelperSvc->toString(primId)<< " twintube: " << m_idHelperSvc->toString(twinId)<<endmsg
359 << " prompthit tdc = " << primHit.tdc() << " twinhit tdc = " << twinHit.tdc()
360 << " tube driftTime = " << primResult<< " second tube driftTime = " << twinResult<<endmsg
361 << " Time difference =" << twin_timedif << " zTwin=" << zTwin<<", errorZ="<<errZTwin);
362 MdtCalibTwinOutput calibResult{primHit,twinHit,primResult, twinResult};
363
364
365 calibResult.setLocZ(zTwin, errZTwin);
366 return calibResult;
367}
368
370 const MuonCalib::MdtRtRelation& rtRelation ) const {
371 if (rtRelation.rt()) {
372 if(driftTime < rtRelation.rt()->tLower() - m_timeWindowLowerBound) {
373 ATH_MSG_VERBOSE( " drift time outside time window "
374 << driftTime << ". Mininum time = "
375 << rtRelation.rt()->tLower() - m_timeWindowLowerBound );
377 } else if (driftTime > rtRelation.rt()->tUpper() + m_timeWindowUpperBound) {
378 ATH_MSG_VERBOSE( " drift time outside time window "
379 << driftTime << ". Maximum time = "
380 << rtRelation.rt()->tUpper() + m_timeWindowUpperBound);
382 }
383 } else {
384 ATH_MSG_WARNING( "No valid rt relation supplied for driftTimeStatus method" );
386 }
388}
389double MdtCalibrationTool::getResolutionFromRt(const EventContext& ctx, const Identifier& moduleID, const double time) const {
390
391 const MuonCalib::MdtFullCalibData* moduleConstants = getCalibConstants(ctx, moduleID);
392 if (!moduleConstants){
393 THROW_EXCEPTION("Failed to retrieve set of calibration constants for "<<m_idHelperSvc->toString(moduleID));
394 }
395 const RtRelationPtr& rtRel{moduleConstants->rtRelation};
396 if (!rtRel) {
397 THROW_EXCEPTION("No rt-relation found for "<<m_idHelperSvc->toString(moduleID));
398 }
399 const double t = std::min(std::max(time, rtRel->rt()->tLower()), rtRel->rt()->tUpper());
400 return rtRel->rtRes()->resolution(t);
401}
#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)
MdtCalibOutput::MdtDriftCircleStatus MdtDriftCircleStatus
MuonCalib::MdtTubeCalibContainer::SingleTubeCalib SingleTubeCalib
MdtCalibrationTool::ToolSettings ToolSettings
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.
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.
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