ATLAS Offline Software
Loading...
Searching...
No Matches
MuonCalib::MultilayerRtDifference Class Reference

#include <MultilayerRtDifference.h>

Collaboration diagram for MuonCalib::MultilayerRtDifference:

Public Member Functions

 MultilayerRtDifference (int min_hits, TDirectory *control_histogram_dir=NULL)
 MultilayerRtDifference (const MultilayerRtDifference &MLRTD)
virtual ~MultilayerRtDifference ()
MultilayerRtDifferenceoperator= (const MultilayerRtDifference &MLRTD)
void Fill (const MdtCalibHitBase &hit, const IRtRelation &rt_relation)
bool DoFit (IRtRelation *rt_relation, const IMdtCalibration::MuonSegVec &seg)
const TF1 * GetFunction () const

Private Attributes

TF1 * m_polfun
MultilayerRtDifference_Histogramsm_histograms
int m_min_number_of_hits

Detailed Description

Definition at line 23 of file MultilayerRtDifference.h.

Constructor & Destructor Documentation

◆ MultilayerRtDifference() [1/2]

MuonCalib::MultilayerRtDifference::MultilayerRtDifference ( int min_hits,
TDirectory * control_histogram_dir = NULL )

Definition at line 147 of file MultilayerRtDifference.cxx.

147 :
148 m_histograms(new MultilayerRtDifference_Histograms(control_histogram_dir)), m_min_number_of_hits(min_hits) {
149 m_polfun = new TF1("polfun", MultilayerRtDifference_ScaleFunction, 4.0, 15.0, 1);
150 }
MultilayerRtDifference_Histograms * m_histograms
Double_t MultilayerRtDifference_ScaleFunction(Double_t *x, Double_t *par)

◆ MultilayerRtDifference() [2/2]

MuonCalib::MultilayerRtDifference::MultilayerRtDifference ( const MultilayerRtDifference & MLRTD)

Definition at line 153 of file MultilayerRtDifference.cxx.

153 : m_min_number_of_hits(MLRTD.m_min_number_of_hits) {
154 m_polfun = new TF1(*MLRTD.m_polfun);
155 m_histograms = new MultilayerRtDifference_Histograms(*MLRTD.m_histograms);
156 }

◆ ~MultilayerRtDifference()

MuonCalib::MultilayerRtDifference::~MultilayerRtDifference ( )
virtual

Definition at line 170 of file MultilayerRtDifference.cxx.

170 {
171 delete m_histograms;
172 delete m_polfun;
173 }

Member Function Documentation

◆ DoFit()

bool MuonCalib::MultilayerRtDifference::DoFit ( IRtRelation * rt_relation,
const IMdtCalibration::MuonSegVec & seg )

Definition at line 183 of file MultilayerRtDifference.cxx.

183 {
184 TProfile *prov_diff = m_histograms->GetProfileDiff(m_min_number_of_hits);
185 if (!prov_diff) {
186 MsgStream log(Athena::getMessageSvc(), "MultilayerRtDifference");
187 log << MSG::WARNING << "MultilayerRtDifference::DoFit: Not enough hits!" << endmsg;
188 return false;
189 }
190 if (prov_diff->Fit("polfun", "Q", "", 4., 15.) != 0) {
191 MsgStream log(Athena::getMessageSvc(), "MultilayerRtDifference");
192 log << MSG::WARNING << "MultilayerRtDifference: Fit of polinomial failed! Not updating scale!" << endmsg;
193 return false;
194 }
195 if (std::abs(m_polfun->GetParameter(0)) < m_polfun->GetParError(0)) {
196 MsgStream log(Athena::getMessageSvc(), "MultilayerRtDifference");
197 log << MSG::WARNING << "MultilayerRtDifference: No Scale update needed! Scale correction: " << m_polfun->GetParameter(0)
198 << " +- " << m_polfun->GetParError(0) << endmsg;
199 return true;
200 }
201 if (!rt_relation) return true;
202 float scale = m_polfun->GetParameter(0);
203 if (std::abs(scale) > 2) scale *= 4;
204 if (rt_relation->hasTmaxDiff()) { scale += rt_relation->GetTmaxDiff(); }
205 rt_relation->SetTmaxDiff(scale);
206 if (seg.empty()) return true;
207 // update input segments
208 for (const auto &segment : seg) {
209 for (const MuonCalibSegment::MdtHitPtr &hit : segment->mdtHOT()) {
210 float old_corr = hit->TemperatureTime();
211 float corr = RtScaleFunction(hit->driftTime(), hit->identify().mdtMultilayer() == 2, *rt_relation);
212 hit->setTemperatureTime(corr);
213 hit->setDriftTime(hit->driftTime() - corr + old_corr);
214 hit->setDriftRadius(rt_relation->radius(hit->driftTime()), hit->sigmaDriftRadius());
215 }
216 }
217 return true;
218 }
#define endmsg
std::shared_ptr< MdtCalibHitBase > MdtHitPtr
typedef for a collection of MdtCalibHitBase s
IMessageSvc * getMessageSvc(bool quiet=false)
float RtScaleFunction(const float t, const bool ml2, const IRtRelation &rtrelation)

◆ Fill()

void MuonCalib::MultilayerRtDifference::Fill ( const MdtCalibHitBase & hit,
const IRtRelation & rt_relation )

Definition at line 175 of file MultilayerRtDifference.cxx.

175 {
176 int ml = hit.identify().mdtMultilayer() - 1;
177 double r_track = std::abs(hit.signedDistanceToTrack());
178 double res = std::abs(hit.driftRadius()) - r_track;
179 double v = rt_relation.driftVelocity(hit.driftTime());
180 m_histograms->GetResHist(ml)->Fill(r_track, res / v);
181 }
std::pair< std::vector< unsigned int >, bool > res

◆ GetFunction()

const TF1 * MuonCalib::MultilayerRtDifference::GetFunction ( ) const
inline

Definition at line 34 of file MultilayerRtDifference.h.

34{ return m_polfun; }

◆ operator=()

MultilayerRtDifference & MuonCalib::MultilayerRtDifference::operator= ( const MultilayerRtDifference & MLRTD)

Definition at line 159 of file MultilayerRtDifference.cxx.

159 {
160 if (this != &MLRTD) {
161 m_min_number_of_hits = MLRTD.m_min_number_of_hits;
162 delete m_polfun;
163 m_polfun = new TF1(*MLRTD.m_polfun);
164 delete m_histograms;
165 m_histograms = new MultilayerRtDifference_Histograms(*MLRTD.m_histograms);
166 }
167 return *this;
168 }

Member Data Documentation

◆ m_histograms

MultilayerRtDifference_Histograms* MuonCalib::MultilayerRtDifference::m_histograms
private

Definition at line 38 of file MultilayerRtDifference.h.

◆ m_min_number_of_hits

int MuonCalib::MultilayerRtDifference::m_min_number_of_hits
private

Definition at line 39 of file MultilayerRtDifference.h.

◆ m_polfun

TF1* MuonCalib::MultilayerRtDifference::m_polfun
private

Definition at line 37 of file MultilayerRtDifference.h.


The documentation for this class was generated from the following files: