257 lines
10 KiB
C++
257 lines
10 KiB
C++
// Geant4 simulation for MuSR
|
|
// AUTHOR: Toni SHIROKA, Paul Scherrer Institut, PSI
|
|
// DATE : 2008-05
|
|
//
|
|
|
|
#include "lem4TabulatedField3D.hh"
|
|
|
|
|
|
lem4TabulatedField3D* lem4TabulatedField3D::pointerToTabulatedField3D=NULL;
|
|
lem4TabulatedField3D* lem4TabulatedField3D::GetInstance() {
|
|
return pointerToTabulatedField3D;
|
|
}
|
|
|
|
|
|
lem4TabulatedField3D::lem4TabulatedField3D( const char* filename, double fieldValue )
|
|
:ffieldValue(fieldValue),invertX(false),invertY(false),invertZ(false)
|
|
{
|
|
pointerToTabulatedField3D=this;
|
|
double lenUnit= meter;
|
|
// double fieldUnit= tesla;
|
|
G4cout << "\n-----------------------------------------------------------"
|
|
<< "\n Magnetic field"
|
|
<< "\n-----------------------------------------------------------"
|
|
<< endl;
|
|
G4cout << " tesla="<<tesla<<G4endl;
|
|
G4cout << " Field set to "<< fieldValue/tesla << " T"<< G4endl;
|
|
G4cout << "\n ---> " "Reading the field grid from " << filename << " ... " << G4endl;
|
|
ifstream file( filename ); // Open the file for reading.
|
|
|
|
// Ignore first blank line
|
|
char buffer[256];
|
|
file.getline(buffer,256);
|
|
|
|
// Read table dimensions
|
|
file >> nx >> ny >> nz; // Note dodgy order
|
|
|
|
G4cout << " [ Number of values x,y,z: "
|
|
<< nx << " " << ny << " " << nz << " ] "
|
|
<< endl;
|
|
|
|
// Set up storage space for table
|
|
xField.resize( nx );
|
|
yField.resize( nx );
|
|
zField.resize( nx );
|
|
int ix, iy, iz;
|
|
for (ix=0; ix<nx; ix++) {
|
|
xField[ix].resize(ny);
|
|
yField[ix].resize(ny);
|
|
zField[ix].resize(ny);
|
|
for (iy=0; iy<ny; iy++) {
|
|
xField[ix][iy].resize(nz);
|
|
yField[ix][iy].resize(nz);
|
|
zField[ix][iy].resize(nz);
|
|
}
|
|
}
|
|
|
|
// Ignore other header information
|
|
// The first line whose second character is '0' is considered to
|
|
// be the last line of the header.
|
|
do {
|
|
file.getline(buffer,256);
|
|
} while ( buffer[1]!='0');
|
|
|
|
// Read in the data
|
|
double xval,yval,zval,bx,by,bz;
|
|
double permeability; // Not used in this example.
|
|
for (ix=0; ix<nx; ix++) {
|
|
for (iy=0; iy<ny; iy++) {
|
|
for (iz=0; iz<nz; iz++) {
|
|
file >> xval >> yval >> zval >> bx >> by >> bz >> permeability;
|
|
if ( ix==0 && iy==0 && iz==0 ) {
|
|
minx = xval * lenUnit;
|
|
miny = yval * lenUnit;
|
|
minz = zval * lenUnit;
|
|
}
|
|
//xField[ix][iy][iz] = bx * fieldUnit * fieldValue;
|
|
//yField[ix][iy][iz] = by * fieldUnit * fieldValue;
|
|
//zField[ix][iy][iz] = bz * fieldUnit * fieldValue;
|
|
// xField[ix][iy][iz] = bx * fieldValue;
|
|
// yField[ix][iy][iz] = by * fieldValue;
|
|
// zField[ix][iy][iz] = bz * fieldValue;
|
|
xField[ix][iy][iz] = bx;
|
|
yField[ix][iy][iz] = by;
|
|
zField[ix][iy][iz] = bz;
|
|
}
|
|
}
|
|
}
|
|
file.close();
|
|
|
|
maxx = xval * lenUnit;
|
|
maxy = yval * lenUnit;
|
|
maxz = zval * lenUnit;
|
|
|
|
G4cout << "\n ---> ... done reading " << endl;
|
|
|
|
// G4cout << " Read values of field from file " << filename << endl;
|
|
G4cout << " ---> assumed the order: x, y, z, Bx, By, Bz "
|
|
<< "\n ---> Min values x,y,z: "
|
|
<< minx/cm << " " << miny/cm << " " << minz/cm << " cm "
|
|
<< "\n ---> Max values x,y,z: "
|
|
<< maxx/cm << " " << maxy/cm << " " << maxz/cm << " cm " << endl;
|
|
|
|
|
|
// Should really check that the limits are not the wrong way around.
|
|
if (maxx < minx) {swap(maxx,minx); invertX = true;}
|
|
if (maxy < miny) {swap(maxy,miny); invertY = true;}
|
|
if (maxz < minz) {swap(maxz,minz); invertZ = true;}
|
|
G4cout << "\nAfter reordering if neccesary"
|
|
<< "\n ---> Min values x,y,z: "
|
|
<< minx/cm << " " << miny/cm << " " << minz/cm << " cm "
|
|
<< " \n ---> Max values x,y,z: "
|
|
<< maxx/cm << " " << maxy/cm << " " << maxz/cm << " cm ";
|
|
|
|
dx = maxx - minx;
|
|
dy = maxy - miny;
|
|
dz = maxz - minz;
|
|
G4cout << "\n ---> Dif values x,y,z (range): "
|
|
<< dx/cm << " " << dy/cm << " " << dz/cm << " cm in z "
|
|
<< "\n-----------------------------------------------------------" << endl;
|
|
}
|
|
|
|
void lem4TabulatedField3D::GetFieldValue(const double point[4],
|
|
double *Bfield ) const
|
|
{
|
|
Bfield[0]=0.; Bfield[1]=0.; Bfield[2]=0.;
|
|
|
|
double x = point[0]+positionOffset[0];
|
|
double y = point[1]+positionOffset[1];
|
|
double z = point[2]+positionOffset[2];
|
|
|
|
|
|
// Check that the point is within the defined region
|
|
if ( x>minx && x<maxx &&
|
|
y>miny && y<maxy &&
|
|
z>minz && z<maxz ) {
|
|
|
|
// Position of given point within region, normalized to the range
|
|
// [0,1]
|
|
double xfraction = (x - minx) / dx;
|
|
double yfraction = (y - miny) / dy;
|
|
double zfraction = (z - minz) / dz;
|
|
|
|
if (invertX) { xfraction = 1 - xfraction;}
|
|
if (invertY) { yfraction = 1 - yfraction;}
|
|
if (invertZ) { zfraction = 1 - zfraction;}
|
|
|
|
// Need addresses of these to pass to modf below.
|
|
// modf uses its second argument as an OUTPUT argument.
|
|
double xdindex, ydindex, zdindex;
|
|
|
|
// Position of the point within the cuboid defined by the
|
|
// nearest surrounding tabulated points
|
|
double xlocal = ( modf(xfraction*(nx-1), &xdindex));
|
|
double ylocal = ( modf(yfraction*(ny-1), &ydindex));
|
|
double zlocal = ( modf(zfraction*(nz-1), &zdindex));
|
|
|
|
// The indices of the nearest tabulated point whose coordinates
|
|
// are all less than those of the given point
|
|
int xindex = static_cast<int>(xdindex);
|
|
int yindex = static_cast<int>(ydindex);
|
|
int zindex = static_cast<int>(zdindex);
|
|
|
|
//cks The following check is necessary - even though xindex and zindex should never be out of range,
|
|
// it may happen (due to some rounding error ?). It is better to leave the check here.
|
|
if ((xindex<0)||(xindex>(nx-2))) {
|
|
std::cout<<"SERIOUS PROBLEM: xindex out of range! xindex="<<xindex<<" x="<<x<<" xfraction="<<xfraction<<std::endl;
|
|
if (xindex<0) xindex=0;
|
|
else xindex=nx-2;
|
|
}
|
|
if ((yindex<0)||(yindex>(ny-2))) {
|
|
std::cout<<"SERIOUS PROBLEM: yindex out of range! yindex="<<yindex<<" y="<<y<<" yfraction="<<yfraction<<std::endl;
|
|
if (yindex<0) yindex=0;
|
|
else yindex=ny-2;
|
|
}
|
|
if ((zindex<0)||(zindex>(nz-2))) {
|
|
std::cout<<"SERIOUS PROBLEM: zindex out of range! zindex="<<zindex<<" z="<<z<<" zfraction="<<zfraction<<std::endl;
|
|
if (zindex<0) zindex=0;
|
|
else zindex=nz-2;
|
|
}
|
|
|
|
#ifdef DEBUG_INTERPOLATING_FIELD
|
|
G4cout << "Local x,y,z: " << xlocal << " " << ylocal << " " << zlocal << endl;
|
|
G4cout << "Index x,y,z: " << xindex << " " << yindex << " " << zindex << endl;
|
|
double valx0z0, mulx0z0, valx1z0, mulx1z0;
|
|
double valx0z1, mulx0z1, valx1z1, mulx1z1;
|
|
valx0z0= table[xindex ][0][zindex]; mulx0z0= (1-xlocal) * (1-zlocal);
|
|
valx1z0= table[xindex+1][0][zindex]; mulx1z0= xlocal * (1-zlocal);
|
|
valx0z1= table[xindex ][0][zindex+1]; mulx0z1= (1-xlocal) * zlocal;
|
|
valx1z1= table[xindex+1][0][zindex+1]; mulx1z1= xlocal * zlocal;
|
|
#endif
|
|
|
|
// Full 3-dimensional version
|
|
Bfield[0] =
|
|
xField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
|
|
xField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
|
|
xField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
|
|
xField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
|
|
xField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
|
|
xField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
|
|
xField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
|
|
xField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
|
|
Bfield[1] =
|
|
yField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
|
|
yField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
|
|
yField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
|
|
yField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
|
|
yField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
|
|
yField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
|
|
yField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
|
|
yField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
|
|
Bfield[2] =
|
|
zField[xindex ][yindex ][zindex ] * (1-xlocal) * (1-ylocal) * (1-zlocal) +
|
|
zField[xindex ][yindex ][zindex+1] * (1-xlocal) * (1-ylocal) * zlocal +
|
|
zField[xindex ][yindex+1][zindex ] * (1-xlocal) * ylocal * (1-zlocal) +
|
|
zField[xindex ][yindex+1][zindex+1] * (1-xlocal) * ylocal * zlocal +
|
|
zField[xindex+1][yindex ][zindex ] * xlocal * (1-ylocal) * (1-zlocal) +
|
|
zField[xindex+1][yindex ][zindex+1] * xlocal * (1-ylocal) * zlocal +
|
|
zField[xindex+1][yindex+1][zindex ] * xlocal * ylocal * (1-zlocal) +
|
|
zField[xindex+1][yindex+1][zindex+1] * xlocal * ylocal * zlocal ;
|
|
|
|
Bfield[0] = Bfield[0] * ffieldValue;
|
|
Bfield[1] = Bfield[1] * ffieldValue;
|
|
Bfield[2] = Bfield[2] * ffieldValue;
|
|
|
|
}
|
|
|
|
if (sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1]+Bfield[2]*Bfield[2])<0.00001*tesla) {
|
|
// Bfield[0] = 0.0;
|
|
// Bfield[1] = 0.0;
|
|
Bfield[2] = 0.00001*tesla;
|
|
}
|
|
// if ((point[2]>-1*cm)&&(point[2]<1*cm)) {
|
|
// G4cout<<"Field at point ("<<point[0]<<","<<point[1]<<","<<point[2]<<","<<point[3]<<" = "
|
|
// << Bfield[0]/tesla<<", "<<Bfield[1]/tesla<<", "<<Bfield[2]/tesla<<G4endl;
|
|
// }
|
|
}
|
|
|
|
G4double lem4TabulatedField3D::GetFieldSetValue() {
|
|
return ffieldValue;
|
|
}
|
|
|
|
void lem4TabulatedField3D::SetFieldValue(double newFieldValue) {
|
|
// // Rescale the magnetic field for a new value of the magnetic field
|
|
// G4double factor = newFieldValue/ffieldValue;
|
|
// for (G4int ix=0; ix<nx; ix++) {
|
|
// for (G4int iy=0; iy<ny; iy++) {
|
|
// for (G4int iz=0; iz<nz; iz++) {
|
|
// xField[ix][iy][iz] = xField[ix][iy][iz] * factor;
|
|
// yField[ix][iy][iz] = yField[ix][iy][iz] * factor;
|
|
// zField[ix][iy][iz] = zField[ix][iy][iz] * factor;
|
|
// }
|
|
// }
|
|
// }
|
|
ffieldValue=newFieldValue;
|
|
G4cout<<"lem4TabulatedField3D.cc: ffieldValue changed to="<< ffieldValue/tesla<<"T "<<G4endl;
|
|
}
|