BmnRoot
Loading...
Searching...
No Matches
BmnAligner.h
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 declaration
9 Alignment algorithm itself
10*/
11
12#ifndef BMNALIGNER_H
13#define BMNALIGNER_H
14
15#include "BmnAlignDefines.h"
16#include "BmnAlignResult.h"
17#include "BmnDataIterator.h"
18#include "BmnDetectorModel.h"
19#include "BmnHitRecord.h"
20#include "BmnMeasureModel.h"
21#include "BmnRamIterator.h"
22#include "BmnRootIterator.h"
24#include "Math/BrentMinimizer1D.h"
25#include "TPad.h"
26#include "TStopwatch.h"
27
28#include <array>
29#include <deque>
30#include <mutex>
31#include <thread> // IWYU pragma: export
32#include <vector>
33
34// different types of hits to be used in test tasks and real tasks
35template<typename HitType = BmnHitRecord>
36class BmnAligner
37{
38 public:
39 BmnAligner( // basic constructor
40 const BmnMeasureModel* measureModel,
41 const BmnDetectorModel* modelDetector,
42 const char* constraintsPath = nullptr,
43 Int_t nThreads = 1);
44 BmnAligner( // main constructor
45 const BmnMeasureModel* measureModel,
46 const BmnDetectorModel* modelDetector,
47 const char* filesListPath,
48 const char* constraintsPath = nullptr,
49 Int_t nThreads = 1);
50 virtual ~BmnAligner();
51
52 // user interface methods
53 Bool_t Align(Bool_t fixEvery = kFALSE);
54 virtual Bool_t PrepareData(const char* tracksPath = nullptr);
55 Bool_t SaveSolution(const char* jsonPath, Bool_t withAlpha = kFALSE);
56 Bool_t LoadSolution(const char* jsonPath, Bool_t withAlpha = kTRUE, Bool_t withMSE = kTRUE);
57 virtual void Draw();
58 // getters
59 Bool_t Initialized() const noexcept { return fInitialized; }
60 Int_t TotalTracks() const noexcept { return fTotalTracks; }
61 Int_t TotalHits() const noexcept { return fTotalHits; }
62 Double_t OmegaScaleFactor() const noexcept { return fOmegaScaleFactor; }
63 BmnAlignResult* GetCurrentResult() const noexcept { return fpResultCurrent; }
64 BmnAlignResult* GetResult(Bool_t last = kTRUE) const noexcept
65 {
66 if (fIterationResults.empty())
67 return nullptr;
68 else if (last)
69 return const_cast<BmnAlignResult*>(&fIterationResults.back());
70 else
71 return const_cast<BmnAlignResult*>(&fIterationResults.front());
72 }
73
74 protected:
75 // main internal methods
76 virtual Bool_t IterateAlignment();
77 Double_t CalculateMSE(const SVectGLT& solution); // does not store result
78 virtual Double_t CalculateMSE(Bool_t withResiduals); // stores result in fResultCurrent
79 Bool_t SolveGlobal(const TMatrixD& mCL, const TVectorD& vBL);
80 virtual void ReportResults() const;
81 // internal service methods
86 inline void AddBlankSolution();
87 inline void StoreZeroSolution(const std::vector<SVectLC>& AlphaZero);
88 inline void MarkTime(TStopwatch& timer);
89 void DoneInTime(TStopwatch& timer, const char* stage);
90
91 typedef std::array<Bool_t, BMN_MODULE_COUNT> DetectorFlags_t;
92 typedef std::array<Int_t, BMN_MODULE_COUNT> DetectorHitCount_t;
93 typedef ROOT::Math::BrentMinimizer1D Minimizer1D_t;
94 // service variables
96 Bool_t fInitialized{kFALSE}, fMSEclcSlow{kFALSE};
97 UInt_t fMSEclcCounter{0u}, fMSEclcTick{1u};
98 Double_t fTotalCpuTime{0.0}, fTotalRealTime{0.0},
99 // calculation state and parameters
103 std::vector<Bool_t>* fpTrackUsable{nullptr}; // TODO: replace with TBits
107 // task definition
110 const char* fpConstraintsPath;
111 std::vector<TVectorD> fConstraints{}; // array of constraint vectors c
112 Int_t fConstraintsCnt{0};
113 // data storage implementation
114 BmnRootIterator* fpRootIterator{nullptr}; // iterator over ROOT files
115 BmnRamIterator<HitType>* fpRamIterator{nullptr}; // iterator over RAM storage
116 BmnDataIterator* fpIteratorMain{nullptr}; // currently used iterator
117 HitType* fpHits{nullptr}; // typed link to raw hits data storage
118 Byte_t* fpHitsData{nullptr}; // raw data storage for hits
119 Int_t* fpFirstHits{nullptr}; // indexes of first hit of each track
120 // iteration results storage
121 BmnAlignResult *fpResultLast{nullptr}, // link to result of previous iteration
122 *fpResultCurrent{nullptr}; // link to result of latest iteration
123 std::deque<BmnAlignResult> fIterationResults; // all results of calculations
124
125 // drawing stuff - works only in ROOT macro mode
126 // TODO: move drawing functionality to some other class (? BmnAlignResult)
127 TPad *fpPad3D{nullptr}, *fpPadG1{nullptr}, *fpPadG2{nullptr}, *fpPadG3{nullptr}, *fpPadG4{nullptr};
128 Double_t fMaxZvalue{0.0};
129
130 // multithreading implementation
131 private:
132 std::vector<BmnDataIterator*> fpIterators; // one iterator per thread
133 std::vector<Int_t> fFirstTracksIdx; // first track to process for each thread
134 std::mutex fMutexGlobal;
135 Int_t fThreadsCnt{1};
136 // thread functions: must have last Int_t argument to support _BMN_ALIGNER_RUN_THREADS macro
137 void ProcessTracks(TVectorD* pvB_dash, TMatrixD* pmC_dash, Int_t threadNum);
138 void CalculateMSEpart(Bool_t withResiduals, Int_t threadNum);
139
140 protected:
141 Int_t GetThreadCount() const noexcept { return fThreadsCnt; }
143};
144
145// macro is used instead of method for simplicity of variable function arguments implementation
146// GetProgressBar() and GetThreadCount() methods are used instead of direct fields access to make
147// the macro available for descendant classes while keeping multithread implementation private
148#define _BMN_ALIGNER_RUN_THREADS(...) \
149 { \
150 std::thread** _THRDS = new std::thread*[GetThreadCount()]; \
151 for (auto i = 0; i < GetThreadCount(); i++) \
152 _THRDS[i] = new std::thread(__VA_ARGS__, i); \
153 for (auto i = 0; i < GetThreadCount(); i++) { \
154 _THRDS[i]->join(); \
155 delete _THRDS[i]; \
156 } \
157 delete[] _THRDS; \
158 }
159
160// inline implementation
161
162template<typename HitType>
163inline BmnAligner<HitType>::BmnAligner(const BmnMeasureModel* measureModel,
164 const BmnDetectorModel* modelDetector,
165 const char* constraintsPath,
166 Int_t nThreads)
167 : fpMeasureModel(measureModel)
168 , fpDetModel(modelDetector)
169 , fpConstraintsPath(constraintsPath)
170 , fThreadsCnt(nThreads)
171{
172 TH1D::AddDirectory(kFALSE); // unbind all historgrams from default file
173 TH1D::StatOverflows(kTRUE); // do not ignore overflows in mean and stddev
174 fBrent.SetNpx(BMN_BRENT_GRID_SIZE); // lower bracketing points number to speed up Brent algorithm
175 fDetectorActive.fill(kTRUE); // parameters of all detectors are active initially
176};
177
178template<typename HitType>
180{
181 fIterationResults.emplace_back(fTotalTracks);
182 fpResultLast = fpResultCurrent;
183 fpResultCurrent = &fIterationResults.back();
184}
185
186template<typename HitType>
187inline void BmnAligner<HitType>::StoreZeroSolution(const std::vector<SVectLC>& AlphaZero)
188{
189 AddBlankSolution();
190 fpResultCurrent->Alpha() = AlphaZero;
191 CalculateMSE(kTRUE);
192}
193
194template<typename HitType>
195inline void BmnAligner<HitType>::MarkTime(TStopwatch& timer)
196{
197 timer.Stop();
198 fTotalCpuTime += timer.CpuTime();
199 fTotalRealTime += timer.RealTime();
200}
201
202#endif // BMNALIGNER_H
#define BMN_BRENT_GRID_SIZE
ROOT::Math::SVector< Double_t, BMN_GLPARAMS_TOTAL > SVectGLT
#define BMN_MODULE_COUNT
HitType
Definition BmnBaseHit.h:15
BmnAligner(const BmnMeasureModel *measureModel, const BmnDetectorModel *modelDetector, const char *filesListPath, const char *constraintsPath=nullptr, Int_t nThreads=1)
virtual ~BmnAligner()
virtual Bool_t IterateAlignment()
DetectorFlags_t fDetectorActive
Definition BmnAligner.h:102
Double_t fTotalCpuTime
Definition BmnAligner.h:98
BmnAlignResult * fpResultCurrent
Definition BmnAligner.h:122
Int_t fConstraintsCnt
Definition BmnAligner.h:112
Bool_t LoadConstraints()
std::vector< TVectorD > fConstraints
Definition BmnAligner.h:111
Double_t fTotalRealTime
Definition BmnAligner.h:98
HitType * fpHits
Definition BmnAligner.h:117
BmnSimpleProgressBar & GetProgressBar() noexcept
Definition BmnAligner.h:142
std::array< Int_t, BMN_MODULE_COUNT > DetectorHitCount_t
Definition BmnAligner.h:92
Int_t fTotalHits
Definition BmnAligner.h:104
void AddBlankSolution()
BmnAlignResult * fpResultLast
Definition BmnAligner.h:121
std::array< Bool_t, BMN_MODULE_COUNT > DetectorFlags_t
Definition BmnAligner.h:91
Bool_t Initialized() const noexcept
Definition BmnAligner.h:59
virtual void Draw()
Int_t fActiveDetCount
Definition BmnAligner.h:101
std::deque< BmnAlignResult > fIterationResults
Definition BmnAligner.h:123
TPad * fpPadG3
Definition BmnAligner.h:127
void StoreZeroSolution(const std::vector< SVectLC > &AlphaZero)
Int_t TotalHits() const noexcept
Definition BmnAligner.h:61
BmnDataIterator * fpIteratorMain
Definition BmnAligner.h:116
ROOT::Math::BrentMinimizer1D Minimizer1D_t
Definition BmnAligner.h:93
BmnSimpleProgressBar fProgressBar
Definition BmnAligner.h:95
Int_t fRegularIters
Definition BmnAligner.h:101
virtual void ReportResults() const
Int_t fIteration
Definition BmnAligner.h:101
UInt_t fMSEclcCounter
Definition BmnAligner.h:97
Int_t fTracksUsableCnt
Definition BmnAligner.h:104
virtual Double_t CalculateMSE(Bool_t withResiduals)
std::vector< Bool_t > * fpTrackUsable
Definition BmnAligner.h:103
Int_t fTotalTracks
Definition BmnAligner.h:104
Double_t fMaxZvalue
Definition BmnAligner.h:128
DetectorHitCount_t fHitsPerDetector
Definition BmnAligner.h:105
Bool_t SaveSolution(const char *jsonPath, Bool_t withAlpha=kFALSE)
BmnRamIterator< HitType > * fpRamIterator
Definition BmnAligner.h:115
Double_t fOmegaScaleFactor
Definition BmnAligner.h:100
TPad * fpPadG4
Definition BmnAligner.h:127
void DrawResiduals()
Int_t TotalTracks() const noexcept
Definition BmnAligner.h:60
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)
BmnAlignResult * GetCurrentResult() const noexcept
Definition BmnAligner.h:63
TPad * fpPad3D
Definition BmnAligner.h:127
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)
void DoneInTime(TStopwatch &timer, const char *stage)
Bool_t fInitialized
Definition BmnAligner.h:96
Double_t OmegaScaleFactor() const noexcept
Definition BmnAligner.h:62
Double_t CalculateMSE(const SVectGLT &solution)
void PrepareDrawing()
Int_t GetThreadCount() const noexcept
Definition BmnAligner.h:141
const BmnDetectorModel * fpDetModel
Definition BmnAligner.h:109
void InitIterators()
const BmnMeasureModel * fpMeasureModel
Definition BmnAligner.h:108
UInt_t fMSEclcTick
Definition BmnAligner.h:97
Minimizer1D_t fBrent
Definition BmnAligner.h:106
const char * fpConstraintsPath
Definition BmnAligner.h:110
BmnAlignResult * GetResult(Bool_t last=kTRUE) const noexcept
Definition BmnAligner.h:64
virtual Bool_t PrepareData(const char *tracksPath=nullptr)
TPad * fpPadG1
Definition BmnAligner.h:127
Bool_t fMSEclcSlow
Definition BmnAligner.h:96
Bool_t SolveGlobal(const TMatrixD &mCL, const TVectorD &vBL)
TPad * fpPadG2
Definition BmnAligner.h:127
void MarkTime(TStopwatch &timer)