Geant4 10.7.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// - 02 Apr. 2008, A.Bagulya: Added SplineInterpolation() and SetSpline()
34// - 19 Jun. 2009, V.Ivanchenko: Removed hidden bin
35// - 15 Mar. 2019 M.Novak: added Value method with the known log-energy value
36// that can avoid the log call in case of log-vectors
37// --------------------------------------------------------------------
38
39#include "G4PhysicsVector.hh"
40#include <iomanip>
41
42// --------------------------------------------------------------
43
45 : useSpline(val)
46{}
47
48// --------------------------------------------------------------
49
51
52// --------------------------------------------------------------
53
55{
56 invdBin = right.invdBin;
57 baseBin = right.baseBin;
59
60 DeleteData();
61 CopyData(right);
62}
63
64// --------------------------------------------------------------
65
67{
68 if(&right == this)
69 {
70 return *this;
71 }
72 invdBin = right.invdBin;
73 baseBin = right.baseBin;
75
76 DeleteData();
77 CopyData(right);
78 return *this;
79}
80
81// --------------------------------------------------------------
82
84{
85 return (this == &right);
86}
87
88// --------------------------------------------------------------
89
91{
92 return (this != &right);
93}
94
95// --------------------------------------------------------------
96
98{
99 useSpline = false;
100 secDerivative.clear();
101}
102
103// --------------------------------------------------------------
104
106{
107 type = vec.type;
108 edgeMin = vec.edgeMin;
109 edgeMax = vec.edgeMax;
111 useSpline = vec.useSpline;
112
113 std::size_t i;
115 for(i = 0; i < numberOfNodes; ++i)
116 {
117 dataVector[i] = (vec.dataVector)[i];
118 }
119 binVector.resize(numberOfNodes);
120 for(i = 0; i < numberOfNodes; ++i)
121 {
122 binVector[i] = (vec.binVector)[i];
123 }
124 if(0 < (vec.secDerivative).size())
125 {
127 for(i = 0; i < numberOfNodes; ++i)
128 {
129 secDerivative[i] = (vec.secDerivative)[i];
130 }
131 }
132}
133
134// --------------------------------------------------------------
135
136G4double G4PhysicsVector::GetLowEdgeEnergy(std::size_t binNumber) const
137{
138 return binVector[binNumber];
139}
140
141// --------------------------------------------------------------
142
143G4bool G4PhysicsVector::Store(std::ofstream& fOut, G4bool ascii) const
144{
145 // Ascii mode
146 if(ascii)
147 {
148 fOut << *this;
149 return true;
150 }
151 // Binary Mode
152
153 // binning
154 fOut.write((char*) (&edgeMin), sizeof edgeMin);
155 fOut.write((char*) (&edgeMax), sizeof edgeMax);
156 fOut.write((char*) (&numberOfNodes), sizeof numberOfNodes);
157
158 // contents
159 std::size_t size = dataVector.size();
160 fOut.write((char*) (&size), sizeof size);
161
162 G4double* value = new G4double[2 * size];
163 for(std::size_t i = 0; i < size; ++i)
164 {
165 value[2 * i] = binVector[i];
166 value[2 * i + 1] = dataVector[i];
167 }
168 fOut.write((char*) (value), 2 * size * (sizeof(G4double)));
169 delete[] value;
170
171 return true;
172}
173
174// --------------------------------------------------------------
175
176G4bool G4PhysicsVector::Retrieve(std::ifstream& fIn, G4bool ascii)
177{
178 // clear properties;
179 dataVector.clear();
180 binVector.clear();
181 secDerivative.clear();
182
183 // retrieve in ascii mode
184 if(ascii)
185 {
186 // binning
187 fIn >> edgeMin >> edgeMax >> numberOfNodes;
188 if(fIn.fail() || numberOfNodes < 2)
189 {
190 return false;
191 }
192 // contents
193 G4int siz = 0;
194 fIn >> siz;
195 if(fIn.fail() || siz != G4int(numberOfNodes))
196 {
197 return false;
198 }
199
200 binVector.reserve(siz);
201 dataVector.reserve(siz);
202 G4double vBin, vData;
203
204 for(G4int i = 0; i < siz; ++i)
205 {
206 vBin = 0.;
207 vData = 0.;
208 fIn >> vBin >> vData;
209 if(fIn.fail())
210 {
211 return false;
212 }
213 binVector.push_back(vBin);
214 dataVector.push_back(vData);
215 }
216
217 // to remove any inconsistency
218 numberOfNodes = siz;
219 edgeMin = binVector[0];
221 return true;
222 }
223
224 // retrieve in binary mode
225 // binning
226 fIn.read((char*) (&edgeMin), sizeof edgeMin);
227 fIn.read((char*) (&edgeMax), sizeof edgeMax);
228 fIn.read((char*) (&numberOfNodes), sizeof numberOfNodes);
229
230 // contents
231 std::size_t size;
232 fIn.read((char*) (&size), sizeof size);
233
234 G4double* value = new G4double[2 * size];
235 fIn.read((char*) (value), 2 * size * (sizeof(G4double)));
236 if(G4int(fIn.gcount()) != G4int(2 * size * (sizeof(G4double))))
237 {
238 delete[] value;
239 return false;
240 }
241
242 binVector.reserve(size);
243 dataVector.reserve(size);
244 for(std::size_t i = 0; i < size; ++i)
245 {
246 binVector.push_back(value[2 * i]);
247 dataVector.push_back(value[2 * i + 1]);
248 }
249 delete[] value;
250
251 // to remove any inconsistency
252 numberOfNodes = size;
253 edgeMin = binVector[0];
255
256 return true;
257}
258
259// --------------------------------------------------------------
260
262{
263 for(std::size_t i = 0; i < numberOfNodes; ++i)
264 {
265 G4cout << binVector[i] / unitE << " " << dataVector[i] / unitV << G4endl;
266 }
267}
268
269// --------------------------------------------------------------------
270
272{
273 std::size_t n = dataVector.size();
274 std::size_t i;
275 for(i = 0; i < n; ++i)
276 {
277 binVector[i] *= factorE;
278 dataVector[i] *= factorV;
279 }
280 secDerivative.clear();
281
282 edgeMin = binVector[0];
283 edgeMax = binVector[n - 1];
284}
285
286// --------------------------------------------------------------
287
289 G4double endPointDerivative)
290// A standard method of computation of second derivatives
291// First derivatives at the first and the last point should be provided
292// See for example W.H. Press et al. "Numerical recipes in C"
293// Cambridge University Press, 1997.
294{
295 if(4 > numberOfNodes) // cannot compute derivatives for less than 4 bins
296 {
298 return;
299 }
300
301 if(!SplinePossible())
302 {
303 return;
304 }
305
306 useSpline = true;
307
308 G4int n = numberOfNodes - 1;
309
310 G4double* u = new G4double[n];
311
312 G4double p, sig, un;
313
314 u[0] = (6.0 / (binVector[1] - binVector[0])) *
315 ((dataVector[1] - dataVector[0]) / (binVector[1] - binVector[0]) -
316 firstPointDerivative);
317
318 secDerivative[0] = -0.5;
319
320 // Decomposition loop for tridiagonal algorithm. secDerivative[i]
321 // and u[i] are used for temporary storage of the decomposed factors.
322
323 for(G4int i = 1; i < n; ++i)
324 {
325 sig =
326 (binVector[i] - binVector[i - 1]) / (binVector[i + 1] - binVector[i - 1]);
327 p = sig * (secDerivative[i - 1]) + 2.0;
328 secDerivative[i] = (sig - 1.0) / p;
329 u[i] =
330 (dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) -
331 (dataVector[i] - dataVector[i - 1]) / (binVector[i] - binVector[i - 1]);
332 u[i] =
333 6.0 * u[i] / (binVector[i + 1] - binVector[i - 1]) - sig * u[i - 1] / p;
334 }
335
336 sig =
337 (binVector[n - 1] - binVector[n - 2]) / (binVector[n] - binVector[n - 2]);
338 p = sig * secDerivative[n - 2] + 2.0;
339 un = (6.0 / (binVector[n] - binVector[n - 1])) *
340 (endPointDerivative - (dataVector[n] - dataVector[n - 1]) /
341 (binVector[n] - binVector[n - 1])) -
342 u[n - 1] / p;
343 secDerivative[n] = un / (secDerivative[n - 1] + 2.0);
344
345 // The back-substitution loop for the triagonal algorithm of solving
346 // a linear system of equations.
347
348 for(G4int k = n - 1; k > 0; --k)
349 {
350 secDerivative[k] *=
351 (secDerivative[k + 1] - u[k] * (binVector[k + 1] - binVector[k - 1]) /
352 (binVector[k + 1] - binVector[k]));
353 }
354 secDerivative[0] = 0.5 * (u[0] - secDerivative[1]);
355
356 delete[] u;
357}
358
359// --------------------------------------------------------------
360
362{
363 // Computation of second derivatives using "Not-a-knot" endpoint conditions
364 // B.I. Kvasov "Methods of shape-preserving spline approximation"
365 // World Scientific, 2000
366
367 if(5 > numberOfNodes) // cannot compute derivatives for less than 4 points
368 {
370 return;
371 }
372
373 if(!SplinePossible())
374 {
375 return;
376 }
377
378 useSpline = true;
379
380 G4int n = numberOfNodes - 1;
381
382 G4double* u = new G4double[n];
383
384 G4double p, sig;
385
386 u[1] = ((dataVector[2] - dataVector[1]) / (binVector[2] - binVector[1]) -
387 (dataVector[1] - dataVector[0]) / (binVector[1] - binVector[0]));
388 u[1] = 6.0 * u[1] * (binVector[2] - binVector[1]) /
389 ((binVector[2] - binVector[0]) * (binVector[2] - binVector[0]));
390
391 // Decomposition loop for tridiagonal algorithm. secDerivative[i]
392 // and u[i] are used for temporary storage of the decomposed factors.
393
394 secDerivative[1] = (2.0 * binVector[1] - binVector[0] - binVector[2]) /
395 (2.0 * binVector[2] - binVector[0] - binVector[1]);
396
397 for(G4int i = 2; i < n - 1; ++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 - 3] + 2.0;
413 u[n - 1] =
414 (dataVector[n] - dataVector[n - 1]) / (binVector[n] - binVector[n - 1]) -
415 (dataVector[n - 1] - dataVector[n - 2]) /
416 (binVector[n - 1] - binVector[n - 2]);
417 u[n - 1] = 6.0 * sig * u[n - 1] / (binVector[n] - binVector[n - 2]) -
418 (2.0 * sig - 1.0) * u[n - 2] / p;
419
420 p = (1.0 + sig) + (2.0 * sig - 1.0) * secDerivative[n - 2];
421 secDerivative[n - 1] = u[n - 1] / p;
422
423 // The back-substitution loop for the triagonal algorithm of solving
424 // a linear system of equations.
425
426 for(G4int k = n - 2; k > 1; --k)
427 {
428 secDerivative[k] *=
429 (secDerivative[k + 1] - u[k] * (binVector[k + 1] - binVector[k - 1]) /
430 (binVector[k + 1] - binVector[k]));
431 }
432 secDerivative[n] =
433 (secDerivative[n - 1] - (1.0 - sig) * secDerivative[n - 2]) / sig;
434 sig = 1.0 - ((binVector[2] - binVector[1]) / (binVector[2] - binVector[0]));
435 secDerivative[1] *= (secDerivative[2] - u[1] / (1.0 - sig));
436 secDerivative[0] = (secDerivative[1] - sig * secDerivative[2]) / (1.0 - sig);
437
438 delete[] u;
439}
440
441// --------------------------------------------------------------
442
444{
445 // A simplified method of computation of second derivatives
446
447 if(3 > numberOfNodes) // cannot compute derivatives for less than 4 bins
448 {
449 useSpline = false;
450 return;
451 }
452
453 if(!SplinePossible())
454 {
455 return;
456 }
457
458 useSpline = true;
459
460 std::size_t n = numberOfNodes - 1;
461
462 for(std::size_t i = 1; i < n; ++i)
463 {
464 secDerivative[i] =
465 3.0 *
466 ((dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) -
467 (dataVector[i] - dataVector[i - 1]) /
468 (binVector[i] - binVector[i - 1])) /
469 (binVector[i + 1] - binVector[i - 1]);
470 }
471 secDerivative[n] = secDerivative[n - 1];
473}
474
475// --------------------------------------------------------------
476
477G4bool G4PhysicsVector::SplinePossible()
478{
479 // Initialise second derivative array. If neighbor energy coincide
480 // or not ordered than spline cannot be applied
481
482 G4bool result = true;
483 for(std::size_t j = 1; j < numberOfNodes; ++j)
484 {
485 if(binVector[j] <= binVector[j - 1])
486 {
487 result = false;
488 useSpline = false;
489 secDerivative.clear();
490 break;
491 }
492 }
493 secDerivative.resize(numberOfNodes, 0.0);
494 return result;
495}
496
497// --------------------------------------------------------------
498
499std::ostream& operator<<(std::ostream& out, const G4PhysicsVector& pv)
500{
501 // binning
502 G4int prec = out.precision();
503 out << std::setprecision(12) << pv.edgeMin << " " << pv.edgeMax << " "
504 << pv.numberOfNodes << G4endl;
505
506 // contents
507 out << pv.dataVector.size() << G4endl;
508 for(std::size_t i = 0; i < pv.dataVector.size(); ++i)
509 {
510 out << pv.binVector[i] << " " << pv.dataVector[i] << G4endl;
511 }
512 out << std::setprecision(prec);
513
514 return out;
515}
516
517//---------------------------------------------------------------
518
519G4double G4PhysicsVector::Value(G4double theEnergy, std::size_t& lastIdx) const
520{
521 G4double y;
522 if(theEnergy <= edgeMin)
523 {
524 lastIdx = 0;
525 y = dataVector[0];
526 }
527 else if(theEnergy >= edgeMax)
528 {
529 lastIdx = numberOfNodes - 1;
530 y = dataVector[lastIdx];
531 }
532 else
533 {
534 lastIdx = FindBin(theEnergy, lastIdx);
535 y = Interpolation(lastIdx, theEnergy);
536 }
537 return y;
538}
539
540//---------------------------------------------------------------
541
543{
544 if(1 >= numberOfNodes)
545 {
546 return 0.0;
547 }
548 G4double y = rand * dataVector[numberOfNodes - 1];
549 std::size_t bin = std::lower_bound(dataVector.begin(), dataVector.end(), y) -
550 dataVector.begin() - 1;
551 bin = std::min(bin, numberOfNodes - 2);
552 G4double res = binVector[bin];
553 G4double del = dataVector[bin + 1] - dataVector[bin];
554 if(del > 0.0)
555 {
556 res += (y - dataVector[bin]) * (binVector[bin + 1] - res) / del;
557 }
558 return res;
559}
560
561//---------------------------------------------------------------
562
564{
566 ed << "Vector type " << type << " length= " << numberOfNodes
567 << " an attempt to put data at index= " << index;
568 G4Exception("G4PhysicsVector::PutValue()", "gl0005", FatalException, ed,
569 "Memory overwritten");
570}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
std::ostream & operator<<(std::ostream &out, const G4PhysicsVector &pv)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4PVDataVector binVector
G4bool operator!=(const G4PhysicsVector &right) const
void ComputeSecondDerivatives(G4double firstPointDerivative, G4double endPointDerivative)
G4PhysicsVectorType type
G4double Value(G4double theEnergy, std::size_t &lastidx) const
std::size_t FindBin(const G4double energy, const std::size_t idx) const
virtual void ScaleVector(G4double factorE, G4double factorV)
G4double GetLowEdgeEnergy(std::size_t binNumber) const
G4bool Store(std::ofstream &fOut, G4bool ascii=false) const
G4double FindLinearEnergy(G4double rand) const
G4bool operator==(const G4PhysicsVector &right) const
G4PVDataVector secDerivative
std::size_t numberOfNodes
virtual ~G4PhysicsVector()
void ComputeSecDerivatives()
G4PhysicsVector(G4bool spline=false)
void FillSecondDerivatives()
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4PVDataVector dataVector
void CopyData(const G4PhysicsVector &vec)
void PrintPutValueError(std::size_t index)
G4PhysicsVector & operator=(const G4PhysicsVector &)
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const