Project
Loading...
Searching...
No Matches
AlignmentHierarchy.cxx
Go to the documentation of this file.
1// Copyright 2019-2026 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
12#include <format>
13#include <fstream>
14#include <sstream>
15#include <fnmatch.h>
16#include <cmath>
17#include <TGeoManager.h>
18#include <TGeoPhysicalNode.h>
19#include <nlohmann/json.hpp>
20
23#include "Framework/Logger.h"
24#include "MathUtils/Utils.h"
25
26namespace o2::its3::align
27{
28
29void HierarchyConstraint::write(std::ostream& os) const
30{
31 os << "!!! " << mName << '\n';
32 os << "Constraint " << mValue << '\n';
33 for (size_t i{0}; i < mLabels.size(); ++i) {
34 os << mLabels[i] << " " << mCoeff[i] << '\n';
35 }
36 os << '\n';
37}
38
39AlignableVolume::AlignableVolume(const char* symName, uint32_t label, uint32_t det, bool sens) : mSymName(symName), mLabel(det, label, sens)
40{
41 init();
42}
43
44AlignableVolume::AlignableVolume(const char* symName, GlobalLabel label) : mSymName(symName), mLabel(label)
45{
46 init();
47}
48
49void AlignableVolume::init()
50{
51 // check if this sym volume actually exists
52 mPNE = gGeoManager->GetAlignableEntry(mSymName.c_str());
53 if (mPNE == nullptr) {
54 LOGP(fatal, "Symbolic volume '{}' has no corresponding alignable entry!", mSymName);
55 }
56 mPN = mPNE->GetPhysicalNode();
57 if (mPN == nullptr) {
58 LOGP(debug, "Adding physical node to {}", mSymName);
59 mPN = gGeoManager->MakePhysicalNode(mPNE->GetTitle());
60 if (mPN == nullptr) {
61 LOGP(fatal, "Failed to make physical node for {}", mSymName);
62 }
63 }
64}
65
67{
68 if (level == 0 && !isRoot()) {
69 LOGP(fatal, "Finalise should be called only from the root node!");
70 }
71 mLevel = level;
72 if (!isLeaf()) {
73 // depth first
74 for (const auto& c : mChildren) {
75 c->finalise(level + 1);
76 }
77 // auto-disable parent RB DOFs if no children are active
78 if (mRigidBody) {
79 int nActiveChildren = 0;
80 for (const auto& c : mChildren) {
81 if (c->isActive()) {
82 ++nActiveChildren;
83 }
84 }
85 if (!nActiveChildren) {
86 for (int iDOF = 0; iDOF < mRigidBody->nDOFs(); ++iDOF) {
87 if (mRigidBody->isFree(iDOF)) {
88 LOGP(warn, "Auto-disabling DOF {} for {} since no active children",
89 mRigidBody->dofName(iDOF), mSymName);
90 mRigidBody->setFree(iDOF, false);
91 }
92 }
93 }
94 }
95 } else {
96 // for sensors we need also to define the transformation from the measurment (TRK) to the local frame (LOC)
97 // need to it with including possible pre-alignment to allow for iterative convergence
98 // (TRK) is defined wrt global z-axis
101 }
102 if (!isRoot()) {
103 // prepare the transformation matrices, e.g. from child frame to parent frame
104 // this is not necessarily just one level transformation
105 TGeoHMatrix mat = *mPN->GetMatrix(); // global matrix (including possible pre-alignment) from this volume to the global frame
106 if (isLeaf()) {
107 mat = mL2G; // for sensor volumes they might have redefined the L2G definition
108 }
109 auto inv = mParent->mPN->GetMatrix()->Inverse(); // global (including possible pre-alignment) from this volume to the global frame
110 mat.MultiplyLeft(inv); // left mult. effectively subtracts the parent transformation which is included in the the childs
111 mL2P = mat; // now this is directly the child to the parent transformation (LOC) (including possible pre-alignment)
112
113 // prepare jacobian from child to parent frame
114 Eigen::Map<const Eigen::Matrix<double, 3, 3, Eigen::RowMajor>> rotL2P(mL2P.GetRotationMatrix());
115 Eigen::Matrix3d rotInv = rotL2P.transpose(); // parent-to-child rotation
116 const double* t = mL2P.GetTranslation(); // child origin in parent frame
117 Eigen::Matrix3d skewT;
118 skewT << 0, -t[2], t[1], t[2], 0, -t[0], -t[1], t[0], 0;
119 mJL2P.setZero();
120 mJL2P.topLeftCorner<3, 3>() = rotInv;
121 mJL2P.topRightCorner<3, 3>() = -rotInv * skewT;
122 mJL2P.bottomRightCorner<3, 3>() = rotInv;
123 mJP2L = mJL2P.inverse();
124 }
125}
126
128{
129 if (isLeaf() || !mRigidBody) {
130 // recurse even if this node has no RB DOFs
131 for (const auto& c : mChildren) {
132 c->writeRigidBodyConstraints(os);
133 }
134 return;
135 }
136
137 for (int iDOF = 0; iDOF < mRigidBody->nDOFs(); ++iDOF) {
138 if (!mRigidBody->isFree(iDOF)) {
139 continue;
140 }
141 double nActiveChildren = 0.;
142 for (const auto& c : mChildren) {
143 if (c->isActive()) {
144 ++nActiveChildren;
145 }
146 }
147 if (nActiveChildren == 0.) {
148 LOGP(fatal, "{} has dof {} active but no active children!", mSymName, mRigidBody->dofName(iDOF));
149 }
150 const double invN = 1.0 / nActiveChildren;
151 HierarchyConstraint con(std::format("DOF {} for {}", mRigidBody->dofName(iDOF), mSymName), 0.0);
152 for (const auto& c : mChildren) {
153 if (!c->mRigidBody) {
154 continue;
155 }
156 for (int jDOF = 0; jDOF < c->mRigidBody->nDOFs(); ++jDOF) {
157 if (!c->mRigidBody->isFree(jDOF)) {
158 continue;
159 }
160 double coeff = invN * c->getJP2L()(iDOF, jDOF);
161 if (std::abs(coeff) > 1e-16f) {
162 con.add(c->getLabel().raw(jDOF), coeff);
163 }
164 }
165 }
166
167 if (con.getSize() > 1) {
168 con.write(os);
169 }
170 }
171 for (const auto& c : mChildren) {
172 c->writeRigidBodyConstraints(os);
173 }
174}
175
176void AlignableVolume::writeParameters(std::ostream& os) const
177{
178 if (isRoot()) {
179 os << "Parameter\n";
180 }
181 if (mRigidBody) {
182 for (int iDOF = 0; iDOF < mRigidBody->nDOFs(); ++iDOF) {
183 os << std::format("{:<10} {:>+15g} {:>+15g} ! {} {} ",
184 mLabel.raw(iDOF), 0.0, (mRigidBody->isFree(iDOF) ? 0.0 : -1.0),
185 (mRigidBody->isFree(iDOF) ? 'V' : 'F'), mRigidBody->dofName(iDOF))
186 << mSymName << '\n';
187 }
188 }
189 if (mCalib) {
190 auto calibLbl = mLabel.asCalib();
191 for (int iDOF = 0; iDOF < mCalib->nDOFs(); ++iDOF) {
192 os << std::format("{:<10} {:>+15g} {:>+15g} ! {} {} ",
193 calibLbl.raw(iDOF), 0.0, (mCalib->isFree(iDOF) ? 0.0 : -1.0),
194 (mCalib->isFree(iDOF) ? 'V' : 'F'), mCalib->dofName(iDOF))
195 << mSymName << '\n';
196 }
197 }
198 for (const auto& c : mChildren) {
199 c->writeParameters(os);
200 }
201}
202
203void AlignableVolume::writeTree(std::ostream& os, int indent) const
204{
205 os << std::string(static_cast<size_t>(indent * 2), ' ') << mSymName << (mLabel.sens() ? " (sens)" : " (pasv)");
206 if (mIsPseudo) {
207 os << " pseudo";
208 } else {
209 int nFreeDofs{0};
210 if (mRigidBody && mRigidBody->nFreeDOFs()) {
211 nFreeDofs += mRigidBody->nFreeDOFs();
212 os << " RB[";
213 for (int i = 0; i < mRigidBody->nDOFs(); ++i) {
214 if (mRigidBody->isFree(i)) {
215 os << " " << mRigidBody->dofName(i) << "(" << mLabel.raw(i) << ")";
216 }
217 }
218 os << " ]";
219 }
220 if (mCalib && mCalib->nFreeDOFs()) {
221 nFreeDofs += mCalib->nFreeDOFs();
222 os << " CAL[";
223 auto calibLbl = mLabel.asCalib();
224 for (int i = 0; i < mCalib->nDOFs(); ++i) {
225 if (mCalib->isFree(i)) {
226 os << " " << mCalib->dofName(i) << "(" << calibLbl.raw(i) << ")";
227 }
228 }
229 os << " ]";
230 }
231 if (!nFreeDofs) {
232 os << " no DOFs";
233 }
234 }
235 os << '\n';
236 for (const auto& c : mChildren) {
237 c->writeTree(os, indent + 2);
238 }
239}
240
241void applyDOFConfig(AlignableVolume* root, const std::string& jsonPath)
242{
243 using json = nlohmann::json;
244 std::ifstream f(jsonPath);
245 if (!f.is_open()) {
246 LOGP(fatal, "Cannot open DOF config file: {}", jsonPath);
247 }
248 auto data = json::parse(f);
249 json rules = data.is_array() ? data : data.value("rules", json::array());
250
251 static const std::map<std::string, int> rbNameToIdx = {
252 {"TX", 0}, {"TY", 1}, {"TZ", 2}, {"RX", 3}, {"RY", 4}, {"RZ", 5}};
253
254 auto matchPattern = [](const std::string& pattern, const std::string& sym) -> bool {
255 if (fnmatch(pattern.c_str(), sym.c_str(), 0) == 0) {
256 return true;
257 }
258 std::string prefixed = "*" + pattern;
259 return fnmatch(prefixed.c_str(), sym.c_str(), 0) == 0;
260 };
261
262 if (data.is_object() && data.contains("defaults")) {
263 json defRule = data["defaults"];
264 defRule["match"] = "*";
265 rules.insert(rules.begin(), defRule);
266 }
267
268 root->traverse([&](AlignableVolume* vol) {
269 const std::string& sym = vol->getSymName();
270 for (const auto& rule : rules) {
271 const auto pattern = rule["match"].get<std::string>();
272 if (!matchPattern(pattern, sym)) {
273 continue;
274 }
275 // rigid body DOFs
276 if (rule.contains("rigidBody")) {
277 const auto& rb = rule["rigidBody"];
278 if (rb.is_string()) {
279 auto s = rb.get<std::string>();
280 if (s == "all" || s == "free") {
281 vol->setRigidBody(std::make_unique<RigidBodyDOFSet>());
282 } else if (s == "fixed") {
283 auto dofSet = std::make_unique<RigidBodyDOFSet>();
284 dofSet->setAllFree(false);
285 vol->setRigidBody(std::move(dofSet));
286 }
287 } else if (rb.is_array()) {
288 auto dofSet = std::make_unique<RigidBodyDOFSet>();
289 dofSet->setAllFree(false);
290 for (const auto& name : rb) {
291 auto it = rbNameToIdx.find(name.get<std::string>());
292 if (it != rbNameToIdx.end()) {
293 dofSet->setFree(it->second, true);
294 }
295 }
296 vol->setRigidBody(std::move(dofSet));
297 } else if (rb.is_object()) {
298 auto dofs = rb.value("dofs", std::string("all"));
299 bool fixed = rb.value("fixed", false);
300 if (dofs == "all") {
301 auto dofSet = std::make_unique<RigidBodyDOFSet>();
302 if (fixed) {
303 dofSet->setAllFree(false);
304 }
305 vol->setRigidBody(std::move(dofSet));
306 } else if (rb["dofs"].is_array()) {
307 auto dofSet = std::make_unique<RigidBodyDOFSet>();
308 dofSet->setAllFree(false);
309 for (const auto& name : rb["dofs"]) {
310 auto it = rbNameToIdx.find(name.get<std::string>());
311 if (it != rbNameToIdx.end()) {
312 dofSet->setFree(it->second, !fixed);
313 }
314 }
315 vol->setRigidBody(std::move(dofSet));
316 }
317 }
318 }
319 // calibration DOFs
320 if (rule.contains("calib")) {
321 const auto& cal = rule["calib"];
322 auto calType = cal.value("type", std::string(""));
323 if (calType == "legendre") {
324 int order = cal.value("order", 3);
325 auto dofSet = std::make_unique<LegendreDOFSet>(order);
326 bool fixed = cal.value("fixed", false);
327 if (fixed) {
328 dofSet->setAllFree(false);
329 }
330 // fix/free individual coefficients by name or index
331 if (cal.contains("free")) {
332 dofSet->setAllFree(false);
333 for (const auto& item : cal["free"]) {
334 if (item.is_number_integer()) {
335 dofSet->setFree(item.get<int>(), true);
336 } else if (item.is_string()) {
337 // match by name e.g. "L(1,0)"
338 for (int k = 0; k < dofSet->nDOFs(); ++k) {
339 if (dofSet->dofName(k) == item.get<std::string>()) {
340 dofSet->setFree(k, true);
341 }
342 }
343 }
344 }
345 }
346 if (cal.contains("fix")) {
347 for (const auto& item : cal["fix"]) {
348 if (item.is_number_integer()) {
349 dofSet->setFree(item.get<int>(), false);
350 } else if (item.is_string()) {
351 for (int k = 0; k < dofSet->nDOFs(); ++k) {
352 if (dofSet->dofName(k) == item.get<std::string>()) {
353 dofSet->setFree(k, false);
354 }
355 }
356 }
357 }
358 }
359 vol->setCalib(std::move(dofSet));
360 }
361 }
362 }
363 });
364}
365
366void writeMillepedeResults(AlignableVolume* root, const std::string& milleResPath, const std::string& outJsonPath, const std::string& injectedJsonPath)
367{
368 using json = nlohmann::json;
369
370 // parse millepede.res: label fittedValue presigma [...]
371 std::ifstream fin(milleResPath);
372 if (!fin.is_open()) {
373 LOGP(fatal, "Cannot open millepede result file: {}", milleResPath);
374 }
375 std::map<uint32_t, double> labelToValue;
376 std::string line;
377 while (std::getline(fin, line)) {
378 if (line.empty() || line[0] == '!' || line[0] == '*') {
379 continue;
380 }
381 if (line.find("Parameter") != std::string::npos) {
382 continue;
383 }
384 std::istringstream iss(line);
385 uint32_t label = 0;
386 double value = NAN, presigma = NAN;
387 if (!(iss >> label >> value >> presigma)) {
388 continue;
389 }
390 if (presigma >= 0.0) { // skip fixed parameters
391 labelToValue[label] = value;
392 }
393 }
394 fin.close();
395 LOGP(info, "Parsed {} not fixed parameters from {}", labelToValue.size(), milleResPath);
396
397 // load injected misalignment if provided (same format as closure test input)
398 // indexed by sensorID
399 std::map<int, std::vector<double>> injRB;
400 std::map<int, std::vector<std::vector<double>>> injMatrix;
401 if (!injectedJsonPath.empty()) {
402 std::ifstream injFile(injectedJsonPath);
403 if (injFile.is_open()) {
404 json injData = json::parse(injFile);
405 for (const auto& item : injData) {
406 int id = item["id"].get<int>();
407 if (item.contains("rigidBody")) {
408 injRB[id] = item["rigidBody"].get<std::vector<double>>();
409 }
410 if (item.contains("matrix")) {
411 injMatrix[id] = item["matrix"].get<std::vector<std::vector<double>>>();
412 }
413 }
414 LOGP(info, "Loaded injected misalignment for {} sensors", injData.size());
415 } else {
416 LOGP(warn, "Cannot open injected misalignment file: {}, writing absolute values", injectedJsonPath);
417 }
418 }
419
420 // collect results per volume that has RB or calib DOFs
421 json output = json::array();
422 root->traverse([&](AlignableVolume* vol) {
423 auto* rb = vol->getRigidBody();
424 auto* cal = vol->getCalib();
425 if ((!rb && !cal) || vol->isPseudo()) {
426 return;
427 }
428 int id = vol->getSensorId();
429 json entry;
430 entry["symName"] = vol->getSymName();
431 entry["id"] = id;
432 bool write = false;
433
434 // rigid body parameters
435 if (rb && rb->nFreeDOFs()) {
436 write = true;
437 json rbArr = json::array();
438 const auto& inj = injRB.contains(id) ? injRB[id] : std::vector<double>{};
439 for (int i = 0; i < rb->nDOFs(); ++i) {
440 uint32_t raw = vol->getLabel().raw(i);
441 auto it = labelToValue.find(raw);
442 double fitted = it != labelToValue.end() ? it->second : 0.0;
443 double ref = i < static_cast<int>(inj.size()) ? inj[i] : 0.0;
444 rbArr.push_back(fitted - ref);
445 }
446 entry["rigidBody"] = rbArr;
447 }
448
449 // calibration (Legendre) parameters
450 if (cal && cal->nFreeDOFs() && cal->type() == DOFSet::Type::Legendre) {
451 write = true;
452 auto* leg = dynamic_cast<const LegendreDOFSet*>(cal);
453 int order = leg->order();
454 auto calibLbl = vol->getLabel().asCalib();
455 const auto& inj = injMatrix.contains(id) ? injMatrix[id] : std::vector<std::vector<double>>{};
456 json matrix = json::array();
457 int idx = 0;
458 for (int i = 0; i <= order; ++i) {
459 json row = json::array();
460 for (int j = 0; j <= i; ++j) {
461 uint32_t raw = calibLbl.raw(idx);
462 auto it = labelToValue.find(raw);
463 double fitted = it != labelToValue.end() ? it->second : 0.0;
464 double ref = (i < static_cast<int>(inj.size()) && j < static_cast<int>(inj[i].size())) ? inj[i][j] : 0.0;
465 row.push_back(fitted - ref);
466 ++idx;
467 }
468 matrix.push_back(row);
469 }
470 entry["matrix"] = matrix;
471 }
472 if (write) {
473 output.push_back(entry);
474 }
475 });
476
477 std::ofstream fout(outJsonPath);
478 if (!fout.is_open()) {
479 LOGP(fatal, "Cannot open output file: {}", outJsonPath);
480 }
481 fout << output.dump(2) << '\n';
482 fout.close();
483 LOGP(info, "Wrote millepede results to {}", outJsonPath);
484}
485
486} // namespace o2::its3::align
General auxilliary methods.
std::ostringstream debug
int32_t i
Definition of the GeometryTGeo class.
void output(const std::map< std::string, ChannelStat > &channels)
Definition rawdump.cxx:197
uint32_t j
Definition RawData.h:0
uint32_t c
Definition RawData.h:2
nlohmann::json json
void traverse(const std::function< void(AlignableVolume *)> &visitor)
void writeParameters(std::ostream &os) const
GlobalLabel getLabel() const noexcept
AlignableVolume * mParent
matrices
void setRigidBody(std::unique_ptr< DOFSet > rb)
void writeTree(std::ostream &os, int indent=0) const
std::string getSymName() const noexcept
void writeRigidBodyConstraints(std::ostream &os) const
AlignableVolume(const AlignableVolume &)=delete
void setCalib(std::unique_ptr< DOFSet > cal)
GlobalLabel asCalib() const noexcept
return a copy of this label with the CALIB bit set (for calibration DOFs on same volume)
constexpr T raw(T dof) const noexcept
produce the raw Millepede label for a given DOF index (rigid body: calib=0 in label)
constexpr bool sens() const noexcept
void add(uint32_t lab, double coeff)
GLuint entry
Definition glcorearb.h:5735
GLsizeiptr size
Definition glcorearb.h:659
GLuint const GLchar * name
Definition glcorearb.h:781
GLdouble f
Definition glcorearb.h:310
GLsizei const GLfloat * value
Definition glcorearb.h:819
GLboolean * data
Definition glcorearb.h:298
GLuint GLsizei const GLchar * label
Definition glcorearb.h:2519
GLint level
Definition glcorearb.h:275
GLuint id
Definition glcorearb.h:650
void writeMillepedeResults(AlignableVolume *root, const std::string &milleResPath, const std::string &outJsonPath, const std::string &injectedJsonPath="")
void applyDOFConfig(AlignableVolume *root, const std::string &jsonPath)
std::vector< int > row
std::array< uint16_t, 5 > pattern