BmnRoot
Loading...
Searching...
No Matches
BmnInnTrackerAlign.cxx
Go to the documentation of this file.
1// @(#)bmnroot/alignment:$Id$
2// Author: Pavel Batyuk <pavel.batyuk@jinr.ru> 2018-11-02
3
5// //
6// BmnInnTrackerAlign //
7// //
8// Interface to get align. corrections for inn. tracker //
9// //
10// //
12
13#include "BmnInnTrackerAlign.h"
14
15#include "TRandom.h"
16#include "UniDetectorParameter.h"
17
18BmnInnTrackerAlign::BmnInnTrackerAlign(Int_t period, Int_t run, TString fileName)
19{
20 fBranchGemCorrs = "BmnGemAlignCorrections";
21 fBranchSilCorrs = "BmnSiliconAlignCorrections";
22
23 // TString alCorrFileName = "";
24
25 if (fileName.Contains("default")) {
26 fFilename = Form("alignment_innTracker_%d.root", UInt_t(gRandom->Integer(UINT32_MAX)));
27 UniDetectorParameter::ReadFile("BM@N", "alignment", period, run, (Char_t*)fFilename.Data());
28 } else
29 fFilename = fileName;
30
31 // Choose appropriate configuration (BM@N or SRC)
32 Bool_t isSRC = (run < 3589) ? kTRUE : kFALSE;
33 TString gPathConfig = gSystem->Getenv("VMCWORKDIR");
34 TString confSi = "SiliconRun" + TString(isSRC ? "SRC" : "") + "Spring2018.xml";
35 TString confGem = "GemRun" + TString(isSRC ? "SRC" : "") + "Spring2018.xml";
36
38 TString gPathSiliconConfig = gPathConfig + "/parameters/silicon/XMLConfigs/";
39 fDetectorSI = new BmnSiliconStationSet(gPathSiliconConfig + confSi);
40
42 TString gPathGemConfig = gPathConfig + "/parameters/gem/XMLConfigs/";
43 fDetectorGEM = new BmnGemStripStationSet(gPathGemConfig + confGem);
44
45 // Get alignment corrections ...
46 TFile* f = nullptr;
47 for (Int_t iOpen = 0; iOpen < 2; iOpen++) {
48 f = new TFile(fFilename.Data());
49 if (iOpen == 0)
50 fCorrsGem = GetGemCorrs(f);
51 else
52 fCorrsSil = GetSiliconCorrs(f);
53 delete f;
54 }
55 remove(fFilename.Data());
56
57 // Get Lorentz corrections ...
58 const Int_t nStat = fDetectorGEM->GetNStations();
59 nStations = fDetectorGEM->GetNStations();
60 fLorCorrs = new Double_t*[nStat];
61 UniDetectorParameter* coeffLorCorrs =
62 UniDetectorParameter::GetDetectorParameter("GEM", "lorentz_shift", period, run);
63 vector<UniValue*> shifts;
64 if (coeffLorCorrs)
65 coeffLorCorrs->GetValue(shifts); //
66 for (Int_t iEle = 0; iEle < nStat; iEle++) {
67 const Int_t nParams = 3; // Parabolic approximation is used
68 fLorCorrs[iEle] = new Double_t[nParams];
69 for (Int_t iParam = 0; iParam < nParams; iParam++)
70 fLorCorrs[iEle][iParam] = (coeffLorCorrs) ? ((LorentzShiftValue*)shifts[iEle])->ls[iParam] : 0.;
71 }
72
73 for (auto shift : shifts)
74 delete shift;
75 delete coeffLorCorrs;
76}
77
79{
80 cout << "GEM corrections: " << endl;
81 for (Int_t iStat = 0; iStat < fDetectorGEM->GetNStations(); iStat++)
82 for (Int_t iMod = 0; iMod < fDetectorGEM->GetGemStation(iStat)->GetNModules(); iMod++) {
83 cout << "Stat# " << iStat << " Mod# " << iMod << " Corrs# ";
84 for (Int_t iCorr = 0; iCorr < 3; iCorr++)
85 cout << fCorrsGem[iStat][iMod][iCorr] << " ";
86 cout << endl;
87 }
88
89 cout << "Lorentz coefficients: " << endl;
90 for (Int_t iStat = 0; iStat < fDetectorGEM->GetNStations(); iStat++)
91 for (Int_t iCoeff = 0; iCoeff < 3; iCoeff++)
92 cout << fLorCorrs[iStat][iCoeff] << " ";
93 cout << endl;
94
95 cout << "SILICON corrections: " << endl;
96 for (Int_t iStat = 0; iStat < fDetectorSI->GetNStations(); iStat++)
97 for (Int_t iMod = 0; iMod < fDetectorSI->GetSiliconStation(iStat)->GetNModules(); iMod++) {
98 cout << "Stat# " << iStat << " Mod# " << iMod << " Corrs# ";
99 for (Int_t iCorr = 0; iCorr < 3; iCorr++)
100 cout << fCorrsSil[iStat][iMod][iCorr] << " ";
101 cout << endl;
102 }
103}
104
106{
107
108 for (Int_t iStat = 0; iStat < fDetectorGEM->GetNStations(); iStat++) {
109 for (Int_t iMod = 0; iMod < fDetectorGEM->GetGemStation(iStat)->GetNModules(); iMod++) {
110 delete[] fCorrsGem[iStat][iMod];
111 }
112 delete[] fCorrsGem[iStat];
113 delete[] fLorCorrs[iStat];
114 }
115 delete[] fCorrsGem;
116 delete[] fLorCorrs;
117
118 for (Int_t iStat = 0; iStat < fDetectorSI->GetNStations(); iStat++) {
119 for (Int_t iMod = 0; iMod < fDetectorSI->GetSiliconStation(iStat)->GetNModules(); iMod++) {
120 delete[] fCorrsSil[iStat][iMod];
121 }
122 delete[] fCorrsSil[iStat];
123 }
124 delete[] fCorrsSil;
125
126 delete fDetectorGEM;
127 delete fDetectorSI;
128}
129
130Double_t*** BmnInnTrackerAlign::GetGemCorrs(TFile* file)
131{
132 const Int_t nStat = fDetectorGEM->GetNStations();
133 const Int_t nParams = 3;
134 Double_t*** corr = new Double_t**[nStat];
135
136 for (Int_t iStat = 0; iStat < nStat; iStat++) {
137 Int_t nModul = fDetectorGEM->GetGemStation(iStat)->GetNModules();
138 corr[iStat] = new Double_t*[nModul];
139
140 for (Int_t iMod = 0; iMod < nModul; iMod++) {
141 corr[iStat][iMod] = new Double_t[nParams];
142
143 for (Int_t iPar = 0; iPar < nParams; iPar++)
144 corr[iStat][iMod][iPar] = 0.;
145 }
146 }
147
148 TTree* t = (TTree*)file->Get("cbmsim");
149 TClonesArray* corrs = nullptr;
150 t->SetBranchAddress(fBranchGemCorrs.Data(), &corrs);
151
152 for (Int_t iEntry = 0; iEntry < t->GetEntries(); iEntry++) {
153 t->GetEntry(iEntry);
154 for (Int_t iCorr = 0; iCorr < corrs->GetEntriesFast(); iCorr++) {
155 BmnGemAlignCorrections* align = (BmnGemAlignCorrections*)corrs->UncheckedAt(iCorr);
156 Int_t iStat = align->GetStation();
157 Int_t iMod = align->GetModule();
158 corr[iStat][iMod][0] = align->GetCorrections().X();
159 corr[iStat][iMod][1] = align->GetCorrections().Y();
160 corr[iStat][iMod][2] = align->GetCorrections().Z();
161 }
162 }
163 delete corrs;
164 t->Delete();
165 return corr;
166}
167
168Double_t*** BmnInnTrackerAlign::GetSiliconCorrs(TFile* file)
169{
170 const Int_t nStat = fDetectorSI->GetNStations();
171 const Int_t nParams = 3;
172 Double_t*** corr = new Double_t**[nStat];
173
174 for (Int_t iStat = 0; iStat < nStat; iStat++) {
175 Int_t nModul = fDetectorSI->GetSiliconStation(iStat)->GetNModules();
176 corr[iStat] = new Double_t*[nModul];
177
178 for (Int_t iMod = 0; iMod < nModul; iMod++) {
179 corr[iStat][iMod] = new Double_t[nParams];
180
181 for (Int_t iPar = 0; iPar < nParams; iPar++)
182 corr[iStat][iMod][iPar] = 0.;
183 }
184 }
185
186 TTree* t = (TTree*)file->Get("cbmsim");
187 TClonesArray* corrs = nullptr;
188 t->SetBranchAddress(fBranchSilCorrs.Data(), &corrs);
189
190 for (Int_t iEntry = 0; iEntry < t->GetEntries(); iEntry++) {
191 t->GetEntry(iEntry);
192 for (Int_t iCorr = 0; iCorr < corrs->GetEntriesFast(); iCorr++) {
193 BmnSiliconAlignCorrections* align = (BmnSiliconAlignCorrections*)corrs->UncheckedAt(iCorr);
194 Int_t iStat = align->GetStation();
195 Int_t iMod = align->GetModule();
196 corr[iStat][iMod][0] = align->GetCorrections().X();
197 corr[iStat][iMod][1] = align->GetCorrections().Y();
198 corr[iStat][iMod][2] = align->GetCorrections().Z();
199 }
200 }
201 delete corrs;
202 t->Delete();
203 return corr;
204}
float f
Definition P4_F32vec4.h:21
BmnGemStripStation * GetGemStation(Int_t station_num)
Double_t *** GetGemCorrs()
Double_t *** GetSiliconCorrs()
BmnSiliconStation * GetSiliconStation(Int_t station_num)
int GetValue(vector< UniValue * > &parameter_value)
get value of detector parameter presented by an array
static UniDetectorParameter * GetDetectorParameter(int value_id)
get detector parameter from the database
static int ReadFile(const char *detector_name, const char *parameter_name, int period_number, int run_number, const char *file_path)
-clang-format