ATLAS Offline Software
Loading...
Searching...
No Matches
T0CalibrationMT.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include <algorithm>
8#include <iostream>
9#include <string>
10#include <vector>
11
13#include "GaudiKernel/MsgStream.h"
23
24#include "TF1.h"
25#include "TFile.h"
26#include "TH1.h"
27#include "TH1F.h"
28#include "TH2F.h"
29
30namespace MuonCalib {
31
33 std::shared_ptr<const T0MTSettings> settings,
34 const std::vector<int> &sort_by,
35 const std::vector<int> &adc_sort_by) :
37 m_settings(settings),
38 m_converged(false),
39 m_name(name),
41 m_sort_by(sort_by),
42 m_adc_sort_by(adc_sort_by) {
43
44 std::string HistoFileName = "T0MT_" + m_name + ".root";
45 m_file = std::make_unique<TFile>(HistoFileName.c_str(), "recreate");
46 m_regiondir = m_file->mkdir(m_name.c_str());
47
48 m_histos.resize(sort_by.size());
49 m_adc_histos.resize(adc_sort_by.size());
50
51 m_tube_ids.resize(sort_by.size());
52 m_adc_tube_ids.resize(adc_sort_by.size());
53 }
54
56 m_file->Write();
57 m_file->Close();
58 }
59
61 for (const MuonCalibSegment::MdtHitPtr& hit : seg.mdtHOT()) {
62 MuonFixedId id = hit->identify();
63 m_nhits_per_tube[id.getIdInt()]++;
64 // get the T0 originally subtracted for this hit
65 int nML = id.mdtMultilayer();
66 int nL = id.mdtTubeLayer();
67 int nT = id.mdtTube();
68 const MdtTubeFitContainer::SingleTubeCalib *stc(nullptr);
69 NtupleStationId sid(id);
70 sid.SetMultilayer(0);
71 std::map<NtupleStationId, MdtTubeFitContainer *>::const_iterator res_it(m_result.find(sid));
72 if (res_it != m_result.end()) {
73 const MdtIdHelper& idHelper{res_it->second->idHelperSvc()->mdtIdHelper()};
74 const Identifier tubeId = idHelper.channelID(id.stationNameString(),
75 id.eta(), id.phi(), nML, nL, nT);
76 stc = res_it->second->getCalib(tubeId);
77 }
78 double oldT0 = 0;
79 if (stc)
80 oldT0 = stc->t0;
81 else {
82 MsgStream log(Athena::getMessageSvc(), "T0CalibrationMT");
83 log << MSG::WARNING << "no Single Tube Calib info found for ML=" << nML << " L=" << nL << " T=" << nT << endmsg;
84 }
85
86 // get histos
87 for (unsigned int i = 0; i < m_sort_by.size(); i++) {
88 T0MTHistos *histos = getHistos(id, i);
89 histos->GetTSpec()->Fill(hit->driftTime() + oldT0 + hit->tubeT0());
90 }
91 for (unsigned int i = 0; i < m_adc_sort_by.size(); i++) {
92 ADCMTHistos *adc_histos = getADCHistos(id, i);
93 adc_histos->GetADCSpec()->Fill(hit->adcCount());
94 }
95 // relative t0 offsets
96 m_rel_tube_t0s[sid].AddHit(*hit);
97 }
98 return true;
99 }
100
102 MsgStream log(Athena::getMessageSvc(), "T0CalibrationMT");
103 std::map<int, MdtTubeFitContainer::SingleTubeFit> full;
104 std::map<int, MdtTubeFitContainer::SingleTubeCalib> st;
105 std::map<int, std::string> fit_by;
106 if (m_settings->FitTime()) {
107 for (unsigned int i = 0; i < m_sort_by.size(); i++) { analyse_tdc(i, full, st, fit_by); }
108 }
109 for (unsigned int i = 0; i < m_adc_sort_by.size(); i++) { analyse_adc(i, full, st); }
110
111 for (auto & it : full) {
112 if (it.first == 0) continue;
113 MuonFixedId fId(it.first);
114 NtupleStationId sid(fId);
115 sid.SetMultilayer(0);
117 MdtTubeFitContainer::SingleTubeCalib &stc(st[it.first]);
118 fi.group_by = fit_by[it.first];
119 int nML = fId.mdtMultilayer();
120 int nL = fId.mdtTubeLayer();
121 int nT = fId.mdtTube();
122 const MdtIdHelper& idHelper{m_result[sid]->idHelperSvc()->mdtIdHelper()};
123 const Identifier tubeId = idHelper.channelID(fId.stationNameString(),
124 fId.eta(), fId.phi(), nML, nL, nT);
125
126 bool setInfo = m_result[sid]->setCalib(stc, tubeId, log);
127 if (!setInfo) {
128 log << MSG::WARNING << "T0CalibrationMT::PROBLEM! could not set SingleTubeCalib info" << endmsg;
129 }
130 fi.n_hits = m_nhits_per_tube[fId.getIdInt()];
132 setInfo = m_result[sid]->setFit(std::move(fi), tubeId, log);
133 if (!setInfo) {
134 MsgStream log(Athena::getMessageSvc(), "T0CalibrationMT");
135 log << MSG::WARNING << "T0CalibrationMT::PROBLEM! could not set SingleTubeFit info" << endmsg;
136 }
137 }
138
140 m_file->Write();
141 return true;
142 }
143
144 bool T0CalibrationMT::analyse_tdc(const int &nr, std::map<int, MdtTubeFitContainer::SingleTubeFit> &full,
145 std::map<int, MdtTubeFitContainer::SingleTubeCalib> &st, std::map<int, std::string> &fit_by_map) {
146 if (m_settings->VerboseLevel() > 1) {
147 MsgStream log(Athena::getMessageSvc(), "T0CalibrationMT");
148 log << MSG::VERBOSE << "T0CalibrationMT::analyse iteration " << m_currentItnum << endmsg;
149 }
150 std::string fit_by("UNKNOWN");
151 switch (m_sort_by[nr]) {
152 case HistogramId::TUBE: fit_by = "TUBE"; break;
153 case HistogramId::LAYER: fit_by = "LAYER"; break;
154 case HistogramId::MULTILAYER: fit_by = "MULTILAYER"; break;
155 case HistogramId::CHAMBER: fit_by = "CHAMBER"; break;
156 case HistogramId::MEZZ_CARD: fit_by = "MEZZ_CARD"; break;
157 }
158
159 for (std::map<HistogramId, std::unique_ptr<T0MTHistos>>::iterator it = m_histos[nr].begin(); it != m_histos[nr].end(); ++it)
160 // loop over m_histos histograms
161 {
162 doTimeFit(it->second.get(), m_tube_ids[nr][it->first], full, st, fit_by_map, fit_by);
163 }
164 return true;
165 }
166
167 bool T0CalibrationMT::analyse_adc(const int &nr, std::map<int, MdtTubeFitContainer::SingleTubeFit> &full,
168 std::map<int, MdtTubeFitContainer::SingleTubeCalib> &st) {
169 for (std::map<HistogramId, std::unique_ptr<ADCMTHistos>>::iterator it = m_adc_histos[nr].begin(); it != m_adc_histos[nr].end();
170 ++it)
171 if (m_settings->FitADC()) { doAdcFit(it->second.get(), m_adc_tube_ids[nr][it->first], full, st); }
172 return true;
173 }
174
176 for (const auto & seg : segs) handleSegment(*seg);
177 analyse();
178 return getResults();
179 }
180
181 void T0CalibrationMT::doTimeFit(T0MTHistos *T0h, const std::set<MuonFixedId> &tube_ids,
182 std::map<int, MdtTubeFitContainer::SingleTubeFit> &fim,
183 std::map<int, MdtTubeFitContainer::SingleTubeCalib> &stcm, std::map<int, std::string> &fit_by_map,
184 const std::string &fit_by) {
185 TDirectory *cwd = gDirectory;
186 m_regiondir->cd();
187 // do t0 fit
189 if (fit_by == "CHAMBER")
191 else if (fit_by == "MULTILAYER")
193 else if (fit_by == "LAYER")
194 theGroup = MdtRelativeTubeT0::LAYER;
195 else if (fit_by == "MEZZ_CARD")
197 if (T0h->FitT0() && m_settings->MinEntriesTime() <= T0h->GetTSpec()->GetEntries()) {
198 const TF1 *fun(T0h->GetT0Function());
199 for (auto tube_id : tube_ids) {
200 if (tube_id.getIdInt() == 0) continue;
201 NtupleStationId stId(tube_id);
202 stId.SetMultilayer(0);
203 double rel_t0(0.0);
204 if (m_settings->T0Settings()->CorrectRelT0s()) rel_t0 = m_rel_tube_t0s[stId].GetRelativeOffset(tube_id, theGroup);
205 MdtTubeFitContainer::SingleTubeFit &fi(fim[tube_id.getIdInt()]);
206 MdtTubeFitContainer::SingleTubeCalib &stc(stcm[tube_id.getIdInt()]);
207 // store parameters
208 fi.statistics = static_cast<int>(T0h->GetTSpec()->GetEntries());
209 fi.chi2Tdc = T0h->T0Chi2();
210 fi.par[0] = fun->GetParameter(T0MTHistos ::T0_PAR_NR_BACK);
211 fi.cov[0] = fun->GetParError(T0MTHistos ::T0_PAR_NR_BACK);
212 fi.par[4] = fun->GetParameter(T0MTHistos ::T0_PAR_NR_T0) + rel_t0;
213 fi.cov[4] = fun->GetParError(T0MTHistos ::T0_PAR_NR_T0);
214 fi.par[6] = fun->GetParameter(T0MTHistos ::T0_PAR_NR_T);
215 fi.cov[6] = fun->GetParError(T0MTHistos ::T0_PAR_NR_T);
216 stc.t0 = fun->GetParameter(T0MTHistos ::T0_PAR_NR_T0) + rel_t0;
217 stc.statusCode = T0h->StatusCode();
218 m_converged = true;
219 fit_by_map[tube_id.getIdInt()] = fit_by;
220 }
221 }
222 if (T0h->FitTmax() && m_settings->MinEntriesTime() <= T0h->GetTSpec()->GetEntries()) {
223 const TF1 *fun(T0h->GetTMaxFunction());
224 // store parameters
225 for (auto tube_id : tube_ids) {
226 if (tube_id.getIdInt() == 0) continue;
227 MdtTubeFitContainer::SingleTubeFit &fi(fim[tube_id.getIdInt()]);
228 // MdtTubeFitContainer::SingleTubeCalib &stc(stcm[*it]);
229 fi.par[5] = fun->GetParameter(T0MTHistos ::TMAX_PAR_NR_TMAX);
230 fi.cov[5] = fun->GetParError(T0MTHistos ::TMAX_PAR_NR_TMAX);
231 fi.chi2TdcEnd = fun->GetChisquare() / fun->GetNDF();
232 }
233 }
234 cwd->cd();
235 }
236
237 void T0CalibrationMT::doAdcFit(ADCMTHistos *T0h, const std::set<MuonFixedId> &tube_ids,
238 std::map<int, MdtTubeFitContainer::SingleTubeFit> &fim,
239 std::map<int, MdtTubeFitContainer::SingleTubeCalib> &stcm) {
240 if (T0h->FitAdc() && m_settings->MinEntriesADC() <= T0h->GetADCSpec()->GetEntries()) {
241 const TF1 *fun(T0h->GetAdcFunction());
242 if (fun == nullptr) return;
243 for (auto tube_id : tube_ids) {
244 if (tube_id.getIdInt() == 0) continue;
245 MdtTubeFitContainer::SingleTubeFit &fi(fim[tube_id.getIdInt()]);
246 MdtTubeFitContainer::SingleTubeCalib &stc(stcm[tube_id.getIdInt()]);
247
248 stc.adcCal = fun->GetParameter(1);
249 for (int i = 0; (i < fun->GetNpar() && i < 4); i++) {
250 fi.adc_par[i] = fun->GetParameter(i);
251 fi.adc_err[i] = fun->GetParError(i);
252 }
253 fi.adc_chi2 = fun->GetChisquare() / fun->GetNDF();
254 }
255 }
256 }
257
258 IMdtCalibration::MdtCalibOutputPtr T0CalibrationMT::getResults() const { return std::make_shared<T0CalibrationOutput>(m_result); }
259
260 ADCMTHistos *T0CalibrationMT::getADCHistos(const MuonFixedId &idtube, unsigned int nr) {
261 HistogramId id;
262 id.Initialize(idtube, m_adc_sort_by[nr]);
263 if (!m_adc_histos[nr][id]) {
264 TDirectory *cwd = gDirectory;
265 m_regiondir->cd();
266 m_adc_histos[nr][id] = std::make_unique<ADCMTHistos>(id.getIdInt(), m_settings.get(), id.HistogramName().c_str());
267 cwd->cd();
268 }
269 m_adc_tube_ids[nr][id].insert(idtube);
270 return m_adc_histos[nr][id].get();
271 }
272
273 T0MTHistos *T0CalibrationMT::getHistos(const MuonFixedId &idtube, unsigned int nr) {
274 HistogramId id;
275 id.Initialize(idtube, m_sort_by[nr]);
276 if (!m_histos[nr][id]) {
277 TDirectory *cwd = gDirectory;
278 m_regiondir->cd();
279 m_histos[nr][id] = std::make_unique<T0MTHistos>(id.getIdInt(), m_settings.get(), id.HistogramName().c_str());
280 cwd->cd();
281 }
282 m_tube_ids[nr][id].insert(idtube);
283 return m_histos[nr][id].get();
284 }
285
287 // This method is called both by the event loop and by the tool.
288 // Only the call from the tool is relevant for this implementation
289 // and should be performed only once.
290
291 if (m_result.size()) return;
292
293 const T0CalibrationOutput *t0Input = dynamic_cast<const T0CalibrationOutput *>(calib_in);
294 if (!calib_in || !t0Input) return;
295 m_result = t0Input->GetMap();
296 for (auto & it : m_result)
297 it.second->setImplementation("T0CalibrationMT");
298 }
299
301
302} // namespace MuonCalib
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define endmsg
TGraphErrors * GetEntries(TH2F *histo)
Identifier channelID(int stationName, int stationEta, int stationPhi, int multilayer, int tubeLayer, int tube) const
Histogram and fitter class for drift time and pulsehight spectra The rising slope is fitted by a ferm...
Definition ADCMTHistos.h:44
TF1 * GetAdcFunction() const
returnd function fitted to adc-spectrum
Definition ADCMTHistos.h:84
TH1F * GetADCSpec()
get adc spectrum
Definition ADCMTHistos.h:67
Identifier class for drift time histograms.
Definition HistogramId.h:24
static const int LAYER
Definition HistogramId.h:45
static const int TUBE
valid values of the sort_by argument of the Initialize function
Definition HistogramId.h:45
static const int MEZZ_CARD
Definition HistogramId.h:45
static const int MULTILAYER
Definition HistogramId.h:45
static const int CHAMBER
Definition HistogramId.h:45
Interface to pass calibration output during calibration.
IMdtCalibration(const std::string &name)
constructor, string used to identify the instance
virtual std::string name() const
returns name (region) of instance
std::shared_ptr< IMdtCalibrationOutput > MdtCalibOutputPtr
std::vector< std::shared_ptr< MuonCalibSegment > > MuonSegVec
A MuonCalibSegment is a reconstructed three dimensional track segment in the MuonSpectrometer.
std::shared_ptr< MdtCalibHitBase > MdtHitPtr
typedef for a collection of MdtCalibHitBase s
const MdtHitVec & mdtHOT() const
retrieve the full set of MdtCalibHitBase s assigned to this segment
Implements fixed identifiers not dependent upon Athena Identifier for internal use in the calibration...
Definition MuonFixedId.h:50
int mdtTubeLayer() const
Mdt specific:
int mdtTube() const
Mdt specific:
unsigned int getIdInt() const
std::string stationNameString() const
int mdtMultilayer() const
Mdt specific:
Station Identifier for sorting calibration data.
void SetMultilayer(const int &ml)
std::unique_ptr< TFile > m_file
pointer to the histogram file
void doTimeFit(T0MTHistos *T0h, const std::set< MuonFixedId > &tube_ids, std::map< int, MdtTubeFitContainer::SingleTubeFit > &fim, std::map< int, MdtTubeFitContainer::SingleTubeCalib > &stcm, std::map< int, std::string > &fit_by_map, const std::string &fit_by)
fit time spectrum
std::vector< std::map< HistogramId, std::unique_ptr< ADCMTHistos > > > m_adc_histos
vector of pointers tube histograms
virtual IMdtCalibration::MdtCalibOutputPtr getResults() const override
T0MTHistos * getHistos(const MuonFixedId &idtube, unsigned int nr)
retrieve pointer for tube idtube histograms
bool analyse_adc(const int &nr, std::map< int, MdtTubeFitContainer::SingleTubeFit > &full, std::map< int, MdtTubeFitContainer::SingleTubeCalib > &st)
std::shared_ptr< const T0MTSettings > m_settings
pointer to the settings
void setInput(const IMdtCalibrationOutput *input) override
unused
void doAdcFit(ADCMTHistos *T0h, const std::set< MuonFixedId > &tube_ids, std::map< int, MdtTubeFitContainer::SingleTubeFit > &fim, std::map< int, MdtTubeFitContainer::SingleTubeCalib > &stcm)
fit adc spectrum
std::vector< std::map< HistogramId, std::set< MuonFixedId > > > m_adc_tube_ids
int m_currentItnum
current iteration (always 1?)
bool m_converged
convergence status
bool handleSegment(MuonCalibSegment &seg)
fill tube spectra
std::vector< std::map< HistogramId, std::set< MuonFixedId > > > m_tube_ids
bool analyse()
extract parameters from spectra
std::map< int, int > m_nhits_per_tube
number of hits per tube
std::map< NtupleStationId, MdtTubeFitContainer * > m_result
tube constants
T0CalibrationMT(const std::string &name, std::shared_ptr< const T0MTSettings > settings, const std::vector< int > &sort_by, const std::vector< int > &adc_sort_by)
constructor
TDirectory * m_regiondir
pointer to the ROOT directory
const std::vector< int > & m_sort_by
std::vector< std::map< HistogramId, std::unique_ptr< T0MTHistos > > > m_histos
vector of pointers tube histograms
virtual IMdtCalibration::MdtCalibOutputPtr analyseSegments(const MuonSegVec &segs) override
new interface function
bool converged() const
return m_converged
bool analyse_tdc(const int &nr, std::map< int, MdtTubeFitContainer::SingleTubeFit > &full, std::map< int, MdtTubeFitContainer::SingleTubeCalib > &st, std::map< int, std::string > &fit_by_map)
std::map< NtupleStationId, MdtRelativeTubeT0 > m_rel_tube_t0s
ADCMTHistos * getADCHistos(const MuonFixedId &idtube, unsigned int nr)
retrieve pointer for tube idtube histograms
const std::vector< int > & m_adc_sort_by
std::string m_name
calibration region name
class for the communication of the results of T0 calibration algorithms
std::map< NtupleStationId, MdtTubeFitContainer * > & GetMap()
Histogram and fitter class for drift time and pulsehight spectra The rising slope is fitted by a ferm...
Definition T0MTHistos.h:40
TH1F * GetTSpec() const
get drift time spectrum
Definition T0MTHistos.h:72
int StatusCode() const
returns status code - the status code applies only to the t0 fit
Definition T0MTHistos.h:102
double T0Chi2() const
returns t0 chi2
Definition T0MTHistos.h:114
bool FitT0()
Perform t0-fit Returns true if fit is successfull.
bool FitTmax()
Performs tmax-fit Returns true if fit is successfull.
const TF1 * GetT0Function() const
returns function fitted to the riding edge of the spectrum
Definition T0MTHistos.h:104
const TF1 * GetTMaxFunction() const
returns function fitted to the riding edge of the spectrum
Definition T0MTHistos.h:110
singleton-like access to IMessageSvc via open function and helper
std::string cwd
Definition listroot.cxx:38
IMessageSvc * getMessageSvc(bool quiet=false)
CscCalcPed - algorithm that finds the Cathode Strip Chamber pedestals from an RDO.
float adcCal
quality flag for the SingleTubeCalib constants: 0 all ok, 1 no hits found, 2 too few hits,...
int n_hits_above_adc_cut
number of hits above adc cut
int statistics
< number of hits used for the fit
float chi2TdcEnd
for MTT0 chi2 of trailing edge fit