Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PurgMagTabulatedField3D.cc
Go to the documentation of this file.
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 // Code developed by:
27 // S.Larsson and J. Generowicz.
28 //
29 // *************************************
30 // * *
31 // * PurgMagTabulatedField3D.cc *
32 // * *
33 // *************************************
34 //
35 // $Id$
36 //
37 
39 #include "G4SystemOfUnits.hh"
40 
41 PurgMagTabulatedField3D::PurgMagTabulatedField3D( const char* filename, double zOffset )
42  :fZoffset(zOffset),invertX(false),invertY(false),invertZ(false)
43 {
44 
45  double lenUnit= meter;
46  double fieldUnit= tesla;
47  G4cout << "\n-----------------------------------------------------------"
48  << "\n Magnetic field"
49  << "\n-----------------------------------------------------------";
50 
51  G4cout << "\n ---> " "Reading the field grid from " << filename << " ... " << endl;
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  << endl;
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  // Ignore other header information
82  // The first line whose second character is '0' is considered to
83  // be the last line of the header.
84  do {
85  file.getline(buffer,256);
86  } while ( buffer[1]!='0');
87 
88  // Read in the data
89  double xval,yval,zval,bx,by,bz;
90  double permeability; // Not used in this example.
91  for (ix=0; ix<nx; ix++) {
92  for (iy=0; iy<ny; iy++) {
93  for (iz=0; iz<nz; iz++) {
94  file >> xval >> yval >> zval >> bx >> by >> bz >> permeability;
95  if ( ix==0 && iy==0 && iz==0 ) {
96  minx = xval * lenUnit;
97  miny = yval * lenUnit;
98  minz = zval * lenUnit;
99  }
100  xField[ix][iy][iz] = bx * fieldUnit;
101  yField[ix][iy][iz] = by * fieldUnit;
102  zField[ix][iy][iz] = bz * fieldUnit;
103  }
104  }
105  }
106  file.close();
107 
108  maxx = xval * lenUnit;
109  maxy = yval * lenUnit;
110  maxz = zval * lenUnit;
111 
112  G4cout << "\n ---> ... done reading " << endl;
113 
114  // G4cout << " Read values of field from file " << filename << endl;
115  G4cout << " ---> assumed the order: x, y, z, Bx, By, Bz "
116  << "\n ---> Min values x,y,z: "
117  << minx/cm << " " << miny/cm << " " << minz/cm << " cm "
118  << "\n ---> Max values x,y,z: "
119  << maxx/cm << " " << maxy/cm << " " << maxz/cm << " cm "
120  << "\n ---> The field will be offset by " << zOffset/cm << " cm " << endl;
121 
122  // Should really check that the limits are not the wrong way around.
123  if (maxx < minx) {swap(maxx,minx); invertX = true;}
124  if (maxy < miny) {swap(maxy,miny); invertY = true;}
125  if (maxz < minz) {swap(maxz,minz); invertZ = true;}
126  G4cout << "\nAfter reordering if neccesary"
127  << "\n ---> Min values x,y,z: "
128  << minx/cm << " " << miny/cm << " " << minz/cm << " cm "
129  << " \n ---> Max values x,y,z: "
130  << maxx/cm << " " << maxy/cm << " " << maxz/cm << " cm ";
131 
132  dx = maxx - minx;
133  dy = maxy - miny;
134  dz = maxz - minz;
135  G4cout << "\n ---> Dif values x,y,z (range): "
136  << dx/cm << " " << dy/cm << " " << dz/cm << " cm in z "
137  << "\n-----------------------------------------------------------" << endl;
138 }
139 
140 void PurgMagTabulatedField3D::GetFieldValue(const double point[4],
141  double *Bfield ) const
142 {
143 
144  double x = point[0];
145  double y = point[1];
146  double z = point[2] + fZoffset;
147 
148  // Check that the point is within the defined region
149  if ( x>=minx && x<=maxx &&
150  y>=miny && y<=maxy &&
151  z>=minz && z<=maxz ) {
152 
153  // Position of given point within region, normalized to the range
154  // [0,1]
155  double xfraction = (x - minx) / dx;
156  double yfraction = (y - miny) / dy;
157  double zfraction = (z - minz) / dz;
158 
159  if (invertX) { xfraction = 1 - xfraction;}
160  if (invertY) { yfraction = 1 - yfraction;}
161  if (invertZ) { zfraction = 1 - zfraction;}
162 
163  // Need addresses of these to pass to modf below.
164  // modf uses its second argument as an OUTPUT argument.
165  double xdindex, ydindex, zdindex;
166 
167  // Position of the point within the cuboid defined by the
168  // nearest surrounding tabulated points
169  double xlocal = ( std::modf(xfraction*(nx-1), &xdindex));
170  double ylocal = ( std::modf(yfraction*(ny-1), &ydindex));
171  double zlocal = ( std::modf(zfraction*(nz-1), &zdindex));
172 
173  // The indices of the nearest tabulated point whose coordinates
174  // are all less than those of the given point
175  int xindex = static_cast<int>(xdindex);
176  int yindex = static_cast<int>(ydindex);
177  int zindex = static_cast<int>(zdindex);
178 
179 
180 #ifdef DEBUG_INTERPOLATING_FIELD
181  G4cout << "Local x,y,z: " << xlocal << " " << ylocal << " " << zlocal << endl;
182  G4cout << "Index x,y,z: " << xindex << " " << yindex << " " << zindex << endl;
183  double valx0z0, mulx0z0, valx1z0, mulx1z0;
184  double valx0z1, mulx0z1, valx1z1, mulx1z1;
185  valx0z0= table[xindex ][0][zindex]; mulx0z0= (1-xlocal) * (1-zlocal);
186  valx1z0= table[xindex+1][0][zindex]; mulx1z0= xlocal * (1-zlocal);
187  valx0z1= table[xindex ][0][zindex+1]; mulx0z1= (1-xlocal) * zlocal;
188  valx1z1= table[xindex+1][0][zindex+1]; mulx1z1= xlocal * zlocal;
189 #endif
190 
191  // Full 3-dimensional version
192  Bfield[0] =
193  xField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
194  xField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
195  xField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
196  xField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
197  xField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
198  xField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
199  xField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
200  xField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
201  Bfield[1] =
202  yField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
203  yField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
204  yField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
205  yField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
206  yField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
207  yField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
208  yField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
209  yField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
210  Bfield[2] =
211  zField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
212  zField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
213  zField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
214  zField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
215  zField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
216  zField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
217  zField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
218  zField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
219 
220  } else {
221  Bfield[0] = 0.0;
222  Bfield[1] = 0.0;
223  Bfield[2] = 0.0;
224  }
225 }
226