CbmRoot
Loading...
Searching...
No Matches
CbmFieldMap.cxx
Go to the documentation of this file.
1/* Copyright (C) 2004-2021 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Mohammad Al-Turany, Denis Bertini [committer], Florian Uhlig, Volker Friese */
4
5// -------------------------------------------------------------------------
6// ----- CbmFieldMap source file -----
7// ----- Created 12/01/04 by M. Al/Turany (CbmField.cxx) -----
8// ----- Redesign 13/02/06 by V. Friese -----
9// -------------------------------------------------------------------------
10#include "CbmFieldMap.h"
11
12#include "CbmFieldMapCreator.h" // for CbmFieldMapCreator
13#include "CbmFieldMapData.h" // for CbmFieldMapData
14#include "CbmFieldPar.h" // for CbmFieldPar
15
16#include <FairField.h> // for FairField
17#include <Logger.h> // for LOG, Logger
18
19#include <TArrayF.h> // for TArrayF
20#include <TFile.h> // for TFile, gFile
21#include <TMath.h> // for Nint
22
23#include <iomanip> // for operator<<, setw
24#include <iostream> // for operator<<, basic_ostream, endl, ost...
25
26#include <assert.h> // for assert
27#include <stdlib.h> // for div, div_t, getenv, exit
28
29using std::cerr;
30using std::cout;
31using std::endl;
32using std::flush;
33using std::right;
34using std::setw;
35using std::showpoint;
36
37
38// ------------- Default constructor ----------------------------------
40 : FairField()
41 , fFileName("")
42 , fScale(1.)
43 , fPosX(0.)
44 , fPosY(0.)
45 , fPosZ(0.)
46 , fXmin(0.)
47 , fXmax(0.)
48 , fXstep(0.)
49 , fYmin(0.)
50 , fYmax(0.)
51 , fYstep(0.)
52 , fZmin(0.)
53 , fZmax(0.)
54 , fZstep(0.)
55 , fNx(0)
56 , fNy(0)
57 , fNz(0)
58 , fBx(nullptr)
59 , fBy(nullptr)
60 , fBz(nullptr)
61 , fBxOrigin(0.)
62 , fByOrigin(0.)
63 , fBzOrigin(0.)
64{
65 // Initilization of arrays is to my knowledge not
66 // possible in member initalization lists
67 for (Int_t i = 0; i < 2; i++) {
68 fHc[i] = 0;
69 for (Int_t j = 0; j < 2; j++) {
70 fHb[i][j] = 0;
71 for (Int_t k = 0; k < 2; k++) {
72 fHa[i][j][k] = 0;
73 }
74 }
75 }
76 // Assign values to data members of base classes
77 // There is no appropriate constructor of the base
78 // class.
79 fName = "";
80 fType = 1;
81}
82// ------------------------------------------------------------------------
83
84
85// ------------- Standard constructor ---------------------------------
86CbmFieldMap::CbmFieldMap(const char* mapName, const char* fileType)
87 : FairField()
88 , fFileName("")
89 , fScale(1.)
90 , fPosX(0.)
91 , fPosY(0.)
92 , fPosZ(0.)
93 , fXmin(0.)
94 , fXmax(0.)
95 , fXstep(0.)
96 , fYmin(0.)
97 , fYmax(0.)
98 , fYstep(0.)
99 , fZmin(0.)
100 , fZmax(0.)
101 , fZstep(0.)
102 , fNx(0)
103 , fNy(0)
104 , fNz(0)
105 , fBx(nullptr)
106 , fBy(nullptr)
107 , fBz(nullptr)
108 , fBxOrigin(0.)
109 , fByOrigin(0.)
110 , fBzOrigin(0.)
111{
112 // Initilization of arrays is to my knowledge not
113 // possible in member initalization lists
114 for (Int_t i = 0; i < 2; i++) {
115 fHc[i] = 0;
116 for (Int_t j = 0; j < 2; j++) {
117 fHb[i][j] = 0;
118 for (Int_t k = 0; k < 2; k++) {
119 fHa[i][j][k] = 0;
120 }
121 }
122 }
123 // Assign values to data members of base classes
124 // There is no appropriate constructor of the base
125 // class.
126 fName = mapName;
127 fTitle = "CbmFieldMap";
128 TString dir = getenv("VMCWORKDIR");
129 fFileName = dir + "/input/" + mapName;
130 if (fileType[0] == 'R') { fFileName += ".root"; }
131 else {
132 fFileName += ".dat";
133 }
134 fType = 1;
135 LOG(info) << "Filename is " << fFileName;
136}
137// ------------------------------------------------------------------------
138
139
140// ------------ Constructor from CbmFieldPar --------------------------
142 : FairField()
143 , fFileName("")
144 , fScale(1.)
145 , fPosX(0.)
146 , fPosY(0.)
147 , fPosZ(0.)
148 , fXmin(0.)
149 , fXmax(0.)
150 , fXstep(0.)
151 , fYmin(0.)
152 , fYmax(0.)
153 , fYstep(0.)
154 , fZmin(0.)
155 , fZmax(0.)
156 , fZstep(0.)
157 , fNx(0)
158 , fNy(0)
159 , fNz(0)
160 , fBx(nullptr)
161 , fBy(nullptr)
162 , fBz(nullptr)
163 , fBxOrigin(0.)
164 , fByOrigin(0.)
165 , fBzOrigin(0.)
166{
167 // Initilization of arrays is to my knowledge not
168 // possible in member initalization lists
169 for (Int_t i = 0; i < 2; i++) {
170 fHc[i] = 0;
171 for (Int_t j = 0; j < 2; j++) {
172 fHb[i][j] = 0;
173 for (Int_t k = 0; k < 2; k++) {
174 fHa[i][j][k] = 0;
175 }
176 }
177 }
178 // Assign values to data members of base classes
179 // There is no appropriate constructor of the base
180 // class.
181 fName = "";
182 fType = 1;
183 if (!fieldPar) { cerr << "-W- CbmFieldConst::CbmFieldMap: empty parameter container!" << endl; }
184 else {
185 fieldPar->MapName(fName);
186 fPosX = fieldPar->GetPositionX();
187 fPosY = fieldPar->GetPositionY();
188 fPosZ = fieldPar->GetPositionZ();
189 fScale = fieldPar->GetScale();
190 TString dir = getenv("VMCWORKDIR");
191 fFileName = dir + "/input/" + fName + ".root";
192 fType = fieldPar->GetType();
193 }
194}
195// ------------------------------------------------------------------------
196
197
198// ------------ Constructor from CbmFieldMapCreator ---------------------
200 : FairField()
201 , fFileName("")
202 , fScale(1.)
203 , fPosX(0.)
204 , fPosY(0.)
205 , fPosZ(0.)
206 , fXmin(0.)
207 , fXmax(0.)
208 , fXstep(0.)
209 , fYmin(0.)
210 , fYmax(0.)
211 , fYstep(0.)
212 , fZmin(0.)
213 , fZmax(0.)
214 , fZstep(0.)
215 , fNx(0)
216 , fNy(0)
217 , fNz(0)
218 , fBx(nullptr)
219 , fBy(nullptr)
220 , fBz(nullptr)
221 , fBxOrigin(0.)
222 , fByOrigin(0.)
223 , fBzOrigin(0.)
224{
225 // Initilization of arrays is to my knowledge not
226 // possible in member initalization lists
227 for (Int_t i = 0; i < 2; i++) {
228 fHc[i] = 0;
229 for (Int_t j = 0; j < 2; j++) {
230 fHb[i][j] = 0;
231 for (Int_t k = 0; k < 2; k++) {
232 fHa[i][j][k] = 0;
233 }
234 }
235 }
236 // Assign values to data members of base classes
237 // There is no appropriate constructor of the base
238 // class.
239 fName = "";
240 fType = 1;
241 if (!creator) { cerr << "-W- CbmFieldMap: no creator given!" << endl; }
242 else {
243 fType = 1;
244 fName = creator->GetMapName();
245 fXmin = creator->GetXmin();
246 fXmax = creator->GetXmax();
247 fYmin = creator->GetYmin();
248 fYmax = creator->GetYmax();
249 fZmin = creator->GetZmin();
250 fZmax = creator->GetZmax();
251 fNx = creator->GetNx();
252 fNy = creator->GetNy();
253 fNz = creator->GetNz();
254 fXstep = (fXmax - fXmin) / Double_t(fNx - 1);
255 fYstep = (fYmax - fYmin) / Double_t(fNy - 1);
256 fZstep = (fZmax - fZmin) / Double_t(fNz - 1);
257 fBx = creator->GetBx();
258 fBy = creator->GetBy();
259 fBz = creator->GetBz();
260 }
261}
262
263// ------------------------------------------------------------------------
264
265
266// ------------ Destructor --------------------------------------------
268{
269 if (fBx) delete fBx;
270 if (fBy) delete fBy;
271 if (fBz) delete fBz;
272}
273// ------------------------------------------------------------------------
274
275
276// ----------- Intialisation ------------------------------------------
278{
279 if (fFileName.EndsWith(".root")) ReadRootFile(fFileName, fName);
280 else if (fFileName.EndsWith(".dat"))
282 else {
283 cerr << "-E- CbmFieldMap::Init: No proper file name defined! (" << fFileName << ")" << endl;
284 Fatal("Init", "No proper file name");
285 }
286 // Fill values needed in the Print() function. This is needed to allow
287 // a constant Print() function.
288 fBxOrigin = GetBx(0., 0., 0.);
289 fByOrigin = GetBy(0., 0., 0.);
290 fBzOrigin = GetBz(0., 0., 0.);
291
292 Print();
293}
294// ------------------------------------------------------------------------
295
296
297// ----- Initialisation from arrays -------------------------------------
298void CbmFieldMap::Init(Int_t nX, Double_t xMin, Double_t xMax, Int_t nY, Double_t yMin, Double_t yMax, Int_t nZ,
299 Double_t zMin, Double_t zMax, TArrayF* bx, TArrayF* by, TArrayF* bz)
300{
301
302 Reset();
303 fNx = nX;
304 fXmin = xMin;
305 fXmax = xMax;
306 fNy = nY;
307 fYmin = yMin;
308 fYmax = yMax;
309 fNz = nZ;
310 fZmin = zMin;
311 fZmax = zMax;
312 Int_t nPoints = nX * nY * nZ;
313 assert(bx->GetSize() == by->GetSize());
314 assert(bz->GetSize() == bx->GetSize());
315 if (bx->GetSize() != nPoints) {
316 LOG(error) << "CbmFieldMap: array size " << bx->GetSize() << " does not match number of grid points.";
317 return;
318 }
319 fBx = new TArrayF(nPoints);
320 fBy = new TArrayF(nPoints);
321 fBz = new TArrayF(nPoints);
322 Int_t index = 0;
323 for (Int_t ix = 0; ix < fNx; ix++) {
324 for (Int_t iy = 0; iy < fNy; iy++) {
325 for (Int_t iz = 0; iz < fNz; iz++) {
326 index = ix * fNy * fNz + iy * fNz + iz;
327 (*fBx)[index] = (*bx)[index];
328 (*fBy)[index] = (*by)[index];
329 (*fBz)[index] = (*bz)[index];
330 } //# iz
331 } //# iy
332 } //# ix
333
334 fXstep = (fXmax - fXmin) / Double_t(fNx - 1);
335 fYstep = (fYmax - fYmin) / Double_t(fNy - 1);
336 fZstep = (fZmax - fZmin) / Double_t(fNz - 1);
337 fBxOrigin = GetBx(0., 0., 0.);
338 fByOrigin = GetBy(0., 0., 0.);
339 fBzOrigin = GetBz(0., 0., 0.);
340
341 LOG(info) << "CbmFieldMap: Initialised from " << nPoints << " grid points";
342 Print();
343}
344// ------------------------------------------------------------------------
345
346
347// ----------- Get x component of the field ---------------------------
348Double_t CbmFieldMap::GetBx(Double_t x, Double_t y, Double_t z)
349{
350
351 Int_t ix = 0;
352 Int_t iy = 0;
353 Int_t iz = 0;
354 Double_t dx = 0.;
355 Double_t dy = 0.;
356 Double_t dz = 0.;
357
358 if (IsInside(x, y, z, ix, iy, iz, dx, dy, dz)) {
359
360 // Get Bx field values at grid cell corners
361 fHa[0][0][0] = fBx->At(ix * fNy * fNz + iy * fNz + iz);
362 fHa[1][0][0] = fBx->At((ix + 1) * fNy * fNz + iy * fNz + iz);
363 fHa[0][1][0] = fBx->At(ix * fNy * fNz + (iy + 1) * fNz + iz);
364 fHa[1][1][0] = fBx->At((ix + 1) * fNy * fNz + (iy + 1) * fNz + iz);
365 fHa[0][0][1] = fBx->At(ix * fNy * fNz + iy * fNz + (iz + 1));
366 fHa[1][0][1] = fBx->At((ix + 1) * fNy * fNz + iy * fNz + (iz + 1));
367 fHa[0][1][1] = fBx->At(ix * fNy * fNz + (iy + 1) * fNz + (iz + 1));
368 fHa[1][1][1] = fBx->At((ix + 1) * fNy * fNz + (iy + 1) * fNz + (iz + 1));
369
370 // Return interpolated field value
371 return Interpolate(dx, dy, dz);
372 }
373
374 return 0.;
375}
376// ------------------------------------------------------------------------
377
378
379// ----------- Get y component of the field ---------------------------
380Double_t CbmFieldMap::GetBy(Double_t x, Double_t y, Double_t z)
381{
382
383 Int_t ix = 0;
384 Int_t iy = 0;
385 Int_t iz = 0;
386 Double_t dx = 0.;
387 Double_t dy = 0.;
388 Double_t dz = 0.;
389
390 if (IsInside(x, y, z, ix, iy, iz, dx, dy, dz)) {
391
392 // Get By field values at grid cell corners
393 fHa[0][0][0] = fBy->At(ix * fNy * fNz + iy * fNz + iz);
394 fHa[1][0][0] = fBy->At((ix + 1) * fNy * fNz + iy * fNz + iz);
395 fHa[0][1][0] = fBy->At(ix * fNy * fNz + (iy + 1) * fNz + iz);
396 fHa[1][1][0] = fBy->At((ix + 1) * fNy * fNz + (iy + 1) * fNz + iz);
397 fHa[0][0][1] = fBy->At(ix * fNy * fNz + iy * fNz + (iz + 1));
398 fHa[1][0][1] = fBy->At((ix + 1) * fNy * fNz + iy * fNz + (iz + 1));
399 fHa[0][1][1] = fBy->At(ix * fNy * fNz + (iy + 1) * fNz + (iz + 1));
400 fHa[1][1][1] = fBy->At((ix + 1) * fNy * fNz + (iy + 1) * fNz + (iz + 1));
401
402 // Return interpolated field value
403 return Interpolate(dx, dy, dz);
404 }
405
406 return 0.;
407}
408// ------------------------------------------------------------------------
409
410
411// ----------- Get z component of the field ---------------------------
412Double_t CbmFieldMap::GetBz(Double_t x, Double_t y, Double_t z)
413{
414
415 Int_t ix = 0;
416 Int_t iy = 0;
417 Int_t iz = 0;
418 Double_t dx = 0.;
419 Double_t dy = 0.;
420 Double_t dz = 0.;
421
422 if (IsInside(x, y, z, ix, iy, iz, dx, dy, dz)) {
423
424 // Get Bz field values at grid cell corners
425 fHa[0][0][0] = fBz->At(ix * fNy * fNz + iy * fNz + iz);
426 fHa[1][0][0] = fBz->At((ix + 1) * fNy * fNz + iy * fNz + iz);
427 fHa[0][1][0] = fBz->At(ix * fNy * fNz + (iy + 1) * fNz + iz);
428 fHa[1][1][0] = fBz->At((ix + 1) * fNy * fNz + (iy + 1) * fNz + iz);
429 fHa[0][0][1] = fBz->At(ix * fNy * fNz + iy * fNz + (iz + 1));
430 fHa[1][0][1] = fBz->At((ix + 1) * fNy * fNz + iy * fNz + (iz + 1));
431 fHa[0][1][1] = fBz->At(ix * fNy * fNz + (iy + 1) * fNz + (iz + 1));
432 fHa[1][1][1] = fBz->At((ix + 1) * fNy * fNz + (iy + 1) * fNz + (iz + 1));
433
434 // Return interpolated field value
435 return Interpolate(dx, dy, dz);
436 }
437
438 return 0.;
439}
440// ------------------------------------------------------------------------
441
442
443// ----------- Check whether a point is inside the map ----------------
444Bool_t CbmFieldMap::IsInside(Double_t x, Double_t y, Double_t z, Int_t& ix, Int_t& iy, Int_t& iz, Double_t& dx,
445 Double_t& dy, Double_t& dz)
446{
447
448 // --- Transform into local coordinate system
449 Double_t xl = x - fPosX;
450 Double_t yl = y - fPosY;
451 Double_t zl = z - fPosZ;
452
453 // --- Check for being outside the map range
454 if (!(xl >= fXmin && xl < fXmax && yl >= fYmin && yl < fYmax && zl >= fZmin && zl < fZmax)) {
455 ix = iy = iz = 0;
456 dx = dy = dz = 0.;
457 return kFALSE;
458 }
459
460 // --- Determine grid cell
461 ix = Int_t((xl - fXmin) / fXstep);
462 iy = Int_t((yl - fYmin) / fYstep);
463 iz = Int_t((zl - fZmin) / fZstep);
464
465
466 // Relative distance from grid point (in units of cell size)
467 dx = (xl - fXmin) / fXstep - Double_t(ix);
468 dy = (yl - fYmin) / fYstep - Double_t(iy);
469 dz = (zl - fZmin) / fZstep - Double_t(iz);
470
471 return kTRUE;
472}
473// ------------------------------------------------------------------------
474
475
476// ---------- Write the map to an ASCII file --------------------------
477void CbmFieldMap::WriteAsciiFile(const char* fileName)
478{
479
480 // Open file
481 LOG(info) << "Writing field map to ASCII file " << fileName;
482 std::ofstream mapFile(fileName);
483 if (!mapFile.is_open()) {
484 LOG(error) << "Could not open file! ";
485 return;
486 }
487
488 // Write field map grid parameters
489 mapFile.precision(6);
490 mapFile << showpoint;
491 if (fType == 1) mapFile << "nosym" << endl;
492 if (fType == 2) mapFile << "sym2" << endl;
493 if (fType == 3) mapFile << "sym3" << endl;
494 mapFile << fXmin << " " << fXmax << " " << fNx << endl;
495 mapFile << fYmin << " " << fYmax << " " << fNy << endl;
496 mapFile << fZmin << " " << fZmax << " " << fNz << endl;
497
498 // Write field values
499 Double_t factor = 10. * fScale; // Takes out scaling and converts kG->T
500 cout << right;
501 Int_t nTot = fNx * fNy * fNz;
502 cout << "-I- CbmFieldMap: " << fNx * fNy * fNz << " entries to write... " << setw(3) << 0 << " % ";
503 Int_t index = 0;
504 div_t modul;
505 Int_t iDiv = TMath::Nint(nTot / 100.);
506 for (Int_t ix = 0; ix < fNx; ix++) {
507 for (Int_t iy = 0; iy < fNy; iy++) {
508 for (Int_t iz = 0; iz < fNz; iz++) {
509 index = ix * fNy * fNz + iy * fNz + iz;
510 modul = div(index, iDiv);
511 if (modul.rem == 0) {
512 Double_t perc = TMath::Nint(100. * index / nTot);
513 cout << "\b\b\b\b\b\b" << setw(3) << perc << " % " << flush;
514 }
515 mapFile << fBx->At(index) / factor << " " << fBy->At(index) / factor << " " << fBz->At(index) / factor << endl;
516 } // z-Loop
517 } // y-Loop
518 } // x-Loop
519 cout << " " << index + 1 << " written" << endl;
520 mapFile.close();
521}
522// ------------------------------------------------------------------------
523
524
525// ------- Write field map to a ROOT file -----------------------------
526void CbmFieldMap::WriteRootFile(const char* fileName, const char* mapName)
527{
528
529 CbmFieldMapData* data = new CbmFieldMapData(mapName, *this);
530
532 TFile* oldFile = gFile;
533 TDirectory* oldDir = gDirectory;
534
535 TFile* file = new TFile(fileName, "RECREATE");
536 data->Write();
537 file->Close();
538
540 gFile = oldFile;
541 gDirectory = oldDir;
542}
543// ------------------------------------------------------------------------
544
545
546// ----- Set the position of the field centre in global coordinates -----
547void CbmFieldMap::SetPosition(Double_t x, Double_t y, Double_t z)
548{
549 fPosX = x;
550 fPosY = y;
551 fPosZ = z;
552}
553// ------------------------------------------------------------------------
554
555
556// --------- Screen output --------------------------------------------
557void CbmFieldMap::Print(Option_t*) const
558{
559 TString type = "Map";
560 if (fType == 2) type = "Map sym2";
561 if (fType == 3) type = "Map sym3";
562 cout << "======================================================" << endl;
563 cout.precision(4);
564 cout << showpoint;
565 cout << "---- " << fTitle << " : " << fName << endl;
566 cout << "----" << endl;
567 cout << "---- Field type : " << type << endl;
568 cout << "----" << endl;
569 cout << "---- Field map grid : " << endl;
570 cout << "---- x = " << setw(4) << fXmin << " to " << setw(4) << fXmax << " cm, " << fNx
571 << " grid points, dx = " << fXstep << " cm" << endl;
572 cout << "---- y = " << setw(4) << fYmin << " to " << setw(4) << fYmax << " cm, " << fNy
573 << " grid points, dy = " << fYstep << " cm" << endl;
574 cout << "---- z = " << setw(4) << fZmin << " to " << setw(4) << fZmax << " cm, " << fNz
575 << " grid points, dz = " << fZstep << " cm" << endl;
576 cout << endl;
577 cout << "---- Field centre position: ( " << setw(6) << fPosX << ", " << setw(6) << fPosY << ", " << setw(6) << fPosZ
578 << ") cm" << endl;
579 cout << "---- Field scaling factor: " << fScale << endl;
580 cout << "----" << endl;
581 cout << "---- Field at origin is ( " << setw(6) << fBxOrigin << ", " << setw(6) << fByOrigin << ", " << setw(6)
582 << fBzOrigin << ") kG" << endl;
583
584 cout << "======================================================" << endl;
585}
586// ------------------------------------------------------------------------
587
588
589// --------- Reset parameters and data (private) ----------------------
591{
592 fPosX = fPosY = fPosZ = 0.;
593 fXmin = fYmin = fZmin = 0.;
594 fXmax = fYmax = fZmax = 0.;
595 fXstep = fYstep = fZstep = 0.;
596 fNx = fNy = fNz = 0;
597 fScale = 1.;
598 if (fBx) {
599 delete fBx;
600 fBx = nullptr;
601 }
602 if (fBy) {
603 delete fBy;
604 fBy = nullptr;
605 }
606 if (fBz) {
607 delete fBz;
608 fBz = nullptr;
609 }
610}
611// ------------------------------------------------------------------------
612
613
614// ----- Read field map from ASCII file (private) ---------------------
615void CbmFieldMap::ReadAsciiFile(const char* fileName)
616{
617
618 Double_t bx = 0., by = 0., bz = 0.;
619
620 // Open file
621 LOG(info) << "CbmFieldMap: Reading field map from ASCII file " << fileName;
622 std::ifstream mapFile(fileName);
623 if (!mapFile.is_open()) {
624 LOG(error) << "CbmFieldMap:ReadAsciiFile: Could not open file!";
625 LOG(fatal) << "CbmFieldMap:ReadAsciiFile: Could not open file!";
626 }
627
628 // Read map type
629 TString type;
630 mapFile >> type;
631 Int_t iType = 0;
632 if (type == "nosym") iType = 1;
633 if (type == "sym2") iType = 2;
634 if (type == "sym3") iType = 3;
635 if (fType != iType) {
636 cout << "-E- CbmFieldMap::ReadAsciiFile: Incompatible map types!" << endl;
637 cout << " Field map is of type " << fType << " but map on file is of type " << iType << endl;
638 Fatal("ReadAsciiFile", "Incompatible map types");
639 }
640
641 // Read grid parameters
642 mapFile >> fXmin >> fXmax >> fNx;
643 mapFile >> fYmin >> fYmax >> fNy;
644 mapFile >> fZmin >> fZmax >> fNz;
645 fXstep = (fXmax - fXmin) / Double_t(fNx - 1);
646 fYstep = (fYmax - fYmin) / Double_t(fNy - 1);
647 fZstep = (fZmax - fZmin) / Double_t(fNz - 1);
648
649 // Create field arrays
650 fBx = new TArrayF(fNx * fNy * fNz);
651 fBy = new TArrayF(fNx * fNy * fNz);
652 fBz = new TArrayF(fNx * fNy * fNz);
653
654 // Read the field values
655 Double_t factor = fScale * 10.; // Factor 10 for T -> kG
656 cout << right;
657 Int_t nTot = fNx * fNy * fNz;
658 cout << "-I- CbmFieldMap: " << nTot << " entries to read... " << setw(3) << 0 << " % ";
659 Int_t index = 0;
660 div_t modul;
661 Int_t iDiv = TMath::Nint(nTot / 100.);
662 for (Int_t ix = 0; ix < fNx; ix++) {
663 for (Int_t iy = 0; iy < fNy; iy++) {
664 for (Int_t iz = 0; iz < fNz; iz++) {
665 if (!mapFile.good())
666 cerr << "-E- CbmFieldMap::ReadAsciiFile: "
667 << "I/O Error at " << ix << " " << iy << " " << iz << endl;
668 index = ix * fNy * fNz + iy * fNz + iz;
669 modul = div(index, iDiv);
670 if (modul.rem == 0) {
671 Double_t perc = TMath::Nint(100. * index / nTot);
672 cout << "\b\b\b\b\b\b" << setw(3) << perc << " % " << flush;
673 }
674 mapFile >> bx >> by >> bz;
675 fBx->AddAt(factor * bx, index);
676 fBy->AddAt(factor * by, index);
677 fBz->AddAt(factor * bz, index);
678 if (mapFile.eof()) {
679 cerr << endl
680 << "-E- CbmFieldMap::ReadAsciiFile: EOF"
681 << " reached at " << ix << " " << iy << " " << iz << endl;
682 mapFile.close();
683 break;
684 }
685 } // z-Loop
686 } // y-Loop0)
687 } // x-Loop
688
689 cout << " " << index + 1 << " read" << endl;
690
691 mapFile.close();
692}
693// ------------------------------------------------------------------------
694
695
696// ------------- Read field map from ROOT file (private) ---------------
697void CbmFieldMap::ReadRootFile(const char* fileName, const char* mapName)
698{
700 TFile* oldFile = gFile;
701 TDirectory* oldDir = gDirectory;
702
703 // Open root file
704 LOG(info) << "CbmFieldMap: Reading field map from ROOT file " << fileName;
705 TFile* file = new TFile(fileName, "READ");
706 if (!(file->IsOpen())) {
707 LOG(error) << "CbmFieldMap:ReadRootFile: Could not open file!";
708 LOG(fatal) << "CbmFieldMap:ReadRootFile: Could not open file!";
709 }
710
711 // Get the field data object
712 CbmFieldMapData* data = nullptr;
713 file->GetObject(mapName, data);
714 if (!data) {
715 cout << "-E- CbmFieldMap::ReadRootFile: data object " << fileName << " not found in file! " << endl;
716 exit(-1);
717 }
718
719 // Get the field parameters
720 SetField(data);
721
722 // Close the root file and delete the data object
723 file->Close();
724 delete data;
725
727 gFile = oldFile;
728 gDirectory = oldDir;
729}
730// ------------------------------------------------------------------------
731
732
733// ------------ Set field parameters and data (private) ----------------
735{
736
737 // Check compatibility
738 if (data->GetType() != fType) {
739 if (!((data->GetType() == 3) && (fType == 5))) // E.Litvinenko
740 {
741 cout << "-E- CbmFieldMap::SetField: Incompatible map types!" << endl;
742 cout << " Field map is of type " << fType << " but map on file is of type " << data->GetType() << endl;
743 Fatal("SetField", "Incompatible map types");
744 }
745 else
746 cout << " CbmFieldMap::SetField: Warning: You are using PosDepScaled "
747 "map (original map type = 3)"
748 << endl;
749 }
750
751
752 fXmin = data->GetXmin();
753 fYmin = data->GetYmin();
754 fZmin = data->GetZmin();
755 fXmax = data->GetXmax();
756 fYmax = data->GetYmax();
757 fZmax = data->GetZmax();
758 fNx = data->GetNx();
759 fNy = data->GetNy();
760 fNz = data->GetNz();
761 fXstep = (fXmax - fXmin) / Double_t(fNx - 1);
762 fYstep = (fYmax - fYmin) / Double_t(fNy - 1);
763 fZstep = (fZmax - fZmin) / Double_t(fNz - 1);
764 if (fBx) delete fBx;
765 if (fBy) delete fBy;
766 if (fBz) delete fBz;
767 fBx = new TArrayF(*(data->GetBx()));
768 fBy = new TArrayF(*(data->GetBy()));
769 fBz = new TArrayF(*(data->GetBz()));
770
771 // Scale and convert from T to kG
772 Double_t factor = fScale * 10.;
773 Int_t index = 0;
774 for (Int_t ix = 0; ix < fNx; ix++) {
775 for (Int_t iy = 0; iy < fNy; iy++) {
776 for (Int_t iz = 0; iz < fNz; iz++) {
777 index = ix * fNy * fNz + iy * fNz + iz;
778 if (fBx) (*fBx)[index] = (*fBx)[index] * factor;
779 if (fBy) (*fBy)[index] = (*fBy)[index] * factor;
780 if (fBz) (*fBz)[index] = (*fBz)[index] * factor;
781 }
782 }
783 }
784}
785// ------------------------------------------------------------------------
786
787
788// ------------ Interpolation in a grid cell (private) -----------------
789Double_t CbmFieldMap::Interpolate(Double_t dx, Double_t dy, Double_t dz)
790{
791
792 // Interpolate in x coordinate
793 fHb[0][0] = fHa[0][0][0] + (fHa[1][0][0] - fHa[0][0][0]) * dx;
794 fHb[1][0] = fHa[0][1][0] + (fHa[1][1][0] - fHa[0][1][0]) * dx;
795 fHb[0][1] = fHa[0][0][1] + (fHa[1][0][1] - fHa[0][0][1]) * dx;
796 fHb[1][1] = fHa[0][1][1] + (fHa[1][1][1] - fHa[0][1][1]) * dx;
797
798 // Interpolate in y coordinate
799 fHc[0] = fHb[0][0] + (fHb[1][0] - fHb[0][0]) * dy;
800 fHc[1] = fHb[0][1] + (fHb[1][1] - fHb[0][1]) * dy;
801
802 // Interpolate in z coordinate
803 return fHc[0] + (fHc[1] - fHc[0]) * dz;
804}
805// ------------------------------------------------------------------------
806
807
ClassImp(CbmConverterManager)
Double_t GetYmax() const
Double_t GetZmax() const
TString GetMapName() const
Double_t GetXmin() const
Double_t GetXmax() const
TArrayF * GetBz() const
Double_t GetYmin() const
Double_t GetZmin() const
TArrayF * GetBy() const
TArrayF * GetBx() const
Double_t GetZmin() const
Double_t GetXmax() const
TArrayF * GetBz() const
Int_t GetNz() const
TArrayF * GetBy() const
Double_t GetYmax() const
Int_t GetNy() const
Int_t GetType() const
Int_t GetNx() const
Double_t GetYmin() const
Double_t GetZmax() const
Double_t GetXmin() const
TArrayF * GetBx() const
Double_t fPosZ
Double_t fBxOrigin
Interpolated field (1-dim)
TArrayF * fBy
Double_t fPosY
TArrayF * GetBz() const
virtual void SetPosition(Double_t x, Double_t y, Double_t z)
Double_t fYmax
Double_t fPosX
Double_t fZmin
Double_t fZmax
virtual void Init()
Double_t fYmin
Double_t fHa[2][2][2]
void WriteRootFile(const char *fileName, const char *mapName)
virtual ~CbmFieldMap()
void ReadAsciiFile(const char *fileName)
virtual void Print(Option_t *="") const
TArrayF * GetBx() const
TArrayF * fBx
Double_t fHb[2][2]
Field at corners of a grid cell.
virtual Bool_t IsInside(Double_t x, Double_t y, Double_t z, Int_t &ix, Int_t &iy, Int_t &iz, Double_t &dx, Double_t &dy, Double_t &dz)
TArrayF * fBz
void SetField(const CbmFieldMapData *data)
Double_t fXmin
TString fFileName
Double_t fYstep
Double_t fZstep
Double_t Interpolate(Double_t dx, Double_t dy, Double_t dz)
Double_t fHc[2]
Interpolated field (2-dim)
void ReadRootFile(const char *fileName, const char *mapName)
void WriteAsciiFile(const char *fileName)
Double_t fXstep
Double_t fBzOrigin
y-component of the field at the origin
Double_t fScale
Double_t fByOrigin
x-component of the field at the origin
TArrayF * GetBy() const
Double_t fXmax
Double_t GetPositionZ() const
Definition CbmFieldPar.h:74
Int_t GetType() const
Definition CbmFieldPar.h:61
void MapName(TString &name)
Definition CbmFieldPar.h:71
Double_t GetPositionY() const
Definition CbmFieldPar.h:73
Double_t GetScale() const
Definition CbmFieldPar.h:75
Double_t GetPositionX() const
Definition CbmFieldPar.h:72