Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysicsVector.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// G4PhysicsVector class implementation
27//
28// Authors:
29// - 02 Dec. 1995, G.Cosmo: Structure created based on object model
30// - 03 Mar. 1996, K.Amako: Implemented the 1st version
31// Revisions:
32// - 11 Nov. 2000, H.Kurashige: Use STL vector for dataVector and binVector
33// --------------------------------------------------------------------
34
35#include "G4PhysicsVector.hh"
36#include <iomanip>
37
38// --------------------------------------------------------------
40 : useSpline(val)
41{}
42
43// --------------------------------------------------------------------
45{
47 if(0 < numberOfNodes)
48 {
49 edgeMin = binVector[0];
51 }
52}
53
54// --------------------------------------------------------------
55G4bool G4PhysicsVector::Store(std::ofstream& fOut, G4bool ascii) const
56{
57 // Ascii mode
58 if(ascii)
59 {
60 fOut << *this;
61 return true;
62 }
63 // Binary Mode
64
65 // binning
66 fOut.write((char*) (&edgeMin), sizeof edgeMin);
67 fOut.write((char*) (&edgeMax), sizeof edgeMax);
68 fOut.write((char*) (&numberOfNodes), sizeof numberOfNodes);
69
70 // contents
71 std::size_t size = dataVector.size();
72 fOut.write((char*) (&size), sizeof size);
73
74 G4double* value = new G4double[2 * size];
75 for(std::size_t i = 0; i < size; ++i)
76 {
77 value[2 * i] = binVector[i];
78 value[2 * i + 1] = dataVector[i];
79 }
80 fOut.write((char*) (value), 2 * size * (sizeof(G4double)));
81 delete[] value;
82
83 return true;
84}
85
86// --------------------------------------------------------------
87G4bool G4PhysicsVector::Retrieve(std::ifstream& fIn, G4bool ascii)
88{
89 // clear properties;
90 dataVector.clear();
91 binVector.clear();
92 secDerivative.clear();
93
94 // retrieve in ascii mode
95 if(ascii)
96 {
97 // binning
98 fIn >> edgeMin >> edgeMax >> numberOfNodes;
99 if(fIn.fail() || numberOfNodes < 2)
100 {
101 return false;
102 }
103 // contents
104 G4int siz = 0;
105 fIn >> siz;
106 if(fIn.fail() || siz != G4int(numberOfNodes))
107 {
108 return false;
109 }
110
111 binVector.reserve(siz);
112 dataVector.reserve(siz);
113 G4double vBin, vData;
114
115 for(G4int i = 0; i < siz; ++i)
116 {
117 vBin = 0.;
118 vData = 0.;
119 fIn >> vBin >> vData;
120 if(fIn.fail())
121 {
122 return false;
123 }
124 binVector.push_back(vBin);
125 dataVector.push_back(vData);
126 }
127 Initialise();
128 return true;
129 }
130
131 // retrieve in binary mode
132 // binning
133 fIn.read((char*) (&edgeMin), sizeof edgeMin);
134 fIn.read((char*) (&edgeMax), sizeof edgeMax);
135 fIn.read((char*) (&numberOfNodes), sizeof numberOfNodes);
136
137 // contents
138 std::size_t size;
139 fIn.read((char*) (&size), sizeof size);
140
141 G4double* value = new G4double[2 * size];
142 fIn.read((char*) (value), 2 * size * (sizeof(G4double)));
143 if(G4int(fIn.gcount()) != G4int(2 * size * (sizeof(G4double))))
144 {
145 delete[] value;
146 return false;
147 }
148
149 binVector.reserve(size);
150 dataVector.reserve(size);
151 for(std::size_t i = 0; i < size; ++i)
152 {
153 binVector.push_back(value[2 * i]);
154 dataVector.push_back(value[2 * i + 1]);
155 }
156 delete[] value;
157
158 Initialise();
159 return true;
160}
161
162// --------------------------------------------------------------
164{
165 for(std::size_t i = 0; i < numberOfNodes; ++i)
166 {
167 G4cout << binVector[i] / unitE << " " << dataVector[i] / unitV
168 << G4endl;
169 }
170}
171
172// --------------------------------------------------------------------
173std::size_t G4PhysicsVector::FindBin(const G4double energy,
174 std::size_t idx) const
175{
176 if(idx + 1 < numberOfNodes &&
177 energy >= binVector[idx] && energy <= binVector[idx])
178 {
179 return idx;
180 }
181 if(energy <= binVector[1])
182 {
183 return 0;
184 }
185 if(energy >= binVector[idxmax])
186 {
187 return idxmax;
188 }
189 return GetBin(energy);
190}
191
192// --------------------------------------------------------------------
194 const G4double factorV)
195{
196 for(std::size_t i = 0; i < numberOfNodes; ++i)
197 {
198 binVector[i] *= factorE;
199 dataVector[i] *= factorV;
200 }
201 Initialise();
202}
203
204// --------------------------------------------------------------------
206 const G4double dir1,
207 const G4double dir2)
208{
209 if(!useSpline) { return; }
210 // cannot compute derivatives for less than 5 points
211 const std::size_t nmin = (stype == G4SplineType::Base) ? 5 : 4;
212 if(nmin > numberOfNodes)
213 {
214 if(0 < verboseLevel)
215 {
216 G4cout << "### G4PhysicsVector: spline cannot be used for "
217 << numberOfNodes << " points - spline disabled"
218 << G4endl;
219 DumpValues();
220 }
221 useSpline = false;
222 return;
223 }
224 // check energies of free vector
226 {
227 for(std::size_t i=0; i<=idxmax; ++i)
228 {
229 if(binVector[i + 1] <= binVector[i])
230 {
231 if(0 < verboseLevel)
232 {
233 G4cout << "### G4PhysicsVector: spline cannot be used, because "
234 << " E[" << i << "]=" << binVector[i]
235 << " >= E[" << i+1 << "]=" << binVector[i + 1]
236 << G4endl;
237 DumpValues();
238 }
239 useSpline = false;
240 return;
241 }
242 }
243 }
244
245 // spline is possible
246 Initialise();
248
249 if(1 < verboseLevel)
250 {
251 G4cout << "### G4PhysicsVector:: FillSecondDerivatives N="
252 << numberOfNodes << G4endl;
253 DumpValues();
254 }
255
256 switch(stype)
257 {
258 case G4SplineType::Base:
259 ComputeSecDerivative1();
260 break;
261
262 case G4SplineType::FixedEdges:
263 ComputeSecDerivative2(dir1, dir2);
264 break;
265
266 default:
267 ComputeSecDerivative0();
268 }
269}
270
271// --------------------------------------------------------------
272void G4PhysicsVector::ComputeSecDerivative0()
273// A simplified method of computation of second derivatives
274{
275 std::size_t n = numberOfNodes - 1;
276
277 for(std::size_t i = 1; i < n; ++i)
278 {
279 secDerivative[i] = 3.0 *
280 ((dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) -
281 (dataVector[i] - dataVector[i - 1]) /
282 (binVector[i] - binVector[i - 1])) /
283 (binVector[i + 1] - binVector[i - 1]);
284 }
287}
288
289// --------------------------------------------------------------
290void G4PhysicsVector::ComputeSecDerivative1()
291// Computation of second derivatives using "Not-a-knot" endpoint conditions
292// B.I. Kvasov "Methods of shape-preserving spline approximation"
293// World Scientific, 2000
294{
295 std::size_t n = numberOfNodes - 1;
296 G4double* u = new G4double[n];
297 G4double p, sig;
298
299 u[1] = ((dataVector[2] - dataVector[1]) / (binVector[2] - binVector[1]) -
300 (dataVector[1] - dataVector[0]) / (binVector[1] - binVector[0]));
301 u[1] = 6.0 * u[1] * (binVector[2] - binVector[1]) /
302 ((binVector[2] - binVector[0]) * (binVector[2] - binVector[0]));
303
304 // Decomposition loop for tridiagonal algorithm. secDerivative[i]
305 // and u[i] are used for temporary storage of the decomposed factors.
306
307 secDerivative[1] = (2.0 * binVector[1] - binVector[0] - binVector[2]) /
308 (2.0 * binVector[2] - binVector[0] - binVector[1]);
309
310 for(std::size_t i = 2; i < n - 1; ++i)
311 {
312 sig =
313 (binVector[i] - binVector[i - 1]) / (binVector[i + 1] - binVector[i - 1]);
314 p = sig * secDerivative[i - 1] + 2.0;
315 secDerivative[i] = (sig - 1.0) / p;
316 u[i] =
317 (dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) -
318 (dataVector[i] - dataVector[i - 1]) / (binVector[i] - binVector[i - 1]);
319 u[i] =
320 (6.0 * u[i] / (binVector[i + 1] - binVector[i - 1])) - sig * u[i - 1] / p;
321 }
322
323 sig =
324 (binVector[n - 1] - binVector[n - 2]) / (binVector[n] - binVector[n - 2]);
325 p = sig * secDerivative[n - 3] + 2.0;
326 u[n - 1] =
327 (dataVector[n] - dataVector[n - 1]) / (binVector[n] - binVector[n - 1]) -
328 (dataVector[n - 1] - dataVector[n - 2]) /
329 (binVector[n - 1] - binVector[n - 2]);
330 u[n - 1] = 6.0 * sig * u[n - 1] / (binVector[n] - binVector[n - 2]) -
331 (2.0 * sig - 1.0) * u[n - 2] / p;
332
333 p = (1.0 + sig) + (2.0 * sig - 1.0) * secDerivative[n - 2];
334 secDerivative[n - 1] = u[n - 1] / p;
335
336 // The back-substitution loop for the triagonal algorithm of solving
337 // a linear system of equations.
338
339 for(std::size_t k = n - 2; k > 1; --k)
340 {
341 secDerivative[k] *=
342 (secDerivative[k + 1] - u[k] * (binVector[k + 1] - binVector[k - 1]) /
343 (binVector[k + 1] - binVector[k]));
344 }
346 (secDerivative[n - 1] - (1.0 - sig) * secDerivative[n - 2]) / sig;
347 sig = 1.0 - ((binVector[2] - binVector[1]) / (binVector[2] - binVector[0]));
348 secDerivative[1] *= (secDerivative[2] - u[1] / (1.0 - sig));
349 secDerivative[0] = (secDerivative[1] - sig * secDerivative[2]) / (1.0 - sig);
350
351 delete[] u;
352}
353
354// --------------------------------------------------------------
355void G4PhysicsVector::ComputeSecDerivative2(G4double firstPointDerivative,
356 G4double endPointDerivative)
357// A standard method of computation of second derivatives
358// First derivatives at the first and the last point should be provided
359// See for example W.H. Press et al. "Numerical recipes in C"
360// Cambridge University Press, 1997.
361{
362 std::size_t n = numberOfNodes - 1;
363 G4double* u = new G4double[n];
364 G4double p, sig, un;
365
366 u[0] = (6.0 / (binVector[1] - binVector[0])) *
367 ((dataVector[1] - dataVector[0]) / (binVector[1] - binVector[0]) -
368 firstPointDerivative);
369
370 secDerivative[0] = -0.5;
371
372 // Decomposition loop for tridiagonal algorithm. secDerivative[i]
373 // and u[i] are used for temporary storage of the decomposed factors.
374
375 for(std::size_t i = 1; i < n; ++i)
376 {
377 sig =
378 (binVector[i] - binVector[i - 1]) / (binVector[i + 1] - binVector[i - 1]);
379 p = sig * (secDerivative[i - 1]) + 2.0;
380 secDerivative[i] = (sig - 1.0) / p;
381 u[i] =
382 (dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) -
383 (dataVector[i] - dataVector[i - 1]) / (binVector[i] - binVector[i - 1]);
384 u[i] =
385 6.0 * u[i] / (binVector[i + 1] - binVector[i - 1]) - sig * u[i - 1] / p;
386 }
387
388 sig =
389 (binVector[n - 1] - binVector[n - 2]) / (binVector[n] - binVector[n - 2]);
390 p = sig * secDerivative[n - 2] + 2.0;
391 un = (6.0 / (binVector[n] - binVector[n - 1])) *
392 (endPointDerivative - (dataVector[n] - dataVector[n - 1]) /
393 (binVector[n] - binVector[n - 1])) -
394 u[n - 1] / p;
395 secDerivative[n] = un / (secDerivative[n - 1] + 2.0);
396
397 // The back-substitution loop for the triagonal algorithm of solving
398 // a linear system of equations.
399
400 for(std::size_t k = n - 1; k > 0; --k)
401 {
402 secDerivative[k] *=
403 (secDerivative[k + 1] - u[k] * (binVector[k + 1] - binVector[k - 1]) /
404 (binVector[k + 1] - binVector[k]));
405 }
406 secDerivative[0] = 0.5 * (u[0] - secDerivative[1]);
407
408 delete[] u;
409}
410
411// --------------------------------------------------------------
412std::ostream& operator<<(std::ostream& out, const G4PhysicsVector& pv)
413{
414 // binning
415 G4long prec = out.precision();
416 out << std::setprecision(12) << pv.edgeMin << " " << pv.edgeMax << " "
417 << pv.numberOfNodes << G4endl;
418
419 // contents
420 out << pv.dataVector.size() << G4endl;
421 for(std::size_t i = 0; i < pv.dataVector.size(); ++i)
422 {
423 out << pv.binVector[i] << " " << pv.dataVector[i] << G4endl;
424 }
425 out.precision(prec);
426
427 return out;
428}
429
430//---------------------------------------------------------------
432{
433 if(0 == numberOfNodes)
434 {
435 return 0.0;
436 }
437 if(1 == numberOfNodes || val <= dataVector[0])
438 {
439 return edgeMin;
440 }
441 if(val >= dataVector[numberOfNodes - 1])
442 {
443 return edgeMax;
444 }
445 std::size_t bin = std::lower_bound(dataVector.cbegin(), dataVector.cend(), val)
446 - dataVector.cbegin() - 1;
447 if(bin > idxmax) { bin = idxmax; }
448 G4double res = binVector[bin];
449 G4double del = dataVector[bin + 1] - dataVector[bin];
450 if(del > 0.0)
451 {
452 res += (val - dataVector[bin]) * (binVector[bin + 1] - res) / del;
453 }
454 return res;
455}
456
457//---------------------------------------------------------------
459 G4double val,
460 const G4String& text)
461{
463 ed << "Vector type: " << type << " length= " << numberOfNodes
464 << "; an attempt to put data at index= " << index
465 << " value= " << val << " in " << text;
466 G4Exception("G4PhysicsVector:", "gl0005",
467 FatalException, ed, "Wrong operation");
468}
469
470//---------------------------------------------------------------
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ T_G4PhysicsFreeVector
std::ostream & operator<<(std::ostream &out, const G4PhysicsVector &pv)
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
std::size_t idxmax
G4PhysicsVectorType type
G4double GetEnergy(const G4double value) const
void PrintPutValueError(std::size_t index, G4double value, const G4String &text)
void ScaleVector(const G4double factorE, const G4double factorV)
G4bool Store(std::ofstream &fOut, G4bool ascii=false) const
std::size_t numberOfNodes
std::vector< G4double > secDerivative
G4PhysicsVector(G4bool spline=false)
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
std::vector< G4double > dataVector
std::vector< G4double > binVector
virtual void Initialise()
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const
std::size_t FindBin(const G4double energy, std::size_t idx) const