Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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) const
 
G4bool Retrieve (std::ifstream &fIn)
 
void SetVerboseLevel (G4int value)
 

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 43 of file G4Physics2DVector.cc.

44  : type(T_G4PhysicsFreeVector),
45  numberOfXNodes(0), numberOfYNodes(0),
46  verboseLevel(0), useBicubic(false)
47 {}
G4Physics2DVector::G4Physics2DVector ( size_t  nx,
size_t  ny 
)
explicit

Definition at line 51 of file G4Physics2DVector.cc.

52  : type(T_G4PhysicsFreeVector),
53  numberOfXNodes(nx), numberOfYNodes(ny),
54  verboseLevel(0), useBicubic(false)
55 {
57 }

Here is the call graph for this function:

G4Physics2DVector::G4Physics2DVector ( const G4Physics2DVector right)

Definition at line 68 of file G4Physics2DVector.cc.

69 {
70  type = right.type;
71 
72  numberOfXNodes = right.numberOfXNodes;
73  numberOfYNodes = right.numberOfYNodes;
74 
75  verboseLevel = right.verboseLevel;
76  useBicubic = right.useBicubic;
77 
78  xVector = right.xVector;
79  yVector = right.yVector;
80 
82  CopyData(right);
83 }
void CopyData(const G4Physics2DVector &vec)

Here is the call graph for this function:

G4Physics2DVector::~G4Physics2DVector ( )

Definition at line 61 of file G4Physics2DVector.cc.

62 {
63  ClearVectors();
64 }

Here is the call graph for this function:

Member Function Documentation

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

Definition at line 189 of file G4Physics2DVector.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Physics2DVector::ClearVectors ( )
protected

Definition at line 122 of file G4Physics2DVector.cc.

123 {
124  for(size_t j=0; j<numberOfYNodes; ++j) {
125  delete value[j];
126  }
127 }
const XML_Char int const XML_Char * value
Definition: expat.h:331

Here is the caller graph for this function:

void G4Physics2DVector::CopyData ( const G4Physics2DVector vec)
protected

Definition at line 131 of file G4Physics2DVector.cc.

132 {
133  for(size_t i=0; i<numberOfXNodes; ++i) {
134  xVector[i] = right.xVector[i];
135  }
136  for(size_t j=0; j<numberOfYNodes; ++j) {
137  yVector[j] = right.yVector[j];
138  G4PV2DDataVector* v0 = right.value[j];
139  for(size_t i=0; i<numberOfXNodes; ++i) {
140  PutValue(i,j,(*v0)[i]);
141  }
142  }
143 }
std::vector< G4double > G4PV2DDataVector
void PutValue(size_t idx, size_t idy, G4double value)

Here is the call graph for this function:

Here is the caller graph for this function:

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 361 of file G4Physics2DVector.cc.

363 {
364  size_t bin;
365  size_t binmax = v.size() - 2;
366 
367  if(z <= v[0]) { bin = 0; }
368  else if(z >= v[binmax]) { bin = binmax; }
369  else {
370  bin = std::lower_bound(v.begin(), v.end(), z) - v.begin() - 1;
371  }
372  return bin;
373 }
size_t G4Physics2DVector::FindBinLocationX ( G4double  x,
size_t  lastidx 
) const
inline

Here is the caller graph for this function:

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

Here is the caller graph for this function:

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

Definition at line 377 of file G4Physics2DVector.cc.

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 }
const XML_Char int const XML_Char * value
Definition: expat.h:331
size_t FindBinLocationY(G4double y, size_t lastidy) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4Physics2DVector::FindLinearX ( G4double  rand,
G4double  y 
) const
inline
size_t G4Physics2DVector::GetLengthX ( ) const
inline

Here is the caller graph for this function:

size_t G4Physics2DVector::GetLengthY ( ) const
inline

Here is the caller graph for this function:

G4PhysicsVectorType G4Physics2DVector::GetType ( ) const
inline
G4double G4Physics2DVector::GetValue ( size_t  idx,
size_t  idy 
) const
inline

Here is the caller graph for this function:

G4double G4Physics2DVector::GetX ( size_t  index) const
inline

Here is the caller graph for this function:

G4double G4Physics2DVector::GetY ( size_t  index) const
inline

Here is the caller graph for this function:

G4Physics2DVector & G4Physics2DVector::operator= ( const G4Physics2DVector right)

Definition at line 87 of file G4Physics2DVector.cc.

88 {
89  if (&right==this) { return *this; }
90  ClearVectors();
91 
92  type = right.type;
93 
94  numberOfXNodes = right.numberOfXNodes;
95  numberOfYNodes = right.numberOfYNodes;
96 
97  verboseLevel = right.verboseLevel;
98  useBicubic = right.useBicubic;
99 
100  PrepareVectors();
101  CopyData(right);
102 
103  return *this;
104 }
void CopyData(const G4Physics2DVector &vec)

Here is the call graph for this function:

void G4Physics2DVector::PrepareVectors ( )
protected

Definition at line 108 of file G4Physics2DVector.cc.

109 {
110  xVector.resize(numberOfXNodes,0.);
111  yVector.resize(numberOfYNodes,0.);
112  value.resize(numberOfYNodes,0);
113  for(size_t j=0; j<numberOfYNodes; ++j) {
115  v->resize(numberOfXNodes,0.);
116  value[j] = v;
117  }
118 }
std::vector< G4double > G4PV2DDataVector
const XML_Char int const XML_Char * value
Definition: expat.h:331

Here is the caller graph for this function:

void G4Physics2DVector::PutValue ( size_t  idx,
size_t  idy,
G4double  value 
)
inline

Here is the caller graph for this function:

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

Definition at line 256 of file G4Physics2DVector.cc.

258 {
259  ClearVectors();
260  numberOfXNodes = vecX.size();
261  numberOfYNodes = vecY.size();
262  PrepareVectors();
263  for(size_t i = 0; i<numberOfXNodes; ++i) {
264  xVector[i] = vecX[i];
265  }
266  for(size_t j = 0; j<numberOfYNodes; ++j) {
267  yVector[j] = vecY[j];
268  }
269 }

Here is the call graph for this function:

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 302 of file G4Physics2DVector.cc.

303 {
304  // initialisation
305  ClearVectors();
306 
307  // binning
308  G4int k;
309  in >> k >> numberOfXNodes >> numberOfYNodes;
310  if (in.fail() || 0 >= numberOfXNodes || 0 >= numberOfYNodes ||
311  numberOfXNodes >= INT_MAX || numberOfYNodes >= INT_MAX) {
312  if( 0 >= numberOfXNodes || numberOfXNodes >= INT_MAX) {
313  numberOfXNodes = 0;
314  }
315  if( 0 >= numberOfYNodes || numberOfYNodes >= INT_MAX) {
316  numberOfYNodes = 0;
317  }
318  return false;
319  }
320  PrepareVectors();
321  type = G4PhysicsVectorType(k);
322 
323  // contents
324  G4double val;
325  for(size_t i = 0; i<numberOfXNodes; ++i) {
326  in >> xVector[i];
327  if (in.fail()) { return false; }
328  }
329  for(size_t j = 0; j<numberOfYNodes; ++j) {
330  in >> yVector[j];
331  if (in.fail()) { return false; }
332  }
333  for(size_t j = 0; j<numberOfYNodes; ++j) {
334  for(size_t i = 0; i<numberOfXNodes; ++i) {
335  in >> val;
336  if (in.fail()) { return false; }
337  PutValue(i, j, val);
338  }
339  }
340  in.close();
341  return true;
342 }
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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Physics2DVector::ScaleVector ( G4double  factor)

Definition at line 347 of file G4Physics2DVector.cc.

348 {
349  G4double val;
350  for(size_t j = 0; j<numberOfYNodes; ++j) {
351  for(size_t i = 0; i<numberOfXNodes; ++i) {
352  val = GetValue(i, j)*factor;
353  PutValue(i, j, val);
354  }
355  }
356 }
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

Here is the call graph for this function:

void G4Physics2DVector::SetBicubicInterpolation ( G4bool  )
inline
void G4Physics2DVector::SetVerboseLevel ( G4int  value)
inline
void G4Physics2DVector::Store ( std::ofstream &  fOut) const

Definition at line 273 of file G4Physics2DVector.cc.

274 {
275  // binning
276  G4int prec = out.precision();
277  out << G4int(type) << " " << numberOfXNodes << " " << numberOfYNodes
278  << G4endl;
279  out << std::setprecision(5);
280 
281  // contents
282  for(size_t i = 0; i<numberOfXNodes-1; ++i) {
283  out << xVector[i] << " ";
284  }
285  out << xVector[numberOfXNodes-1] << G4endl;
286  for(size_t j = 0; j<numberOfYNodes-1; ++j) {
287  out << yVector[j] << " ";
288  }
289  out << yVector[numberOfYNodes-1] << G4endl;
290  for(size_t j = 0; j<numberOfYNodes; ++j) {
291  for(size_t i = 0; i<numberOfXNodes-1; ++i) {
292  out << GetValue(i, j) << " ";
293  }
294  out << GetValue(numberOfXNodes-1,j) << G4endl;
295  }
296  out.precision(prec);
297  out.close();
298 }
G4double GetValue(size_t idx, size_t idy) const
int G4int
Definition: G4Types.hh:78
static const double prec
Definition: RanecuEngine.cc:58
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4double G4Physics2DVector::Value ( G4double  x,
G4double  y,
size_t &  lastidx,
size_t &  lastidy 
) const

Definition at line 147 of file G4Physics2DVector.cc.

149 {
150  G4double x = xx;
151  G4double y = yy;
152 
153  // no interpolation outside the table
154  if(x < xVector[0]) {
155  x = xVector[0];
156  } else if(x > xVector[numberOfXNodes - 1]) {
157  x = xVector[numberOfXNodes - 1];
158  }
159  if(y < yVector[0]) {
160  y = yVector[0];
161  } else if(y > yVector[numberOfYNodes - 1]) {
162  y = yVector[numberOfYNodes - 1];
163  }
164 
165  // find bins
166  idx = FindBinLocationX(x, idx);
167  idy = FindBinLocationY(y, idy);
168 
169  // interpolate
170  if(useBicubic) {
171  return BicubicInterpolation(x, y, idx, idy);
172  } else {
173  G4double x1 = xVector[idx];
174  G4double x2 = xVector[idx+1];
175  G4double y1 = yVector[idy];
176  G4double y2 = yVector[idy+1];
177  G4double v11= GetValue(idx, idy);
178  G4double v12= GetValue(idx+1, idy);
179  G4double v21= GetValue(idx, idy+1);
180  G4double v22= GetValue(idx+1, idy+1);
181  return ((y2 - y)*(v11*(x2 - x) + v12*(x - x1)) +
182  ((y - y1)*(v21*(x2 - x) + v22*(x - x1))))/((x2 - x1)*(y2 - y1));
183  }
184 }
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

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Physics2DVector::Value ( G4double  x,
G4double  y 
) const

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