ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
HadrontherapyMagneticField3D.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file HadrontherapyMagneticField3D.cc
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // Hadrontherapy advanced example for Geant4
27 // See more at: https://twiki.cern.ch/twiki/bin/view/Geant4/AdvancedExamplesHadrontherapy
28 
30 #include "G4SystemOfUnits.hh"
31 #include "G4AutoLock.hh"
32 
33 namespace{ G4Mutex MyHadrontherapyLock=G4MUTEX_INITIALIZER; }
34 
35 using namespace std;
36 
38  :fXoffset(xOffset),invertX(false),invertY(false),invertZ(false)
39 {
40  //The format file is: X Y Z Ex Ey Ez
41 
42  double lenUnit= meter;
43  double fieldUnit= tesla;
44  G4cout << "\n-----------------------------------------------------------"
45  << "\n Magnetic field"
46  << "\n-----------------------------------------------------------";
47 
48 
49  G4cout << "\n ---> " "Reading the field grid from " << filename << " ... " << G4endl;
50  G4AutoLock lock(&MyHadrontherapyLock);
51 
52  ifstream file( filename ); // Open the file for reading.
53 
54  // Ignore first blank line
55  char buffer[256];
56  file.getline(buffer,256);
57 
58  // Read table dimensions
59  file >> nx >> ny >> nz; // Note dodgy order
60 
61  G4cout << " [ Number of values x,y,z: "
62  << nx << " " << ny << " " << nz << " ] "
63  << G4endl;
64 
65  // Set up storage space for table
66  xField.resize( nx );
67  yField.resize( nx );
68  zField.resize( nx );
69  int ix, iy, iz;
70  for (ix=0; ix<nx; ix++) {
71  xField[ix].resize(ny);
72  yField[ix].resize(ny);
73  zField[ix].resize(ny);
74  for (iy=0; iy<ny; iy++) {
75  xField[ix][iy].resize(nz);
76  yField[ix][iy].resize(nz);
77  zField[ix][iy].resize(nz);
78  }
79  }
80 
81  // Read in the data
82  G4double xval=0.;
83  G4double yval=0.;
84  G4double zval=0.;
85  G4double bx=0.;
86  G4double by=0.;
87  G4double bz=0.;
88  for (ix=0; ix<nx; ix++) {
89  for (iy=0; iy<ny; iy++) {
90  for (iz=0; iz<nz; iz++) {
91  file >> xval >> yval >> zval >> bx >> by >> bz ;
92  if ( ix==0 && iy==0 && iz==0 ) {
93  minx = xval * lenUnit;
94  miny = yval * lenUnit;
95  minz = zval * lenUnit;
96  }
97  xField[ix][iy][iz] = bx * fieldUnit;
98  yField[ix][iy][iz] = by * fieldUnit;
99  zField[ix][iy][iz] = bz * fieldUnit;
100  }
101  }
102  }
103  file.close();
104 
105  lock.unlock();
106 
107  maxx = xval * lenUnit;
108  maxy = yval * lenUnit;
109  maxz = zval * lenUnit;
110 
111  G4cout << "\n ---> ... done reading " << G4endl;
112 
113  // G4cout << " Read values of field from file " << filename << G4endl;
114  G4cout << " ---> assumed the order: x, y, z, Bx, By, Bz "
115  << "\n ---> Min values x,y,z: "
116  << minx/cm << " " << miny/cm << " " << minz/cm << " cm "
117  << "\n ---> Max values x,y,z: "
118  << maxx/cm << " " << maxy/cm << " " << maxz/cm << " cm "
119  << "\n ---> The field will be offset by " << xOffset/cm << " cm " << G4endl;
120 
121  // Should really check that the limits are not the wrong way around.
122  if (maxx < minx) {swap(maxx,minx); invertX = true;}
123  if (maxy < miny) {swap(maxy,miny); invertY = true;}
124  if (maxz < minz) {swap(maxz,minz); invertZ = true;}
125  G4cout << "\nAfter reordering if neccesary"
126  << "\n ---> Min values x,y,z: "
127  << minx/cm << " " << miny/cm << " " << minz/cm << " cm "
128  << " \n ---> Max values x,y,z: "
129  << maxx/cm << " " << maxy/cm << " " << maxz/cm << " cm ";
130 
131  dx = maxx - minx;
132  dy = maxy - miny;
133  dz = maxz - minz;
134  G4cout << "\n ---> Dif values x,y,z (range): "
135  << dx/cm << " " << dy/cm << " " << dz/cm << " cm in z "
136  << "\n-----------------------------------------------------------" << G4endl;
137 }
138 
140  double *Bfield ) const
141 {
142  double x = point[0]+ fXoffset;
143  double y = point[1];
144  double z = point[2];
145 
146  // Position of given point within region, normalized to the range
147  // [0,1]
148  double xfraction = (x - minx) / dx;
149  double yfraction = (y - miny) / dy;
150  double zfraction = (z - minz) / dz;
151 
152  if (invertX) { xfraction = 1 - xfraction;}
153  if (invertY) { yfraction = 1 - yfraction;}
154  if (invertZ) { zfraction = 1 - zfraction;}
155 
156  // Need addresses of these to pass to modf below.
157  // modf uses its second argument as an OUTPUT argument.
158  double xdindex, ydindex, zdindex;
159 
160  // Position of the point within the cuboid defined by the
161  // nearest surrounding tabulated points
162  double xlocal = ( std::modf(xfraction*(nx-1), &xdindex));
163  double ylocal = ( std::modf(yfraction*(ny-1), &ydindex));
164  double zlocal = ( std::modf(zfraction*(nz-1), &zdindex));
165 
166  // The indices of the nearest tabulated point whose coordinates
167  // are all less than those of the given point
168  int xindex = static_cast<int>(std::floor(xdindex));
169  int yindex = static_cast<int>(std::floor(ydindex));
170  int zindex = static_cast<int>(std::floor(zdindex));
171 
172  // Check that the point is within the defined region
173  if ((xindex < 0) || (xindex >= nx - 1) ||
174  (yindex < 0) || (yindex >= ny - 1) ||
175  (zindex < 0) || (zindex >= nz - 1))
176  {
177  Bfield[0] = 0.0;
178  Bfield[1] = 0.0;
179  Bfield[2] = 0.0;
180  }
181  else
182  {
183 
184 #ifdef DEBUG_INTERPOLATING_FIELD
185  G4cout << "Local x,y,z: " << xlocal << " " << ylocal << " " << zlocal << G4endl;
186  G4cout << "Index x,y,z: " << xindex << " " << yindex << " " << zindex << G4endl;
187  double valx0z0, mulx0z0, valx1z0, mulx1z0;
188  double valx0z1, mulx0z1, valx1z1, mulx1z1;
189  valx0z0= table[xindex ][0][zindex]; mulx0z0= (1-xlocal) * (1-zlocal);
190  valx1z0= table[xindex+1][0][zindex]; mulx1z0= xlocal * (1-zlocal);
191  valx0z1= table[xindex ][0][zindex+1]; mulx0z1= (1-xlocal) * zlocal;
192  valx1z1= table[xindex+1][0][zindex+1]; mulx1z1= xlocal * zlocal;
193 #endif
194 
195  // Full 3-dimensional version
196  Bfield[0] =
197  xField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
198  xField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
199  xField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
200  xField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
201  xField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
202  xField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
203  xField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
204  xField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
205 
206  Bfield[1] =
207  yField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
208  yField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
209  yField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
210  yField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
211  yField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
212  yField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
213  yField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
214  yField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
215 
216  Bfield[2] =
217  zField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
218  zField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
219  zField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
220  zField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
221  zField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
222  zField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
223  zField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
224  zField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
225  }
226 
227 //In order to obtain the output file with the magnetic components read from a particle passing in the magnetic field
228 /* std::ofstream MagneticField("MagneticField.out", std::ios::app);
229  MagneticField<< Bfield[0] << '\t' << " "
230  << Bfield[1] << '\t' << " "
231  << Bfield[2] << '\t' << " "
232  << point[0] << '\t' << " "
233  << point[1] << '\t' << " "
234  << point[2] << '\t' << " "
235  << std::endl;*/
236 
237 }