Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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)
 
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 GetX (size_t index) const
 
G4double GetY (size_t index) const
 
G4double GetValue (size_t idx, size_t idy) 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)
 
G4double GetLastX () const
 
G4double GetLastY () const
 
G4double GetLastValue () const
 
size_t GetLastBinX () const
 
size_t GetLastBinY () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

void PrepareVectors ()
 
void ClearVectors ()
 
void CopyData (const G4Physics2DVector &vec)
 
void ComputeValue (G4double x, G4double y)
 
void BicubicInterpolation (size_t idx, size_t idy)
 
size_t FindBinLocation (G4double z, const G4PV2DDataVector &)
 
void FindBin (G4double z, const G4PV2DDataVector &, size_t &lastidx)
 
void FindBinLocationX (G4double x)
 
void FindBinLocationY (G4double y)
 

Detailed Description

Definition at line 60 of file G4Physics2DVector.hh.

Constructor & Destructor Documentation

◆ G4Physics2DVector() [1/3]

G4Physics2DVector::G4Physics2DVector ( )

Definition at line 47 of file G4Physics2DVector.cc.

49 numberOfXNodes(0), numberOfYNodes(0),
50 verboseLevel(0), useBicubic(false)
51{
52 cache = new G4Physics2DVectorCache();
53}
@ T_G4PhysicsFreeVector

◆ G4Physics2DVector() [2/3]

G4Physics2DVector::G4Physics2DVector ( size_t  nx,
size_t  ny 
)

Definition at line 57 of file G4Physics2DVector.cc.

59 numberOfXNodes(nx), numberOfYNodes(ny),
60 verboseLevel(0), useBicubic(false)
61{
62 cache = new G4Physics2DVectorCache();
64}

◆ G4Physics2DVector() [3/3]

G4Physics2DVector::G4Physics2DVector ( const G4Physics2DVector right)

Definition at line 76 of file G4Physics2DVector.cc.

77{
78 type = right.type;
79
80 numberOfXNodes = right.numberOfXNodes;
81 numberOfYNodes = right.numberOfYNodes;
82
83 verboseLevel = right.verboseLevel;
84 useBicubic = right.useBicubic;
85
86 xVector = right.xVector;
87 yVector = right.yVector;
88
89 cache = new G4Physics2DVectorCache();
91 CopyData(right);
92}
void CopyData(const G4Physics2DVector &vec)

◆ ~G4Physics2DVector()

G4Physics2DVector::~G4Physics2DVector ( )

Definition at line 68 of file G4Physics2DVector.cc.

69{
70 delete cache;
72}

Member Function Documentation

◆ BicubicInterpolation()

void G4Physics2DVector::BicubicInterpolation ( size_t  idx,
size_t  idy 
)
protected

Definition at line 206 of file G4Physics2DVector.cc.

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

Referenced by ComputeValue().

◆ ClearVectors()

void G4Physics2DVector::ClearVectors ( )
protected

Definition at line 132 of file G4Physics2DVector.cc.

133{
134 for(size_t j=0; j<numberOfYNodes; ++j) {
135 delete value[j];
136 }
137}

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

◆ ComputeValue()

void G4Physics2DVector::ComputeValue ( G4double  x,
G4double  y 
)
protected

Definition at line 157 of file G4Physics2DVector.cc.

158{
159 if(xx != cache->lastBinX) {
160 if(xx <= xVector[0]) {
161 cache->lastX = xVector[0];
162 cache->lastBinX = 0;
163 } else if(xx >= xVector[numberOfXNodes-1]) {
164 cache->lastX = xVector[numberOfXNodes-1];
165 cache->lastBinX = numberOfXNodes-2;
166 } else {
167 cache->lastX = xx;
169 }
170 }
171 if(yy != cache->lastBinY) {
172 if(yy <= yVector[0]) {
173 cache->lastY = yVector[0];
174 cache->lastBinY = 0;
175 } else if(yy >= yVector[numberOfYNodes-1]) {
176 cache->lastY = yVector[numberOfYNodes-1];
177 cache->lastBinY = numberOfYNodes-2;
178 } else {
179 cache->lastY = yy;
181 }
182 }
183 size_t idx = cache->lastBinX;
184 size_t idy = cache->lastBinY;
185 if(useBicubic) {
186 BicubicInterpolation(idx, idy);
187 } else {
188 G4double x1 = xVector[idx];
189 G4double x2 = xVector[idx+1];
190 G4double y1 = yVector[idy];
191 G4double y2 = yVector[idy+1];
192 G4double x = cache->lastX;
193 G4double y = cache->lastY;
194 G4double v11= GetValue(idx, idy);
195 G4double v12= GetValue(idx+1, idy);
196 G4double v21= GetValue(idx, idy+1);
197 G4double v22= GetValue(idx+1, idy+1);
198 cache->lastValue =
199 ((y2 - y)*(v11*(x2 - x) + v12*(x - x1)) +
200 ((y - y1)*(v21*(x2 - x) + v22*(x - x1))))/((x2 - x1)*(y2 - y1));
201 }
202}
void FindBinLocationY(G4double y)
void BicubicInterpolation(size_t idx, size_t idy)
void FindBinLocationX(G4double x)

◆ CopyData()

void G4Physics2DVector::CopyData ( const G4Physics2DVector vec)
protected

Definition at line 141 of file G4Physics2DVector.cc.

142{
143 for(size_t i=0; i<numberOfXNodes; ++i) {
144 xVector[i] = right.xVector[i];
145 }
146 for(size_t j=0; j<numberOfYNodes; ++j) {
147 yVector[j] = right.yVector[j];
148 G4PV2DDataVector* v0 = right.value[j];
149 for(size_t i=0; i<numberOfXNodes; ++i) {
150 PutValue(i,j,(*v0)[i]);
151 }
152 }
153}
std::vector< G4double > G4PV2DDataVector
void PutValue(size_t idx, size_t idy, G4double value)

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

◆ FindBin()

void G4Physics2DVector::FindBin ( G4double  z,
const G4PV2DDataVector ,
size_t &  lastidx 
)
inlineprotected

◆ FindBinLocation()

size_t G4Physics2DVector::FindBinLocation ( G4double  z,
const G4PV2DDataVector v 
)
protected

Definition at line 373 of file G4Physics2DVector.cc.

375{
376 size_t lowerBound = 0;
377 size_t upperBound = v.size() - 2;
378
379 while (lowerBound <= upperBound)
380 {
381 size_t midBin = (lowerBound + upperBound)/2;
382 if( z < v[midBin] ) { upperBound = midBin-1; }
383 else { lowerBound = midBin+1; }
384 }
385
386 return upperBound;
387}

◆ FindBinLocationX()

void G4Physics2DVector::FindBinLocationX ( G4double  x)
inlineprotected

Referenced by ComputeValue().

◆ FindBinLocationY()

void G4Physics2DVector::FindBinLocationY ( G4double  y)
inlineprotected

Referenced by ComputeValue().

◆ GetLastBinX()

size_t G4Physics2DVector::GetLastBinX ( ) const
inline

◆ GetLastBinY()

size_t G4Physics2DVector::GetLastBinY ( ) const
inline

◆ GetLastValue()

G4double G4Physics2DVector::GetLastValue ( ) const
inline

◆ GetLastX()

G4double G4Physics2DVector::GetLastX ( ) const
inline

◆ GetLastY()

G4double G4Physics2DVector::GetLastY ( ) const
inline

◆ GetLengthX()

size_t G4Physics2DVector::GetLengthX ( ) const
inline

◆ GetLengthY()

size_t G4Physics2DVector::GetLengthY ( ) const
inline

◆ GetType()

G4PhysicsVectorType G4Physics2DVector::GetType ( ) const
inline

◆ GetValue()

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

◆ GetVerboseLevel()

G4int G4Physics2DVector::GetVerboseLevel ( ) const
inline

◆ GetX()

G4double G4Physics2DVector::GetX ( size_t  index) const
inline

◆ GetY()

G4double G4Physics2DVector::GetY ( size_t  index) const
inline

◆ operator=()

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

Definition at line 96 of file G4Physics2DVector.cc.

97{
98 if (&right==this) { return *this; }
100
101 type = right.type;
102
103 numberOfXNodes = right.numberOfXNodes;
104 numberOfYNodes = right.numberOfYNodes;
105
106 verboseLevel = right.verboseLevel;
107 useBicubic = right.useBicubic;
108
109 cache->Clear();
111 CopyData(right);
112
113 return *this;
114}

◆ PrepareVectors()

void G4Physics2DVector::PrepareVectors ( )
protected

Definition at line 118 of file G4Physics2DVector.cc.

119{
120 xVector.resize(numberOfXNodes,0.);
121 yVector.resize(numberOfYNodes,0.);
122 value.resize(numberOfYNodes,0);
123 for(size_t j=0; j<numberOfYNodes; ++j) {
125 v->resize(numberOfXNodes,0.);
126 value[j] = v;
127 }
128}

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

◆ PutValue()

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

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

◆ PutVectors()

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

Definition at line 274 of file G4Physics2DVector.cc.

276{
277 ClearVectors();
278 numberOfXNodes = vecX.size();
279 numberOfYNodes = vecY.size();
281 if(!cache) { cache = new G4Physics2DVectorCache(); }
282 cache->Clear();
283 for(size_t i = 0; i<numberOfXNodes; ++i) {
284 xVector[i] = vecX[i];
285 }
286 for(size_t j = 0; j<numberOfYNodes; ++j) {
287 yVector[j] = vecY[j];
288 }
289}

◆ PutX()

void G4Physics2DVector::PutX ( size_t  idx,
G4double  value 
)
inline

◆ PutY()

void G4Physics2DVector::PutY ( size_t  idy,
G4double  value 
)
inline

◆ Retrieve()

G4bool G4Physics2DVector::Retrieve ( std::ifstream &  fIn)

Definition at line 322 of file G4Physics2DVector.cc.

323{
324 // initialisation
325 cache->Clear();
326 ClearVectors();
327
328 // binning
329 G4int k;
330 in >> k >> numberOfXNodes >> numberOfYNodes;
331 if (in.fail()) { return false; }
333 type = G4PhysicsVectorType(k);
334
335 // contents
336 G4double val;
337 for(size_t i = 0; i<numberOfXNodes; ++i) {
338 in >> xVector[i];
339 if (in.fail()) { return false; }
340 }
341 for(size_t j = 0; j<numberOfYNodes; ++j) {
342 in >> yVector[j];
343 if (in.fail()) { return false; }
344 }
345 for(size_t j = 0; j<numberOfYNodes; ++j) {
346 for(size_t i = 0; i<numberOfXNodes; ++i) {
347 in >> val;
348 if (in.fail()) { return false; }
349 PutValue(i, j, val);
350 }
351 }
352 in.close();
353 return true;
354}
G4PhysicsVectorType
int G4int
Definition: G4Types.hh:66

◆ ScaleVector()

void G4Physics2DVector::ScaleVector ( G4double  factor)

Definition at line 359 of file G4Physics2DVector.cc.

360{
361 G4double val;
362 for(size_t j = 0; j<numberOfYNodes; ++j) {
363 for(size_t i = 0; i<numberOfXNodes; ++i) {
364 val = GetValue(i, j)*factor;
365 PutValue(i, j, val);
366 }
367 }
368}

◆ SetBicubicInterpolation()

void G4Physics2DVector::SetBicubicInterpolation ( G4bool  )
inline

◆ SetVerboseLevel()

void G4Physics2DVector::SetVerboseLevel ( G4int  value)
inline

◆ Store()

void G4Physics2DVector::Store ( std::ofstream &  fOut)

Definition at line 293 of file G4Physics2DVector.cc.

294{
295 // binning
296 G4int prec = out.precision();
297 out << G4int(type) << " " << numberOfXNodes << " " << numberOfYNodes
298 << G4endl;
299 out << std::setprecision(5);
300
301 // contents
302 for(size_t i = 0; i<numberOfXNodes-1; ++i) {
303 out << xVector[i] << " ";
304 }
305 out << xVector[numberOfXNodes-1] << G4endl;
306 for(size_t j = 0; j<numberOfYNodes-1; ++j) {
307 out << yVector[j] << " ";
308 }
309 out << yVector[numberOfYNodes-1] << G4endl;
310 for(size_t j = 0; j<numberOfYNodes; ++j) {
311 for(size_t i = 0; i<numberOfXNodes-1; ++i) {
312 out << GetValue(i, j) << " ";
313 }
314 out << GetValue(numberOfXNodes-1,j) << G4endl;
315 }
316 out.precision(prec);
317 out.close();
318}
#define G4endl
Definition: G4ios.hh:52

◆ Value()


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