Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
HadrontherapyMagneticField3D Class Reference

#include <HadrontherapyMagneticField3D.hh>

Inheritance diagram for HadrontherapyMagneticField3D:
Collaboration diagram for HadrontherapyMagneticField3D:

Public Member Functions

 HadrontherapyMagneticField3D (const char *filename, double xOffset)
 
void GetFieldValue (const double Point[4], double *Bfield) const
 
- Public Member Functions inherited from G4MagneticField
 G4MagneticField ()
 
virtual ~G4MagneticField ()
 
 G4MagneticField (const G4MagneticField &r)
 
G4MagneticFieldoperator= (const G4MagneticField &p)
 
G4bool DoesFieldChangeEnergy () const
 
- Public Member Functions inherited from G4ElectroMagneticField
 G4ElectroMagneticField ()
 
virtual ~G4ElectroMagneticField ()
 
 G4ElectroMagneticField (const G4ElectroMagneticField &r)
 
G4ElectroMagneticFieldoperator= (const G4ElectroMagneticField &p)
 
- Public Member Functions inherited from G4Field
 G4Field (G4bool gravityOn=false)
 
 G4Field (const G4Field &)
 
virtual ~G4Field ()
 
G4Fieldoperator= (const G4Field &p)
 
G4bool IsGravityActive () const
 
void SetGravityActive (G4bool OnOffFlag)
 
virtual G4FieldClone () const
 

Detailed Description

Definition at line 40 of file HadrontherapyMagneticField3D.hh.

Constructor & Destructor Documentation

HadrontherapyMagneticField3D::HadrontherapyMagneticField3D ( const char *  filename,
double  xOffset 
)

Definition at line 35 of file HadrontherapyMagneticField3D.cc.

36  :fXoffset(xOffset),invertX(false),invertY(false),invertZ(false)
37 {
38  //The format file is: X Y Z Ex Ey Ez
39 
40  double lenUnit= meter;
41  double fieldUnit= tesla;
42  G4cout << "\n-----------------------------------------------------------"
43  << "\n Magnetic field"
44  << "\n-----------------------------------------------------------";
45 
46 
47  G4cout << "\n ---> " "Reading the field grid from " << filename << " ... " << endl;
48  G4AutoLock lock(&MyHadrontherapyLock);
49 
50  ifstream file( filename ); // Open the file for reading.
51 
52  // Ignore first blank line
53  char buffer[256];
54  file.getline(buffer,256);
55 
56  // Read table dimensions
57  file >> nx >> ny >> nz; // Note dodgy order
58 
59  G4cout << " [ Number of values x,y,z: "
60  << nx << " " << ny << " " << nz << " ] "
61  << endl;
62 
63  // Set up storage space for table
64  xField.resize( nx );
65  yField.resize( nx );
66  zField.resize( nx );
67  int ix, iy, iz;
68  for (ix=0; ix<nx; ix++) {
69  xField[ix].resize(ny);
70  yField[ix].resize(ny);
71  zField[ix].resize(ny);
72  for (iy=0; iy<ny; iy++) {
73  xField[ix][iy].resize(nz);
74  yField[ix][iy].resize(nz);
75  zField[ix][iy].resize(nz);
76  }
77  }
78 
79  // Read in the data
80  G4double xval=0.;
81  G4double yval=0.;
82  G4double zval=0.;
83  G4double bx=0.;
84  G4double by=0.;
85  G4double bz=0.;
86  for (ix=0; ix<nx; ix++) {
87  for (iy=0; iy<ny; iy++) {
88  for (iz=0; iz<nz; iz++) {
89  file >> xval >> yval >> zval >> bx >> by >> bz ;
90  if ( ix==0 && iy==0 && iz==0 ) {
91  minx = xval * lenUnit;
92  miny = yval * lenUnit;
93  minz = zval * lenUnit;
94  }
95  xField[ix][iy][iz] = bx * fieldUnit;
96  yField[ix][iy][iz] = by * fieldUnit;
97  zField[ix][iy][iz] = bz * fieldUnit;
98  }
99  }
100  }
101  file.close();
102 
103  lock.unlock();
104 
105  maxx = xval * lenUnit;
106  maxy = yval * lenUnit;
107  maxz = zval * lenUnit;
108 
109  G4cout << "\n ---> ... done reading " << endl;
110 
111  // G4cout << " Read values of field from file " << filename << endl;
112  G4cout << " ---> assumed the order: x, y, z, Bx, By, Bz "
113  << "\n ---> Min values x,y,z: "
114  << minx/cm << " " << miny/cm << " " << minz/cm << " cm "
115  << "\n ---> Max values x,y,z: "
116  << maxx/cm << " " << maxy/cm << " " << maxz/cm << " cm "
117  << "\n ---> The field will be offset by " << xOffset/cm << " cm " << endl;
118 
119  // Should really check that the limits are not the wrong way around.
120  if (maxx < minx) {swap(maxx,minx); invertX = true;}
121  if (maxy < miny) {swap(maxy,miny); invertY = true;}
122  if (maxz < minz) {swap(maxz,minz); invertZ = true;}
123  G4cout << "\nAfter reordering if neccesary"
124  << "\n ---> Min values x,y,z: "
125  << minx/cm << " " << miny/cm << " " << minz/cm << " cm "
126  << " \n ---> Max values x,y,z: "
127  << maxx/cm << " " << maxy/cm << " " << maxz/cm << " cm ";
128 
129  dx = maxx - minx;
130  dy = maxy - miny;
131  dz = maxz - minz;
132  G4cout << "\n ---> Dif values x,y,z (range): "
133  << dx/cm << " " << dy/cm << " " << dz/cm << " cm in z "
134  << "\n-----------------------------------------------------------" << endl;
135 }
static constexpr double tesla
Definition: G4SIunits.hh:268
#define buffer
Definition: xmlparse.cc:628
static constexpr double meter
Definition: G4SIunits.hh:82
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Member Function Documentation

void HadrontherapyMagneticField3D::GetFieldValue ( const double  Point[4],
double *  Bfield 
) const
virtual

Implements G4MagneticField.

Definition at line 137 of file HadrontherapyMagneticField3D.cc.

139 {
140  double x = point[0]+ fXoffset;
141  double y = point[1];
142  double z = point[2];
143 
144  // Position of given point within region, normalized to the range
145  // [0,1]
146  double xfraction = (x - minx) / dx;
147  double yfraction = (y - miny) / dy;
148  double zfraction = (z - minz) / dz;
149 
150  if (invertX) { xfraction = 1 - xfraction;}
151  if (invertY) { yfraction = 1 - yfraction;}
152  if (invertZ) { zfraction = 1 - zfraction;}
153 
154  // Need addresses of these to pass to modf below.
155  // modf uses its second argument as an OUTPUT argument.
156  double xdindex, ydindex, zdindex;
157 
158  // Position of the point within the cuboid defined by the
159  // nearest surrounding tabulated points
160  double xlocal = ( std::modf(xfraction*(nx-1), &xdindex));
161  double ylocal = ( std::modf(yfraction*(ny-1), &ydindex));
162  double zlocal = ( std::modf(zfraction*(nz-1), &zdindex));
163 
164  // The indices of the nearest tabulated point whose coordinates
165  // are all less than those of the given point
166  int xindex = static_cast<int>(std::floor(xdindex));
167  int yindex = static_cast<int>(std::floor(ydindex));
168  int zindex = static_cast<int>(std::floor(zdindex));
169 
170  // Check that the point is within the defined region
171  if ((xindex < 0) || (xindex >= nx - 1) ||
172  (yindex < 0) || (yindex >= ny - 1) ||
173  (zindex < 0) || (zindex >= nz - 1))
174  {
175  Bfield[0] = 0.0;
176  Bfield[1] = 0.0;
177  Bfield[2] = 0.0;
178  }
179  else
180  {
181 
182 #ifdef DEBUG_INTERPOLATING_FIELD
183  G4cout << "Local x,y,z: " << xlocal << " " << ylocal << " " << zlocal << endl;
184  G4cout << "Index x,y,z: " << xindex << " " << yindex << " " << zindex << endl;
185  double valx0z0, mulx0z0, valx1z0, mulx1z0;
186  double valx0z1, mulx0z1, valx1z1, mulx1z1;
187  valx0z0= table[xindex ][0][zindex]; mulx0z0= (1-xlocal) * (1-zlocal);
188  valx1z0= table[xindex+1][0][zindex]; mulx1z0= xlocal * (1-zlocal);
189  valx0z1= table[xindex ][0][zindex+1]; mulx0z1= (1-xlocal) * zlocal;
190  valx1z1= table[xindex+1][0][zindex+1]; mulx1z1= xlocal * zlocal;
191 #endif
192 
193  // Full 3-dimensional version
194  Bfield[0] =
195  xField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
196  xField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
197  xField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
198  xField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
199  xField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
200  xField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
201  xField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
202  xField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
203 
204  Bfield[1] =
205  yField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
206  yField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
207  yField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
208  yField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
209  yField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
210  yField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
211  yField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
212  yField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
213 
214  Bfield[2] =
215  zField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
216  zField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
217  zField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
218  zField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
219  zField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
220  zField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
221  zField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
222  zField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
223  }
224 
225 //In order to obtain the output file with the magnetic components read from a particle passing in the magnetic field
226 /* std::ofstream MagneticField("MagneticField.out", std::ios::app);
227  MagneticField<< Bfield[0] << '\t' << " "
228  << Bfield[1] << '\t' << " "
229  << Bfield[2] << '\t' << " "
230  << point[0] << '\t' << " "
231  << point[1] << '\t' << " "
232  << point[2] << '\t' << " "
233  << G4endl;*/
234 
235 }
tuple x
Definition: test.py:50
G4GLOB_DLL std::ostream G4cout
tuple z
Definition: test.py:28

The documentation for this class was generated from the following files: