Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4Physics2DVector Class Reference

#include <G4Physics2DVector.hh>

Public Member Functions

 G4Physics2DVector ()
 
 G4Physics2DVector (size_t nx, size_t ny)
 
 G4Physics2DVector (const G4Physics2DVector &)
 
G4Physics2DVectoroperator= (const G4Physics2DVector &)
 
 ~G4Physics2DVector ()
 
G4double Value (G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
 
G4double Value (G4double x, G4double y) const
 
void PutX (size_t idx, G4double value)
 
void PutY (size_t idy, G4double value)
 
void PutValue (size_t idx, size_t idy, G4double value)
 
void PutVectors (const std::vector< G4double > &vecX, const std::vector< G4double > &vecY)
 
void ScaleVector (G4double factor)
 
G4double FindLinearX (G4double rand, G4double y, size_t &lastidy) const
 
G4double FindLinearX (G4double rand, G4double y) const
 
G4double GetX (size_t index) const
 
G4double GetY (size_t index) const
 
G4double GetValue (size_t idx, size_t idy) const
 
size_t FindBinLocationX (G4double x, size_t lastidx) const
 
size_t FindBinLocationY (G4double y, size_t lastidy) const
 
size_t GetLengthX () const
 
size_t GetLengthY () const
 
G4PhysicsVectorType GetType () const
 
void SetBicubicInterpolation (G4bool)
 
void Store (std::ofstream &fOut)
 
G4bool Retrieve (std::ifstream &fIn)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

void PrepareVectors ()
 
void ClearVectors ()
 
void CopyData (const G4Physics2DVector &vec)
 
G4double BicubicInterpolation (G4double x, G4double y, size_t idx, size_t idy) const
 
size_t FindBinLocation (G4double z, const G4PV2DDataVector &) const
 
size_t FindBin (G4double z, const G4PV2DDataVector &, size_t idz, size_t idzmax) const
 

Detailed Description

Definition at line 62 of file G4Physics2DVector.hh.

Constructor & Destructor Documentation

G4Physics2DVector::G4Physics2DVector ( )

Definition at line 47 of file G4Physics2DVector.cc.

48  : type(T_G4PhysicsFreeVector),
49  numberOfXNodes(0), numberOfYNodes(0),
50  verboseLevel(0), useBicubic(false)
51 {}
G4Physics2DVector::G4Physics2DVector ( size_t  nx,
size_t  ny 
)

Definition at line 55 of file G4Physics2DVector.cc.

References PrepareVectors().

56  : type(T_G4PhysicsFreeVector),
57  numberOfXNodes(nx), numberOfYNodes(ny),
58  verboseLevel(0), useBicubic(false)
59 {
61 }
G4Physics2DVector::G4Physics2DVector ( const G4Physics2DVector right)

Definition at line 72 of file G4Physics2DVector.cc.

References CopyData(), and PrepareVectors().

73 {
74  type = right.type;
75 
76  numberOfXNodes = right.numberOfXNodes;
77  numberOfYNodes = right.numberOfYNodes;
78 
79  verboseLevel = right.verboseLevel;
80  useBicubic = right.useBicubic;
81 
82  xVector = right.xVector;
83  yVector = right.yVector;
84 
86  CopyData(right);
87 }
void CopyData(const G4Physics2DVector &vec)
G4Physics2DVector::~G4Physics2DVector ( )

Definition at line 65 of file G4Physics2DVector.cc.

References ClearVectors().

66 {
67  ClearVectors();
68 }

Member Function Documentation

G4double G4Physics2DVector::BicubicInterpolation ( G4double  x,
G4double  y,
size_t  idx,
size_t  idy 
) const
protected

Definition at line 193 of file G4Physics2DVector.cc.

References GetValue().

Referenced by Value().

195 {
196  // Bicubic interpolation according to
197  // 1. H.M. Antia, "Numerical Methods for Scientists and Engineers",
198  // MGH, 1991.
199  // 2. W.H. Press et al., "Numerical recipes. The Art of Scientific
200  // Computing", Cambridge University Press, 2007.
201  G4double x1 = xVector[idx];
202  G4double x2 = xVector[idx+1];
203  G4double y1 = yVector[idy];
204  G4double y2 = yVector[idy+1];
205  G4double f1 = GetValue(idx, idy);
206  G4double f2 = GetValue(idx+1, idy);
207  G4double f3 = GetValue(idx+1, idy+1);
208  G4double f4 = GetValue(idx, idy+1);
209 
210  G4double dx = x2 - x1;
211  G4double dy = y2 - y1;
212 
213  G4double h1 = (x - x1)/dx;
214  G4double h2 = (y - y1)/dy;
215 
216  G4double h12 = h1*h1;
217  G4double h13 = h12*h1;
218  G4double h22 = h2*h2;
219  G4double h23 = h22*h2;
220 
221  // Three derivatives at each of four points (1-4) defining the
222  // subregion are computed by numerical centered differencing from
223  // the functional values already tabulated on the grid.
224 
225  G4double f1x = DerivativeX(idx, idy, dx);
226  G4double f2x = DerivativeX(idx+1, idy, dx);
227  G4double f3x = DerivativeX(idx+1, idy+1, dx);
228  G4double f4x = DerivativeX(idx, idy+1, dx);
229 
230  G4double f1y = DerivativeY(idx, idy, dy);
231  G4double f2y = DerivativeY(idx+1, idy, dy);
232  G4double f3y = DerivativeY(idx+1, idy+1, dy);
233  G4double f4y = DerivativeY(idx, idy+1, dy);
234 
235  G4double dxy = dx*dy;
236  G4double f1xy = DerivativeXY(idx, idy, dxy);
237  G4double f2xy = DerivativeXY(idx+1, idy, dxy);
238  G4double f3xy = DerivativeXY(idx+1, idy+1, dxy);
239  G4double f4xy = DerivativeXY(idx, idy+1, dxy);
240 
241  return
242  f1 + f1y*h2 + (3*(f4-f1) - 2*f1y - f4y)*h22 + (2*(f1 - f4) + f1y + f4y)*h23
243  + f1x*h1 + f1xy*h1*h2 +(3*(f4x - f1x) - 2*f1xy - f4xy)*h1*h22
244  + (2*(f1x - f4x) + f1xy + f4xy)*h1*h23
245  + (3*(f2 - f1) - 2*f1x - f2x)*h12 + (3*f2y - 3*f1y - 2*f1xy - f2xy)*h12*h2
246  + (9*(f1 - f2 + f3 - f4) + 6*f1x + 3*f2x - 3*f3x - 6*f4x + 6*f1y - 6*f2y
247  - 3*f3y + 3*f4y + 4*f1xy + 2*f2xy + f3xy + 2*f4xy)*h12*h22
248  + (6*(-f1 + f2 - f3 + f4) - 4*f1x - 2*f2x + 2*f3x + 4*f4x - 3*f1y
249  + 3*f2y + 3*f3y - 3*f4y - 2*f1xy - f2xy - f3xy - 2*f4xy)*h12*h23
250  + (2*(f1 - f2) + f1x + f2x)*h13 + (2*(f1y - f2y) + f1xy + f2xy)*h13*h2
251  + (6*(-f1 + f2 -f3 + f4) + 3*(-f1x - f2x + f3x + f4x) - 4*f1y
252  + 4*f2y + 2*f3y - 2*f4y - 2*f1xy - 2*f2xy - f3xy - f4xy)*h13*h22
253  + (4*(f1 - f2 + f3 - f4) + 2*(f1x + f2x - f3x - f4x)
254  + 2*(f1y - f2y - f3y + f4y) + f1xy + f2xy + f3xy + f4xy)*h13*h23;
255 }
G4double GetValue(size_t idx, size_t idy) const
double G4double
Definition: G4Types.hh:76
void G4Physics2DVector::ClearVectors ( )
protected

Definition at line 126 of file G4Physics2DVector.cc.

Referenced by operator=(), PutVectors(), Retrieve(), and ~G4Physics2DVector().

127 {
128  for(size_t j=0; j<numberOfYNodes; ++j) {
129  delete value[j];
130  }
131 }
const XML_Char int const XML_Char * value
void G4Physics2DVector::CopyData ( const G4Physics2DVector vec)
protected

Definition at line 135 of file G4Physics2DVector.cc.

References PutValue().

Referenced by G4Physics2DVector(), and operator=().

136 {
137  for(size_t i=0; i<numberOfXNodes; ++i) {
138  xVector[i] = right.xVector[i];
139  }
140  for(size_t j=0; j<numberOfYNodes; ++j) {
141  yVector[j] = right.yVector[j];
142  G4PV2DDataVector* v0 = right.value[j];
143  for(size_t i=0; i<numberOfXNodes; ++i) {
144  PutValue(i,j,(*v0)[i]);
145  }
146  }
147 }
std::vector< G4double > G4PV2DDataVector
void PutValue(size_t idx, size_t idy, G4double value)
size_t G4Physics2DVector::FindBin ( G4double  z,
const G4PV2DDataVector ,
size_t  idz,
size_t  idzmax 
) const
inlineprotected
size_t G4Physics2DVector::FindBinLocation ( G4double  z,
const G4PV2DDataVector v 
) const
protected

Definition at line 359 of file G4Physics2DVector.cc.

361 {
362  size_t lowerBound = 0;
363  size_t upperBound = v.size() - 2;
364 
365  while (lowerBound <= upperBound)
366  {
367  size_t midBin = (lowerBound + upperBound)/2;
368  if( z < v[midBin] ) { upperBound = midBin-1; }
369  else { lowerBound = midBin+1; }
370  }
371 
372  return upperBound;
373 }
G4double z
Definition: TRTMaterials.hh:39
size_t G4Physics2DVector::FindBinLocationX ( G4double  x,
size_t  lastidx 
) const
inline

Referenced by Value().

size_t G4Physics2DVector::FindBinLocationY ( G4double  y,
size_t  lastidy 
) const
inline

Referenced by FindLinearX(), and Value().

G4double G4Physics2DVector::FindLinearX ( G4double  rand,
G4double  y,
size_t &  lastidy 
) const

Definition at line 377 of file G4Physics2DVector.cc.

References FindBinLocationY().

379 {
380  G4double y = yy;
381 
382  // no interpolation outside the table
383  if(y < yVector[0]) {
384  y = yVector[0];
385  } else if(y > yVector[numberOfYNodes - 1]) {
386  y = yVector[numberOfYNodes - 1];
387  }
388 
389  // find bins
390  idy = FindBinLocationY(y, idy);
391 
392  G4double x1 = InterpolateLinearX(*(value[idy]), rand);
393  G4double x2 = InterpolateLinearX(*(value[idy+1]), rand);
394  G4double res = x1;
395  G4double del = yVector[idy+1] - yVector[idy];
396  if(del != 0.0) {
397  res += (x2 - x1)*(y - yVector[idy])/del;
398  }
399  return res;
400 }
size_t FindBinLocationY(G4double y, size_t lastidy) const
const XML_Char int const XML_Char * value
double G4double
Definition: G4Types.hh:76
G4double G4Physics2DVector::FindLinearX ( G4double  rand,
G4double  y 
) const
inline
size_t G4Physics2DVector::GetLengthX ( ) const
inline
size_t G4Physics2DVector::GetLengthY ( ) const
inline
G4PhysicsVectorType G4Physics2DVector::GetType ( ) const
inline
G4double G4Physics2DVector::GetValue ( size_t  idx,
size_t  idy 
) const
inline
G4int G4Physics2DVector::GetVerboseLevel ( ) const
inline
G4double G4Physics2DVector::GetX ( size_t  index) const
inline
G4double G4Physics2DVector::GetY ( size_t  index) const
inline
G4Physics2DVector & G4Physics2DVector::operator= ( const G4Physics2DVector right)

Definition at line 91 of file G4Physics2DVector.cc.

References ClearVectors(), CopyData(), and PrepareVectors().

92 {
93  if (&right==this) { return *this; }
94  ClearVectors();
95 
96  type = right.type;
97 
98  numberOfXNodes = right.numberOfXNodes;
99  numberOfYNodes = right.numberOfYNodes;
100 
101  verboseLevel = right.verboseLevel;
102  useBicubic = right.useBicubic;
103 
104  PrepareVectors();
105  CopyData(right);
106 
107  return *this;
108 }
void CopyData(const G4Physics2DVector &vec)
void G4Physics2DVector::PrepareVectors ( )
protected

Definition at line 112 of file G4Physics2DVector.cc.

References test::v.

Referenced by G4Physics2DVector(), operator=(), PutVectors(), and Retrieve().

113 {
114  xVector.resize(numberOfXNodes,0.);
115  yVector.resize(numberOfYNodes,0.);
116  value.resize(numberOfYNodes,0);
117  for(size_t j=0; j<numberOfYNodes; ++j) {
119  v->resize(numberOfXNodes,0.);
120  value[j] = v;
121  }
122 }
std::vector< G4double > G4PV2DDataVector
const XML_Char int const XML_Char * value
void G4Physics2DVector::PutValue ( size_t  idx,
size_t  idy,
G4double  value 
)
inline

Referenced by CopyData(), Retrieve(), and ScaleVector().

void G4Physics2DVector::PutVectors ( const std::vector< G4double > &  vecX,
const std::vector< G4double > &  vecY 
)

Definition at line 260 of file G4Physics2DVector.cc.

References ClearVectors(), and PrepareVectors().

262 {
263  ClearVectors();
264  numberOfXNodes = vecX.size();
265  numberOfYNodes = vecY.size();
266  PrepareVectors();
267  for(size_t i = 0; i<numberOfXNodes; ++i) {
268  xVector[i] = vecX[i];
269  }
270  for(size_t j = 0; j<numberOfYNodes; ++j) {
271  yVector[j] = vecY[j];
272  }
273 }
void G4Physics2DVector::PutX ( size_t  idx,
G4double  value 
)
inline
void G4Physics2DVector::PutY ( size_t  idy,
G4double  value 
)
inline
G4bool G4Physics2DVector::Retrieve ( std::ifstream &  fIn)

Definition at line 306 of file G4Physics2DVector.cc.

References ClearVectors(), INT_MAX, PrepareVectors(), and PutValue().

Referenced by G4OpticalSurface::ReadDichroicFile().

307 {
308  // initialisation
309  ClearVectors();
310 
311  // binning
312  G4int k;
313  in >> k >> numberOfXNodes >> numberOfYNodes;
314  if (in.fail() || 0 >= numberOfXNodes || 0 >= numberOfYNodes ||
315  numberOfXNodes >= INT_MAX || numberOfYNodes >= INT_MAX) {
316  return false;
317  }
318  PrepareVectors();
319  type = G4PhysicsVectorType(k);
320 
321  // contents
322  G4double val;
323  for(size_t i = 0; i<numberOfXNodes; ++i) {
324  in >> xVector[i];
325  if (in.fail()) { return false; }
326  }
327  for(size_t j = 0; j<numberOfYNodes; ++j) {
328  in >> yVector[j];
329  if (in.fail()) { return false; }
330  }
331  for(size_t j = 0; j<numberOfYNodes; ++j) {
332  for(size_t i = 0; i<numberOfXNodes; ++i) {
333  in >> val;
334  if (in.fail()) { return false; }
335  PutValue(i, j, val);
336  }
337  }
338  in.close();
339  return true;
340 }
G4PhysicsVectorType
int G4int
Definition: G4Types.hh:78
void PutValue(size_t idx, size_t idy, G4double value)
#define INT_MAX
Definition: templates.hh:111
double G4double
Definition: G4Types.hh:76
void G4Physics2DVector::ScaleVector ( G4double  factor)

Definition at line 345 of file G4Physics2DVector.cc.

References GetValue(), and PutValue().

346 {
347  G4double val;
348  for(size_t j = 0; j<numberOfYNodes; ++j) {
349  for(size_t i = 0; i<numberOfXNodes; ++i) {
350  val = GetValue(i, j)*factor;
351  PutValue(i, j, val);
352  }
353  }
354 }
G4double GetValue(size_t idx, size_t idy) const
void PutValue(size_t idx, size_t idy, G4double value)
double G4double
Definition: G4Types.hh:76
void G4Physics2DVector::SetBicubicInterpolation ( G4bool  )
inline
void G4Physics2DVector::SetVerboseLevel ( G4int  value)
inline
void G4Physics2DVector::Store ( std::ofstream &  fOut)

Definition at line 277 of file G4Physics2DVector.cc.

References G4endl, and GetValue().

278 {
279  // binning
280  G4int prec = out.precision();
281  out << G4int(type) << " " << numberOfXNodes << " " << numberOfYNodes
282  << G4endl;
283  out << std::setprecision(5);
284 
285  // contents
286  for(size_t i = 0; i<numberOfXNodes-1; ++i) {
287  out << xVector[i] << " ";
288  }
289  out << xVector[numberOfXNodes-1] << G4endl;
290  for(size_t j = 0; j<numberOfYNodes-1; ++j) {
291  out << yVector[j] << " ";
292  }
293  out << yVector[numberOfYNodes-1] << G4endl;
294  for(size_t j = 0; j<numberOfYNodes; ++j) {
295  for(size_t i = 0; i<numberOfXNodes-1; ++i) {
296  out << GetValue(i, j) << " ";
297  }
298  out << GetValue(numberOfXNodes-1,j) << G4endl;
299  }
300  out.precision(prec);
301  out.close();
302 }
G4double GetValue(size_t idx, size_t idy) const
int G4int
Definition: G4Types.hh:78
#define G4endl
Definition: G4ios.hh:61
G4double G4Physics2DVector::Value ( G4double  x,
G4double  y,
size_t &  lastidx,
size_t &  lastidy 
) const

Definition at line 151 of file G4Physics2DVector.cc.

References BicubicInterpolation(), FindBinLocationX(), FindBinLocationY(), GetValue(), and test::x.

Referenced by G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom(), G4SeltzerBergerModel::ComputeDXSectionPerAtom(), G4SeltzerBergerModel::SampleSecondaries(), and G4LivermoreBremsstrahlungModel::SampleSecondaries().

153 {
154  G4double x = xx;
155  G4double y = yy;
156 
157  // no interpolation outside the table
158  if(x < xVector[0]) {
159  x = xVector[0];
160  } else if(x > xVector[numberOfXNodes - 1]) {
161  x = xVector[numberOfXNodes - 1];
162  }
163  if(y < yVector[0]) {
164  y = yVector[0];
165  } else if(y > yVector[numberOfYNodes - 1]) {
166  y = yVector[numberOfYNodes - 1];
167  }
168 
169  // find bins
170  idx = FindBinLocationX(x, idx);
171  idy = FindBinLocationY(y, idy);
172 
173  // interpolate
174  if(useBicubic) {
175  return BicubicInterpolation(x, y, idx, idy);
176  } else {
177  G4double x1 = xVector[idx];
178  G4double x2 = xVector[idx+1];
179  G4double y1 = yVector[idy];
180  G4double y2 = yVector[idy+1];
181  G4double v11= GetValue(idx, idy);
182  G4double v12= GetValue(idx+1, idy);
183  G4double v21= GetValue(idx, idy+1);
184  G4double v22= GetValue(idx+1, idy+1);
185  return ((y2 - y)*(v11*(x2 - x) + v12*(x - x1)) +
186  ((y - y1)*(v21*(x2 - x) + v22*(x - x1))))/((x2 - x1)*(y2 - y1));
187  }
188 }
G4double GetValue(size_t idx, size_t idy) const
size_t FindBinLocationX(G4double x, size_t lastidx) const
G4double BicubicInterpolation(G4double x, G4double y, size_t idx, size_t idy) const
size_t FindBinLocationY(G4double y, size_t lastidy) const
double G4double
Definition: G4Types.hh:76
G4double G4Physics2DVector::Value ( G4double  x,
G4double  y 
) const

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