BmnRoot
Loading...
Searching...
No Matches
BmnAlignerApply.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 BmnAlignerApply class implementation
9*/
10
11#include "BmnAlignerApply.h"
12
13#include "BmnGemStripStationSet.h"
14#include "BmnSiliconStationSet.h"
15#include "DstEventHeader.h"
16#include "FairRootManager.h"
17#include "FairRunAna.h"
18#include "TObjString.h"
19#include "TSystem.h"
20#include "fairlogger/Logger.h"
21
22#include <iostream>
23
25{
26 FairRootManager* ioman = FairRootManager::Instance();
27 if (!ioman) {
28 LOG(error) << " <!> BmnAligner::Init() : FairRootManager not found";
29 return kERROR;
30 }
31 fpStsHitsArray = (TClonesArray*)ioman->GetObject("StsHit");
32 if (!fpStsHitsArray) {
33 LOG(error) << " <!> BmnAligner::Init() : Failed to get StsHit array";
34 return kERROR;
35 }
36 try {
37 if (!LoadSolution(fAlignFilePath, kFALSE, kFALSE)) {
38 LOG(error) << " <!> BmnAlignerApply::Init() : Failed to load alignment solution from file "
39 << fAlignFilePath;
40 return kERROR;
41 }
42 } catch (std::exception& e) {
43 LOG(error) << " <!> BmnAlignerApply::Init() : " << e.what();
44 return kERROR;
45 } catch (...) {
46 LOG(error) << " <!> BmnAlignerApply::Init() : Unknown exception";
47 return kERROR;
48 }
49 if (!fLorentzCorrectionApply)
50 return kSUCCESS;
51
52 /*
53 * Lorentz corrections section-----------------------------------------------------------------
54 *
55 * The code below is directly taken from BmnToCbmHitConverter::ApplyAlignment(CbmStsHit*).
56 * It implements the Lorentz correction of hit positions and was not revised, except few
57 * lines marked with tag PLV. Indentation is the kept as in the original code.
58 */
59
60 FairField* mField = FairRunAna::Instance()->GetField();
61 fieldFlag = (TMath::Abs(mField->GetBy(0, 0, 135)) > 0.5) ? kTRUE : kFALSE;
62 TString gPathConfig = gSystem->Getenv("VMCWORKDIR");
63 // PLV : begin
64 // fGemConfigFile and fSilConfigFile initialization moved from ctor
65 // original code:
66 // , fGemConfigFile("GemRun8.xml")
67 // , fSilConfigFile("SiliconRun8.xml")
68 // modified code:
69 TString fGemConfigFile("GemRun8.xml"), fSilConfigFile("SiliconRun8.xml");
70 // BmnGemStripStationSet and BmnSiliconStationSet moved from class fields to local variables
71 // original code:
72 // GemStationSet = new BmnGemStripStationSet(gPathConfig + "/parameters/gem/XMLConfigs/" + fGemConfigFile);
73 // SilStationSet = new BmnSiliconStationSet(gPathConfig + "/parameters/silicon/XMLConfigs/" + fSilConfigFile);
74 // modified code:
75 BmnGemStripStationSet* GemStationSet =
76 new BmnGemStripStationSet(gPathConfig + "/parameters/gem/XMLConfigs/" + fGemConfigFile);
77 BmnSiliconStationSet* SilStationSet =
78 new BmnSiliconStationSet(gPathConfig + "/parameters/silicon/XMLConfigs/" + fSilConfigFile);
79 // PLV : end
80
81 static int irun = -1;
82 // PLV : begin
83 // zStat is not used in calculations
84 // static Double_t zStat[19] = {0};
85 // PLV : end
86
87 if (irun <= 0)
88 irun = ((DstEventHeader*)FairRootManager::Instance()->GetObject("DstEventHeader."))->GetRunId(); // AZ-010723
89
90 // Lorentz shifts
91 // AZ-010723 if (hLorCor[0][0] == nullptr) {
92 if (irun <= 0)
93 // PLV : begin
94 // all the section moved from method ApplyAlignment() - return type is different
95 // original code:
96 // return;
97 // modified code:
98 return kERROR;
99 // PLV : end
100 else if (hLorCorX[0][0] == nullptr) {
101 // Station Z-positions
102 int nsta = SilStationSet->GetNStations();
103 // PLV : begin
104 // zStat is not used in calculations
105 // for (int ista = 0; ista < nsta; ++ista)
106 // zStat[ista] = SilStationSet->GetZStationPosition(ista);
107 // PLV : end
108 int nsta1 = GemStationSet->GetNStations();
109 // PLV : begin
110 // zStat is not used in calculations
111 // for (int ista = 0; ista < nsta1; ++ista)
112 // zStat[ista + nsta] = GemStationSet->GetZStationPosition(ista);
113 // PLV : end
114 // PLV : begin
115 // non-necessary console output removed
116 // if (fVerbose > 0)
117 // for (int j = 0; j < nsta + nsta1; ++j)
118 // cout << " Stat: " << j << " " << zStat[j] << endl;
119 // PLV : end
120
121 map<int, vector<Double_t>> lorCorX[2], lorCorY[2];
122 int nx, ny;
123
124 for (int igem = 0; igem < 2; ++igem) {
125 ReadCorrections(irun, igem, lorCorX, lorCorY, nx, ny);
126 Double_t dx = 160.0 / nx, dy = 45.0 / ny;
127 // Double_t dx = 20.0, dy = 15.0;
128
129 for (int iupd = 0; iupd < 2; ++iupd) {
130 int istb = (igem == 0) ? 0 : nsta;
131 int iste = (igem == 0) ? nsta : nsta + nsta1;
132 if (iste > (int)lorCorX[0].size())
133 --iste; // exclude last small GEM
134
135 for (int ist = istb; ist < iste; ++ist) {
136 TString hName = "hLorCorX_";
137 hName += iupd;
138 hName += "_";
139 hName += ist;
140 // hLorCor[iupd][ist] = new TH2D(hName, "", 8, -80, 80, 3, 0, 45);
141 int nxfsd = 0;
142 if (igem)
143 hLorCorX[iupd][ist] = new TH2D(hName, "", nx, -80, 80, ny, 0, 45);
144 else {
145 nxfsd = lorCorX[0][ist + 1].size();
146 hLorCorX[iupd][ist] = new TH2D(hName, "", nxfsd, 0, nxfsd, ny, 0, ny);
147 }
148
149 hName = "hLorCorY_";
150 hName += iupd;
151 hName += "_";
152 hName += ist;
153 if (igem)
154 hLorCorY[iupd][ist] = new TH2D(hName, "", nx, -80, 80, ny, 0, 45);
155 else
156 hLorCorY[iupd][ist] = new TH2D(hName, "", nxfsd, 0, nxfsd, ny, 0, ny);
157
158 if (igem) {
159 // for (int j = 0; j < nx*ny; ++j) {
160 for (int j = 0; j < 24; ++j) {
161 // int ireg1 = j - 13;
162 // if (ireg1 < 0) ++ireg1;
163 int ireg1 = j - 12;
164 int ireg = TMath::Abs(ireg1);
165 if (ireg1 < 0)
166 --ireg;
167 int iy2 = ireg / 4;
168 int ix2 = ireg % 4;
169 Double_t x = (ix2 + 0.5) * dx;
170 if (ireg1 < 0)
171 x = -x;
172 Double_t y = (iy2 + 0.5) * dy;
173 hLorCorX[iupd][ist]->Fill(x, y, lorCorX[iupd][ist + 1][j]); // X-corrections
174 if (lorCorY[iupd].find(ist + 1) != lorCorY[iupd].end())
175 hLorCorY[iupd][ist]->Fill(x, y, lorCorY[iupd][ist + 1][j]); // Y-corrections
176 } // for (int j = 0;
177 } else {
178 for (int j = 0; j < nxfsd; ++j) {
179 hLorCorX[iupd][ist]->Fill(j + 0.5, 0.5, lorCorX[iupd][ist + 1][j]); // X-corrections
180 hLorCorY[iupd][ist]->Fill(j + 0.5, 0.5, lorCorY[iupd][ist + 1][j]); // Y-corrections
181 }
182 } // if (igem)
183 } // for (int ist = istb;
184 } // for (int iupd = 0;
185 } // for (int igem = 0;
186 } // else if (hLorCorX[0][0] == nullptr)
187
188 /*
189 * End of Lorentz corrections section ---------------------------------------------------------
190 */
191 return kSUCCESS;
192}
193
194void BmnAlignerApply::Exec(Option_t* opt)
195{
196 for (Int_t i = 0, n = fpStsHitsArray->GetEntriesFast(); i < n; i++) {
197 CbmStsHit* hit = static_cast<CbmStsHit*>(fpStsHitsArray->At(i));
198 SVectGL& align = fpResultCurrent->A(fpDetModel->IdxFromHit(hit->GetDetectorID()));
199 // keep alignment application code inside BmnMeasureModel
200 // to prevent inconsistency in sign selection and rotation sequence
201 fpMeasureModel->AlignHit(*hit, align);
202 if (!fLorentzCorrectionApply)
203 continue;
204
205 /*
206 * Lorentz corrections section-----------------------------------------------------------------
207 *
208 * The code below is directly taken from BmnToCbmHitConverter::ApplyAlignment(CbmStsHit*).
209 * It implements the Lorentz correction of hit positions and was not revised, except few
210 * lines marked with tag PLV. Indentation is the kept as in the original code.
211 */
212
213 int ista = hit->GetStationNr() - 1;
214 int isec = hit->GetSectorNr() - 1;
215 if (fieldFlag && hit->GetSystemId() == kGEM) {
216 // Lorentz corrections
217 int iupd = isec / 4;
218 Double_t xh = hit->GetX(), yh = TMath::Abs(hit->GetY()), xh1 = xh, yh1 = yh;
219 TH2D* h2 = hLorCorX[iupd][ista];
220 if (h2 == nullptr)
221 // PLV : begin
222 // in case of non-initialized h2 we must proceed to the next hit
223 // original code:
224 // return;
225 // modified code:
226 continue;
227 // PLV : end
228 if (xh1 >= h2->GetXaxis()->GetXmax())
229 xh1 = h2->GetXaxis()->GetXmax() - 0.1;
230 else if (xh1 < h2->GetXaxis()->GetXmin())
231 xh1 = h2->GetXaxis()->GetXmin() + 0.1;
232 if (yh1 >= h2->GetYaxis()->GetXmax())
233 yh1 = h2->GetYaxis()->GetXmax() - 0.1;
234 else if (yh1 < h2->GetYaxis()->GetXmin())
235 yh1 = h2->GetYaxis()->GetXmin() + 0.1;
236 // if (yh < h2->GetYaxis()->GetXmin() || yh >= h2->GetYaxis()->GetXmax()) cout << " yy " << isec << " " <<
237 // yh << endl;
238 Double_t dx = h2->Interpolate(xh1, yh1);
239 // dx = 0; //AZ-060723 - no corrections
240 hit->SetX(xh - dx * 0.65);
241 h2 = hLorCorY[iupd][ista];
242 Double_t dy = h2->Interpolate(xh1, yh1);
243 if (iupd)
244 dy = -dy;
245 // cout << iupd << " " << ista << " " << dy << endl;
246 Double_t ynew = yh - dy * 0.65;
247 if (iupd)
248 ynew = -ynew;
249 hit->SetY(ynew);
250 } else if (hit->GetSystemId() != kGEM) {
251 int iupd = isec / hLorCorX[0][ista]->GetNbinsX();
252 int ireg = isec % hLorCorX[0][ista]->GetNbinsX();
253 TH2D* h2 = hLorCorX[iupd][ista];
254 Double_t dx = h2->GetBinContent(ireg + 1, 1);
255 hit->SetX(hit->GetX() - dx * 0.65);
256 h2 = hLorCorY[iupd][ista];
257 Double_t dy = h2->GetBinContent(ireg + 1, 1);
258 hit->SetY(hit->GetY() - dy * 0.65);
259 }
260
261 /*
262 * End of Lorentz corrections section ---------------------------------------------------------
263 */
264
265 } // end of loop over hits
266}
267
269{
270 // we could do hLorCorX and hLorCorY cleanup here but
271 // it would interfere with ROOT logic of working with histograms
272}
273
274/*
275 * Lorentz corrections section-----------------------------------------------------------------
276 *
277 * The code below is directly taken from BmnToCbmHitConverter::ApplyAlignment(CbmStsHit*).
278 * It implements the Lorentz correction of hit positions and was not revised, except few
279 * lines marked with tag PLV. Indentation is the kept as in the original code.
280 */
281
282// PLV : begin
283// method ReadCorrections() header modified for BmnAlignerApply class
284// original code:
285// void BmnToCbmHitConverter::ReadCorrections(int irun,
286// int igem,
287// std::map<int, std::vector<Double_t>>* lorCorX,
288// std::map<int, std::vector<Double_t>>* lorCorY,
289// int& nx,
290// int& ny)
291// modified code:
292void BmnAlignerApply::ReadCorrections(int irun,
293 int igem,
294 std::map<int, std::vector<Double_t>>* lorCorX,
295 std::map<int, std::vector<Double_t>>* lorCorY,
296 int& nx,
297 int& ny)
298// PLV : end
299{
300 // Read file with "effective" Lorentz-corrections
301
302 ifstream fin;
303 TString fileName, chline;
304
305 fileName = gSystem->ExpandPathName("$VMCWORKDIR");
306 // AZ-030825 fileName += "/input/";
307 fileName += "/parameters/reco/"; // AZ-030825
308 if (igem)
309 fileName += "LorentzGemCorr.txt";
310 else
311 fileName += "LorentzFsdCorr.txt";
312 fin.open(fileName.Data());
313 if (!fin.is_open()) { // AZ-040825
314 cout << " !!! Lorentz correction file " << fileName << " not found !!!" << endl;
315 exit(1);
316 }
317
318 map<int, std::size_t> runsMap;
319 TObjArray* tokens = nullptr;
320
321 while (1) {
322 // while (!fin.eof()) {
323 size_t pos = fin.tellg();
324 chline.ReadLine(fin);
325 if (fin.eof())
326 break;
327 if (chline.Contains("Run")) {
328 // new run
329 // cout << chline << endl;
330 tokens = chline.Tokenize(" ");
331 TString run = ((TObjString*)tokens->UncheckedAt(1))->String();
332 run.Atoi();
333 // runsMap[run.Atoi()] = fin.tellg();
334 runsMap[run.Atoi()] = pos;
335 // cout << run.Atoi() << " " << pos << endl;
336 }
337 if (tokens) {
338 tokens->Delete();
339 delete tokens;
340 tokens = nullptr;
341 }
342 }
343 // if (chline[0] == '*') continue; // comment
344 // cout << runsMap.size() << endl;
345
346 auto iter = runsMap.upper_bound(irun);
347 if (iter != runsMap.begin())
348 --iter;
349 // PLV : begin
350 // non-necessary console output removed
351 // cout << " ---> Using corrections for run " << iter->first << " " << iter->second << " " << igem << endl;
352 // PLV : end
353
354 fin.clear();
355 fin.seekg(iter->second);
356 chline.ReadLine(fin);
357 // cout << chline << endl;
358
359 // Read detectors
360 while (1) {
361 chline.ReadLine(fin);
362 if (fin.eof())
363 break;
364 if (chline.Contains("Station")) {
365 // cout << chline << endl;
366 tokens = chline.Tokenize(" ");
367 TString stat = ((TObjString*)tokens->UncheckedAt(1))->String();
368 int ista = stat.Atoi();
369 int iupd = 0; // half-station
370 if (chline.Contains("lower"))
371 iupd = 1;
372 nx = ((TObjString*)tokens->UncheckedAt(3))->String().Atoi();
373 ny = ((TObjString*)tokens->UncheckedAt(4))->String().Atoi();
374 int ixdata = 0;
375 if (chline.Contains("X-"))
376 ixdata = 1;
377 // cout << ista << " " << iupd << " " << nx << " " << ny << endl;
378 if (nx * ny < 1)
379 cout << " !!! Wrong input line !!! " << chline << " " << nx << " " << ny << endl;
380 if (tokens) {
381 tokens->Delete();
382 delete tokens;
383 tokens = nullptr;
384 }
385 chline.ReadLine(fin);
386 tokens = chline.Tokenize(", ");
387 int nreg = tokens->GetEntriesFast();
388 if (nreg != nx * ny) {
389 cout << " !!! Wrong input line !!! " << chline << " " << nreg << " " << nx * ny << endl;
390 exit(0);
391 }
392 auto lorCor = vector<Double_t>(nreg);
393 // cout << tokens->GetEntriesFast() << " " << ((TObjString*)tokens->UncheckedAt(0))->String() << endl;
394 for (int itok = 0; itok < nreg; ++itok) {
395 Double_t corr = ((TObjString*)tokens->UncheckedAt(itok))->String().Atof();
396 // lorCor[iupd][ista][itok] = corr;
397 lorCor[itok] = corr;
398 }
399 if (ixdata)
400 lorCorX[iupd][ista] = lorCor;
401 else
402 lorCorY[iupd][ista] = lorCor;
403 if (tokens) {
404 tokens->Delete();
405 delete tokens;
406 tokens = nullptr;
407 }
408 } else if (chline.Contains("*")) {
409 continue; // comment
410 } else
411 break;
412 } // while (1)
413 fin.close();
414}
415
416/*
417 * End of Lorentz corrections section ---------------------------------------------------------
418 */
ROOT::Math::SVector< Double_t, BMN_GLOBAL_PARAMS_PD > SVectGL
int i
Definition P4_F32vec4.h:22
@ kGEM
SVectGL & A(Int_t index) noexcept
virtual void Exec(Option_t *opt)
virtual void Finish()
virtual InitStatus Init()
BmnAlignResult * fpResultCurrent
Definition BmnAligner.h:122
Bool_t LoadSolution(const char *jsonPath, Bool_t withAlpha=kTRUE, Bool_t withMSE=kTRUE)
const BmnDetectorModel * fpDetModel
Definition BmnAligner.h:109
const BmnMeasureModel * fpMeasureModel
Definition BmnAligner.h:108
Int_t IdxFromHit(Int_t hitID) const
void AlignHit(FairHit &hit, SVectGL &alignValues) const
virtual Int_t GetStationNr() const
Definition CbmStsHit.h:66
Int_t GetSectorNr() const
Definition CbmStsHit.h:68
Int_t GetSystemId() const
Definition CbmStsHit.h:64
-clang-format