FairRoot/PandaRoot
read_jost_fieldmaps.cxx
Go to the documentation of this file.
1 // Panda fieldmaps created with tosca by jost luehning are read
2 // and transformed to panda compatible fieldmaps
3 
4 #include<TString.h>
5 #include<TStopwatch.h>
6 #include<TROOT.h>
7 #include<FairModule.h>
8 #include<PndCave.h>
9 #include<PndMagnet.h>
10 #include<PndPipe.h>
11 #include<PndMvdDetector.h>
12 #include<FairRunSim.h>
13 #include<PndLmdDetector.h>
14 #include<FairPrimaryGenerator.h>
15 #include<FairBoxGenerator.h>
16 #include<PndMultiField.h>
17 #include<PndTransMap.h>
18 #include<PndDipoleMap.h>
19 #include<PndSolenoidMap.h>
20 #include<FairParRootFileIo.h>
21 #include<PndMultiFieldPar.h>
22 #include<map>
23 #include<vector>
24 #include<fstream>
25 #include<iostream>
26 #include<iomanip>
27 #include<sstream>
28 #include<TH1.h>
29 #include<TTree.h>
30 #include<TFile.h>
31 #include<TCanvas.h>
32 #include<TVector3.h>
33 #include "TPolyLine3D.h"
34 #include "TView.h"
35 #include "TCanvas.h"
36 #include "TRandom.h"
37 #include<TApplication.h>
38 #include "TGLViewer.h"
39 #include <TGeoManager.h>
40 #include "TGeoMaterial.h"
41 #include "TGeoMedium.h"
42 #include "TGeoVolume.h"
43 #include "TString.h"
44 #include "TSystem.h"
45 #include "FairGeoLoader.h"
46 #include "FairGeoInterface.h"
47 #include "FairGeoMedia.h"
48 #include "FairGeoBuilder.h"
49 
50 using namespace std;
51 
52 // for testing purposes the default
53 // field map name is modified
54 // Moreover in case of setting a string of
55 // length > 0
56 // files are stored in the $VMCWORKDIR/input/ folder
57 const string name_modification = "";
58 
59 vector<string> fieldmapnames;
60 
61 TTree* gtree_values(NULL);
62 
63 double gBx, gBy, gBz, gDBx, gDBy, gDBz;
64 
65 // read fem calculated fieldmap and transform it
66 // fieldtype 1: no symetry
67 // fieldtype 2: Solenoidfieled
68 // fieldtype 3: Dipolefield
69 // fieldtype 4: Transition region
70 // filename_out is the pandaroot filename without the extenstion
71 // in contrast to filename_in having the extension
72 void transform_jost_fieldmap(string filename_in, string filename_out, int fieldtype){
73  cout << " reading file " << filename_in << endl;
74  ifstream jostmapFile(filename_in.c_str());
75  if (!jostmapFile.is_open()) {
76  cout << "ReadAsciiFile: Could not open file! " << endl;
77  return;
78  }
79  // coordinate min max ranges
80  double x_min = -1.e9;
81  double x_max = -1.e9;
82  double y_min = -1.e9;
83  double y_max = -1.e9;
84  double z_min = -1.e9;
85  double z_max = -1.e9;
86  // number of entries
87  int n_x = 0;
88  int n_y = 0;
89  int n_z = 0;
90  // spacing (maximum allowed precision is 1 mm)
91  int d_x = 0;
92  int d_y = 0;
93  int d_z = 0;
94  // store the grid axes
95  map <double, int> x_counts;
96  map <double, int> y_counts;
97  map <double, int> z_counts;
98  map <double, int>::iterator it_counts, it_counts_prev;
99  map <string, TVector3> field_in;
100  map <string, TVector3>::iterator it_field_in;
101  // first pass: analyzing parameters
102  while (1)
103  {
104  string line;
105  getline (jostmapFile , line);
106  if (!jostmapFile.good()) break;
107  // go through the line
108  bool isvalidline = true;
109  // I'm expecting 6 entries x,y,z and Bx, By, Bz;
110  char value[100];
111  value[0] = '\0';
112  unsigned int ivalchar = 0;
113  double values[6];
114  unsigned int ival = 0;
115  stringstream key;
116  for (unsigned int ichar = 0; ichar < line.size(); ichar++){
117  if (line[ichar] == '#'){
118  isvalidline = false;
119  break;
120  }
121  //if (line[ichar] == '\n' || line[ichar] == '\0'){
122  //ival = 0;
123  // break;
124  //}
125  if (line[ichar] == ' ' || line[ichar] == '\n' || line[ichar] == '\0' || ichar == line.size()-1){
126  if (value[0] != '\0' && ivalchar != 0){
127  if (ival > 5) {
128  cout << " hmm, more than 6 values per line? -> Fatal " << endl;
129  break;
130  }
131  value[ivalchar] = '\n';
132  values[ival] = atof(value);
133  value[0] = '\0';
134  ivalchar = 0;
135  ival++;
136  }
137  if (line[ichar] == '\n' || line[ichar] == '\0'){
138  ival = 0;
139  break;
140  } else continue;
141  }
142  value[ivalchar] = line[ichar];
143  ivalchar++;
144  }
145  if (ival != 6 && ival > 0 && isvalidline) {
146  cout << " not enough values in this line!!! -> Fatal " << ival << endl;
147  cout << line << endl;
148  continue;
149  } else {
150  if (isvalidline){
151  // analyze the values
152  if (x_min == -1.e9 || values[0] < x_min) x_min = values[0];
153  if (x_max == -1.e9 || values[0] > x_max) x_max = values[0];
154  if (y_min == -1.e9 || values[1] < y_min) y_min = values[1];
155  if (y_max == -1.e9 || values[1] > y_max) y_max = values[1];
156  if (z_min == -1.e9 || values[2] < z_min) z_min = values[2];
157  if (z_max == -1.e9 || values[2] > z_max) z_max = values[2];
158  x_counts[values[0]]++;
159  y_counts[values[1]]++;
160  z_counts[values[2]]++;
161  // maximum resolution in the grid spacing will be 1 mm
162  key << ((int)floor(values[0])) << "_" << ((int)floor(values[1])) << "_" << ((int)floor(values[2]));
163  TVector3 fieldval(values[3], values[4], values[5]);
164  field_in[key.str()] = fieldval;
165  }
166  }
167  }
168  jostmapFile.close();
169 
170  cout << " finished reading " << endl;
171 
172  n_x = x_counts.size();
173  n_y = y_counts.size();
174  n_z = z_counts.size();
175  d_x = (int)floor((x_max - x_min)/(n_x-1));
176  d_y = (int)floor((y_max - y_min)/(n_y-1));
177  d_z = (int)floor((z_max - z_min)/(n_z-1));
178 
179  cout << "grid size is " << n_x << " X " << n_y << " X " << n_z << endl;
180  cout << "grid spacing is [mm] " << d_x << " " << d_y << " " << d_z << endl;
181  cout << x_min << " < x [mm] < " << x_max << endl;
182  cout << y_min << " < y [mm] < " << y_max << endl;
183  cout << z_min << " < z [mm] < " << z_max << endl;
184 
185  cout << "checking consistency " << endl;
186 
187  // get spacing and check grid in the file
188  for (it_counts = x_counts.begin(); it_counts != x_counts.end(); it_counts++){
189  if (it_counts == x_counts.begin()) {
190  it_counts_prev = it_counts;
191  continue;
192  }
193  double spacing = it_counts->first - it_counts_prev->first;
194  if (fabs(spacing - d_x) > 1e-5)
195  cout << " warning: grid spacing does not match for dx: " << spacing << " vs " << d_x << endl;
196  if (it_counts->second != it_counts_prev->second)
197  cout << " warning: grid is not fully determined in x! " << endl;
198  it_counts_prev = it_counts;
199  }
200  for (it_counts = y_counts.begin(); it_counts != y_counts.end(); it_counts++){
201  if (it_counts == y_counts.begin()) {
202  it_counts_prev = it_counts;
203  continue;
204  }
205  double spacing = it_counts->first - it_counts_prev->first;
206  if (fabs(spacing - d_y) > 1e-5)
207  cout << " warning: grid spacing does not match for dy: " << spacing << " vs " << d_y << endl;
208  if (it_counts->second != it_counts_prev->second)
209  cout << " warning: grid is not fully determined in y! " << endl;
210  it_counts_prev = it_counts;
211  }
212  for (it_counts = z_counts.begin(); it_counts != z_counts.end(); it_counts++){
213  if (it_counts == z_counts.begin()) {
214  it_counts_prev = it_counts;
215  continue;
216  }
217  double spacing = it_counts->first - it_counts_prev->first;
218  if (fabs(spacing - d_z) > 1e-5)
219  cout << " warning: grid spacing does not match for dz: " << spacing << " vs " << d_z << endl;
220  if (it_counts->second != it_counts_prev->second)
221  cout << " warning: grid is not fully determined in z! " << endl;
222  it_counts_prev = it_counts;
223  }
224  cout << " done " << endl;
225  // now we can try to transform the maps
226  // units are assumed to be gaus and mm
227  // pandaroot needs kGaus and cm
228 
229  cout << "writing " << filename_out << " files" << endl;
230  ofstream mapFile((filename_out+".dat").c_str());
231  if (!mapFile.is_open()) {
232  cout << "WriteAsciiFile: Could not open file! " << endl;
233  return;
234  }
235 
236  // Write field map grid parameters
237  mapFile.precision(6);
238  mapFile << showpoint;
239  if (fieldtype == 1)
240  mapFile << "nosym" << endl;
241  if (fieldtype == 2)
242  mapFile << "Solenoid" << endl;
243  if (fieldtype == 3)
244  mapFile << "Dipole" << endl;
245  if (fieldtype == 4)
246  mapFile << "Trans" << endl;
247  //if (map_to_corr->GetUnit() == 10.0)
248  // mapFile << "T" << endl;
249  //if (map_to_corr->GetUnit() == 0.001)
250  double unit = 0.001;
251  mapFile << "G" << endl; // units are assumed in Gauss
252  //if (map_to_corr->GetUnit() == 1.0)
253  // mapFile << "kG" << endl;
254 
255  mapFile << x_min/10. << " " << x_max/10. << " "
256  << n_x << endl;
257  mapFile << y_min/10. << " " << y_max/10. << " "
258  << n_y << endl;
259  mapFile << z_min/10. << " " << z_max/10. << " "
260  << n_z << endl;
261 
262  // Write field values
263  double scale = 1.;
264  Double_t factor = unit * scale;
265  cout << right;
266  Int_t nTot = n_x * n_y * n_z;
267  cout << "-I- PndFieldMap: " << n_x * n_y * n_z << " entries to write... "
268  << setw(3) << 0 << " % ";
269  Int_t index = 0;
270  div_t modul;
271  Int_t iDiv = TMath::Nint(nTot / 100.);
272 
273  for (int x = 0; x < n_x; x++) {
274  for (int y = 0; y < n_y; y++) {
275  for (int z = 0; z < n_z; z++) {
276  index = x * n_y * n_z + y
277  * n_z + z;
278  //double Bx1 = (*(map_to_corr->GetBx()))[index];
279  //double Bx2 = (*(map_to_interpol->GetBx()))[index];
280  //double By1 = (*(map_to_corr->GetBy()))[index];
281  //double By2 = (*(map_to_interpol->GetBy()))[index];
282  //double Bz1 = (*(map_to_corr->GetBz()))[index];
283  //double Bz2 = (*(map_to_interpol->GetBz()))[index];
284  stringstream key;
285  key << (int)floor(x_min + x * d_x) << "_";
286  key << (int)floor(y_min + y * d_y) << "_";
287  key << (int)floor(z_min + z * d_z);
288  it_field_in = field_in.find(key.str());
289  if (it_field_in != field_in.end()){
290  double new_Bx = it_field_in->second.X();
291  double new_By = it_field_in->second.Y();
292  double new_Bz = it_field_in->second.Z();
293 
294  //cout << " changing \t" << Bx1 << '\t' << By1 << '\t' << Bz1 << '\n' ;
295  //cout << " to \t" << new_Bx << '\t' << new_By << '\t' << new_Bz << '\n' << endl;
296  modul = div(index, iDiv);
297  if (modul.rem == 0) {
298  Double_t perc = TMath::Nint(100. * index / nTot);
299  cout << "\b\b\b\b\b\b" << setw(3) << perc << " % " << flush;
300  }
301  //mapFile << new_Bx / factor << " " << new_By
302  // / factor << " " << new_Bz / factor << endl;
303  mapFile << new_Bx << " " << new_By << " " << new_Bz << endl;
304  } else {
305  cout << " fatal: no value to key " << key << " found!" << endl;
306  }
307  }
308  }
309  }
310 
311  // for (Int_t ix = 0; ix < fNx; ix++) {
312  // for (Int_t iy = 0; iy < fNy; iy++) {
313  // for (Int_t iz = 0; iz < fNz; iz++) {
314  // index = ix * fNy * fNz + iy * fNz + iz;
315  // } // z-Loop
316  // } // y-Loop
317  // } // x-Loop
318  cout << " " << index + 1 << " written" << endl;
319  mapFile.close();
320 
321  // test
322  //map_to_corr->WriteAsciiFile((fileName+".dat").c_str());
323 
324  // It is indeed annoying that I did not find any accessors to the field values to
325  // modify. So I have to write an ascii file, read it and save it as a root file.
326 
327  cout << " moving ascii file " << filename_out << " to the PANDAROOT input folder " << endl;
328  // move it to the input directory of PANDAROOT
329  string dir = getenv("VMCWORKDIR");
330  string moveto = dir + "/input/" + filename_out + ".dat";
331  system(("mv "+filename_out+".dat"+" "+moveto).c_str());
332  cout << " converting ASCII files to ROOT files " << endl;
333  PndFieldMap* fieldmap_converted = NULL;
334  //gfield = 0;
335  if (fieldtype == 2) {
336  fieldmap_converted = new PndSolenoidMap(filename_out.c_str(), "A");
337  //gfield = 3;
338  }
339  if (fieldtype == 3) {
340  fieldmap_converted = new PndDipoleMap(filename_out.c_str(), "A");
341  //gfield = 3;
342  }
343  if (fieldtype == 4) {
344  fieldmap_converted = new PndTransMap (filename_out.c_str(), "A");
345  //gfield = 4;
346  }
347  if (fieldmap_converted){
348  fieldmap_converted->Init();
349  stringstream _filename;
350  _filename << filename_out;// << "_" <<
351  fieldmap_converted->WriteRootFile((_filename.str()+".root").c_str(),_filename.str().c_str());
352  cout << " converted " << n_x << " x " << n_y << " x " << n_z << endl;
353  cout << " to " << fieldmap_converted->GetNx() << " x " << fieldmap_converted->GetNy() << " x " << fieldmap_converted->GetNz() << endl;
354  cout << " values " << endl;
355  /************* feature of displaying results ************************/
356  if (1){
357  static TApplication* myapp = NULL;
358  static TCanvas *canvas_fieldlines = NULL;
359  if (!myapp){
360  myapp = new TApplication("myapp",0,0);
361  /*
362  TGeoManager *geom = new TGeoManager("geom", "geo manager");
363  TGeoMaterial *something = new TGeoMaterial("something", 9.01, 4, 1.848);
364  something->SetUniqueID(1);
365  TGeoMedium *somemedium = new TGeoMedium("something", 3, 5, 0, 1, 10, 2, 0.1e11, 0.2, 0.2e-3, 0.2e-1);
366  TGeoVolume *top_vol = gGeoManager->MakeBox("cave", somemedium, 5, 5, 350);
367  top_vol->SetVisibility(0);
368  TGeoVolume *beam_pipe_vol = gGeoManager->MakeTube("beampipe",somemedium,2.0,2.1,5);
369  beam_pipe_vol->SetVisibility(0);
370  beam_pipe_vol->SetLineColor(5);
371  top_vol->AddNode(beam_pipe_vol, 1);
372  geom->SetTopVolume(top_vol);
373  gGeoManager->CloseGeometry();
374  canvas_fieldlines = new TCanvas("canvas_fieldlines", "field lines", 600, 600);
375  top_vol->Draw("ogl");
376  cout << " drawing done " << endl;
377  */
378  gROOT->Macro("$VMCWORKDIR/gconfig/rootlogon.C");
379  TString vmcWorkdir = gSystem->Getenv("VMCWORKDIR");
380  // materials and media
381  FairGeoLoader* geoLoad = new FairGeoLoader("TGeo", "FairGeoLoader");
382  FairGeoInterface* geoFace = geoLoad->getGeoInterface();
383  geoFace->setMediaFile(vmcWorkdir + "/geometry/media_pnd.geo");
384  geoFace->readMedia();
385  geoFace->print();
386  FairGeoMedia* geoMedia = geoFace->getMedia();
387  FairGeoMedium *FairMediumAir = geoMedia->getMedium("air");
388  FairGeoBuilder* geoBuild = geoLoad->getGeoBuilder();
389  Int_t nmed=geoBuild->createMedium(FairMediumAir);
390  TGeoVolume *top_vol = gGeoManager->MakeBox("cave", gGeoManager->GetMedium("air"), 400, 400, 1000);
391  top_vol->SetVisibility(0);
392  TFile* solenoid_file = new TFile("$VMCWORKDIR/geometry/FullSolenoid.root" ,"OPEN");
393  TGeoVolume * solenoid = (TGeoVolume*) solenoid_file->Get("topNode");
394  if (!solenoid){
395  cout << " 'Error: no solenoid file found to display " << endl;
396  } else {
397  solenoid->SetLineColor(12);
398  solenoid->SetTransparency(50);
399  top_vol->AddNode(solenoid, 1);
400  }
401  TFile* beampipe_file = new TFile("$VMCWORKDIR/geometry/beampipe_201309.root" ,"OPEN");
402  TGeoVolume * beampipe = (TGeoVolume*) beampipe_file->Get("pipeassembly");
403  if (!beampipe){
404  cout << " 'Error: no beam pipe file found to display " << endl;
405  } else {
406  beampipe->SetTransparency(50);
407  top_vol->AddNode(beampipe, 1);
408  }
409  gGeoManager->SetTopVolume(top_vol);
410  top_vol->Draw("ogl");
411  }
412  //static TGLViewer *view(NULL);
413  //static TView *view(NULL);
414  //if (!view) {
415  //view = TView::CreateView(1);
416  //view = (TGLViewer *)gPad->GetViewer3D();
417  //view->Draw("ogl");
418  //view->SetRange(x_min/10., y_min/10., z_min/10., x_max/10., y_max/10., z_max/10.);
419  //view->SetRange(-231, -231, -172, 231, 231, 343);
420  //view->SetRange(-5, -5, -172, 5, 5, 343);
421  //}
422  // store the end points of the previous visualization
423  static vector<double>xstart_prev;
424  static vector<double>ystart_prev;
425  static vector<double>zstart_prev;
426  static vector<int> line_color;
427  // initilize with a starting grid
428  if (xstart_prev.size()==0){
429  //const Int_t n = 100;
430  //create a raster for one plane to start the field lines from
431  for (double ix = -2.; ix < 2.; ix += .2){
432  for (double iy = -2.; iy < 2.; iy += .2){
433  double r2 = (ix*ix + iy*iy);
434  if (r2 > 4) continue;
435  xstart_prev.push_back(ix);
436  ystart_prev.push_back(iy);
437  zstart_prev.push_back(-50.);
438  int color = 12 + floor(7./4.*r2);
439  line_color.push_back(color);
440  }
441  }
442  }
443  // continue drawing from discontinued points
444  for (unsigned int istartpoint = 0; istartpoint < xstart_prev.size(); istartpoint++){
445  int npoints = 1000;
446  double stepsize = d_x/10.;
447  if (d_y/10. < stepsize) stepsize = d_y/10.;
448  if (d_z/10. < stepsize) stepsize = d_z/10.;
449  vector<double> xs;
450  vector<double> ys;
451  vector<double> zs;
452  xs.push_back(xstart_prev[istartpoint]);
453  ys.push_back(ystart_prev[istartpoint]);
454  zs.push_back(zstart_prev[istartpoint]);
455  //cout << endl;
456  for (int ipoint = 1; ipoint < npoints; ipoint++){
457  // get the previous position
458  double pval[3] = {xs[ipoint-1],ys[ipoint-1],zs[ipoint-1]};
459  // check the validity of the field map range
460  //cout << x_min/10. << " " << pval[0] << " " << x_max/10. << endl;
461  //cout << y_min/10. << " " << pval[1] << " " << y_max/10. << endl;
462  //cout << z_min/10. << " " << pval[2] << " " << z_max/10. << endl << endl;
463  if (!(x_min/10. < pval[0] && pval[0] < x_max/10.)||
464  !(y_min/10. < pval[1] && pval[1] < y_max/10.)||
465  !(z_min/10. < pval[2] && pval[2] < z_max/10.)){
466  //if (xs.size()>2)
467  //cout << " break due to out of range " << endl;
468  xstart_prev[istartpoint] = pval[0];
469  ystart_prev[istartpoint] = pval[1];
470  zstart_prev[istartpoint] = pval[2];
471  break;
472  }
473  // get the corresponding field value
474  // what is in fact the direction of the field line
475  double Bval[3];
476  fieldmap_converted->GetBxyz(pval, Bval);
477  // scale to the stepsize and add to the position
478  double fieldstrength = sqrt(Bval[0]*Bval[0]+Bval[1]*Bval[1]+Bval[2]*Bval[2]);
479  pval[0] = pval[0] + stepsize/fieldstrength*Bval[0];
480  pval[1] = pval[1] + stepsize/fieldstrength*Bval[1];
481  pval[2] = pval[2] + stepsize/fieldstrength*Bval[2];
482  // store the next point
483  xs.push_back(pval[0]);
484  ys.push_back(pval[1]);
485  zs.push_back(pval[2]);
486  }
487  // create a line out of the points
488  if (xs.size() > 2){
489  //cout << " creating polyline " << endl;
490  TPolyLine3D* fieldline = new TPolyLine3D(xs.size(), &xs[0], &ys[0], &zs[0]);
491  fieldline->SetLineColor(line_color[istartpoint]);
492  fieldline->Draw("same");
493  // mirror it
494  for (unsigned int ipoint = 0; ipoint < xs.size(); ipoint++){
495  xs[ipoint] = -xs[ipoint];
496  }
497  fieldline = new TPolyLine3D(xs.size(), &xs[0], &ys[0], &zs[0]);
498  fieldline->SetLineColor(line_color[istartpoint]);
499  fieldline->Draw("same");
500  for (unsigned int ipoint = 0; ipoint < xs.size(); ipoint++){
501  ys[ipoint] = -ys[ipoint];
502  }
503  fieldline = new TPolyLine3D(xs.size(), &xs[0], &ys[0], &zs[0]);
504  fieldline->SetLineColor(line_color[istartpoint]);
505  fieldline->Draw("same");
506  for (unsigned int ipoint = 0; ipoint < xs.size(); ipoint++){
507  xs[ipoint] = -xs[ipoint];
508  }
509  fieldline = new TPolyLine3D(xs.size(), &xs[0], &ys[0], &zs[0]);
510  fieldline->SetLineColor(line_color[istartpoint]);
511  fieldline->Draw("same");
512  }
513  }
514  static int counter(0);
515  counter++;
516  if (counter == 5){
517  myapp->Run();
518  }
519  }
520  /************* end of displaying results ****************************/
521  }
522 }
523 
525  //Load basic libraries
526  //gROOT->Macro("$VMCWORKDIR/gconfig/rootlogon.C");
527 
528  transform_jost_fieldmap("s1301_z-172_hsc.dat", "SolenoidMap1", 2);
529  transform_jost_fieldmap("s1301_z-40_hsc.dat", "SolenoidMap2", 2);
530  transform_jost_fieldmap("s1301_z180_hsc.dat", "SolenoidMap3", 2);
531  transform_jost_fieldmap("s1301_z248_hsc.dat", "SolenoidMap4", 2);
532  transform_jost_fieldmap("p1301_z283_hsc_0150.dat", "TransMap.0150", 4);
533  transform_jost_fieldmap("p1301_z283_hsc_0406.dat", "TransMap.0406", 4);
534 
535 }
536 
537 int main(int argc, char **argv) {
539  return 0;
540 }
541 
int main(int argc, char **argv)
virtual void Init()
double gBx
FairGeoLoader * geoLoad
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void WriteRootFile(const char *fileName, const char *mapName)
Int_t GetNy() const
Definition: PndFieldMap.h:98
double gBy
TGeoManager * gGeoManager
FairGeoMedium * FairMediumAir
const string name_modification
int counter
Definition: ZeeAnalysis.C:59
FairGeoMedia * geoMedia
Double_t
double gBz
TTree * gtree_values(NULL)
Int_t GetNx() const
Definition: PndFieldMap.h:97
Double_t z
double gDBz
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
Double_t x
double gDBx
void read_jost_fieldmaps()
double gDBy
FairGeoInterface * geoFace
Double_t y
Int_t GetNz() const
Definition: PndFieldMap.h:99
int Nint(float x)
Definition: PndCAMath.h:117
vector< string > fieldmapnames
double r2
FairGeoBuilder * geoBuild
void transform_jost_fieldmap(string filename_in, string filename_out, int fieldtype)