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