ATLAS Offline Software
Loading...
Searching...
No Matches
MdtAnalyticRtCalibAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
5
8
12
14#include "TGraph.h"
15#include "TGraphErrors.h"
16#include "TH1F.h"
17#include "TCanvas.h"
18#include "TLegend.h"
19#include "TROOT.h"
20
21#include <format>
22
23using namespace MuonCalib;
24using namespace Acts;
25namespace MuonCalibR4{
26
28 ATH_CHECK(m_idHelperSvc.retrieve());
29 ATH_CHECK(m_readKey.initialize());
30 ATH_CHECK(m_writeKey.initialize());
31 if (m_saveDiagnostic) {
32 ATH_CHECK(m_visualSvc.retrieve());
33 m_clientToken.canvasLimit = -1;
34 m_clientToken.preFixName = "AnalyticMdtCalib";
35 ATH_CHECK(m_visualSvc->registerClient(m_clientToken));
36 }
37 return StatusCode::SUCCESS;
38 }
39 StatusCode MdtAnalyticRtCalibAlg::execute(const EventContext& ctx) const {
40 SG::WriteCondHandle writeHandle{m_writeKey, ctx};
41 if(writeHandle.isValid()) {
42 ATH_MSG_DEBUG("CondHandle " << writeHandle.fullKey() << " is already valid.");
43 return StatusCode::SUCCESS;
44 }
46 SG::ReadCondHandle readHandle{m_readKey, ctx};
47 ATH_CHECK(readHandle.isValid());
48 writeHandle.addDependency(readHandle);
50 auto writeCdo = std::make_unique<MdtCalibDataContainer>(m_idHelperSvc.get(), readHandle->granularity());
52 const MdtIdHelper& idHelper{m_idHelperSvc->mdtIdHelper()};
53 std::map<RtRelationPtr, RtRelationPtr, CalibParamSorter> translatedRts{CalibParamSorter{std::pow(0.1, m_precCutOff)}};
54
55 RtRelationPtr dummyFillerRt{};
56 MdtFullCalibData::CorrectionPtr dummyFillerCorr{};
57 std::unordered_set<Identifier> missedElements{};
58 for (auto itr = idHelper.detectorElement_begin(); itr != idHelper.detectorElement_end(); ++itr){
59 const Identifier& detId{*itr};
60 if (!readHandle->hasDataForChannel(detId, msgStream())){
61 ATH_MSG_VERBOSE("There's no calibration data available for "<<m_idHelperSvc->toStringDetEl(detId));
62 missedElements.insert(detId);
63 continue;
64 }
65 const MdtFullCalibData* copyMe = readHandle->getCalibData(detId, msgStream());
66
67 if (copyMe->corrections && !writeCdo->storeData(detId, copyMe->corrections, msgStream())) {
68 return StatusCode::FAILURE;
69 }
71 const std::vector<Identifier> tubes = tubeIds(detId);
72 if (std::ranges::none_of(tubes, [copyMe](const Identifier& tubeId){
73 return !copyMe->tubeCalib->getCalib(tubeId);
74 })) {
76 if (!writeCdo->storeData(detId, copyMe->tubeCalib, msgStream())) {
77 return StatusCode::FAILURE;
78 }
79 } else {
80 ATH_MSG_VERBOSE("Tube calibration constants invalid for: "<<m_idHelperSvc->toStringDetEl(detId));
81 missedElements.insert(detId);
82 }
83
84 RtRelationPtr& translated = translatedRts[copyMe->rtRelation];
85 if (!translated) {
86 translated = translateRt(ctx, detId, *copyMe->rtRelation);
87 }
88 if (!dummyFillerRt) {
89 dummyFillerRt = translated;
90 dummyFillerCorr = copyMe->corrections;
91 }
92 if (!writeCdo->storeData(detId, translated, msgStream())) {
93 ATH_MSG_FATAL("Failed to store rt relation for "<<m_idHelperSvc->toStringDetEl(detId));
94 return StatusCode::FAILURE;
95 }
96 }
97 if (!m_fillMissingCh) {
98 missedElements.clear();
99 }
100 for (const Identifier& detId: missedElements) {
101 ATH_MSG_WARNING("Initial container did not contain any constants for "<<m_idHelperSvc->toString(detId)
102 <<". Insert instead the first translated rt-relation for those channels");
103 const bool fillRt = !writeCdo->hasDataForChannel(detId, msgStream()) ||
104 !writeCdo->getCalibData(detId, msgStream())->rtRelation;
105 if(fillRt) {
106 if (!writeCdo->storeData(detId, dummyFillerRt, msgStream())) {
107 ATH_MSG_FATAL("Failed to store rt relation for "<<m_idHelperSvc->toStringDetEl(detId));
108 return StatusCode::FAILURE;
109 }
110 if (dummyFillerCorr && !writeCdo->storeData(detId, dummyFillerCorr, msgStream())) {
111 return StatusCode::FAILURE;
112 }
113 }
114 if (writeCdo->getCalibData(detId, msgStream())->tubeCalib) {
115 continue;
116 }
117 MdtCalibDataContainer::TubeContainerPtr dummyT0cont = std::make_unique<MdtTubeCalibContainer>(m_idHelperSvc.get(), detId);
119 dummyT0Calib.t0 = m_missingT0;
120 for (const Identifier& tubeId : tubeIds(detId)) {
121 if (!dummyT0cont->setCalib(dummyT0Calib, tubeId, msgStream())) {
122 return StatusCode::FAILURE;
123 }
124 }
125 if (!writeCdo->storeData(detId, std::move(dummyT0cont), msgStream())) {
126 return StatusCode::FAILURE;
127 }
128 }
129 ATH_CHECK(writeHandle.record(std::move(writeCdo)));
130 return StatusCode::SUCCESS;
131 }
132 std::vector<Identifier> MdtAnalyticRtCalibAlg::tubeIds(const Identifier& chId) const {
133 std::vector<Identifier> tubes{};
134 const MdtIdHelper& idHelper{m_idHelperSvc->mdtIdHelper()};
135 for (int ml = 1; ml <= idHelper.multilayerMax(chId); ++ml) {
136 const Identifier detId = idHelper.multilayerID(chId, ml);
137 for (int layer = 1; layer <= idHelper.tubeLayerMax(detId); ++layer) {
138 for (int tube = 1; tube <= idHelper.tubeMax(detId); ++tube) {
139 tubes.emplace_back(idHelper.channelID(detId, ml, layer, tube));
140 }
141 }
142 }
143 return tubes;
144 }
145
147 MdtAnalyticRtCalibAlg::translateRt(const EventContext& ctx,
148 const Identifier& detId,
149 const MdtRtRelation& inRel) const {
150
151 const std::vector<SamplePoint> rtPoints = fetchDataPoints(*inRel.rt(), *inRel.rtRes());
152 IRtRelationPtr rt{};
153 ITrRelationPtr tr{};
154 unsigned int order{0};
155 while (order++ <= m_maxOrder && !rt) {
156 ATH_MSG_DEBUG("Attempt to fit rt relation for "<<m_idHelperSvc->toStringDetEl(detId)
157 <<" using a polynomial of order "<<order);
158 switch(m_polyTypeRt) {
159 case toUnderlying(PolyType::ChebyChev):
160 rt = RtFromPoints::getRtChebyshev(rtPoints, order);
161 break;
162 case toUnderlying(PolyType::Legendre):
163 rt = RtFromPoints::getRtLegendre(rtPoints, order);
164 break;
165 case toUnderlying(PolyType::Simple):
166 rt = RtFromPoints::getRtSimplePoly(rtPoints, order);
167 break;
168 default:
169 break;
170 }
171 const unsigned rtNdoF = rtPoints.size() - rt->nDoF();
172 const double rtChi2 = calculateChi2(rtPoints, *rt) / rtNdoF;
173 ATH_MSG_VERBOSE("Fit has a chi2 of "<<rtChi2<<". Cut off at "<<m_chiCutOff);
174 if (rtChi2 > m_chiCutOff) {
175 rt.reset();
176 }
177 }
178 if (!rt) {
179 ATH_MSG_ERROR("Failed to fit rt relation for "<<m_idHelperSvc->toStringDetEl(detId));
180 return nullptr;
181 }
182 const std::vector<SamplePoint> trPoints = swapCoordinates(rtPoints, *rt);
183 order = 0;
184 while (order++ <= m_maxOrder && !tr) {
185 ATH_MSG_DEBUG("Now continue with tr fit for "<<m_idHelperSvc->toStringDetEl(detId)
186 <<" using a polynomial of order "<<order);
187 switch(m_polyTypeTr) {
188 case toUnderlying(PolyType::Legendre):
189 tr = RtFromPoints::getTrLegendre(trPoints, order);
190 break;
191 case toUnderlying(PolyType::ChebyChev):
192 tr = RtFromPoints::getTrChebyshev(trPoints, order);
193 break;
194 case toUnderlying(PolyType::Simple):
195 tr = RtFromPoints::getTrSimplePoly(trPoints, order);
196 break;
197 default:
198 break;
199 }
200 const unsigned trNdoF = rtPoints.size() - tr->nDoF();
201 const double trChi2 = calculateChi2(trPoints, *tr) / trNdoF;
203 ATH_MSG_VERBOSE("T-R fit resulted in a chi2/nDoF for tr: "<<trChi2<<" ("<<trNdoF<<")");
204 if (trChi2 > m_chiCutOff) {
205 tr.reset();
206 }
207 }
208 if (!tr) {
209 ATH_MSG_FATAL("Failed to fit tr relation for "<<m_idHelperSvc->toStringDetEl(detId));
210 return nullptr;
211 }
212 auto rtReso = inRel.smartReso();
213 if (m_fitRtReso) {
214 std::unique_ptr<IRtResolution> fittedReso{};
215 std::vector<SamplePoint> rtResoPoints = fetchResolution(rtPoints, m_relUncReso);
216 order = 0;
217 while (order++ <= m_maxOrderReso && !fittedReso) {
218 ATH_MSG_VERBOSE("Finally fit the resolution function "<<m_idHelperSvc->toStringDetEl(detId)
219 <<" using a polynomial of order "<<order);
220 fittedReso = RtFromPoints::getResoChebyshev(rtPoints, rt, m_relUncReso, order);
221 const unsigned nDoF = rtResoPoints.size() - fittedReso->nDoF();
222 const double chi2 = calculateChi2(rtResoPoints, *fittedReso) / nDoF;
223 ATH_MSG_VERBOSE("Fit has a chi2 of "<<chi2<<". Cut off at "<<m_chiCutOff);
224 if (chi2 > m_chiCutOff) {
225 if (order == m_maxOrderReso) {
226 drawResoFunc(ctx, detId, rtResoPoints, *fittedReso);
227 }
228 fittedReso.reset();
229 }
230 }
231
232 if (fittedReso) {
233 drawResoFunc(ctx, detId, rtResoPoints, *fittedReso);
234 rtReso = std::move(fittedReso);
235 } else {
236 ATH_MSG_WARNING("Despite of having a "<<m_maxOrderReso
237 <<" no rt resolution function could be fitted for "<<m_idHelperSvc->toStringDetEl(detId));
238 }
239 }
240 auto finalRt = std::make_unique<MdtRtRelation>(std::move(rt), rtReso, std::move(tr));
241 drawRt(ctx, detId, rtPoints, *finalRt);
242 return finalRt;
243 }
244 void MdtAnalyticRtCalibAlg::drawRt(const EventContext& ctx,
245 const Identifier& detId,
246 const std::vector<SamplePoint>& rtPoints,
247 const MdtRtRelation& inRel) const {
248 if (!m_saveDiagnostic) {
249 return;
250 }
251 const std::string chName = std::format("{:}{:d}{:}{:d}M{:1d}",
252 m_idHelperSvc->stationNameString(detId),
253 std::abs(m_idHelperSvc->stationEta(detId)),
254 m_idHelperSvc->stationEta(detId)> 0? 'A' : 'C',
255 m_idHelperSvc->stationPhi(detId),
256 m_idHelperSvc->mdtIdHelper().multilayer(detId));
257 auto canvas = m_visualSvc->prepareCanvas(ctx, m_clientToken, std::format("RtRelation_{:}", chName) );
258 canvas->setAxisTitles("drift time [ns]", "drift radius [mm]");
259 auto dataGraph = std::make_unique<TGraphErrors>();
260 for (const SamplePoint& dataPoint : rtPoints) {
261 dataGraph->SetPoint(dataGraph->GetN(), dataPoint.x1(), dataPoint.x2());
262 dataGraph->SetPointError(dataGraph->GetN()-1, 0., dataPoint.error());
263 canvas->expandPad(dataPoint.x1(), dataPoint.x2());
264 }
265
267 auto rtGraph = std::make_unique<TGraph>();
268 auto trGraph = std::make_unique<TGraph>();
269
270 double driftTime{inRel.rt()->tLower()};
271 while (driftTime <= inRel.rt()->tUpper()) {
272 const double radius = inRel.rt()->radius(driftTime);
273 const double backRadius = std::clamp(radius,inRel.tr()->minRadius(), inRel.tr()->maxRadius());
274 const double backTime = inRel.tr()->driftTime(backRadius).value_or(-666.);
275 rtGraph->SetPoint(rtGraph->GetN(), driftTime, radius);
276 trGraph->SetPoint(trGraph->GetN(), backTime, radius);
277 canvas->expandPad(driftTime, radius);
278 canvas->expandPad(backTime, radius);
279 driftTime+=1.;
280 }
281
282 const unsigned rtNDoF= rtPoints.size() - inRel.rt()->nDoF();
283 const unsigned trNDoF= rtPoints.size() - inRel.tr()->nDoF();
284 const double rtChi2 = calculateChi2(rtPoints, *inRel.rt())/ rtNDoF;
285 const double trChi2 = calculateChi2(swapCoordinates(rtPoints, *inRel.rt()), *inRel.tr())/ trNDoF;
286
287 auto legend = std::make_unique<TLegend>(0.2,0.7,0.6, 0.9,chName.c_str());
288 legend->SetFillStyle(0);
289 legend->SetBorderSize(0);
290 legend->AddEntry(rtGraph.get(),std::format("{:}, order: {:d}, #chi^{{2}}: {:.3f}({:d})",
291 inRel.rt()->name(), inRel.rt()->nDoF(), rtChi2, rtNDoF).c_str(),"L");
292
293 legend->AddEntry(trGraph.get(),std::format("{:}, order: {:d}, #chi^{{2}}: {:.3f}({:d})",
294 inRel.tr()->name(), inRel.tr()->nDoF(), trChi2, trNDoF).c_str(),"L");
295
296 dataGraph->SetMarkerSize(0);
297 rtGraph->SetLineColor(kRed);
298 trGraph->SetLineColor(kBlue);
299
300 canvas->add(std::move(dataGraph), "P");
301 canvas->add(std::move(rtGraph), "C");
302 canvas->add(std::move(trGraph), "C");
303 canvas->add(std::move(legend));
304 }
305
306 void MdtAnalyticRtCalibAlg::drawResoFunc(const EventContext& ctx,
307 const Identifier& detId,
308 const std::vector<SamplePoint>& resoPoints,
309 const IRtResolution& inReso) const{
310
311 const std::string chName = std::format("{:}{:d}{:}{:d}M{:1d}",
312 m_idHelperSvc->stationNameString(detId),
313 std::abs(m_idHelperSvc->stationEta(detId)),
314 m_idHelperSvc->stationEta(detId)> 0? 'A' : 'C',
315 m_idHelperSvc->stationPhi(detId),
316 m_idHelperSvc->mdtIdHelper().multilayer(detId));
317
318 auto canvas = m_visualSvc->prepareCanvas(ctx, m_clientToken, std::format("RtReso_{:}", chName) );
319 canvas->setAxisTitles("drift time [ns]", "#sigma(r_{drift}) [mm]");
320 auto dataGraph = std::make_unique<TGraphErrors>();
321 for (const SamplePoint& dataPoint : resoPoints) {
322 dataGraph->SetPoint(dataGraph->GetN(), dataPoint.x1(), dataPoint.x2());
323 dataGraph->SetPointError(dataGraph->GetN()-1, 0., dataPoint.error());
324 canvas->expandPad(dataPoint.x1(), dataPoint.x2());
325 }
326 auto resoGraph = std::make_unique<TGraph>();
327 const auto [tLow, tHigh] = interval(resoPoints);
328 const auto [resoLow, resoHigh] = minMax(resoPoints);
329 canvas->expandPad(tLow, resoLow);
330 canvas->expandPad(tHigh, resoHigh);
331 double driftTime{tLow};
332 while (driftTime <= tHigh) {
333 const double evalReso = inReso.resolution(driftTime);
334 resoGraph->SetPoint(resoGraph->GetN(), driftTime, evalReso);
335 driftTime+=0.5;
336 }
337 canvas->add(std::move(dataGraph), "P");
338 canvas->add(std::move(resoGraph), "C");
339
340 const double chi2 = calculateChi2(resoPoints, inReso) / (resoPoints.size() - inReso.nDoF());
341 canvas->add(MuonValR4::drawLabel(std::format("{:}, {:}, order: {:1d} #chi^{{2}}: {:.3f}", chName, inReso.name(),
342 inReso.nDoF(), chi2), 0.2, 0.8));
343 }
344}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
void reset(std::unique_ptr< DerivType > newObj)
Overload the pointer.
Identifier multilayerID(const Identifier &channeldID) const
static int tubeLayerMax()
static int multilayerMax()
int tubeMax() const
Identifier channelID(int stationName, int stationEta, int stationPhi, int multilayer, int tubeLayer, int tube) const
Gaudi::Property< bool > m_fillMissingCh
At the end of the translation, it's checked whether all channels have been assigned.
Gaudi::Property< unsigned > m_precCutOff
Precision cut-off to treat two incoming rt-relations as equivalent.
virtual StatusCode execute(const EventContext &ctx) const override final
Gaudi::Property< float > m_missingT0
Default t0 constant to use, in case there's.
Gaudi::Property< unsigned > m_maxOrder
Maximum order of the polynomial in use.
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
void drawRt(const EventContext &ctx, const Identifier &detId, const std::vector< MuonCalib::SamplePoint > &rtPoints, const MuonCalib::MdtRtRelation &inRel) const
SG::WriteCondHandleKey< MuonCalib::MdtCalibDataContainer > m_writeKey
Gaudi::Property< bool > m_fitRtReso
Toggle whether the resolution shall be also converted into a polynomial.
Gaudi::Property< double > m_relUncReso
Assignment of the relative uncertainty on each resolution data point.
MuonCalib::MdtFullCalibData::RtRelationPtr RtRelationPtr
void drawResoFunc(const EventContext &ctx, const Identifier &detId, const std::vector< MuonCalib::SamplePoint > &resoPoints, const MuonCalib::IRtResolution &inReso) const
ServiceHandle< MuonValR4::IRootVisualizationService > m_visualSvc
Service handle of the visualization service.
Gaudi::Property< bool > m_saveDiagnostic
Save diagnostic histograms.
Gaudi::Property< unsigned > m_maxOrderReso
Maximal order to use for the resolution.
Gaudi::Property< int > m_polyTypeRt
Toggle the polynomial for the Rt-relation: ChebyChev or Legendre.
Gaudi::Property< int > m_polyTypeTr
Toggle the polynomial for the Rt-relation: ChebyChev or Legendre.
MuonValR4::IRootVisualizationService::ClientToken m_clientToken
Token to be presented to the visualization service.
Gaudi::Property< float > m_chiCutOff
Stop incrementing the order once the chi2CutOff is reached.
virtual StatusCode initialize() override final
RtRelationPtr translateRt(const EventContext &ctx, const Identifier &detId, const MuonCalib::MdtRtRelation &inRel) const
Translates the rt / tr & resolution relation from a look-up table into the requested polynomial type.
std::vector< Identifier > tubeIds(const Identifier &chId) const
SG::ReadCondHandleKey< MuonCalib::MdtCalibDataContainer > m_readKey
virtual std::string name() const =0
virtual unsigned nDoF() const =0
Returns the number of degrees of freedom of the relation function.
virtual double tLower() const =0
Returns the lower time covered by the r-t.
virtual double radius(double t) const =0
returns drift radius for a given time
virtual double tUpper() const =0
Returns the upper time covered by the r-t.
Generic interface to retrieve the resolution on the drift radius as a function of the drift time.
virtual unsigned nDoF() const =0
Returns the number of degrees of freedom of the relation function.
virtual double resolution(double t, double bgRate=0.0) const =0
returns resolution for a give time and background rate
virtual double maxRadius() const =0
Returns the maximum drift-radius.
virtual unsigned nDoF() const =0
Returns the number of degrees of freedom of the tr relation.
virtual std::optional< double > driftTime(const double r) const =0
Interface method for fetching the drift-time from the radius Returns a nullopt if the time is out of ...
virtual double minRadius() const =0
Returns the minimum drift-radius.
MdtFullCalibData::TubeContainerPtr TubeContainerPtr
class which holds calibration constants per rt-region
const ITrRelation * tr() const
t(r) relationship
const IRtResolutionPtr & smartReso() const
const IRtResolution * rtRes() const
resolution
const IRtRelation * rt() const
rt relation
const SingleTubeCalib * getCalib(const Identifier &tubeId) const
return calibration constants of a single tube
bool setCalib(SingleTubeCalib val, const Identifier &tubeId, MsgStream &msg)
set the calibration constants of a single tube
static std::unique_ptr< ITrRelation > getTrLegendre(const std::vector< SamplePoint > &dataPoints, const unsigned order)
Converts a list of t(r) data points into a t(r) relation expressed as a series of legendre polynomial...
static std::unique_ptr< IRtRelation > getRtChebyshev(const std::vector< SamplePoint > &dataPoints, const unsigned order)
Converts a list of r-t data points into a r(t) relation expressed as a series of chebychev polynomial...
static std::unique_ptr< ITrRelation > getTrChebyshev(const std::vector< SamplePoint > &dataPoints, const unsigned order)
Converts a list of r-t data points into a t(r) relation expressed as a series of chebychev polynomial...
static std::unique_ptr< IRtRelation > getRtSimplePoly(const std::vector< SamplePoint > &dataPoints, const unsigned order)
Converts a list of r(t) data points into a r(t) relation expressed as a series of elementary monomoni...
static std::unique_ptr< IRtResolution > getResoChebyshev(const std::vector< SamplePoint > &dataPoints, const unsigned order)
Converts a list of reso - t into a reso(t) relation expressed as a series of chebychev polynomials.
static std::unique_ptr< IRtRelation > getRtLegendre(const std::vector< SamplePoint > &dataPoints, const unsigned order)
Converts a list of r-t data points into a r(t) relation expressed as a series of legendre polynomials...
static std::unique_ptr< ITrRelation > getTrSimplePoly(const std::vector< SamplePoint > &dataPoints, const unsigned order)
Converts a list of t(r) data points into a t(r) relation expressed as a series of elementary monomoni...
This class provides a sample point for the BaseFunctionFitter.
Definition SamplePoint.h:15
const_id_iterator detectorElement_begin() const
Iterators over full set of ids.
const_id_iterator detectorElement_end() const
void addDependency(const EventIDRange &range)
StatusCode record(const EventIDRange &range, T *t)
record handle, with explicit range DEPRECATED
const DataObjID & fullKey() const
double chi2(TH1 *h0, TH1 *h1)
CscCalcPed - algorithm that finds the Cathode Strip Chamber pedestals from an RDO.
std::pair< double, double > interval(const std::vector< SamplePoint > &points)
Returns the interval covered by the sample points.
std::pair< double, double > minMax(const std::vector< SamplePoint > &points)
Returns the minimum & maximum values covered by the sample points.
GeoModel::TransientConstSharedPtr< IRtRelation > IRtRelationPtr
Definition IRtRelation.h:17
GeoModel::TransientConstSharedPtr< ITrRelation > ITrRelationPtr
Definition ITrRelation.h:16
std::vector< SamplePoint > swapCoordinates(const std::vector< SamplePoint > &points, const IRtRelation &rtRel)
Creates a new vector of samples points with x1 exchanged by x2 and vice-versa.
double calculateChi2(const std::vector< SamplePoint > &dataPoints, const IRtRelation &rtRel)
Returns the chi2 of the rt-relation w.r.t.
std::vector< SamplePoint > fetchResolution(const std::vector< SamplePoint > &points, const double uncert)
Creates a new vector of sample points where the x2 is assigned to the uncertainty and the uncertainty...
std::vector< SamplePoint > fetchDataPoints(const IRtRelation &rtRel, const double relUnc)
Constructs a list of sample points from the rt-relation.
std::unique_ptr< TLatex > drawLabel(const std::string &text, const double xPos, const double yPos, const unsigned int fontSize=18)
Create a TLatex label,.
const std::string & chName(ChIndex index)
convert ChIndex into a string
Helper struct to group Mdt calibration constants which are equivalent within the target precision.
class which holds the full set of calibration constants for a given tube
GeoModel::TransientConstSharedPtr< MdtCorFuncSet > CorrectionPtr