BmnRoot
Loading...
Searching...
No Matches
BmnFieldMap.cxx
Go to the documentation of this file.
1// ---------------------------------------------------------------------------------
2// ----- BmnFieldMap header file -------------
3// ----- Created 03/02/2015 by P. Batyuk -------------
4// ---------------------------------------------------------------------------------
5
6#include "BmnFieldMap.h"
7
8#include "BmnFieldPoint.h"
9#include "TFile.h"
10#include "TMath.h"
11#include "TSystem.h"
12#include "TTree.h"
13
14#include <iomanip>
15#include <iostream>
16using namespace std;
17
19 : FairField()
20 , fScale(1.)
21 , fPosX(0.)
22 , fPosY(0.)
23 , fPosZ(0.)
24 , fXmin(0.)
25 , fXmax(0.)
26 , fXstep(0.)
27 , fYmin(0.)
28 , fYmax(0.)
29 , fYstep(0.)
30 , fZmin(0.)
31 , fZmax(0.)
32 , fZstep(0.)
33 , fNx(0)
34 , fNy(0)
35 , fNz(0)
36 , fBx(nullptr)
37 , fBy(nullptr)
38 , fBz(nullptr)
39 , fDebugInfo(kFALSE)
40 , fIsOff(kFALSE)
41{
42 // Initilization of arrays is to my knowledge not
43 // possible in member initalization lists
44 for (Int_t i = 0; i < 2; i++) {
45 fHc[i] = 0;
46 for (Int_t j = 0; j < 2; j++) {
47 fHb[i][j] = 0;
48 for (Int_t k = 0; k < 2; k++) {
49 fHa[i][j][k] = 0;
50 }
51 }
52 }
53
54 // fName stores file name or path (using $VMCWORKDIR var) with magnetic field map
55 fName = "";
56 // fType: 1- nosym magnetic field, 2 - sym2, 3 - sym3
57 fType = 1;
58}
59
60BmnFieldMap::BmnFieldMap(const char* mapFileName)
61 : BmnFieldMap()
62{
63 fName = mapFileName;
64 fType = 1;
65}
66
68 : BmnFieldMap()
69{
70 if (!fieldPar)
71 LOG(error) << "BmnFieldMap::BmnFieldMap: empty parameter container";
72 else {
73 fPosX = fieldPar->GetPositionX();
74 fPosY = fieldPar->GetPositionY();
75 fPosZ = fieldPar->GetPositionZ();
76 fScale = fieldPar->GetScale();
77 fIsOff = fieldPar->IsFieldOff();
78 fieldPar->MapName(fName);
79 fType = fieldPar->GetType();
80 }
81}
82
84{
85 if (fBx)
86 delete fBx;
87 if (fBy)
88 delete fBy;
89 if (fBz)
90 delete fBz;
91}
92
94{
95 if (!fName.Contains("/")) {
96 TString localDbPath = gSystem->Getenv("DBL_FILE_PATH");
97 if (localDbPath == "")
98 fName = "$VMCWORKDIR/input/" + fName;
99 else
100 fName = "$DBL_FILE_PATH/storage/" + fName;
101 }
102
103 if (fName.Contains("extrap")) {
104 LOG(info) << "Extrapolated magnetic field map is used";
105 gSystem->ExpandPathName(fName);
107 } else {
108 LOG(info) << "Basic magnetic field map is used";
109 gSystem->ExpandPathName(fName);
111 }
115}
116
117void BmnFieldMap::WriteAsciiFile(const char* asciiFileName)
118{
119 TString strFileName(asciiFileName);
120 gSystem->ExpandPathName(strFileName);
121 // Open file
122 LOG(info) << "BmnFieldMap: Writing field map to ASCII file " << strFileName.Data() << endl;
123 ofstream mapFile(strFileName.Data());
124 if (!mapFile.is_open()) {
125 LOG(error) << "BmnFieldMap:ReadAsciiFile: Could not open file" << endl;
126 return;
127 }
128
129 // Write field map grid parameters
130 mapFile.precision(6);
131 mapFile << showpoint;
132 if (fType == 1)
133 mapFile << "nosym" << endl;
134 if (fType == 2)
135 mapFile << "sym2" << endl;
136 if (fType == 3)
137 mapFile << "sym3" << endl;
138 mapFile << fXmin << " " << fXmax << " " << fNx << endl;
139 mapFile << fYmin << " " << fYmax << " " << fNy << endl;
140 mapFile << fZmin << " " << fZmax << " " << fNz << endl;
141
142 // Write field values
143 Double_t factor = 10. * fScale; // Takes out scaling and converts kG->T
144 cout << right;
145 Int_t nTot = fNx * fNy * fNz;
146 LOG(info) << "BmnFieldMap: " << fNx * fNy * fNz << " entries to write... " << setw(3) << 0 << " % ";
147 Int_t index = 0;
148 div_t modul;
149 Int_t iDiv = TMath::Nint(nTot / 100.);
150 for (Int_t ix = 0; ix < fNx; ix++) {
151 for (Int_t iy = 0; iy < fNy; iy++) {
152 for (Int_t iz = 0; iz < fNz; iz++) {
153 index = ix * fNy * fNz + iy * fNz + iz;
154 modul = div(index, iDiv);
155 if (modul.rem == 0) {
156 Double_t perc = TMath::Nint(100. * index / nTot);
157 cout << "\b\b\b\b\b\b" << setw(3) << perc << " % " << flush;
158 }
159 mapFile << fBx->At(index) / factor << " " << fBy->At(index) / factor << " " << fBz->At(index) / factor
160 << endl;
161 } // z-Loop
162 } // y-Loop
163 } // x-Loop
164 cout << " " << index + 1 << " written" << endl;
165 mapFile.close();
166}
167
168void BmnFieldMap::WriteRootFile(const char* rootFileName, const char* mapFileName)
169{
170 BmnFieldMapData* data = new BmnFieldMapData(mapFileName, *this);
171 TFile* oldFile = gFile;
172 TString strFileName(rootFileName);
173 gSystem->ExpandPathName(strFileName);
174 TFile* file = TFile::Open(strFileName, "RECREATE");
175 data->Write();
176 delete file;
177 if (oldFile)
178 oldFile->cd();
179}
180
181void BmnFieldMap::SetPosition(Double_t x, Double_t y, Double_t z)
182{
183 fPosX = x;
184 fPosY = y;
185 fPosZ = z;
189}
190
191void BmnFieldMap::Print(Option_t*) const
192{
193 TString type = "Map";
194 if (fType == 2)
195 type = "Map sym2";
196 if (fType == 3)
197 type = "Map sym3";
198 cout << "======================================================" << endl;
199 cout.precision(4);
200 cout << showpoint;
201 cout << "---- " << fTitle << " : " << fName << endl;
202 cout << "----" << endl;
203 cout << "---- Field type : " << type << endl;
204 cout << "----" << endl;
205 cout << "---- Field map grid : " << endl;
206 cout << "---- x = " << setw(4) << fXmin << " to " << setw(4) << fXmax << " cm, " << fNx
207 << " grid points, dx = " << fXstep << " cm" << endl;
208 cout << "---- y = " << setw(4) << fYmin << " to " << setw(4) << fYmax << " cm, " << fNy
209 << " grid points, dy = " << fYstep << " cm" << endl;
210 cout << "---- z = " << setw(4) << fZmin << " to " << setw(4) << fZmax << " cm, " << fNz
211 << " grid points, dz = " << fZstep << " cm" << endl;
212 cout << endl;
213 cout << "---- Target position: ( " << setw(6) << fPosX << ", " << setw(6) << fPosY << ", " << setw(6) << fPosZ
214 << ") cm" << endl;
215 cout << "---- Field scaling factor: " << fScale << endl;
216 if (fIsOff)
217 cout << "---- Field is off ----" << endl;
218 cout << "----" << endl;
219 cout << "---- Field at Target Position ( " << setw(6) << fPosBx << ", " << setw(6) << fPosBy << ", " << setw(6)
220 << fPosBz << ") kG" << endl;
221 cout << "======================================================" << endl;
222}
223
225{
226 fPosX = fPosY = fPosZ = 0.;
227 fPosBx = fPosBy = fPosBz = 0.;
228 fXmin = fYmin = fZmin = 0.;
229 fXmax = fYmax = fZmax = 0.;
230 fXstep = fYstep = fZstep = 0.;
231 fNx = fNy = fNz = 0;
232 fScale = 1.;
233 fIsOff = kFALSE;
234 if (fBx) {
235 delete fBx;
236 fBx = NULL;
237 }
238 if (fBy) {
239 delete fBy;
240 fBy = NULL;
241 }
242 if (fBz) {
243 delete fBz;
244 fBz = NULL;
245 }
246}
247
248void BmnFieldMap::ReadAsciiFile(const char* asciiFileName)
249{
250 Double_t bx = 0., by = 0., bz = 0.;
251
252 TString strFileName(asciiFileName);
253 gSystem->ExpandPathName(strFileName);
254 // Open file
255 LOG(info) << "BmnFieldMap: Reading field map from ASCII file " << strFileName.Data();
256 ifstream mapFile(strFileName.Data());
257 if (!mapFile.is_open()) {
258 LOG(error) << "BmnFieldMap:ReadAsciiFile: Could not open file!";
259 LOG(fatal) << "BmnFieldMap:ReadAsciiFile: Could not open file!";
260 }
261
262 // Read map type
263 TString type;
264 mapFile >> type;
265 Int_t iType = 0;
266 if (type == "nosym")
267 iType = 1;
268 if (type == "sym2")
269 iType = 2;
270 if (type == "sym3")
271 iType = 3;
272 if (fType != iType) {
273 cout << "-E- BmnFieldMap::ReadAsciiFile: Incompatible map types!" << endl;
274 cout << " Field map is of type " << fType << " but map on file is of type " << iType << endl;
275 Fatal("ReadAsciiFile", "Incompatible map types");
276 }
277
278 // Read grid parameters
279 mapFile >> fXmin >> fXmax >> fNx;
280 mapFile >> fYmin >> fYmax >> fNy;
281 mapFile >> fZmin >> fZmax >> fNz;
282 fXstep = (fXmax - fXmin) / Double_t(fNx - 1);
283 fYstep = (fYmax - fYmin) / Double_t(fNy - 1);
284 fZstep = (fZmax - fZmin) / Double_t(fNz - 1);
285
286 // Create field arrays
287 fBx = new TArrayF(fNx * fNy * fNz);
288 fBy = new TArrayF(fNx * fNy * fNz);
289 fBz = new TArrayF(fNx * fNy * fNz);
290
291 // Read the field values
292 Double_t factor = fScale * 10.; // Factor 10 for T -> kG
293 cout << right;
294 Int_t nTot = fNx * fNy * fNz;
295 // cout << "-I- BmnFieldMap: " << nTot << " entries to read... " << setw(3) << 0 << " % ";
296 Int_t index = 0;
297 div_t modul;
298 Int_t iDiv = TMath::Nint(nTot / 100.);
299 for (Int_t ix = 0; ix < fNx; ix++) {
300 for (Int_t iy = 0; iy < fNy; iy++) {
301 for (Int_t iz = 0; iz < fNz; iz++) {
302 if (!mapFile.good())
303 cerr << "-E- BmnFieldMap::ReadAsciiFile: "
304 << "I/O Error at " << ix << " " << iy << " " << iz << endl;
305 index = ix * fNy * fNz + iy * fNz + iz;
306 modul = div(index, iDiv);
307 // if (modul.rem == 0) {
308 // Double_t perc = TMath::Nint(100. * index / nTot);
309 // cout << "\b\b\b\b\b\b" << setw(3) << perc << " % " << flush;
310 //}
311 mapFile >> bx >> by >> bz;
312 fBx->AddAt(factor * bx, index);
313 fBy->AddAt(factor * by, index);
314 fBz->AddAt(factor * bz, index);
315 if (mapFile.eof()) {
316 cerr << endl
317 << "-E- BmnFieldMap::ReadAsciiFile: EOF"
318 << " reached at " << ix << " " << iy << " " << iz << endl;
319 mapFile.close();
320 break;
321 }
322 } // z-Loop
323 } // y-Loop0)
324 } // x-Loop
325
326 cout << " " << index + 1 << " read" << endl;
327
328 fXmin += fPosX;
329 fXmax += fPosX;
330
331 fYmin += fPosY;
332 fYmax += fPosY;
333
334 fZmin += fPosZ;
335 fZmax += fPosZ;
336
337 mapFile.close();
338}
339
340void BmnFieldMap::ReadRootFileNewFormat(const char* rootFileName)
341{
342 TString strFileName(rootFileName);
343 gSystem->ExpandPathName(strFileName);
344 // Open root file
345 LOG(info) << "BmnFieldMap: Reading field map from ROOT file " << strFileName.Data();
346 TFile* file = TFile::Open(strFileName, "READ");
347 if ((file == nullptr) || (!file->IsOpen()))
348 LOG(fatal) << "BmnFieldMap:ReadRootFileNewFormat: Could not open the file with a field map";
349
350 BmnFieldPoint point;
351 TTree* tree = (TTree*)file->Get("Field_map");
352 tree->SetBranchAddress("field_points", &point);
353
354 fNx = 126;
355 fXmin = -155.78;
356 fXmax = 144.82;
357 fXstep = (fXmax - fXmin) / (fNx - 1);
358
359 fNy = 47;
360 fYmin = -37.67;
361 fYmax = 54.33;
362 fYstep = (fYmax - fYmin) / (fNy - 1);
363
364 fNz = 241;
365 fZmin = -162.16;
366 fZmax = 439.04;
367 fZstep = (fZmax - fZmin) / (fNz - 1);
368
369 Int_t nTot = fNx * fNy * fNz;
370 fBx = new TArrayF(nTot);
371 fBy = new TArrayF(nTot);
372 fBz = new TArrayF(nTot);
373
374 Int_t index = 0;
375
376 for (Int_t i = 0; i < tree->GetEntries(); i++) {
377 tree->GetEntry(i);
378 // В дереве, сформированном из CSV файлов, данные записаны в одном порядке, а BmnFieldMap они хранятся в другом.
379 // Поэтому сначала из общего сквозного индекса получаем отдельные индексы по каждому измерению.
380 // При этом индекс по x инвертируем, так как направление оси в Bmn и в рамке с замерами поля противоположные.
381 Int_t iz = i % fNz;
382 Int_t ix = fNx - ((i / fNz) % fNx) - 1;
383 Int_t iy = ((i / fNz) / fNx) % fNy;
384 // После этого формируем итоговый сквозной индекс для карты поля в формате BmnFieldMap
385 index = ix * fNy * fNz + iy * fNz + iz;
386 // И пишем по этому индексу компоненты, предварительно скалируя их.
387 fBx->AddAt(fScale * point.GetBx(), index);
388 fBy->AddAt(fScale * point.GetBy(), index);
389 fBz->AddAt(fScale * point.GetBz(), index);
390 }
391
392 delete file;
393}
394
395void BmnFieldMap::ReadRootFileNewFormatExtrap(const char* rootFileName)
396{
397 TString strFileName(rootFileName);
398 gSystem->ExpandPathName(strFileName);
399 // Open root file
400 LOG(info) << "BmnFieldMap: Reading field map from ROOT file " << strFileName.Data();
401 TFile* file = TFile::Open(strFileName, "READ");
402 if ((file == nullptr) || (!file->IsOpen()))
403 LOG(fatal) << "BmnFieldMap:ReadRootFileNewFormatExtrap: Could not open the file with a field map";
404
405 BmnFieldPoint point;
406 TTree* tree = (TTree*)file->Get("Field_map");
407 tree->SetBranchAddress("field_points", &point);
408
409 fNx = 206;
410 fXmin = -251.972;
411 fXmax = 241.012;
412 fXstep = (fXmax - fXmin) / (fNx - 1);
413
414 fNy = 87;
415 fYmin = -77.67;
416 fYmax = 94.33;
417 fYstep = (fYmax - fYmin) / (fNy - 1);
418
419 fNz = 341;
420 fZmin = -287.41;
421 fZmax = 564.29;
422 fZstep = (fZmax - fZmin) / (fNz - 1);
423
424 Int_t nTot = fNx * fNy * fNz;
425 fBx = new TArrayF(nTot);
426 fBy = new TArrayF(nTot);
427 fBz = new TArrayF(nTot);
428
429 Int_t index = 0;
430
431 for (Int_t i = 0; i < tree->GetEntries(); i++) {
432 tree->GetEntry(i);
433 // В дереве, сформированном из CSV файлов, данные записаны в одном порядке, а BmnFieldMap они хранятся в другом.
434 // Поэтому сначала из общего сквозного индекса получаем отдельные индексы по каждому измерению.
435 // При этом индекс по x инвертируем, так как направление оси в Bmn и в рамке с замерами поля противоположные.
436 Int_t iz = i % fNz;
437 Int_t ix = fNx - ((i / fNz) % fNx) - 1;
438 Int_t iy = ((i / fNz) / fNx) % fNy;
439 // После этого формируем итоговый сквозной индекс для карты поля в формате BmnFieldMap
440 index = ix * fNy * fNz + iy * fNz + iz;
441 // И пишем по этому индексу компоненты, предварительно скалируя их.
442 fBx->AddAt(fScale * point.GetBx(), index);
443 fBy->AddAt(fScale * point.GetBy(), index);
444 fBz->AddAt(fScale * point.GetBz(), index);
445 }
446
447 delete file;
448}
449
450void BmnFieldMap::ReadRootFile(const char* rootFileName)
451{
452 TString* type = 0;
453 coordinate_info_t X, Y, Z;
454 vector<Double_t>* read_Field = NULL;
455
456 TString strFileName(rootFileName);
457 gSystem->ExpandPathName(strFileName);
458 // Open root file
459 LOG(info) << "BmnFieldMap: Reading field map from ROOT file " << strFileName.Data();
460 TFile* file = TFile::Open(strFileName, "READ");
461 if ((file == nullptr) || (!file->IsOpen()))
462 LOG(fatal) << "BmnFieldMap:ReadRootFile: Could not open the file with a field map";
463
464 TTree* tree = (TTree*)file->Get("Main_info");
465 TTree* t = (TTree*)file->Get("Field_map");
466 tree->SetBranchAddress("Field_type", &type);
467 tree->SetBranchAddress("Main_info_X", &X);
468 tree->SetBranchAddress("Main_info_Y", &Y);
469 tree->SetBranchAddress("Main_info_Z", &Z);
470 t->SetBranchAddress("field_map", &read_Field);
471
472 tree->GetEntry(0);
473 Int_t iType = 0;
474
475 if (*type == "nosym")
476 iType = 1;
477 if (*type == "sym2")
478 iType = 2;
479 if (*type == "sym3")
480 iType = 3;
481 if (fType != iType) {
482 LOG(error) << "BmnFieldMap::ReadRootFile: Incompatible map types!";
483 LOG(info) << " Field map is of type " << fType << " but map on file is of type " << iType;
484 LOG(fatal) << "ReadRootFile: Incompatible map types";
485 }
486
487 fNx = X.N;
488 fXmin = X.min + fPosX;
489 fXmax = X.max + fPosX;
490 fXstep = X.step;
491
492 fNy = Y.N;
493 fYmin = Y.min + fPosY;
494 fYmax = Y.max + fPosY;
495 fYstep = Y.step;
496
497 fNz = Z.N;
498 fZmin = Z.min + fPosZ;
499 fZmax = Z.max + fPosZ;
500 fZstep = Z.step;
501
502 fBx = new TArrayF(fNx * fNy * fNz);
503 fBy = new TArrayF(fNx * fNy * fNz);
504 fBz = new TArrayF(fNx * fNy * fNz);
505
506 Double_t factor = fScale;
507 Int_t nTot = fNx * fNy * fNz;
508 Int_t index = 0;
509 div_t modul;
510 Int_t iDiv = TMath::Nint(nTot / 100.);
511
512 for (Int_t ix = 0; ix < fNx; ix++) {
513 for (Int_t iy = 0; iy < fNy; iy++) {
514 for (Int_t iz = 0; iz < fNz; iz++) {
515 modul = div(index, iDiv);
516 // if (modul.rem == 0) {
517 // Double_t perc = TMath::Nint(100. * index / nTot);
518 // }
519 index = ix * fNy * fNz + iy * fNz + iz;
520 t->GetEntry(index);
521 vector<Double_t> v = (*read_Field);
522 fBx->AddAt(factor * v[0], index);
523 fBy->AddAt(factor * v[1], index);
524 fBz->AddAt(factor * v[2], index);
525 } // z-Loop
526 } // y-Loop
527 } // x-Loop
528
529 t->Delete();
530 tree->Delete();
531 delete type;
532 delete read_Field;
533 delete file;
534}
535
537{
538 // Check compatibility
539 if (data->GetType() != fType) {
540 if (!((data->GetType() == 3) && (fType == 5))) // E.Litvinenko
541 {
542 LOG(error) << "BmnFieldMap::SetField: Incompatible map types!";
543 LOG(info) << " Field map is of type " << fType << " but map on file is of type " << data->GetType();
544 LOG(fatal) << "SetField: Incompatible map types";
545 } else
546 LOG(info) << " BmnFieldMap::SetField: Warning: You are using PosDepScaled map (original map type = 3)";
547 }
548
549 fXmin = data->GetXmin();
550 fYmin = data->GetYmin();
551 fZmin = data->GetZmin();
552 fXmax = data->GetXmax();
553 fYmax = data->GetYmax();
554 fZmax = data->GetZmax();
555 fNx = data->GetNx();
556 fNy = data->GetNy();
557 fNz = data->GetNz();
558 fXstep = (fXmax - fXmin) / Double_t(fNx - 1);
559 fYstep = (fYmax - fYmin) / Double_t(fNy - 1);
560 fZstep = (fZmax - fZmin) / Double_t(fNz - 1);
561 if (fBx)
562 delete fBx;
563 if (fBy)
564 delete fBy;
565 if (fBz)
566 delete fBz;
567 fBx = new TArrayF(*(data->GetBx()));
568 fBy = new TArrayF(*(data->GetBy()));
569 fBz = new TArrayF(*(data->GetBz()));
570
571 // Scale and convert from T to kG
572 Double_t factor = fScale * 10.;
573 Int_t index = 0;
574 for (Int_t ix = 0; ix < fNx; ix++) {
575 for (Int_t iy = 0; iy < fNy; iy++) {
576 for (Int_t iz = 0; iz < fNz; iz++) {
577 index = ix * fNy * fNz + iy * fNz + iz;
578 if (fBx)
579 (*fBx)[index] = (*fBx)[index] * factor;
580 if (fBy)
581 (*fBy)[index] = (*fBy)[index] * factor;
582 if (fBz)
583 (*fBz)[index] = (*fBz)[index] * factor;
584 }
585 }
586 }
587}
588
589Double_t BmnFieldMap::Interpolate(Double_t dx, Double_t dy, Double_t dz)
590{
591 // Interpolate in x coordinate
592 fHb[0][0] = fHa[0][0][0] + (fHa[1][0][0] - fHa[0][0][0]) * dx;
593 fHb[1][0] = fHa[0][1][0] + (fHa[1][1][0] - fHa[0][1][0]) * dx;
594 fHb[0][1] = fHa[0][0][1] + (fHa[1][0][1] - fHa[0][0][1]) * dx;
595 fHb[1][1] = fHa[0][1][1] + (fHa[1][1][1] - fHa[0][1][1]) * dx;
596
597 // Interpolate in y coordinate
598 fHc[0] = fHb[0][0] + (fHb[1][0] - fHb[0][0]) * dy;
599 fHc[1] = fHb[0][1] + (fHb[1][1] - fHb[0][1]) * dy;
600
601 // Interpolate in z coordinate
602 return fHc[0] + (fHc[1] - fHc[0]) * dz;
603}
__m128 v
Definition P4_F32vec4.h:1
int i
Definition P4_F32vec4.h:22
TArrayF * GetBy() const
Double_t GetYmax() const
Double_t GetZmin() const
Double_t GetXmin() const
TArrayF * GetBz() const
Int_t GetNz() const
Int_t GetType() const
Double_t GetXmax() const
Int_t GetNy() const
TArrayF * GetBx() const
Int_t GetNx() const
Double_t GetYmin() const
Double_t GetZmax() const
void ReadRootFileNewFormat(const char *rootFileName)
Double_t fHa[2][2][2]
TArrayF * fBy
void SetField(const BmnFieldMapData *data)
Bool_t fIsOff
Double_t fPosZ
Double_t fPosBx
TArrayF * GetBx() const
Definition BmnFieldMap.h:96
void ReadRootFileNewFormatExtrap(const char *rootFileName)
void ReadRootFile(const char *rootFileName)
Double_t fPosY
TArrayF * GetBy() const
Definition BmnFieldMap.h:97
TArrayF * fBz
Double_t Interpolate(Double_t dx, Double_t dy, Double_t dz)
Double_t fScale
Double_t fPosX
Double_t fHc[2]
Interpolated field (2-dim)
virtual ~BmnFieldMap()
void SetPosition(Double_t x, Double_t y, Double_t z)
Double_t fZmin
TArrayF * GetBz() const
Definition BmnFieldMap.h:98
Double_t fYmax
Double_t fHb[2][2]
Field at corners of a grid cell.
void WriteRootFile(const char *rootFileName, const char *mapFileName)
virtual void Init()
Double_t fPosBy
Double_t fXstep
Double_t fZstep
Double_t fPosBz
Double_t fXmin
TArrayF * fBx
void ReadAsciiFile(const char *asciiFileName)
Double_t fZmax
Double_t fYmin
Double_t fYstep
void WriteAsciiFile(const char *asciiFileName)
Double_t fXmax
virtual void Print(Option_t *) const
Double_t GetPositionZ() const
Definition BmnFieldPar.h:69
Bool_t IsFieldOff() const
Definition BmnFieldPar.h:71
Int_t GetType() const
Definition BmnFieldPar.h:56
Double_t GetScale() const
Definition BmnFieldPar.h:70
Double_t GetPositionY() const
Definition BmnFieldPar.h:68
Double_t GetPositionX() const
Definition BmnFieldPar.h:67
void MapName(TString &name)
Definition BmnFieldPar.h:66
Double_t GetBy()
Double_t GetBz()
Double_t GetBx()
STL namespace.