BmnRoot
Loading...
Searching...
No Matches
BmnAligner.cxx
Go to the documentation of this file.
1/*
2 BM@N alignment routine
3 BM@N experiment at NICA complex, JINR, 2025
4
5 Department: Math & Soft Group of HEP lab
6 Author: Igor Polev, polev@jinr.ru
7
8 BmnAligner class implementation
9*/
10
11#undef BMN_TEST_ALGORITHM // special test mode
12
13#include "BmnAligner.h"
14
15#include "BmnAlignDefines.h"
16#include "BmnDataIterator.h"
17#include "Math/Functor.h"
18#include "TCanvas.h"
19#include "TDecompLU.h"
20#include "TDecompSVD.h"
21#include "TError.h"
22#include "TPolyLine3D.h"
23#include "TRandom.h"
24#include "fairlogger/Logger.h"
25
26#include <fstream>
27#include <iomanip>
28#include <iostream>
29#include <nlohmann/json.hpp>
30
31#define BMN_MSE_SLOW_THRESHOLD 2.0 // time threshold in seconds to show MSE calculation progress bar
32
33template<typename HitType>
35 const BmnDetectorModel* modelDetector,
36 const char* filesListPath,
37 const char* constraintsPath,
38 Int_t nThreads)
39 : BmnAligner(measureModel, modelDetector, constraintsPath, nThreads)
40{
41 TStopwatch timer;
42 LOG(state) << " > BmnAligner : initialization started";
43
44 // initialize ROOT file iterator
45 LOG(info) << " ROOT files list: " << filesListPath;
48 LOG(error) << "<!> ROOT iterator initialization failed";
49 return;
50 }
56 // all track are usable by default
57 // TODO : replace std::vector with TBits
59 fpTrackUsable = new std::vector<Bool_t>(fTracksUsableCnt, kTRUE);
60
61 // RAM mode:
62 // If RAM requiered to store all hits is reasonable, we will read
63 // all hits from files and itarate them in RAM during calculations.
64 // In this case we need special iterator.
65 auto ram_size = fTotalHits * sizeof(HitType);
66 if (ram_size <= BMN_MAX_HIT_RAM_SIZE) {
67 fpFirstHits = new Int_t[fTotalTracks + 1]; // one extra element for iterator convenience
68 fpHitsData = new Byte_t[ram_size];
69 fpHits = reinterpret_cast<HitType*>(fpHitsData);
72 };
73 LOG(info) << Form(" %.1f GB of RAM required", ram_size / 1073741824.0);
74 LOG_IF(info, fpIteratorMain == fpRamIterator) << " RAM storage selected";
75 LOG_IF(info, fpIteratorMain == fpRootIterator) << " file storage selected";
77 // prepare set of iterators for parallelization by data
80 fInitialized = kTRUE;
81 DoneInTime(timer, "initialization");
84template<typename HitType>
86{
87 if (fpTrackUsable)
88 delete fpTrackUsable;
89 if (fpRootIterator)
90 delete fpRootIterator;
91 if (fpRamIterator)
92 delete fpRamIterator;
93 if (fpHitsData)
94 delete[] fpHitsData;
95 if (fpFirstHits)
96 delete[] fpFirstHits;
97 if (fThreadsCnt > 1)
98 for (auto pItr : fpIterators)
99 delete pItr;
100}
101
102// usage of option tracksPath is not implemented in this class
103template<typename HitType>
104Bool_t BmnAligner<HitType>::PrepareData(const char* tracksPath)
105{
106 TStopwatch timer;
107 LOG(state) << " > BmnAligner : data preparation started";
108 if (!fInitialized) {
109 LOG(error) << " <!> BmnAligner is not initialized";
110 return kFALSE;
111 }
112 LOG_IF(info, fpIteratorMain == fpRamIterator) << " loading data from ROOT files into RAM...";
113 LOG_IF(info, fpIteratorMain == fpRootIterator) << " parsing ROOT data files...";
114 fProgressBar.Init();
115
116 std::vector<SVectLC> arrAlpha(fTotalTracks); // zero approximation of linear tracks
117
118 Int_t detID, iTrack{0}, iHit{0}; // global hit index
119 Double_t Dx{0.0}, Dy{0.0};
120 fpRootIterator->ResetAll();
121 while (!fpRootIterator->EndOfTracks()) {
122 // in RAM mode we store index of first hits in track for memory optimization purposes
123 if (fpIteratorMain == fpRamIterator)
124 fpFirstHits[iTrack] = iHit;
125 // init linear regression terms
126 Double_t xSum{0.0}, ySum{0.0}, zSum{0.0}, z2Sum{0.0}, zxSum{0.0}, zySum{0.0};
127
128 // iterate hits in track
129 while (!fpRootIterator->EndOfHits()) {
130 detID = fpRootIterator->HitDetectorID();
131 // check if detector element encoding is valid
132 if (fpDetModel->UnknownID(detID)) {
133 const BmnRootIterator::Hit_t* hit{fpRootIterator->GetHit()};
134 LOG(error) << " <!> Unknown detector module encoding:"
135 << "\n iterator::HitDetectorID() = " << detID
136 << "\n FairHit::GetDetectorID() = " << hit->GetDetectorID()
137 << "\n CbmStsHit::GetSystemId() = " << hit->GetSystemId()
138 << "\n CbmStsHit::GetStationNr() = " << hit->GetStationNr()
139 << "\n CbmStsHit::GetSectorNr() = " << hit->GetSectorNr()
140 << "\n CbmStsHit::GetSensorNr() = " << hit->GetSensorNr();
141 return kFALSE;
142 }
143 // count hits for each detector element
144 fHitsPerDetector[fpDetModel->Idx(detID)]++;
145
146 // inplace hit copy
147 if (fpIteratorMain == fpRamIterator)
148 new (&fpHits[iHit]) HitType(fpRootIterator->GetHit(), *fpDetModel);
149
150 Double_t x{fpRootIterator->HitX()}, y{fpRootIterator->HitY()}, z{fpRootIterator->HitZ()};
151 // calculate linear regression terms
152 xSum += x;
153 ySum += y;
154 zSum += z;
155 z2Sum += z * z;
156 zxSum += z * x;
157 zySum += z * y;
158 // find max Z value (for visualization only)
159 if (z > fMaxZvalue)
160 fMaxZvalue = z;
161 // calculate average hit position errors
162 Dx += fpRootIterator->GetHit()->GetDx();
163 Dy += fpRootIterator->GetHit()->GetDy();
164
165 iHit++;
166 fpRootIterator->NextHit();
167 } // loop by hits
168
169 // calculate zero approximation of tracks as linear regression
170 Double_t hitCnt{static_cast<Double_t>(fpRootIterator->HitsInTrack())},
171 factor{1.0 / (hitCnt * z2Sum - zSum * zSum)};
172 arrAlpha[iTrack] = SVectLC((xSum * z2Sum - zSum * zxSum) * factor, (hitCnt * zxSum - zSum * xSum) * factor,
173 (ySum * z2Sum - zSum * zySum) * factor, (hitCnt * zySum - zSum * ySum) * factor);
174
175 iTrack++;
176 fpRootIterator->NextTrack();
177 } // loop by tracks
178 if (fpIteratorMain == fpRamIterator)
179 fpFirstHits[iTrack] = iHit; // last syntetic element for iterator convenience
180 StoreZeroSolution(arrAlpha);
181
182 // MSE scale factors calculation
183 Dx /= iHit;
184 Dy /= iHit;
185 fOmegaScaleFactor = (Dx + Dy) * 0.5;
186 fOmegaScaleFactor *= fOmegaScaleFactor;
187
188 // evaluate MSE calculation time
189 TStopwatch timerMSE;
190 timerMSE.Stop();
191 timerMSE.Reset();
192 SVectGLT point;
193 for (Int_t i = 0; i < BMN_GLOBAL_PARAMS_PD; i++) {
194 BMN_RND.RndmArray(BMN_GLPARAMS_TOTAL, point.Array());
195 timerMSE.Start(kFALSE);
196 CalculateMSE(point);
197 timerMSE.Stop();
198 }
199 fMSEclcSlow = timerMSE.RealTime() / BMN_GLOBAL_PARAMS_PD > BMN_MSE_SLOW_THRESHOLD;
200 if (!fMSEclcSlow)
201 fMSEclcTick = static_cast<UInt_t>(BMN_MSE_SLOW_THRESHOLD / timerMSE.RealTime() * BMN_GLOBAL_PARAMS_PD);
202
203 fProgressBar.Clear();
204 if (fpIteratorMain == fpRamIterator) {
205 LOG(info) << " " << iTrack << " tracks loaded";
206 LOG(info) << " " << iHit << " hits loaded";
207 }
208 LOG(info) << " average hit position errors: Dx = " << Dx << ", Dy = " << Dy;
209 DoneInTime(timer, "data preparation");
210 return LoadConstraints();
211}
212
213template<typename HitType>
215{
216 TStopwatch timer;
217 LOG(state) << " > BmnAligner : loading constraints...";
218 if (!fpConstraintsPath) {
219 LOG(warning) << " constraints file in not provided, proceeding without constraints";
220 return kTRUE;
221 }
222 std::ifstream ifs(fpConstraintsPath);
223 if (!ifs.good()) {
224 LOG(error) << " <!> LoadConstraints: failed to open JSON file " << *fpConstraintsPath;
225 return kFALSE;
226 }
227 nlohmann::json jsonData;
228 try {
229 ifs >> jsonData;
230 } catch (const std::exception& e) {
231 LOG(error) << " <!> LoadConstraints: failed to parse JSON file " << *fpConstraintsPath;
232 LOG(error) << e.what();
233 return kFALSE;
234 }
235 if (!jsonData.is_array()) {
236 LOG(error) << " <!> LoadConstraints: top-level JSON is not an array in " << *fpConstraintsPath;
237 return kFALSE;
238 }
239 fConstraintsCnt = jsonData.size();
240 if (fConstraintsCnt <= 0) {
241 LOG(error) << " <!> LoadConstraints: empty constraints array in " << *fpConstraintsPath;
242 return kFALSE;
243 }
244
245 fConstraints.resize(fConstraintsCnt);
246 for (Int_t i = 0; i < fConstraintsCnt; i++) {
247 const auto& jsonLine = jsonData[i];
248 if (!jsonLine.is_array()) {
249 LOG(error) << " <!> LoadConstraints: constraint " << i << " is not an array in " << fpConstraintsPath;
250 return kFALSE;
251 }
252 Int_t lineSize = jsonLine.size();
253 if (lineSize > BMN_GLPARAMS_TOTAL) {
254 LOG(warning) << " constraint " << i << " size exceeds global parameters count, truncating";
255 lineSize = BMN_GLPARAMS_TOTAL;
256 } else if (lineSize < BMN_GLPARAMS_TOTAL)
257 LOG(warning) << " constraint " << i << " size is less then global parameters count, padding with zeroes";
258 fConstraints[i].ResizeTo(BMN_GLPARAMS_TOTAL);
259 for (Int_t k = 0; k < lineSize; k++)
260 fConstraints[i](k) = jsonLine[k].get<Double_t>();
261 }
262
263 LOG(info) << " constraints loaded: " << fConstraintsCnt;
264 DoneInTime(timer, "loading constraints");
265 return kTRUE;
266}
267
268template<typename HitType>
269Bool_t BmnAligner<HitType>::LoadSolution(const char* jsonPath, Bool_t withAlpha, Bool_t withMSE)
270{
271 if (!jsonPath)
272 return kFALSE;
273 TStopwatch timer;
274 LOG(state) << " > BmnAligner : loading solution from file " << jsonPath;
275
276 std::ifstream ifs(jsonPath);
277 if (!ifs.good()) {
278 LOG(error) << " <!> LoadSolution: failed to open JSON file";
279 return kFALSE;
280 }
281 using namespace nlohmann;
282 json js;
283 Int_t elementsLoaded{0};
284 try {
285 ifs >> js;
286
287 // check if solution is suitable
288 json jsAlignmentType = BMN_CORRECTIONS;
289 if (js.contains("CorrectionValues")) {
290 if (js["CorrectionValues"] != jsAlignmentType) {
291 LOG(error) << " <!> LoadSolution: wrong alignment type in JSON file\n"
292 << " expected: " << jsAlignmentType << ", in file: " << js["CorrectionValues"];
293 return kFALSE;
294 }
295 } else
296 LOG(warning) << " CorrectionValues key is missing: assuming alignment type is correct";
297 if (js.contains("DetectorElements")) {
298 if (static_cast<Int_t>(js["DetectorElements"]) != BMN_MODULE_COUNT) {
299 LOG(error) << " <!> LoadSolution: wrong detector configuration in JSON file\n"
300 << " modules expected: " << BMN_MODULE_COUNT << ", in file: " << js["DetectorElements"];
301 return kFALSE;
302 }
303 } else
304 LOG(warning) << " DetectorElements key is missing: assuming detector configuration is correct";
305 if (!js.contains("CorrectionsPerDetector") || !js["CorrectionsPerDetector"].is_array()) {
306 LOG(error) << " LoadSolution: CorrectionsPerDetector key is missing or not an array";
307 return kFALSE;
308 }
309
310 // load solution for each stored detector element
311 Int_t elementsCount{static_cast<Int_t>(js["CorrectionsPerDetector"].size())};
312 for (Int_t i = 0; i < elementsCount; i++) {
313 json& element = js["CorrectionsPerDetector"][i];
314 if (!element.contains("detector_id")) {
315 LOG(warning) << " ATTENTION! DetectorID is missing for some element, skipping.";
316 continue;
317 }
318 Int_t detID = element["detector_id"].get<Int_t>();
319 if (fpDetModel->UnknownID(detID)) {
320 LOG(warning) << " ATTENTION! Unregistered detector element ID: " << detID << ", skipping.";
321 continue;
322 }
323 if (!element.contains("values") || !element["values"].is_array()
324 || element["values"].size() != BMN_GLOBAL_PARAMS_PD)
325 {
326 LOG(warning) << " ATTENTION! Correction values for the element " << detID
327 << " are invalid, skipping.";
328 continue;
329 }
330 if (!elementsLoaded)
331 AddBlankSolution();
332 Int_t detIdx = fpDetModel->Idx(detID);
333 json& values = element["values"];
334 SVectGL& A = fpResultCurrent->A(detIdx);
335 for (Int_t j = 0; j < BMN_GLOBAL_PARAMS_PD; j++)
336 A(j) = values[j].get<Double_t>();
337 if (element.contains("isFixed")) {
338 Bool_t isFixed = element["isFixed"].get<Bool_t>();
339 if (fDetectorActive[detIdx] == isFixed) {
340 fDetectorActive[detIdx] = not isFixed;
341 if (isFixed)
342 fActiveDetCount--;
343 else
344 fActiveDetCount++;
345 }
346 }
347 elementsLoaded++;
348 }
349 fpResultCurrent->ConcatenateA();
350 LOG(info) << " alignment corrections loaded for " << elementsLoaded << " detector elements";
351
352 // TODO: load tracks approximation form JSON file if available
353 if (withAlpha) {
354 if (fpResultLast) {
355 LOG(info) << " loading tracks approximation from last knownsolution";
356 fpResultCurrent->Alpha() = fpResultLast->Alpha();
357 } else
358 LOG(warning) << " ATTENTION! Local parameters approximation not found.";
359 }
360 } catch (std::exception& e) {
361 LOG(error) << " <!> LoadSolution: failed to parse JSON file " << jsonPath;
362 LOG(error) << e.what();
363 return kFALSE;
364 }
365 if (withMSE)
366 CalculateMSE(kTRUE);
367
368 DoneInTime(timer, "solution loaded");
369 return kTRUE;
370}
371
372template<typename HitType>
373Bool_t BmnAligner<HitType>::SaveSolution(const char* jsonPath, Bool_t withAlpha)
374{
375 TStopwatch timer;
376 LOG(state) << " > BmnAligner : saving solution to file " << jsonPath;
377 std::ofstream ofs(jsonPath);
378 if (!ofs.good()) {
379 LOG(error) << " <!> SaveSolution: failed to open output file";
380 return kFALSE;
381 }
382
383 BmnAlignResult& solution{*GetCurrentResult()};
384 using namespace nlohmann;
385 json js; // all data to be saved
386
387 // alignment solution
388 js["CorrectionValues"] = BMN_CORRECTIONS;
389 js["DetectorElements"] = BMN_MODULE_COUNT;
390 js["CorrectionsPerDetector"] = json::array();
391 for (Int_t i = 0; i < BMN_MODULE_COUNT; i++) {
392 json element;
393 element["detector_id"] = fpDetModel->EncodedID(i);
394 element["values"] = json::array();
395 for (Int_t j = 0; j < BMN_GLOBAL_PARAMS_PD; j++)
396 element["values"].push_back(solution.A(i)(j));
397 js["CorrectionsPerDetector"].push_back(element);
398 }
399
400 // tracks approximation
401 if (withAlpha) {
402 js["MinHitsPerTrack"] = BMN_MIN_HITS_PER_TRACK;
403 js["TracksCount"] = fTotalTracks;
404 js["TracksParameters"] = json::array();
405 for (Int_t i = 0; i < fTotalTracks; i++) {
406 json element = json::array();
407 for (Int_t j = 0; j < BMN_LOCAL_PARAMS_PT; j++)
408 element.push_back(solution.Alpha(i)(j));
409 js["TracksParameters"].push_back(element);
410 }
411 }
412
413 // write data
414 try {
415 ofs << std::setw(4) << js << std::endl;
416 } catch (std::exception& e) {
417 LOG(error) << " <!> SaveSolution: failed to write JSON data to file\n" << e.what();
418 return kFALSE;
419 }
420
421 DoneInTime(timer, "solution save");
422 return kTRUE;
423}
424
425/*
426 * Main alignment loop.
427 * Due to non-linearity of dependence between measurements and alignment parameters,
428 * iterative alignment is required.
429 * fixEvery - after main alignmenet procedure fix corrections for detector elements
430 * one by one and repeat alignment (seems to be useless)
431 */
432template<typename HitType>
433Bool_t BmnAligner<HitType>::Align(Bool_t fixEvery)
434{
435 TStopwatch timer;
436 LOG(state) << " > BmnAligner : iterative alignment started";
437
438 // check detector elements for data sufficiency
439 Int_t minHits{static_cast<Int_t>(BMN_MIN_HITS_PORTION * fTotalTracks / fpDetModel->MaxModulesInStation())};
440 for (Int_t i = 0; i < BMN_MODULE_COUNT; i++) {
441 if (!fDetectorActive[i])
442 continue;
443 if (fHitsPerDetector[i] < minHits) {
444 LOG(warning) << Form(" detector element %d has only %d hits of %d required, it will be fixed",
445 fpDetModel->EncodedID(i), fHitsPerDetector[i], minHits);
446 fDetectorActive[i] = kFALSE;
447 fActiveDetCount--;
448 }
449 }
450
451 // iterative detector elements deactivation - optional feature
452 while (fActiveDetCount >= BMN_MIN_ACTIVE_DETECTORS) {
453 LOG(info) << " " << fActiveDetCount << " detector elements active of " << BMN_MODULE_COUNT << " total";
454
455 // main alignment loop
456 Int_t regularRetries{BMN_MAX_REGULAR_RETRY}, retriesLast;
457 Double_t progressAbs{0.0}, progressRel{BMN_MIN_REL_PROGRESS};
458 while (std::abs(progressRel) >= BMN_MIN_REL_PROGRESS || regularRetries > 0)
459 // regularRetries is never negative, but it is more safe - in case of future code changes
460 {
461 if (progressAbs < 0.0) { // no progress at all
462 if (!fRegularIters) // regularization as a last resort
463 fRegularIters = BMN_MIN_REGULAR_ITERS;
464 else {
465 LOG(error) << " <!> BmnAligner::Align : algorithm is diverging";
466 return kFALSE;
467 }
468 } else if (progressRel < BMN_MIN_REL_PROGRESS && !fRegularIters) // progress is poor
469 fRegularIters = BMN_MIN_REGULAR_ITERS; // but we have regularization to try
470
471 // perform alignment iteration
472 if (!IterateAlignment()) {
473 LOG(error) << " <!> BmnAligner::Align : iteration failed";
474 return kFALSE;
475 }
476 progressAbs = fpResultLast->GetValueMSE() - fpResultCurrent->GetValueMSE();
477 progressRel = progressAbs / fpResultLast->GetValueMSE();
478 LOG(info) << Form(" = %+.3f%% relative progress, MSE change : %.3e -> %.3e\n", progressRel * 100.0,
479 fpResultLast->GetValueMSE(), fpResultCurrent->GetValueMSE());
480
481 retriesLast = regularRetries;
482 if (fRegularIters) {
483 // decide about regularization
484 if (progressRel >= BMN_REL_PROGRESS_GOOD)
485 fRegularIters = BMN_MIN_REGULAR_ITERS; // keep going!
486 else
487 fRegularIters--; // attention! be careful to keep it non-negative
488 if (!fRegularIters && regularRetries)
489 regularRetries--;
490 }
491 LOG_IF(info, retriesLast != regularRetries) << " " << regularRetries << " regularization retries left";
492 } // main alignment loop
493
494 // fix best aligned detector element (optional)
495 if (fixEvery) {
496 BmnAlignResult::IdxValuePair_t MSElist{fpResultCurrent->GetSortedMSE()};
497 for (const auto& [idx, mse] : MSElist)
498 if (fDetectorActive[idx]) {
499 fDetectorActive[idx] = kFALSE;
500 fActiveDetCount--;
501 break;
502 }
503 } else
504 break;
505 } // detector elements deactivation loop
506
507 LOG(info) << " " << fIterationResults.size() - 1 << " iterations until algorithm convergence criteria reached";
508 DoneInTime(timer, "iterative alignment");
509 ReportResults();
510 return kTRUE;
511}
512
513/*
514 * One iteration of alignment. Consists of two parts:
515 * 1. parallel processing of tracks
516 * 2. final calculation of (global) alignment parameters
517 *
518 * Parallelized (intermediate) calculations implemented in ProcessTracks()
519 * method below, here only parallelization control is performed. Thus,
520 * calulations code is splitted between two methods: ProcessTracks() and
521 * IterateAlignment(). This code is marked by noticable comments.
522 *
523 * Calculations in ProcessTracks() performed according to:
524 * Volker Blobel, Claus Kleinwort
525 * "A New Method for the High-Precision Alignment of Track Detectors"
526 * http://arxiv.org/abs/hep-ex/0208021v1
527 * also widely known as Millepede method.
528 * This paper is referred below as _VBCK_. Variable naming is also corresponds to _VBCK_.
529 *
530 * Attention! In _VBCK_ (12) term z_k refers to track model which is incorrect
531 * in non-linear case. In general case it should refer to a _difference_ between
532 * measured value and track model. More accurate definition is given in documentation
533 * on Millepede-II software package published here:
534 * https://www.desy.de/~kleinwrt/MP2/doc/html/draftman_page.html
535 */
536template<typename HitType>
538{
539 TStopwatch timer;
540 fIteration++;
541 LOG(state) << " > #" << fIteration + 1 << " iteration started";
542 if (!fInitialized) {
543 LOG(error) << " <!> BmnAligner is not initialized";
544 return kFALSE;
545 }
546 AddBlankSolution();
547
548 // prepare storage for parallel calculations
549 const Int_t ActiveParamsCnt{fActiveDetCount * BMN_GLOBAL_PARAMS_PD};
550 TVectorD vB_dash(ActiveParamsCnt); // b' from _VBCK_ (17)
551 TMatrixD mC_dash(ActiveParamsCnt, ActiveParamsCnt); // C' from _VBCK_ (17)
552
553 // run Millepede-like dimension reduction calculations in parallel
554 LOG(info) << " " << fThreadsCnt << " threads started to perform dimension reduction (Millepede)";
555 fProgressBar.Init();
556 TStopwatch timerMille;
557 _BMN_ALIGNER_RUN_THREADS(&BmnAligner::ProcessTracks, this, &vB_dash, &mC_dash)
558 LOG(info) << Form(" %.1fs (REAL) / %.1fs (CPU) took to complete dimension reduction", timerMille.RealTime(),
559 timerMille.CpuTime());
560 fProgressBar.Clear();
561
562 /*
563 * ---------------------------------------------------
564 * Final calculations of global parameters: _VBCK_(16)
565 * ---------------------------------------------------
566 */
567
568 /*
569 * Prepare system of linear equations for global parameters
570 * --------------------------------------------------------
571 * 1. make a copy of vB_dash and mC_dash (original matrices will be used too)
572 * 2. append with user defined constraints taking into account active detectors
573 * 3. remove zero constraints if any (because of inactive detectors)
574 * 4. iteratively solve linear systems with different regularization factors
575 * 5. check results of iterative regularization
576 * 6. select optimal regularization factor if possible
577 */
578
579 TStopwatch timerGlobal;
580 fMSEclcCounter = 0u;
581
582 TVectorD vBL(vB_dash);
583 TMatrixD mCL(mC_dash);
584 if (fConstraintsCnt) {
585 // add user defined constraints to the system
586 vBL.ResizeTo(ActiveParamsCnt + fConstraintsCnt);
587 mCL.ResizeTo(ActiveParamsCnt + fConstraintsCnt, ActiveParamsCnt + fConstraintsCnt);
588 Double_t weight;
589 std::vector<Bool_t> constActive(fConstraintsCnt, false);
590 for (Int_t i = 0, k = 0; i < BMN_GLPARAMS_TOTAL; i++) {
591 if (!fDetectorActive[i / BMN_GLOBAL_PARAMS_PD])
592 continue;
593 for (Int_t j = 0; j < fConstraintsCnt; j++) {
594 weight = fConstraints[j](i);
595 if (weight) {
596 mCL(k, ActiveParamsCnt + j) = weight;
597 mCL(ActiveParamsCnt + j, k) = weight;
598 constActive[j] = kTRUE;
599 }
600 }
601 k++;
602 }
603 // remove zero constraints if any
604 TMatrixD row(1, mCL.GetNcols());
605 Int_t dstIdx{ActiveParamsCnt - 1}, // one less then first constraint
606 maxIdx{ActiveParamsCnt + fConstraintsCnt - 1};
607 for (Int_t srcIdx = ActiveParamsCnt, i = 0; srcIdx <= maxIdx; srcIdx++, i++) {
608 if (!constActive[i])
609 continue;
610 if (++dstIdx == srcIdx)
611 continue;
612 mCL.GetSub(srcIdx, srcIdx, 0, maxIdx, row);
613 mCL.SetSub(dstIdx, 0, row);
614 mCL.SetSub(0, dstIdx, row.T()); // mCL is symmetric
615 }
616 dstIdx++; // total number of active rows in mCL
617 mCL.ResizeTo(dstIdx, dstIdx);
618 vBL.ResizeTo(dstIdx);
619 LOG(info) << " " << dstIdx - ActiveParamsCnt << " active parameter constraints applied";
620 }
621 // TODO: remove inactivated constraints from fConstraints
622
623 if (!fRegularIters) { // no regularization
624 LOG(info) << " solving non-regularized system for global parameters...";
625 if (!SolveGlobal(mCL, vBL)) {
626 LOG(warning) << " iteration failed to improve MSE";
627 *fpResultCurrent = *fpResultLast;
628 }
629 } else { // Levenberg-Marquardt (LM) regularization
630
631 LOG(info) << " starting Levenberg-Marquardt regularization...";
632 // find crude range for LM factor in which to perform 1D optimization
633 Double_t valueMSEbest{fpResultLast->GetValueMSE()}, LMfactorBest{0.0}, LMfactor{1.0};
634 Int_t LMiterarion{0};
635 Bool_t LMworking{kTRUE};
636 while ((LMworking || LMiterarion < BMN_LM_ITERATIONS_MIN) && LMiterarion < BMN_LM_ITERATIONS_MAX) {
637 if (LMfactor > 1.0)
638 // apply normalization factor to diagonal elements
639 for (Int_t i = 0; i < ActiveParamsCnt; i++)
640 mCL(i, i) = mC_dash(i, i) * LMfactor;
641
642 if (SolveGlobal(mCL, vBL)) {
643 if (LMfactor == 1.0)
644 // at first check if non-regularized solution is good enough
645 if (1.0 - fpResultCurrent->GetValueMSE() / valueMSEbest >= BMN_REL_PROGRESS_GOOD) {
646 // cancel regularization mode but keep retry count
647 fRegularIters = 0;
648 break;
649 }
650
651 // evaluate results of successful iteration
652 if (fpResultCurrent->GetValueMSE() < valueMSEbest) {
653 valueMSEbest = fpResultCurrent->GetValueMSE();
654 LMfactorBest = LMfactor;
655 } else
656 LMworking = kFALSE; // stop if MSE is not improved
657 } else
658 LOG(warning) << Form(" LM regularization failed (F=%.3e)", LMfactor);
659
660 LMfactor *= BMN_LM_MULT_STEP;
661 LMiterarion++;
662 } // end of Levenberg-Marquardt regularization loop
663
664 // regularization factor optimization function
665 ROOT::Math::Functor1D fnLM([&mCL, &mC_dash, &vBL, this, ActiveParamsCnt](Double_t factor) {
666 TMatrixD mCLcopy{mCL};
667 for (Int_t i = 0; i < ActiveParamsCnt; i++)
668 mCLcopy(i, i) = mC_dash(i, i) * factor;
669 if (!SolveGlobal(mCLcopy, vBL)) // expected to allways succeed
670 throw std::runtime_error(Form("Levenberg-Marquardt failed (F=%.3e)", factor));
671 return fpResultCurrent->GetValueMSE();
672 });
673 // classification of Levenberg-Marquardt results
674 if (!fRegularIters) {
675 CalculateMSE(kTRUE);
676 LOG(info) << Form(
677 " regularization cancelled : sufficient %+.1f%% of progress achieved using non-regularized method",
678 100.0 * (fpResultCurrent->GetValueMSE() / fpResultLast->GetValueMSE() - 1.0));
679
680 } else if (LMfactorBest < 1.0) {
681 LOG(warning) << " iteration failed to improve MSE";
682 *fpResultCurrent = *fpResultLast;
683
684 } else if (LMfactorBest == 1.0 || LMworking) { // it works like XOR here
685 // LMfactorBest has boundary values
686 if (LMworking)
687 LOG(info) // max possible LMfactor
688 << Form(" F = %.3e : LM regularization approached gradient descent", LMfactorBest);
689 else
690 LOG(info) // min possible LMfactor (1.0)
691 << " F = 1.0 : LM regularization did not improve results";
692 fnLM(LMfactorBest);
693 CalculateMSE(kTRUE);
694
695 } else { // LMfactorBest > 1.0 and LMworking is false
696 // search for optimal LMfactor between LMfactorBest and LMfactorBest * BMN_LM_MULT_STEP
697
698 ROOT::Math::BrentMinimizer1D brent;
699 brent.SetNpx(BMN_BRENT_GRID_SIZE);
700 brent.SetFunction(fnLM, LMfactorBest, LMfactorBest * BMN_LM_MULT_STEP);
701 Bool_t success{brent.Minimize(static_cast<Int_t>(0.1 / BMN_1D_STEP_TOLERANCE), BMN_1D_STEP_TOLERANCE * 10.0,
703 if (success)
704 LOG(info) << Form(" F = %.3e : optimal LM factor value", brent.XMinimum());
705 else
706 LOG(warning) << Form(" F = %.3e : rough LM factor value used, 1D optimization failed",
707 LMfactorBest);
708 fnLM(success ? brent.XMinimum() : LMfactorBest);
709 CalculateMSE(kTRUE);
710
711 } // end if (LM results classification)
712 } // end if (regularization)
713 LOG(info) << Form(" %d MSE value calculations performed during iteration", fMSEclcCounter);
714 LOG(info) << Form(" %.1fs (REAL) / %.1fs (CPU) took to complete parameters estimation", timerGlobal.RealTime(),
715 timerGlobal.CpuTime());
716
717 // don't call MarkTime() here - it will be done in Align()
718 LOG(state) << Form(" < #%d iteration took %.1fs (REAL) / %.1fs (CPU)", fIteration + 1, timer.RealTime(),
719 timer.CpuTime());
720 return kTRUE;
721}
722
723template<typename HitType>
724Bool_t BmnAligner<HitType>::SolveGlobal(const TMatrixD& mCL, const TVectorD& vBL)
725{
726 /*
727 * Solve system of linear equations for global parameters
728 * ------------------------------------------------------
729 * 1. get direct solution using LU or SVD decomposition
730 * 2. store Lagrange multipliers and SVD singular values (for analysis if required)
731 * 3. append solution with zeroes for inactive detectors
732 * 4. find optimal step size in proposed direction
733 * 5. calculate new approximation of global parameters
734 * 6. in case of any failure proceed to next regularization factor
735 */
736 // TODO: not necessary to solve full system with large LMfactor
737 const Int_t ActiveParamsCnt{fActiveDetCount * BMN_GLOBAL_PARAMS_PD};
738
739 // Suppress ROOT error messages: LU decomposition expected to fail sometimes
740 auto oldErrorLevel{gErrorIgnoreLevel};
741 gErrorIgnoreLevel = kBreak;
742 TDecompLU mLU{mCL};
743 Bool_t success{mLU.Decompose()};
744 TDecompSVD mSVD{mCL}; // we need SVD anyway to evaluate singularity of system
745 mSVD.Decompose();
746 gErrorIgnoreLevel = oldErrorLevel; // Restore original error level
747 fpResultCurrent->SVDSig().ResizeTo(mSVD.GetSig());
748 fpResultCurrent->SVDSig() = mSVD.GetSig();
749
750 TVectorD vAL{success ? mLU.Invert(success) * vBL : mSVD.Solve(vBL, success)};
751 if (!success) {
752 LOG(warning) << " <!> global matrix inversion unsuccessful";
753 return kFALSE;
754 }
755
756 // separate Lagrange multipliers
757 Int_t activeConstCnt{vAL.GetNrows() - ActiveParamsCnt};
758 if (activeConstCnt) {
759 fpResultCurrent->L().ResizeTo(activeConstCnt);
760 vAL.GetSub(ActiveParamsCnt, vAL.GetNrows() - 1, fpResultCurrent->L());
761 }
762 // expand vAL to full vector of global parameters
763 Double_t *dst{fpResultCurrent->A().begin()}, *src{vAL.GetMatrixArray()};
764 for (Int_t i = 0; i < BMN_MODULE_COUNT; i++)
765 if (fDetectorActive[i])
766 for (Int_t j = 0; j < BMN_GLOBAL_PARAMS_PD; j++)
767 *dst++ = *src++;
768 else
770
771 if (fRegularIters) { // find optimal step in proposed direction
772 // 1D function to minimize
773 ROOT::Math::Functor1D fnMSE([this](Double_t step) {
774 SVectGLT point{fpResultLast->A() + fpResultCurrent->A() * step};
775 return CalculateMSE(point);
776 });
777 // search for initial segment containing the minimum
778 Double_t valueMSE, stepMin{0.0}, stepMax{1.0};
779 valueMSE = fnMSE(stepMax);
780 while (valueMSE < fpResultLast->GetValueMSE()) {
781 stepMin = stepMax;
782 stepMax *= 2.0;
783 valueMSE = fnMSE(stepMax);
784 }
785 // find minimum using Brent method
786 constexpr static const Int_t iterations{static_cast<Int_t>(1.0 / BMN_1D_STEP_TOLERANCE)};
787 fBrent.SetFunction(fnMSE, stepMin, stepMax);
788 success = fBrent.Minimize(iterations, BMN_1D_STEP_TOLERANCE, BMN_1D_STEP_TOLERANCE * 0.1);
789 if (!success) {
790 LOG(warning) << " Brent method failed to calculate step size";
791 return kFALSE;
792 }
793
794 // select new approximation of global parameters
795 fpResultCurrent->A() *= fBrent.XMinimum();
796 } // end if (regularization)
797
798 fpResultCurrent->A() += fpResultLast->A();
799 fpResultCurrent->SeparateA();
800 CalculateMSE(kFALSE);
801
802 return kTRUE;
803}
804
805/*
806 * Millepede-like dimension reduction algorithm.
807 * For detailes see:
808 * http://arxiv.org/abs/hep-ex/0208021v1
809 * https://www.desy.de/~kleinwrt/MP2/doc/html/draftman_page.html
810 *
811 * Each thread processes separate subset of all tracks defined by iterator.
812 * Synchronization is necessary at result storage stage only.
813 */
814template<typename HitType>
815void BmnAligner<HitType>::ProcessTracks(TVectorD* pvB_dash, TMatrixD* pmC_dash, Int_t threadNum)
816{
817 // prepare vector of zeroes to be used for variables initialization
818 constexpr static const Int_t SMatrMaxSize{
819 std::max(std::max(static_cast<Int_t>(SMatrGL::kSize), static_cast<Int_t>(SMatrLC::kSize)),
820 static_cast<Int_t>(SMatrGP::kSize))};
821 static const std::vector<Double_t> zeroes(SMatrMaxSize, 0.0);
822 static const auto zStart = zeroes.begin();
823 const Int_t ActiveParamsCnt{fActiveDetCount * BMN_GLOBAL_PARAMS_PD};
824
825 // ROOT matrices and vectors are initialized with zeroes during construction
826 // accelerated ROOT vectors and matrices
827 SVectGL vGradG,
828 vB[BMN_MODULE_COUNT]; // get memory for full set of detectors,
829 SMatrGL mC[BMN_MODULE_COUNT]; // not for active only, for loops optimization
831 SMatrLC mGamma;
832 SVectLC vBeta, vDelta, vGradL;
833 // standard non-accelerated aliases to variables of accelerated types
834 TVectorD vaDelta(BMN_LOCAL_PARAMS_PT);
835 TMatrixD maGamma(BMN_LOCAL_PARAMS_PT, BMN_LOCAL_PARAMS_PT);
836 vaDelta.Use(BMN_LOCAL_PARAMS_PT, vDelta.Array());
837 maGamma.Use(BMN_LOCAL_PARAMS_PT, BMN_LOCAL_PARAMS_PT, mGamma.Array());
838 // non-accelerated ROOT matrices and vectors due to runtime-defined sizes
839 TMatrixD mGfull(ActiveParamsCnt, BMN_LOCAL_PARAMS_PT), mCfull(ActiveParamsCnt, ActiveParamsCnt),
840 mC_part(ActiveParamsCnt, ActiveParamsCnt);
841 TVectorD vBfull(ActiveParamsCnt), vB_part(ActiveParamsCnt);
842
843 Double_t z, omega;
844
845 // iterate tracks
846 BmnDataIterator* pItr{fpIterators[threadNum]};
847 Int_t iTrack{fFirstTracksIdx[threadNum]};
848 pItr->ResetAll();
849 while (!pItr->EndOfTracks()) {
850
851 // skip unusable tracks
852 if (!fpTrackUsable->at(iTrack)) {
853 iTrack++;
854 pItr->NextTrack();
855 continue;
856 }
857
858 /*
859 * ---------------------------------------------------
860 * Local parameters section
861 * ---------------------------------------------------
862 */
863
864 // reset cumulative variables to zero
865 vBeta.SetElements(zStart, SVectLC::kSize);
866 mGamma.SetElements(zStart, SMatrLC::kSize);
867
868 // get previous approximation of local parameters for this track
869 SVectLC& lc{fpResultLast->Alpha(iTrack)};
870 // iterate hits
871 pItr->ResetHits();
872 while (!pItr->EndOfHits()) {
873 // TODO: cache gradient calculations - it does not depend on x and y
874
875 // get previous approximation of global parameters for this hit (detector)
876 const SVectGL& gl{fpResultLast->A(fpDetModel->Idx(pItr->HitDetectorID()))};
877
878 // notice: the code below is independent of corrections sign selection
879 z = pItr->HitZ();
880
881 omega = fOmegaScaleFactor * pItr->HitWx();
882 vGradL = fpMeasureModel->GradXlocal(z, gl, lc);
883 mGamma += omega * ROOT::Math::TensorProd(vGradL, vGradL);
884 vBeta += (omega * (pItr->HitX() - fpMeasureModel->X(z, gl, lc))) * vGradL;
885
886 omega = fOmegaScaleFactor * pItr->HitWy();
887 vGradL = fpMeasureModel->GradYlocal(z, gl, lc);
888 mGamma += omega * ROOT::Math::TensorProd(vGradL, vGradL);
889 vBeta += (omega * (pItr->HitY() - fpMeasureModel->Y(z, gl, lc))) * vGradL;
890
891 pItr->NextHit();
892 } // loop by hits
893
894 // local parameters alpha : _VBCK_(11)
895 // pay attantion that inverse of Gamma is stored in mGamma
896 // TODO: check if Cholesky inversion is aplicable
897 if (!mGamma.Invert())
898 throw std::runtime_error("matrix inversion failed while calculating local parameters");
899 vDelta = mGamma * vBeta;
900 // store updated local parameters (no sycn is needed between threads)
901 fpResultCurrent->Alpha(iTrack) = fpResultLast->Alpha(iTrack) + vDelta;
902
903 /*
904 * ---------------------------------------------------
905 * Global parameters section
906 * ---------------------------------------------------
907 */
908
909 // TODO: put global section together with local in same loop - can't use updated local parameter values
910
911 // reset cumulative variables to zero
912 for (Int_t i = 0; i < BMN_MODULE_COUNT; i++) {
913 if (!fDetectorActive[i])
914 continue; // working only with active detectors
915 vB[i].SetElements(zStart, SVectGL::kSize);
916 mC[i].SetElements(zStart, SMatrGL::kSize);
917 mG[i].SetElements(zStart, SMatrGP::kSize);
918 }
919
920 // iterate hits
921 pItr->ResetHits();
922 while (!pItr->EndOfHits()) {
923 const Int_t iDet{fpDetModel->Idx(pItr->HitDetectorID())};
924 if (!fDetectorActive[iDet]) { // working only with active detectors
925 pItr->NextHit();
926 continue;
927 }
928
929 // We use global parameters approximation from previous iteration,
930 // although we have new, more precise values, at this moment.
931 // It is intended: usage of new values is not consistent with
932 // mathematical model of MILLEPEDE method of dimemsionality reduction.
933 const SVectGL& gl{fpResultLast->A(iDet)};
934
935 // notice: the code below is independent of corrections sign selection
936 z = pItr->HitZ();
937
938 omega = fOmegaScaleFactor * pItr->HitWx();
939 vGradL = fpMeasureModel->GradXlocal(z, gl, lc);
940 vGradG = fpMeasureModel->GradXglobal(z, gl, lc);
941 mG[iDet] += omega * ROOT::Math::TensorProd(vGradG, vGradL);
942 mC[iDet] += omega * ROOT::Math::TensorProd(vGradG, vGradG);
943 vB[iDet] += (omega * (pItr->HitX() - fpMeasureModel->X(z, gl, lc))) * vGradG;
944
945 omega = fOmegaScaleFactor * pItr->HitWy();
946 vGradL = fpMeasureModel->GradYlocal(z, gl, lc);
947 vGradG = fpMeasureModel->GradYglobal(z, gl, lc);
948 mG[iDet] += omega * ROOT::Math::TensorProd(vGradG, vGradL);
949 mC[iDet] += omega * ROOT::Math::TensorProd(vGradG, vGradG);
950 vB[iDet] += (omega * (pItr->HitY() - fpMeasureModel->Y(z, gl, lc))) * vGradG;
951
952 pItr->NextHit();
953 } // loop by hits
954
955 // collecting C_dash, b_dash : _VBCK_(17)
956 TVectorD vaB(BMN_GLOBAL_PARAMS_PD);
958 for (Int_t i = 0, j = 0; i < BMN_MODULE_COUNT; i++) {
959 if (!fDetectorActive[i])
960 continue; // working only with active detectors
961 vaB.Use(BMN_GLOBAL_PARAMS_PD, vB[i].Array());
962 maG.Use(BMN_GLOBAL_PARAMS_PD, BMN_LOCAL_PARAMS_PT, mG[i].Array());
963 maC.Use(BMN_GLOBAL_PARAMS_PD, BMN_GLOBAL_PARAMS_PD, mC[i].Array());
964 vBfull.SetSub(j, vaB);
965 mGfull.SetSub(j, 0, maG);
966 mCfull.SetSub(j, j, maC);
967 j += BMN_GLOBAL_PARAMS_PD; // destination index
968 }
969 TMatrixD mGfullT(TMatrixD::kTransposed, mGfull);
970 mC_part += mCfull - mGfull * maGamma * mGfullT;
971 vB_part += vBfull - mGfull * vaDelta;
972
973 iTrack++;
974 pItr->NextTrack();
975 } // loop by tracks
976
977 // store final result
978 std::lock_guard<std::mutex> lock(fMutexGlobal);
979 *pmC_dash += mC_part;
980 *pvB_dash += vB_part;
981}
982
983template<typename HitType>
985{
986 // to calculate MSE at arbitrary point we substitute this point in place of
987 // final solution to keep cumbersome multithreaded MSE calculation scheme intact
988 BmnAlignResult *pBackupSolution{fpResultCurrent}, tempSolution(*fpResultCurrent);
989 tempSolution.A() = solution;
990 tempSolution.SeparateA();
991 fpResultCurrent = &tempSolution;
992 Double_t valueMSE{CalculateMSE(kFALSE)};
993 fpResultCurrent = pBackupSolution;
994 return valueMSE;
995}
996
997template<typename HitType>
998Double_t BmnAligner<HitType>::CalculateMSE(Bool_t withResiduals)
999{
1000 fMSEclcCounter++; // MSE calculation is expensive because it requires to traverse all hits
1001 fpResultCurrent->ResetValueMSE();
1002 if (withResiduals) {
1003 fpResultCurrent->ResidualsX().Reset();
1004 fpResultCurrent->ResidualsY().Reset();
1005 }
1006 if (fMSEclcSlow)
1007 fProgressBar.Init();
1008 _BMN_ALIGNER_RUN_THREADS(&BmnAligner::CalculateMSEpart, this, withResiduals)
1009 if (fMSEclcSlow)
1010 fProgressBar.Clear();
1011 else if (0 == fMSEclcCounter % fMSEclcTick)
1012 std::cout << " " << fMSEclcCounter << " MSE calculations\r" << std::flush;
1013 return fpResultCurrent->GetValueMSE();
1014}
1015
1016template<typename HitType>
1017void BmnAligner<HitType>::CalculateMSEpart(Bool_t withResiduals, Int_t threadNum)
1018{
1019 Double_t delta, omega;
1020 std::array<Double_t, BMN_MODULE_COUNT> arrMSE;
1021 TH1D *pResX, *pResY;
1022
1023 arrMSE.fill(0.0);
1024 if (withResiduals) {
1025 pResX = new TH1D(Form("residualsX_%lu", fpResultCurrent->GetHistCounter()), "Partial residuals by X",
1027 pResY = new TH1D(Form("residualsY_%lu", fpResultCurrent->GetHistCounter()), "Partial residuals by Y",
1029 }
1030
1031 Int_t iTrack;
1032 BmnDataIterator* pItr;
1033 if (threadNum == -1) { // single-thread computation
1034 iTrack = 0;
1035 pItr = fpIteratorMain;
1036 } else { // multi-thread computation
1037 iTrack = fFirstTracksIdx[threadNum];
1038 pItr = fpIterators[threadNum];
1039 }
1040 pItr->ResetAll();
1041 while (!pItr->EndOfTracks()) {
1042 if (!fpTrackUsable->at(iTrack)) {
1043 iTrack++;
1044 pItr->NextTrack();
1045 continue;
1046 }
1047 // notice: the code below is independent of corrections sign selection
1048 const SVectLC& lc{fpResultCurrent->Alpha(iTrack)};
1049 while (!pItr->EndOfHits()) {
1050 Int_t detectorIdx{fpDetModel->Idx(pItr->HitDetectorID())};
1051 const SVectGL& gl{fpResultCurrent->A(detectorIdx)};
1052
1053 omega = fOmegaScaleFactor * pItr->HitWx();
1054 delta = pItr->HitX() - fpMeasureModel->X(pItr->HitZ(), gl, lc);
1055 arrMSE[detectorIdx] += omega * delta * delta;
1056 if (withResiduals && TMath::Abs(delta) > BMN_ALMOST_ZERO)
1057 // ignore almost precise zero residuals which occure
1058 // by construction of local parameters zero approximation
1059 pResX->Fill(delta);
1060
1061 omega = fOmegaScaleFactor * pItr->HitWy();
1062 delta = pItr->HitY() - fpMeasureModel->Y(pItr->HitZ(), gl, lc);
1063 arrMSE[detectorIdx] += omega * delta * delta;
1064 if (withResiduals && TMath::Abs(delta) > BMN_ALMOST_ZERO)
1065 // ignore almost precise zero residuals which occure
1066 // by construction of local parameters zero approximation
1067 pResY->Fill(delta);
1068
1069 pItr->NextHit();
1070 }
1071 iTrack++;
1072 pItr->NextTrack();
1073 }
1074
1075 // store results
1076 std::lock_guard<std::mutex> lock(fMutexGlobal);
1077 for (Int_t iDet = 0; iDet < BMN_MODULE_COUNT; iDet++)
1078 fpResultCurrent->AddValueMSE(iDet, arrMSE[iDet]);
1079 if (withResiduals) {
1080 fpResultCurrent->ResidualsX().Add(pResX);
1081 fpResultCurrent->ResidualsY().Add(pResY);
1082 delete pResX;
1083 delete pResY;
1084 }
1085}
1086
1087template<typename HitType>
1089{
1090 BmnAlignResult &before{*GetResult(kFALSE)}, &after{*GetResult(kTRUE)};
1091 Int_t idxLast{after.SVDSig().GetNrows() - 1};
1092 static const char* drawLine{"================================================================================"};
1093 LOG(info) << "";
1094 LOG(info) << drawLine;
1095 LOG(info) << Form("%.2e -> %.2e (%+.1f%%) : MSE value change", before.GetValueMSE(), after.GetValueMSE(),
1096 100.0 * (after.GetValueMSE() - before.GetValueMSE()) / before.GetValueMSE());
1097 LOG(info) << Form(
1098 "%.2e -> %.2e (%+.1f%%) : mean error in X change", before.ResidualsX().GetMean(), after.ResidualsX().GetMean(),
1099 100.0 * (after.ResidualsX().GetMean() - before.ResidualsX().GetMean()) / before.ResidualsX().GetMean());
1100 LOG(info) << Form(
1101 "%.2e -> %.2e (%+.1f%%) : mean error in Y change", before.ResidualsY().GetMean(), after.ResidualsY().GetMean(),
1102 100.0 * (after.ResidualsY().GetMean() - before.ResidualsY().GetMean()) / before.ResidualsY().GetMean());
1103 LOG(info) << Form("%.2e -> %.2e (%+.1f%%) : error stddev in X change", before.ResidualsX().GetStdDev(),
1104 after.ResidualsX().GetStdDev(),
1105 100.0 * (after.ResidualsX().GetStdDev() - before.ResidualsX().GetStdDev())
1106 / before.ResidualsX().GetStdDev());
1107 LOG(info) << Form("%.2e -> %.2e (%+.1f%%) : error stddev in Y change", before.ResidualsY().GetStdDev(),
1108 after.ResidualsY().GetStdDev(),
1109 100.0 * (after.ResidualsY().GetStdDev() - before.ResidualsY().GetStdDev())
1110 / before.ResidualsY().GetStdDev());
1111 LOG(info) << Form("SVD minimum singular value : %.2e", after.SVDSig()(idxLast));
1112 LOG(info) << Form("SVD condition number : %.2e", after.SVDSig()(0) / after.SVDSig()(idxLast));
1113 LOG(info) << drawLine;
1114 Int_t cpuTime{static_cast<Int_t>(fTotalCpuTime)}, realTime{static_cast<Int_t>(fTotalRealTime)};
1115 LOG(info) << Form("total CPU time spent .......... %dh %02dm %02ds", cpuTime / 3600, (cpuTime / 60) % 60,
1116 cpuTime % 60);
1117 LOG(info) << Form("total real time spent ......... %dh %02dm %02ds", realTime / 3600, (realTime / 60) % 60,
1118 realTime % 60);
1119 LOG(info) << Form("multithreading boost factor ... x%.1f", fTotalCpuTime / fTotalRealTime);
1120}
1121
1122template<typename HitType>
1124{
1125 // check number of requested threads
1126 Int_t max_threads = fTotalTracks / BMN_MIN_TRACKS_PER_THREAD;
1127 if (!max_threads)
1128 max_threads = 1;
1129 if (max_threads < fThreadsCnt) {
1130 LOG(warning) << " too many threads requested for the size of dataset, using " << max_threads
1131 << " threads instead of " << fThreadsCnt;
1132 fThreadsCnt = max_threads;
1133 }
1134 if (fpIteratorMain == fpRootIterator && fThreadsCnt > 1) {
1135 LOG(warning) << " parallelization implemented for RAM mode only";
1136 fThreadsCnt = 1;
1137 }
1138
1139 if (fThreadsCnt > 1) {
1140 // create separate RamIterator for each thread
1141 fFirstTracksIdx.resize(fThreadsCnt);
1142 fpIterators.resize(fThreadsCnt);
1143 Int_t tracksPerThread = fTotalTracks / fThreadsCnt;
1144 for (Int_t i = 0, iTrack = 0; i < fThreadsCnt; i++, iTrack += tracksPerThread) {
1145 fFirstTracksIdx[i] = iTrack;
1146 fpIterators[i] = new BmnRamIterator(fpRamIterator->Slice(
1147 iTrack,
1148 i == fThreadsCnt - 1 // last thread case
1149 ? fTotalTracks - tracksPerThread * i // the rest of tracks for the last thread
1150 : tracksPerThread));
1151 }
1152 } else {
1153 fFirstTracksIdx.resize(1);
1154 fpIterators.resize(1);
1155 fFirstTracksIdx[0] = 0;
1156 fpIterators[0] = fpIteratorMain;
1157 }
1158}
1159
1160template<typename HitType>
1161void BmnAligner<HitType>::DoneInTime(TStopwatch& timer, const char* stage)
1162{
1163 MarkTime(timer);
1164 LOG(state) << Form(" < BmnAligner : %s completed in %.1fs (REAL) / %.1fs (CPU)", stage, timer.RealTime(),
1165 timer.CpuTime());
1166}
1167
1168// visualization works in ROOT macro mode only
1169template<typename HitType>
1171{
1172 // prepare canvas and pads
1173 static const UInt_t ww{1900}, wh{1010};
1174 TCanvas* pCnv = new TCanvas( // RAM is not released in macro mode
1175 "canvas", "BM@N alignment");
1176 pCnv->SetRealAspectRatio();
1177 pCnv->SetFixedAspectRatio();
1178 pCnv->SetCanvasSize(ww, wh);
1179 pCnv->SetWindowSize(ww + 20, wh + 70);
1180
1181 // divide canvas into pads
1182 // pay attention: dynamic RAM is not released in ROOT macro mode
1183 fpPad3D = new TPad("pad3D", "3D", 0.0, 0.34, 0.6, 1.0);
1184 fpPadG1 = new TPad("padG1", "G1", 0.0, 0.0, 0.6, 0.34);
1185 fpPadG2 = new TPad("padG2", "G2", 0.6, 0.67, 1.0, 1.0);
1186 fpPadG3 = new TPad("padG3", "G3", 0.6, 0.34, 1.0, 0.67);
1187 fpPadG4 = new TPad("padG4", "G4", 0.6, 0.0, 1.0, 0.34);
1188 fpPad3D->Draw();
1189 fpPadG1->Draw();
1190 fpPadG2->Draw();
1191 fpPadG3->Draw();
1192 fpPadG4->Draw();
1193}
1194
1195template<typename HitType>
1197{
1198 // Draw residuals histograms
1199 TH1D* pHist;
1200
1201 fpPadG3->cd();
1202 fpPadG3->SetTitle("Residuals by Y");
1203 // fpPadG3->SetLogy();
1204 pHist = &GetResult(kFALSE)->ResidualsY();
1205 pHist->SetFillColor(kCyan);
1206 pHist->SetLineColor(kCyan);
1207 pHist->GetYaxis()->SetRangeUser(0.0, GetResult(kTRUE)->ResidualsY().GetMaximum() * 1.1);
1208 pHist->Draw();
1209 pHist = &GetResult(kTRUE)->ResidualsY();
1210 pHist->SetLineColor(kRed);
1211 pHist->Draw("CSAME");
1212
1213 fpPadG2->cd();
1214 fpPadG2->SetTitle("Residuals by X");
1215 // fpPadG2->SetLogy();
1216 pHist = &GetResult(kFALSE)->ResidualsX();
1217 pHist->SetFillColor(kTeal);
1218 pHist->SetLineColor(kTeal);
1219 pHist->GetYaxis()->SetRangeUser(0.0, GetResult(kTRUE)->ResidualsX().GetMaximum() * 1.1);
1220 pHist->Draw();
1221 pHist = &GetResult(kTRUE)->ResidualsX();
1222 pHist->SetLineColor(kRed);
1223 pHist->Draw("CSAME");
1224}
1225
1226template<typename HitType>
1228{
1229 PrepareDrawing();
1230
1231 // draw some tracks in 3D pad
1232
1233 fpPad3D->cd();
1234 fpPad3D->SetTitle("Some random tracks");
1235 const Double_t maxZdraw{fMaxZvalue * 1.1};
1236 const Double_t trackProbability{static_cast<Double_t>(BMN_DRAW_TRACKS_CNT)
1237 / static_cast<Double_t>(fTracksUsableCnt)};
1238 TRandom trueRnd(0);
1239 TPolyLine3D *pTrackLines = new TPolyLine3D[BMN_DRAW_TRACKS_CNT],
1240 *pApproxLines = new TPolyLine3D[BMN_DRAW_TRACKS_CNT];
1241
1242 fProgressBar.Init();
1243 Int_t iTrack{0}, iTrackDrawn{0};
1244 fpIteratorMain->ResetAll();
1245 while (iTrackDrawn < BMN_DRAW_TRACKS_CNT && !fpIteratorMain->EndOfTracks()) {
1246 // skip tracks according to probability
1247 if (!fpTrackUsable->at(iTrack) || trueRnd.Uniform() > trackProbability) {
1248 fpIteratorMain->NextTrack();
1249 iTrack++;
1250 continue;
1251 }
1252
1253 // draw track itself
1254 Int_t iHit{0};
1255 while (!fpIteratorMain->EndOfHits()) {
1256 // x and z are swapped for visual convenience
1257 pTrackLines[iTrackDrawn].SetPoint(iHit, fpIteratorMain->HitZ(), fpIteratorMain->HitY(),
1258 fpIteratorMain->HitX());
1259 fpIteratorMain->NextHit();
1260 iHit++;
1261 } // loop by hits
1262 pTrackLines[iTrackDrawn].SetLineColor(kBlue);
1263 pTrackLines[iTrackDrawn].Draw();
1264
1265 // draw liner approximation of the track
1266 const SVectLC& Alpha{GetResult()->Alpha(iTrack)};
1267 pApproxLines[iTrackDrawn].SetPoint(0, 0.0, Alpha(iY0), Alpha(iX0));
1268 pApproxLines[iTrackDrawn].SetPoint(1, maxZdraw, Alpha(iY0) + maxZdraw * Alpha(iY1),
1269 Alpha(iX0) + maxZdraw * Alpha(iX1));
1270 pApproxLines[iTrackDrawn].SetLineColor(kRed);
1271 pApproxLines[iTrackDrawn].SetLineStyle(kDotted);
1272 pApproxLines[iTrackDrawn].Draw();
1273
1274 iTrackDrawn++;
1275 fpIteratorMain->NextTrack();
1276 iTrack++;
1277 } // loop by tracks
1278 fProgressBar.Clear();
1279
1280 DrawResiduals();
1281}
1282
1283// Instantiate the template
1284template class BmnAligner<>;
1285
1286#include "BmnATestHit.h"
1287template class BmnAligner<BmnATestHit>;
#define BMN_MIN_REL_PROGRESS
#define BMN_CORRECTIONS
ROOT::Math::SMatrix< Double_t, BMN_GLOBAL_PARAMS_PD, BMN_GLOBAL_PARAMS_PD > SMatrGL
#define BMN_BRENT_GRID_SIZE
#define BMN_RESIDUALS_BINS
@ iY1
@ iY0
@ iX1
@ iX0
ROOT::Math::SVector< Double_t, BMN_GLPARAMS_TOTAL > SVectGLT
#define BMN_LM_ITERATIONS_MAX
ROOT::Math::SMatrix< Double_t, BMN_GLOBAL_PARAMS_PD, BMN_LOCAL_PARAMS_PT > SMatrGP
#define BMN_DRAW_TRACKS_CNT
ROOT::Math::SMatrix< Double_t, BMN_LOCAL_PARAMS_PT, BMN_LOCAL_PARAMS_PT > SMatrLC
#define BMN_MIN_HITS_PORTION
ROOT::Math::SVector< Double_t, BMN_GLOBAL_PARAMS_PD > SVectGL
#define BMN_MIN_TRACKS_PER_THREAD
#define BMN_1D_STEP_TOLERANCE
#define BMN_RESIDUALS_RANGE_X
#define BMN_REL_PROGRESS_GOOD
#define BMN_MIN_HITS_PER_TRACK
#define BMN_LM_MULT_STEP
#define BMN_MIN_ACTIVE_DETECTORS
#define BMN_RESIDUALS_RANGE_Y
ROOT::Math::SVector< Double_t, BMN_LOCAL_PARAMS_PT > SVectLC
#define BMN_MODULE_COUNT
#define BMN_MIN_REGULAR_ITERS
#define BMN_LM_ITERATIONS_MIN
#define BMN_MAX_HIT_RAM_SIZE
#define BMN_ALMOST_ZERO
#define BMN_MAX_REGULAR_RETRY
#define BMN_LOCAL_PARAMS_PT
#define BMN_GLOBAL_PARAMS_PD
#define BMN_MSE_SLOW_THRESHOLD
HitType
Definition BmnBaseHit.h:15
const Float_t delta
Distance between GEM-stations.
Definition BmnMwpcHit.cxx:8
int i
Definition P4_F32vec4.h:22
#define _BMN_ALIGNER_RUN_THREADS(...)
Definition BmnAligner.h:148
TVectorD & SVDSig() noexcept
std::array< std::pair< Int_t, Double_t >, BMN_MODULE_COUNT > IdxValuePair_t
Bool_t LoadConstraints()
HitType * fpHits
Definition BmnAligner.h:117
Int_t fTotalHits
Definition BmnAligner.h:104
BmnDataIterator * fpIteratorMain
Definition BmnAligner.h:116
virtual ~BmnAligner()
BmnSimpleProgressBar fProgressBar
Definition BmnAligner.h:95
virtual Bool_t IterateAlignment()
Int_t fTracksUsableCnt
Definition BmnAligner.h:104
std::vector< Bool_t > * fpTrackUsable
Definition BmnAligner.h:103
Int_t fTotalTracks
Definition BmnAligner.h:104
virtual Bool_t PrepareData(const char *tracksPath=nullptr)
Bool_t SaveSolution(const char *jsonPath, Bool_t withAlpha=kFALSE)
BmnRamIterator< HitType > * fpRamIterator
Definition BmnAligner.h:115
void DrawResiduals()
Int_t * fpFirstHits
Definition BmnAligner.h:119
BmnRootIterator * fpRootIterator
Definition BmnAligner.h:114
Bool_t LoadSolution(const char *jsonPath, Bool_t withAlpha=kTRUE, Bool_t withMSE=kTRUE)
Byte_t * fpHitsData
Definition BmnAligner.h:118
Bool_t Align(Bool_t fixEvery=kFALSE)
BmnAligner(const BmnMeasureModel *measureModel, const BmnDetectorModel *modelDetector, const char *constraintsPath=nullptr, Int_t nThreads=1)
Definition BmnAligner.h:163
void DoneInTime(TStopwatch &timer, const char *stage)
Bool_t fInitialized
Definition BmnAligner.h:96
Double_t CalculateMSE(const SVectGLT &solution)
void PrepareDrawing()
const BmnDetectorModel * fpDetModel
Definition BmnAligner.h:109
void InitIterators()
virtual void ReportResults() const
virtual void Draw()
Bool_t SolveGlobal(const TMatrixD &mCL, const TVectorD &vBL)
virtual Double_t HitX() const =0
virtual Double_t HitY() const =0
virtual Double_t HitZ() const =0
virtual Double_t HitWx() const =0
virtual void ResetAll()=0
virtual void NextTrack()=0
Bool_t EndOfTracks() const noexcept
virtual Int_t HitDetectorID() const =0
virtual void NextHit()=0
virtual Double_t HitWy() const =0
Bool_t EndOfHits() const noexcept
Bool_t Initialized() const noexcept
Int_t GetTotalTracks() const noexcept
Int_t GetTotalHits() const noexcept
Int_t GetSystemId() const
Definition CbmStsHit.h:64
a class to store JSON values
Definition json.hpp:17282
bool contains(KeyT &&key) const
check the existence of an element in a JSON object
Definition json.hpp:19640
size_type size() const noexcept
returns the number of elements
Definition json.hpp:19842
constexpr bool is_array() const noexcept
return whether value is an array
Definition json.hpp:18518
void push_back(basic_json &&val)
add an object to an array
Definition json.hpp:19986
auto get() const noexcept(noexcept(std::declval< const basic_json_t & >().template get_impl< ValueType >(detail::priority_tag< 4 > {}))) -> decltype(std::declval< const basic_json_t & >().template get_impl< ValueType >(detail::priority_tag< 4 > {}))
get a (pointer) value (explicit)
Definition json.hpp:18898
namespace for Niels Lohmann
Definition json.hpp:98