CbmRoot
Loading...
Searching...
No Matches
CbmFieldMapDistorted.cxx
Go to the documentation of this file.
1/* Copyright (C) 2008-2021 Justus-Liebig-Universitaet Giessen, Giessen
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Elena Litvinenko, Elena Lebedeva [committer], Florian Uhlig */
4
5// -------------------------------------------------------------------------
6// ----- CbmFieldMapDistorted source file -----
7// ----- Created 18/01/08 by E.Litvinenko -----
8// -------------------------------------------------------------------------
10
11#include "CbmFieldConst.h" // for CbmFieldConst
12#include "CbmFieldMapSym1.h" // for CbmFieldMapSym1
13#include "CbmFieldMapSym2.h" // for CbmFieldMapSym2
14#include "CbmFieldMapSym3.h" // for CbmFieldMapSym3
15#include "CbmFieldPar.h" // for CbmFieldPar, kTypeDistorted
16
17#include <FairField.h> // for FairField
18
19#include <TFile.h> // for TFile, gFile
20#include <TFormula.h> // for TFormula
21#include <TH1.h> // for TH1D
22#include <TSystem.h> // for TSystem, gSystem, kFileExists
23
24#include <iostream> // for operator<<, endl, ostream, basic_ostream
25
26#include <string.h> // for strlen
27
28using std::cerr;
29using std::cout;
30using std::endl;
31
32// ------------- Default constructor ----------------------------------
34 : CbmFieldMap()
35 , fParentField(nullptr)
36 , fTypeOfParent(1)
37 , fDistortionFilename("")
38 , fBxDistortionFormulaMult(nullptr)
39 , fBxDistortionFormulaAdd(nullptr)
40 , fByDistortionFormulaMult(nullptr)
41 , fByDistortionFormulaAdd(nullptr)
42 , fBzDistortionFormulaMult(nullptr)
43 , fBzDistortionFormulaAdd(nullptr)
44{
45 fType = kTypeDistorted;
46}
47// ------------------------------------------------------------------------
48
49
50// ------------- Standard constructor (with FieldMap Parent Field ) ---------------------------------
51CbmFieldMapDistorted::CbmFieldMapDistorted(const char* mapName, const char* pfDistortionFilename,
52 const char* parentName, const char* fileTypeParent, Int_t pfTypeOfParent)
53 : CbmFieldMap()
54 , fParentField(nullptr)
55 , fTypeOfParent(pfTypeOfParent)
56 , fDistortionFilename(pfDistortionFilename)
57 , fBxDistortionFormulaMult(nullptr)
58 , fBxDistortionFormulaAdd(nullptr)
59 , fByDistortionFormulaMult(nullptr)
60 , fByDistortionFormulaAdd(nullptr)
61 , fBzDistortionFormulaMult(nullptr)
62 , fBzDistortionFormulaAdd(nullptr)
63{
64 fName = mapName;
65 fType = kTypeDistorted;
66
67 switch (pfTypeOfParent) {
68 case 3: fParentField = new CbmFieldMapSym3(parentName, fileTypeParent); break;
69 case 2: fParentField = new CbmFieldMapSym2(parentName, fileTypeParent); break;
70 case 5: fParentField = new CbmFieldMapSym1(parentName, fileTypeParent); break;
71 default: fParentField = new CbmFieldMap(parentName, fileTypeParent); break;
72 }
73}
74// ------------------------------------------------------------------------
75
76// ------------- Constructor (with Constant Parent Field ) -------------------------
77CbmFieldMapDistorted::CbmFieldMapDistorted(Double_t xMin, Double_t xMax, Double_t yMin, Double_t yMax, Double_t zMin,
78 Double_t zMax, Double_t bX, Double_t bY, Double_t bZ, const char* mapName,
79 const char* pfDistortionFilename, const char* parentname)
80 : CbmFieldMap()
81 , fParentField(new CbmFieldConst(parentname, xMin, xMax, yMin, yMax, zMin, zMax, bX, bY, bZ))
82 , fTypeOfParent(0)
83 , fDistortionFilename(pfDistortionFilename)
84 , fBxDistortionFormulaMult(nullptr)
85 , fBxDistortionFormulaAdd(nullptr)
86 , fByDistortionFormulaMult(nullptr)
87 , fByDistortionFormulaAdd(nullptr)
88 , fBzDistortionFormulaMult(nullptr)
89 , fBzDistortionFormulaAdd(nullptr)
90{
91 fName = mapName;
92 fType = kTypeDistorted;
93}
94
95
96// ------------- Constructor from CbmFieldPar -------------------------
98 : CbmFieldMap()
99 , fParentField(nullptr)
100 , fTypeOfParent(0)
101 , fDistortionFilename("")
102 , fBxDistortionFormulaMult(nullptr)
103 , fBxDistortionFormulaAdd(nullptr)
104 , fByDistortionFormulaMult(nullptr)
105 , fByDistortionFormulaAdd(nullptr)
106 , fBzDistortionFormulaMult(nullptr)
107 , fBzDistortionFormulaAdd(nullptr)
108{
109 if (fieldPar) {
110 fieldPar->MapName(fName);
111 fType = fieldPar->GetType();
112
113 if (fType == kTypeDistorted) { // normal case of distorted field map parameters as input
114
115 fTypeOfParent = fieldPar->GetTypeOfParent();
116
117 TString parentName;
118 fieldPar->GetParentName(parentName);
119
120 switch (fTypeOfParent) {
121 case 3: fParentField = new CbmFieldMapSym3(parentName.Data(), "R"); break;
122 case 2: fParentField = new CbmFieldMapSym2(parentName.Data(), "R"); break;
123 case 5: fParentField = new CbmFieldMapSym1(parentName.Data(), "R"); break;
124 case 0:
126 parentName.Data(), fieldPar->GetXmin(), fieldPar->GetXmax(), fieldPar->GetYmin(), fieldPar->GetYmax(),
127 fieldPar->GetZmin(), fieldPar->GetZmax(), fieldPar->GetBx(), fieldPar->GetBy(), fieldPar->GetBz());
128 break;
129 default: fParentField = new CbmFieldMap(parentName.Data(), "R"); break;
130 }
131
132 if ((fParentField) && (fTypeOfParent)) { // for field map parent
134 ->SetPosition(fieldPar->GetPositionX(), fieldPar->GetPositionY(), fieldPar->GetPositionZ());
135 ((CbmFieldMap*) fParentField)->SetScale(fieldPar->GetScale());
136 }
137
139 }
140 else { // try to create a distorted field map from the parameters of constant field or normal field map as input
141 switch (fType) {
142 case 3: fParentField = new CbmFieldMapSym3(fieldPar); break;
143 case 2: fParentField = new CbmFieldMapSym2(fieldPar); break;
144 case 5: fParentField = new CbmFieldMapSym1(fieldPar); break;
145 case 0: fParentField = new CbmFieldConst(fieldPar); break;
146 default: fParentField = new CbmFieldMap(fieldPar); break;
147 }
148 fType = kTypeDistorted;
150 }
151 }
152}
153// ------------------------------------------------------------------------
154
155
156// ------------ Destructor --------------------------------------------
158// ------------------------------------------------------------------------
159
160
161// ----------- Intialisation ------------------------------------------
163{
164 fParentField->Init();
165 fParentField->Print("");
168
169 // Fill values needed in the Print() function. This is needed to allow
170 // a constant Print() function.
171 Double_t po[3], b[3];
172 po[0] = 0;
173 po[1] = 0;
174 po[2] = 0;
175 GetFieldValue(po, b);
176 fBxOrigin = b[0];
177 fByOrigin = b[1];
178 fBzOrigin = b[2];
179}
180// ------------------------------------------------------------------------
181
182
183// ------------- Set from parent CbmField -------------------------
185{
186 fTypeOfParent = 0;
187 fParentField = field;
188
189 if (field) {
190 fTypeOfParent = field->GetType();
191
192 CbmFieldConst* field0 = 0;
193 CbmFieldMap* field1 = 0;
194
195 switch (fTypeOfParent) {
196
197 case 0:
198 field0 = (CbmFieldConst*) field;
199 fXmin = field0->GetXmin();
200 fXmax = field0->GetXmax();
201 fYmin = field0->GetYmin();
202 fZmax = field0->GetYmax();
203 fZmin = field0->GetZmin();
204 fZmax = field0->GetZmax();
205 break;
206
207 default:
208 field1 = (CbmFieldMap*) field;
209 fXmin = field1->GetXmin();
210 fXmax = field1->GetXmax();
211 fYmin = field1->GetYmin();
212 fZmax = field1->GetYmax();
213 fZmin = field1->GetZmin();
214 fZmax = field1->GetZmax();
215
216 fNx = field1->GetNx();
217 fNy = field1->GetNy();
218 fNz = field1->GetNz();
219
220 fXstep = field1->GetXstep();
221 fYstep = field1->GetYstep();
222 fZstep = field1->GetZstep();
223
224 fScale = field1->GetScale();
225 fPosX = field1->GetPositionX();
226 fPosY = field1->GetPositionY();
227 fPosZ = field1->GetPositionZ();
228 break;
229 }
230 }
231}
232// ------------------------------------------------------------------------
233
234// ----------- Read Distortion Formulas from Distortion File ------------------------------------------
236{
238 TFile* oldFile = gFile;
239 TDirectory* oldDir = gDirectory;
240
241 if (filename) {
242 if (strlen(filename)) { fDistortionFilename = filename; }
243 }
244 if (fDistortionFilename.Data()) {
245 if (strlen(fDistortionFilename.Data())) {
246 if (gSystem->AccessPathName(fDistortionFilename + ".root", kFileExists)) {
247 cerr << "CbmFieldMapDistorted::ReadDistortionInformation Warning: file " << (fDistortionFilename.Data())
248 << " not exists yet !!!" << endl;
249 }
250 else {
251 TFile* f = new TFile(fDistortionFilename + ".root");
252 if (f) {
253 fBxDistortionFormulaMult = f->Get<TFormula>("BxDistortionFormulaMult");
254 LOG_IF(fatal, !fBxDistortionFormulaMult) << "fBxDistortionFormulaMult not found in file " << f->GetName();
255 fBxDistortionFormulaAdd = f->Get<TFormula>("BxDistortionFormulaAdd");
256 LOG_IF(fatal, !fBxDistortionFormulaAdd) << "fBxDistortionFormulaAdd not found in file " << f->GetName();
257 fByDistortionFormulaMult = f->Get<TFormula>("ByDistortionFormulaMult");
258 LOG_IF(fatal, !fByDistortionFormulaMult) << "fByDistortionFormulaMult not found in file " << f->GetName();
259 fByDistortionFormulaAdd = f->Get<TFormula>("ByDistortionFormulaAdd");
260 LOG_IF(fatal, !fByDistortionFormulaAdd) << "fByDistortionFormulaAdd not found in file " << f->GetName();
261 fBzDistortionFormulaMult = f->Get<TFormula>("BzDistortionFormulaMult");
262 LOG_IF(fatal, !fBzDistortionFormulaMult) << "fBzDistortionFormulaMult not found in file " << f->GetName();
263 fBzDistortionFormulaAdd = f->Get<TFormula>("BzDistortionFormulaAdd");
264 LOG_IF(fatal, !fBzDistortionFormulaAdd) << "fBzDistortionFormulaAdd not found in file " << f->GetName();
265 f->Close();
266 }
267 else {
268 cerr << "CbmFieldMapDistorted::ReadDistortionInformation ERROR: file " << (fDistortionFilename.Data())
269 << " can not be read !!!" << endl;
270 }
271 }
272 }
273 }
275 gFile = oldFile;
276 gDirectory = oldDir;
277}
278
279// ----------- Write Distortion Formulas to Distortion File ------------------------------------------
281{
283 TFile* oldFile = gFile;
284 TDirectory* oldDir = gDirectory;
285
286 if (filename) {
287 if (strlen(filename)) { fDistortionFilename = filename; }
288 }
289
290 if (fDistortionFilename.Data()) {
291 if (strlen(fDistortionFilename.Data())) {
292 TFile* f = new TFile(fDistortionFilename + ".root", "RECREATE");
293 if (f) {
294
295 if (fBxDistortionFormulaMult) fBxDistortionFormulaMult->Write("BxDistortionFormulaMult");
296 if (fBxDistortionFormulaAdd) fBxDistortionFormulaAdd->Write("BxDistortionFormulaAdd");
297
298 if (fByDistortionFormulaMult) fByDistortionFormulaMult->Write("ByDistortionFormulaMult");
299 if (fByDistortionFormulaAdd) fByDistortionFormulaAdd->Write("ByDistortionFormulaAdd");
300
301 if (fBzDistortionFormulaMult) fBzDistortionFormulaMult->Write("BzDistortionFormulaMult");
302 if (fBzDistortionFormulaAdd) fBzDistortionFormulaAdd->Write("BzDistortionFormulaAdd");
303
304 f->Write();
305 f->Close();
306 }
307 }
308 }
310 gFile = oldFile;
311 gDirectory = oldDir;
312}
313
314// ---------------Getter and Setter for Distortion Formulas------------------------------
315TFormula* CbmFieldMapDistorted::GetDistortionFormula(const char* component_option, const char* action_option)
316{
317 // component_option: "x","y","z"; action_option: "m","a"
318
319 TString co = component_option;
320 TString ao = action_option;
321
322 if (co.Contains("y") && ao.Contains("m")) return fByDistortionFormulaMult;
323 if (co.Contains("y") && ao.Contains("a")) return fByDistortionFormulaAdd;
324
325 if (co.Contains("x") && ao.Contains("m")) return fBxDistortionFormulaMult;
326 if (co.Contains("x") && ao.Contains("a")) return fBxDistortionFormulaAdd;
327
328 if (co.Contains("z") && ao.Contains("m")) return fBzDistortionFormulaMult;
329 if (co.Contains("z") && ao.Contains("a")) return fBzDistortionFormulaAdd;
330
331 return 0;
332}
333
334Bool_t CbmFieldMapDistorted::SetDistortionFormula(TFormula* parDistortionFormula, const char* component_option,
335 const char* action_option)
336{
337 TString co = component_option;
338 TString ao = action_option;
339 Int_t counter = 0;
340
341 if (co.Contains("y") && ao.Contains("m")) {
342 fByDistortionFormulaMult = parDistortionFormula;
343 counter++;
344 }
345 if (co.Contains("y") && ao.Contains("a")) {
346 fByDistortionFormulaAdd = parDistortionFormula;
347 counter++;
348 }
349
350 if (co.Contains("x") && ao.Contains("m")) {
351 fBxDistortionFormulaMult = parDistortionFormula;
352 counter++;
353 }
354 if (co.Contains("x") && ao.Contains("a")) {
355 fBxDistortionFormulaAdd = parDistortionFormula;
356 counter++;
357 }
358
359 if (co.Contains("z") && ao.Contains("m")) {
360 fBzDistortionFormulaMult = parDistortionFormula;
361 counter++;
362 }
363 if (co.Contains("z") && ao.Contains("a")) {
364 fBzDistortionFormulaAdd = parDistortionFormula;
365 counter++;
366 }
367 return (counter > 0);
368}
369
370Bool_t CbmFieldMapDistorted::SetDistortionFormula(const char* parDistortionFormulaText, const char* component_option,
371 const char* action_option)
372{
373 TString co = component_option;
374 TString ao = action_option;
375 Int_t counter = 0;
376
377 if (co.Contains("y") && ao.Contains("m")) {
378 fByDistortionFormulaMult = new TFormula("ByDistortionFormulaMult", parDistortionFormulaText);
379 counter++;
380 }
381 if (co.Contains("y") && ao.Contains("a")) {
382 fByDistortionFormulaAdd = new TFormula("ByDistortionFormulaAdd", parDistortionFormulaText);
383 counter++;
384 }
385
386 if (co.Contains("x") && ao.Contains("m")) {
387 fBxDistortionFormulaMult = new TFormula("BxDistortionFormulaMult", parDistortionFormulaText);
388 counter++;
389 }
390 if (co.Contains("x") && ao.Contains("a")) {
391 fBxDistortionFormulaAdd = new TFormula("BxDistortionFormulaAdd", parDistortionFormulaText);
392 counter++;
393 }
394
395 if (co.Contains("z") && ao.Contains("m")) {
396 fBzDistortionFormulaMult = new TFormula("BzDistortionFormulaMult", parDistortionFormulaText);
397 counter++;
398 }
399 if (co.Contains("z") && ao.Contains("a")) {
400 fBzDistortionFormulaAdd = new TFormula("BzDistortionFormulaAdd", parDistortionFormulaText);
401 counter++;
402 }
403 return (counter > 0);
404}
405// ------------------------------------------------------------------------
406
407// ----------- Get x component of the field ---------------------------
408Double_t CbmFieldMapDistorted::GetBx(Double_t x, Double_t y, Double_t z)
409{
410 Double_t bx = fParentField->GetBx(x, y, z);
411
412 if (fBxDistortionFormulaMult) { bx *= fBxDistortionFormulaMult->Eval(x, y, z); }
413
414 if (fBxDistortionFormulaAdd) { bx += fBxDistortionFormulaAdd->Eval(x, y, z); }
415
416 return bx;
417}
418
419// ----------- Get y component of the field ---------------------------
420Double_t CbmFieldMapDistorted::GetBy(Double_t x, Double_t y, Double_t z)
421{
422 Double_t by = fParentField->GetBy(x, y, z);
423
424 if (fByDistortionFormulaMult) { by *= fByDistortionFormulaMult->Eval(x, y, z); }
425 if (fByDistortionFormulaAdd) { by += fByDistortionFormulaAdd->Eval(x, y, z); }
426
427 return by;
428}
429
430// ----------- Get z component of the field ---------------------------
431Double_t CbmFieldMapDistorted::GetBz(Double_t x, Double_t y, Double_t z)
432{
433 Double_t bz = fParentField->GetBz(x, y, z);
434
435 if (fBzDistortionFormulaMult) { bz *= fBzDistortionFormulaMult->Eval(x, y, z); }
436 if (fBzDistortionFormulaAdd) { bz += fBzDistortionFormulaAdd->Eval(x, y, z); }
437
438 return bz;
439}
440
441void CbmFieldMapDistorted::Print(Option_t*) const
442{
443 cout << "=============================================================" << endl;
444 cout << "---- " << fTitle << " : " << fName << endl;
445
446 cout << "== Parent Field: ==" << endl;
447 fParentField->Print("");
448 cout << "==============================" << endl;
449
450 cout << "== Distortion Information File : ==" << endl;
451 cout << fDistortionFilename << endl;
452 cout << "==============================" << endl;
453
454 cout << "== Bx Distortion Formula Mult : ==" << endl;
455 if (fBxDistortionFormulaMult) cout << fBxDistortionFormulaMult->GetExpFormula("p") << endl;
456 cout << "== Bx Distortion Formula Add : ==" << endl;
457 if (fBxDistortionFormulaAdd) cout << fBxDistortionFormulaAdd->GetExpFormula("p") << endl;
458 cout << "==============================" << endl;
459
460 cout << "== By Distortion Formula Mult : ==" << endl;
461 if (fByDistortionFormulaMult) cout << fByDistortionFormulaMult->GetExpFormula("p") << endl;
462 cout << "== By Distortion Formula Add : ==" << endl;
463 if (fByDistortionFormulaAdd) cout << fByDistortionFormulaAdd->GetExpFormula("p") << endl;
464 cout << "==============================" << endl;
465
466 cout << "== Bz Distortion Formula Mult : ==" << endl;
467 if (fBzDistortionFormulaMult) cout << fBzDistortionFormulaMult->GetExpFormula("p") << endl;
468 cout << "== Bz Distortion Formula Add : ==" << endl;
469 if (fBzDistortionFormulaAdd) cout << fBzDistortionFormulaAdd->GetExpFormula("p") << endl;
470 cout << "==============================" << endl;
471
472
473 cout << "---- Distorted field at origin is ( " << fBxOrigin << ", " << fByOrigin << ", " << fBzOrigin << ") kG"
474 << endl;
475
476 cout << "=============================================================" << endl;
477}
478// ------------------------------------------------------------------------
479
480
481// ----- Set the position of the field centre in global coordinates -----
482void CbmFieldMapDistorted::SetPosition(Double_t x, Double_t y, Double_t z)
483{
484 fPosX = x;
485 fPosY = y;
486 fPosZ = z;
487 if (fTypeOfParent) ((CbmFieldMap*) fParentField)->SetPosition(x, y, z);
488}
489// ------------------------------------------------------------------------
490
491// ---------------------Set a global field scaling factor------------------
493{
494 fScale = factor;
495 if (fTypeOfParent) ((CbmFieldMap*) fParentField)->SetScale(factor);
496}
497// ------------------------------------------------------------------------
498
499// ---------------------Plot distorted and parend field (By component)-----
500void CbmFieldMapDistorted::PlotBy(Int_t n, Double_t zmin, Double_t zmax)
501{
502 TH1D *h, *hp;
503 Int_t n0 = n;
504 if (!n) n0 = 1;
505 Double_t dz = (zmax - zmin) / n0;
506
507 h = new TH1D("hField", fName, n, zmin, zmax);
508 for (Int_t i = 0; i < n0; i++)
509 h->SetBinContent(i + 1, GetBy(0., 0., zmin + i * dz));
510 h->SetLineWidth(2);
511 h->Draw();
512
513 if (fParentField) {
514 hp = new TH1D("hFieldParent", fParentField->GetName(), n, zmin, zmax);
515 for (Int_t i = 0; i < n0; i++)
516 hp->SetBinContent(i + 1, fParentField->GetBy(0., 0., zmin + i * dz));
517 hp->SetLineColor(2);
518 hp->Draw("sames");
519 }
520}
521// ------------------------------------------------------------------------
const int kTypeDistorted
Definition CbmFieldPar.h:33
Double_t GetXmax() const
Double_t GetZmax() const
Double_t GetZmin() const
Double_t GetYmin() const
Double_t GetXmin() const
Double_t GetYmax() const
TFormula * fByDistortionFormulaMult
getter/setter options: ("x","m"), ("x","a")
virtual void SetFromParent(FairField *field)
virtual void SetPosition(Double_t x, Double_t y, Double_t z)
virtual TFormula * GetDistortionFormula(const char *component_option="y", const char *action_option="m")
void ReadDistortionInformation(const char *filename=0)
void PlotBy(Int_t n=250, Double_t zmin=-50, Double_t zmax=450)
void WriteDistortionInformation(const char *filename=0)
virtual void SetScale(Double_t factor)
virtual Bool_t SetDistortionFormula(TFormula *parDistortionFormula, const char *component_option="y", const char *action_option="m")
TFormula * fBzDistortionFormulaMult
getter/setter options: ("y","m"), ("y","a")
virtual void Print(Option_t *="") const
Double_t fPosZ
Double_t fBxOrigin
Interpolated field (1-dim)
Int_t GetNy() const
Double_t GetZmax() const
Double_t fPosY
TArrayF * GetBz() const
Double_t GetXstep() const
Double_t fPosX
Double_t fZmin
Double_t fZmax
Double_t fYmin
Double_t GetYmin() const
TArrayF * GetBx() const
Double_t GetPositionY() const
Double_t fXmin
Double_t GetYmax() const
Double_t GetYstep() const
Double_t GetPositionZ() const
Double_t fYstep
Double_t fZstep
Double_t GetXmin() const
Int_t GetNz() const
Double_t GetPositionX() const
Int_t GetNx() const
Double_t fXstep
Double_t fBzOrigin
y-component of the field at the origin
Double_t GetXmax() const
Double_t GetZmin() const
Double_t GetZstep() const
Double_t fScale
Double_t fByOrigin
x-component of the field at the origin
Double_t GetScale() const
TArrayF * GetBy() const
Double_t fXmax
Double_t GetXmin() const
Definition CbmFieldPar.h:62
Double_t GetBy() const
Definition CbmFieldPar.h:69
Double_t GetXmax() const
Definition CbmFieldPar.h:63
Double_t GetPositionZ() const
Definition CbmFieldPar.h:74
Int_t GetType() const
Definition CbmFieldPar.h:61
Int_t GetTypeOfParent() const
Definition CbmFieldPar.h:79
Double_t GetBx() const
Definition CbmFieldPar.h:68
Double_t GetYmin() const
Definition CbmFieldPar.h:64
void GetParentName(TString &parentname)
Definition CbmFieldPar.h:78
Double_t GetZmax() const
Definition CbmFieldPar.h:67
void GetDistortionFilename(TString &filename)
Definition CbmFieldPar.h:77
void MapName(TString &name)
Definition CbmFieldPar.h:71
Double_t GetPositionY() const
Definition CbmFieldPar.h:73
Double_t GetZmin() const
Definition CbmFieldPar.h:66
Double_t GetBz() const
Definition CbmFieldPar.h:70
Double_t GetScale() const
Definition CbmFieldPar.h:75
Double_t GetYmax() const
Definition CbmFieldPar.h:65
Double_t GetPositionX() const
Definition CbmFieldPar.h:72
Data class with information on a STS local track.