10 #include <top/modules/TOPAligner/TOPAlignerModule.h>
11 #include <top/geometry/TOPGeometryPar.h>
12 #include <top/reconstruction_cpp/TOPTrack.h>
13 #include <framework/logging/Logger.h>
44 setDescription(
"Alignment of TOP");
48 addParam(
"targetModule", m_targetMid,
49 "Module to be aligned. Must be 1 <= Mid <= 16.", 1);
50 addParam(
"maxFails", m_maxFails,
51 "Maximum number of consecutive failed iterations before resetting the procedure", 100);
52 addParam(
"sample", m_sample,
53 "sample type: one of dimuon, bhabha or cosmics", std::string(
"dimuon"));
54 addParam(
"minMomentum", m_minMomentum,
55 "minimal track momentum if sample is cosmics", 1.0);
56 addParam(
"deltaEcms", m_deltaEcms,
57 "c.m.s energy window (half size) if sample is dimuon or bhabha", 0.1);
58 addParam(
"dr", m_dr,
"cut on POCA in r", 2.0);
59 addParam(
"dz", m_dz,
"cut on POCA in abs(z)", 4.0);
60 addParam(
"minZ", m_minZ,
61 "minimal local z of extrapolated hit", -130.0);
62 addParam(
"maxZ", m_maxZ,
63 "maximal local z of extrapolated hit", 130.0);
64 addParam(
"outFileName", m_outFileName,
65 "Root output file name containing alignment results",
66 std::string(
"TopAlignPars.root"));
67 addParam(
"stepPosition", m_stepPosition,
"step size for translations", 1.0);
68 addParam(
"stepAngle", m_stepAngle,
"step size for rotations", 0.01);
69 addParam(
"stepTime", m_stepTime,
"step size for t0", 0.05);
71 for (
const auto& parName : m_align.getParameterNames()) names += parName +
", ";
74 addParam(
"parInit", m_parInit,
75 "initial parameter values in the order [" + names +
"]. "
76 "If list is too short, the missing ones are set to 0.", m_parInit);
77 auto parFixed = m_parFixed;
78 addParam(
"parFixed", m_parFixed,
"list of names of parameters to be fixed. "
79 "Valid names are: " + names, parFixed);
83 void TOPAlignerModule::initialize()
87 const auto* geo = TOPGeometryPar::Instance()->getGeometry();
88 if (not geo->isModuleIDValid(m_targetMid)) {
89 B2ERROR(
"Target module ID = " << m_targetMid <<
" is invalid.");
94 if (not(m_sample ==
"dimuon" or m_sample ==
"bhabha" or m_sample ==
"cosmics")) {
95 B2ERROR(
"Invalid sample type '" << m_sample <<
"'");
97 if (m_sample ==
"bhabha") m_chargedStable = Const::electron;
102 m_selector.setMinMomentum(m_minMomentum);
103 m_selector.setDeltaEcms(m_deltaEcms);
104 m_selector.setCutOnPOCA(m_dr, m_dz);
105 m_selector.setCutOnLocalZ(m_minZ, m_maxZ);
109 m_align.setModuleID(m_targetMid);
110 m_align.setSteps(m_stepPosition, m_stepAngle, m_stepTime);
111 m_align.setParameters(m_parInit);
112 for (
const auto& parName : m_parFixed) {
113 m_align.fixParameter(parName);
118 m_digits.isRequired();
119 m_tracks.isRequired();
120 m_extHits.isRequired();
121 m_recBunch.isOptional();
125 m_file = TFile::Open(m_outFileName.c_str(),
"RECREATE");
126 if (m_file->IsZombie()) {
127 B2FATAL(
"Couldn't open file '" << m_outFileName <<
"' for writing!");
133 m_alignTree =
new TTree(
"alignTree",
"TOP alignment results");
134 m_alignTree->Branch(
"ModuleId", &m_targetMid);
135 m_alignTree->Branch(
"iter", &m_iter);
136 m_alignTree->Branch(
"ntrk", &m_ntrk);
137 m_alignTree->Branch(
"errorCode", &m_errorCode);
138 m_alignTree->Branch(
"iterPars", &m_vAlignPars);
139 m_alignTree->Branch(
"iterParsErr", &m_vAlignParsErr);
140 m_alignTree->Branch(
"valid", &m_valid);
141 m_alignTree->Branch(
"numPhot", &m_numPhot);
142 m_alignTree->Branch(
"x", &m_x);
143 m_alignTree->Branch(
"y", &m_y);
144 m_alignTree->Branch(
"z", &m_z);
145 m_alignTree->Branch(
"p", &m_p);
146 m_alignTree->Branch(
"theta", &m_theta);
147 m_alignTree->Branch(
"phi", &m_phi);
148 m_alignTree->Branch(
"r_poca", &m_pocaR);
149 m_alignTree->Branch(
"z_poca", &m_pocaZ);
150 m_alignTree->Branch(
"x_poca", &m_pocaX);
151 m_alignTree->Branch(
"y_poca", &m_pocaY);
152 m_alignTree->Branch(
"Ecms", &m_cmsE);
153 m_alignTree->Branch(
"charge", &m_charge);
154 m_alignTree->Branch(
"PDG", &m_PDG);
158 void TOPAlignerModule::event()
165 if (m_recBunch.isValid()) {
166 if (not m_recBunch->isReconstructed())
return;
171 for (
const auto& track : m_tracks) {
175 if (not trk.
isValid())
continue;
181 if (not m_selector.isSelected(trk))
continue;
184 int err = m_align.iterate(trk, m_chargedStable);
190 }
else if (m_countFails <= m_maxFails) {
193 B2INFO(
"Reached maximum allowed number of failed iterations. "
194 "Resetting TOPalign object");
200 m_vAlignPars = m_align.getParameters();
201 m_vAlignParsErr = m_align.getErrors();
202 m_ntrk = m_align.getNumUsedTracks();
204 m_valid = m_align.isValid();
205 m_numPhot = m_align.getNumOfPhotons();
208 const auto& localPosition = m_selector.getLocalPosition();
209 m_x = localPosition.X();
210 m_y = localPosition.Y();
211 m_z = localPosition.Z();
212 const auto& localMomentum = m_selector.getLocalMomentum();
213 m_p = localMomentum.Mag();
214 m_theta = localMomentum.Theta();
215 m_phi = localMomentum.Phi();
216 const auto& pocaPosition = m_selector.getPOCAPosition();
217 m_pocaR = pocaPosition.Perp();
218 m_pocaZ = pocaPosition.Z();
219 m_pocaX = pocaPosition.X();
220 m_pocaY = pocaPosition.Y();
221 m_cmsE = m_selector.getCMSEnergy();
229 TString resMsg =
"M= ";
230 resMsg += m_align.getModuleID();
234 resMsg += m_errorCode;
236 resMsg += m_align.isValid();
237 for (
auto par : m_vAlignPars) {
248 void TOPAlignerModule::terminate()
252 m_alignTree->Write();
254 TH1F valid(
"valid",
"status valid", 16, 0.5, 16.5);
255 valid.SetXTitle(
"slot ID");
256 valid.SetBinContent(m_targetMid, m_valid);
259 TH1F ntrk(
"ntrk",
"number of tracks", 16, 0.5, 16.5);
260 ntrk.SetXTitle(
"slot ID");
261 ntrk.SetBinContent(m_targetMid, m_ntrk);
264 std::string name, title;
265 name =
"results_slot" + to_string(m_targetMid);
266 title =
"alignment parameters, slot " + to_string(m_targetMid);
267 int npar = m_align.getParams().size();
268 TH1F h0(name.c_str(), title.c_str(), npar, 0, npar);
269 const auto& par = m_align.getParams();
270 const auto& err = m_align.getErrors();
271 for (
int i = 0; i < npar; i++) {
272 h0.SetBinContent(i + 1, par[i]);
273 h0.SetBinError(i + 1, err[i]);
277 name =
"errMatrix_slot" + to_string(m_targetMid);
278 title =
"error matrix, slot " + to_string(m_targetMid);
279 TH2F h1(name.c_str(), title.c_str(), npar, 0, npar, npar, 0, npar);
280 const auto& errMatrix = m_align.getErrorMatrix();
281 for (
int i = 0; i < npar; i++) {
282 for (
int k = 0; k < npar; k++) {
283 h1.SetBinContent(i + 1, k + 1, errMatrix[i][k]);
288 name =
"corMatrix_slot" + to_string(m_targetMid);
289 title =
"correlation matrix, slot " + to_string(m_targetMid);
290 TH2F h2(name.c_str(), title.c_str(), npar, 0, npar, npar, 0, npar);
291 std::vector<double> diag;
292 for (
int i = 0; i < npar; i++) {
293 double d = errMatrix[i][i];
294 if (d != 0) d = 1.0 / sqrt(d);
297 for (
int i = 0; i < npar; i++) {
298 for (
int k = 0; k < npar; k++) {
299 h2.SetBinContent(i + 1, k + 1, diag[i] * diag[k] * errMatrix[i][k]);
307 B2RESULT(
"TOPAligner: slot = " << m_targetMid <<
", status = successful, "
308 <<
"iterations = " << m_iter <<
", tracks used = " << m_ntrk);
310 B2RESULT(
"TOPAligner: slot = " << m_targetMid <<
", status = failed, "
311 <<
"error code = " << m_errorCode
312 <<
", iterations = " << m_iter <<
", tracks used = " << m_ntrk);
Reconstructed track at TOP.
int getPDGCode() const
Returns PDG code of associated MCParticle (returns 0 if none)
double getCharge() const
Returns charge.
bool isValid() const
Checks if track is successfully constructed.
int getModuleID() const
Returns slot ID.
Utility for the track selection - used in various calibration modules.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.