Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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//
27// $Id$
28//
29//
30// --------------------------------------------------------------
31// GEANT 4 class implementation file
32//
33// G4PhysicsVector.cc
34//
35// History:
36// 02 Dec. 1995, G.Cosmo : Structure created based on object model
37// 03 Mar. 1996, K.Amako : Implemented the 1st version
38// 01 Jul. 1996, K.Amako : Hidden bin from the user introduced
39// 12 Nov. 1998, K.Amako : A bug in GetVectorLength() fixed
40// 11 Nov. 2000, H.Kurashige : use STL vector for dataVector and binVector
41// 18 Jan. 2001, H.Kurashige : removed ptrNextTable
42// 09 Mar. 2001, H.Kurashige : added G4PhysicsVector type
43// 05 Sep. 2008, V.Ivanchenko : added protections for zero-length vector
44// 11 May 2009, A.Bagulya : added new implementation of methods
45// ComputeSecondDerivatives - first derivatives at edge points
46// should be provided by a user
47// FillSecondDerivatives - default computation base on "not-a-knot"
48// algorithm
49// 19 Jun. 2009, V.Ivanchenko : removed hidden bin
50// 17 Nov. 2009, H.Kurashige : use pointer for DataVector
51// 04 May 2010 H.Kurashige : use G4PhyscisVectorCache
52// 28 May 2010 H.Kurashige : Stop using pointers to G4PVDataVector
53// 16 Aug. 2011 H.Kurashige : Add dBin, baseBin and verboseLevel
54// --------------------------------------------------------------
55
56#include "G4PhysicsVector.hh"
57#include <iomanip>
58
60
61// --------------------------------------------------------------
62
64 : type(T_G4PhysicsVector),
65 edgeMin(0.), edgeMax(0.), numberOfNodes(0),
66 useSpline(spline),
67 dBin(0.), baseBin(0.),
68 verboseLevel(0)
69{
71}
72
73// --------------------------------------------------------------
74
76{
77 delete cache; cache =0;
78}
79
80// --------------------------------------------------------------
81
83{
85 dBin = right.dBin;
86 baseBin = right.baseBin;
88
89 DeleteData();
90 CopyData(right);
91}
92
93// --------------------------------------------------------------
94
96{
97 if (&right==this) { return *this; }
98 dBin = right.dBin;
99 baseBin = right.baseBin;
101
102 DeleteData();
103 CopyData(right);
104 return *this;
105}
106
107// --------------------------------------------------------------
108
110{
111 return (this == &right);
112}
113
114// --------------------------------------------------------------
115
117{
118 return (this != &right);
119}
120
121// --------------------------------------------------------------
122
124{
125 secDerivative.clear();
126}
127
128// --------------------------------------------------------------
129
131{
132 type = vec.type;
133 edgeMin = vec.edgeMin;
134 edgeMax = vec.edgeMax;
138 cache->lastBin = vec.GetLastBin();
139 useSpline = vec.useSpline;
140
141 size_t i;
142 dataVector.clear();
143 for(i=0; i<(vec.dataVector).size(); i++){
144 dataVector.push_back( (vec.dataVector)[i] );
145 }
146 binVector.clear();
147 for(i=0; i<(vec.binVector).size(); i++){
148 binVector.push_back( (vec.binVector)[i] );
149 }
150 secDerivative.clear();
151 for(i=0; i<(vec.secDerivative).size(); i++){
152 secDerivative.push_back( (vec.secDerivative)[i] );
153 }
154}
155
156// --------------------------------------------------------------
157
159{
160 return binVector[binNumber];
161}
162
163// --------------------------------------------------------------
164
165G4bool G4PhysicsVector::Store(std::ofstream& fOut, G4bool ascii)
166{
167 // Ascii mode
168 if (ascii)
169 {
170 fOut << *this;
171 return true;
172 }
173 // Binary Mode
174
175 // binning
176 fOut.write((char*)(&edgeMin), sizeof edgeMin);
177 fOut.write((char*)(&edgeMax), sizeof edgeMax);
178 fOut.write((char*)(&numberOfNodes), sizeof numberOfNodes);
179
180 // contents
181 size_t size = dataVector.size();
182 fOut.write((char*)(&size), sizeof size);
183
184 G4double* value = new G4double[2*size];
185 for(size_t i = 0; i < size; ++i)
186 {
187 value[2*i] = binVector[i];
188 value[2*i+1]= dataVector[i];
189 }
190 fOut.write((char*)(value), 2*size*(sizeof (G4double)));
191 delete [] value;
192
193 return true;
194}
195
196// --------------------------------------------------------------
197
198G4bool G4PhysicsVector::Retrieve(std::ifstream& fIn, G4bool ascii)
199{
200 // clear properties;
202 cache->lastValue =0.;
203 cache->lastBin =0;
204 dataVector.clear();
205 binVector.clear();
206 secDerivative.clear();
207
208 // retrieve in ascii mode
209 if (ascii){
210 // binning
211 fIn >> edgeMin >> edgeMax >> numberOfNodes;
212 if (fIn.fail()) { return false; }
213 // contents
214 G4int siz=0;
215 fIn >> siz;
216 if (fIn.fail()) { return false; }
217 if (siz<=0)
218 {
219#ifdef G4VERBOSE
220 G4cerr << "G4PhysicsVector::Retrieve():";
221 G4cerr << " Invalid vector size: " << siz << G4endl;
222#endif
223 return false;
224 }
225
226 binVector.reserve(siz);
227 dataVector.reserve(siz);
228 G4double vBin, vData;
229
230 for(G4int i = 0; i < siz ; i++)
231 {
232 vBin = 0.;
233 vData= 0.;
234 fIn >> vBin >> vData;
235 if (fIn.fail()) { return false; }
236 binVector.push_back(vBin);
237 dataVector.push_back(vData);
238 }
239
240 // to remove any inconsistency
241 numberOfNodes = siz;
242 edgeMin = binVector[0];
244 return true ;
245 }
246
247 // retrieve in binary mode
248 // binning
249 fIn.read((char*)(&edgeMin), sizeof edgeMin);
250 fIn.read((char*)(&edgeMax), sizeof edgeMax);
251 fIn.read((char*)(&numberOfNodes), sizeof numberOfNodes );
252
253 // contents
254 size_t size;
255 fIn.read((char*)(&size), sizeof size);
256
257 G4double* value = new G4double[2*size];
258 fIn.read((char*)(value), 2*size*(sizeof(G4double)) );
259 if (G4int(fIn.gcount()) != G4int(2*size*(sizeof(G4double))) )
260 {
261 delete [] value;
262 return false;
263 }
264
265 binVector.reserve(size);
266 dataVector.reserve(size);
267 for(size_t i = 0; i < size; ++i)
268 {
269 binVector.push_back(value[2*i]);
270 dataVector.push_back(value[2*i+1]);
271 }
272 delete [] value;
273
274 // to remove any inconsistency
275 numberOfNodes = size;
276 edgeMin = binVector[0];
278
279 return true;
280}
281
282// --------------------------------------------------------------
283
284void
286{
287 size_t n = dataVector.size();
288 size_t i;
289 if(n > 0) {
290 for(i=0; i<n; ++i) {
291 binVector[i] *= factorE;
292 dataVector[i] *= factorV;
293 }
294 }
295 // n = secDerivative.size();
296 // if(n > 0) { for(i=0; i<n; ++i) { secDerivative[i] *= factorV; } }
297 secDerivative.clear();
298
299 edgeMin *= factorE;
300 edgeMax *= factorE;
301 cache->lastEnergy = factorE*(cache->lastEnergy);
302 cache->lastValue = factorV*(cache->lastValue);
303}
304
305// --------------------------------------------------------------
306
307void
309 G4double endPointDerivative)
310 // A standard method of computation of second derivatives
311 // First derivatives at the first and the last point should be provided
312 // See for example W.H. Press et al. "Numerical recipes in C"
313 // Cambridge University Press, 1997.
314{
315 if(4 > numberOfNodes) // cannot compute derivatives for less than 4 bins
316 {
318 return;
319 }
320
321 if(!SplinePossible()) { return; }
322
323 G4int n = numberOfNodes-1;
324
325 G4double* u = new G4double [n];
326
327 G4double p, sig, un;
328
329 u[0] = (6.0/(binVector[1]-binVector[0]))
330 * ((dataVector[1]-dataVector[0])/(binVector[1]-binVector[0])
331 - firstPointDerivative);
332
333 secDerivative[0] = - 0.5;
334
335 // Decomposition loop for tridiagonal algorithm. secDerivative[i]
336 // and u[i] are used for temporary storage of the decomposed factors.
337
338 for(G4int i=1; i<n; ++i)
339 {
340 sig = (binVector[i]-binVector[i-1]) / (binVector[i+1]-binVector[i-1]);
341 p = sig*(secDerivative[i-1]) + 2.0;
342 secDerivative[i] = (sig - 1.0)/p;
343 u[i] = (dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i])
344 - (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]);
345 u[i] = 6.0*u[i]/(binVector[i+1]-binVector[i-1]) - sig*u[i-1]/p;
346 }
347
348 sig = (binVector[n-1]-binVector[n-2]) / (binVector[n]-binVector[n-2]);
349 p = sig*secDerivative[n-2] + 2.0;
350 un = (6.0/(binVector[n]-binVector[n-1]))
351 *(endPointDerivative -
352 (dataVector[n]-dataVector[n-1])/(binVector[n]-binVector[n-1])) - u[n-1]/p;
353 secDerivative[n] = un/(secDerivative[n-1] + 2.0);
354
355 // The back-substitution loop for the triagonal algorithm of solving
356 // a linear system of equations.
357
358 for(G4int k=n-1; k>0; --k)
359 {
360 secDerivative[k] *=
361 (secDerivative[k+1] -
362 u[k]*(binVector[k+1]-binVector[k-1])/(binVector[k+1]-binVector[k]));
363 }
364 secDerivative[0] = 0.5*(u[0] - secDerivative[1]);
365
366 delete [] u;
367}
368
369// --------------------------------------------------------------
370
372 // Computation of second derivatives using "Not-a-knot" endpoint conditions
373 // B.I. Kvasov "Methods of shape-preserving spline approximation"
374 // World Scientific, 2000
375{
376 if(5 > numberOfNodes) // cannot compute derivatives for less than 4 points
377 {
379 return;
380 }
381
382 if(!SplinePossible()) { return; }
383
384 G4int n = numberOfNodes-1;
385
386 G4double* u = new G4double [n];
387
388 G4double p, sig;
389
390 u[1] = ((dataVector[2]-dataVector[1])/(binVector[2]-binVector[1]) -
391 (dataVector[1]-dataVector[0])/(binVector[1]-binVector[0]));
392 u[1] = 6.0*u[1]*(binVector[2]-binVector[1])
393 / ((binVector[2]-binVector[0])*(binVector[2]-binVector[0]));
394
395 // Decomposition loop for tridiagonal algorithm. secDerivative[i]
396 // and u[i] are used for temporary storage of the decomposed factors.
397
398 secDerivative[1] = (2.0*binVector[1]-binVector[0]-binVector[2])
399 / (2.0*binVector[2]-binVector[0]-binVector[1]);
400
401 for(G4int i=2; i<n-1; ++i)
402 {
403 sig = (binVector[i]-binVector[i-1]) / (binVector[i+1]-binVector[i-1]);
404 p = sig*secDerivative[i-1] + 2.0;
405 secDerivative[i] = (sig - 1.0)/p;
406 u[i] = (dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i])
407 - (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]);
408 u[i] = (6.0*u[i]/(binVector[i+1]-binVector[i-1])) - sig*u[i-1]/p;
409 }
410
411 sig = (binVector[n-1]-binVector[n-2]) / (binVector[n]-binVector[n-2]);
412 p = sig*secDerivative[n-3] + 2.0;
413 u[n-1] = (dataVector[n]-dataVector[n-1])/(binVector[n]-binVector[n-1])
414 - (dataVector[n-1]-dataVector[n-2])/(binVector[n-1]-binVector[n-2]);
415 u[n-1] = 6.0*sig*u[n-1]/(binVector[n]-binVector[n-2])
416 - (2.0*sig - 1.0)*u[n-2]/p;
417
418 p = (1.0+sig) + (2.0*sig-1.0)*secDerivative[n-2];
419 secDerivative[n-1] = u[n-1]/p;
420
421 // The back-substitution loop for the triagonal algorithm of solving
422 // a linear system of equations.
423
424 for(G4int k=n-2; k>1; --k)
425 {
426 secDerivative[k] *=
427 (secDerivative[k+1] -
428 u[k]*(binVector[k+1]-binVector[k-1])/(binVector[k+1]-binVector[k]));
429 }
430 secDerivative[n] = (secDerivative[n-1] - (1.0-sig)*secDerivative[n-2])/sig;
431 sig = 1.0 - ((binVector[2]-binVector[1])/(binVector[2]-binVector[0]));
432 secDerivative[1] *= (secDerivative[2] - u[1]/(1.0-sig));
433 secDerivative[0] = (secDerivative[1] - sig*secDerivative[2])/(1.0-sig);
434
435 delete [] u;
436}
437
438// --------------------------------------------------------------
439
440void
442 // A simplified method of computation of second derivatives
443{
444 if(!SplinePossible()) { return; }
445
446 if(3 > numberOfNodes) // cannot compute derivatives for less than 4 bins
447 {
448 useSpline = false;
449 return;
450 }
451
452 size_t n = numberOfNodes-1;
453
454 for(size_t i=1; i<n; ++i)
455 {
456 secDerivative[i] =
457 3.0*((dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i]) -
458 (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]))
459 /(binVector[i+1]-binVector[i-1]);
460 }
463}
464
465// --------------------------------------------------------------
466
467G4bool G4PhysicsVector::SplinePossible()
468 // Initialise second derivative array. If neighbor energy coincide
469 // or not ordered than spline cannot be applied
470{
471 secDerivative.clear();
472 if(!useSpline) { return useSpline; }
474 for(size_t j=0; j<numberOfNodes; ++j)
475 {
476 secDerivative.push_back(0.0);
477 if(j > 0)
478 {
479 if(binVector[j]-binVector[j-1] <= 0.) { useSpline = false; }
480 }
481 }
482 return useSpline;
483}
484
485// --------------------------------------------------------------
486
487std::ostream& operator<<(std::ostream& out, const G4PhysicsVector& pv)
488{
489 // binning
490 out << std::setprecision(12) << pv.edgeMin << " "
491 << pv.edgeMax << " " << pv.numberOfNodes << G4endl;
492
493 // contents
494 out << pv.dataVector.size() << G4endl;
495 for(size_t i = 0; i < pv.dataVector.size(); i++)
496 {
497 out << pv.binVector[i] << " " << pv.dataVector[i] << G4endl;
498 }
499 out << std::setprecision(6);
500
501 return out;
502}
503
504//---------------------------------------------------------------
505
506void G4PhysicsVector::ComputeValue(G4double theEnergy)
507{
508 // Use cache for speed up - check if the value 'theEnergy' lies
509 // between the last energy and low edge of of the
510 // bin of last call, then the last bin location is used.
511
512 if( theEnergy < cache->lastEnergy
513 && theEnergy >= binVector[cache->lastBin]) {
514 cache->lastEnergy = theEnergy;
515 Interpolation(cache->lastBin);
516
517 } else if( theEnergy <= edgeMin ) {
518 cache->lastBin = 0;
521
522 } else if( theEnergy >= edgeMax ) {
526
527 } else {
528 cache->lastBin = FindBinLocation(theEnergy);
529 cache->lastEnergy = theEnergy;
530 Interpolation(cache->lastBin);
531 }
532}
533
@ T_G4PhysicsVector
G4Allocator< G4PhysicsVector > aPVAllocator
std::ostream & operator<<(std::ostream &out, const G4PhysicsVector &pv)
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4PVDataVector binVector
void ComputeSecondDerivatives(G4double firstPointDerivative, G4double endPointDerivative)
G4PhysicsVectorType type
size_t GetLastBin() const
virtual size_t FindBinLocation(G4double theEnergy) const =0
virtual void ScaleVector(G4double factorE, G4double factorV)
G4PhysicsVectorCache * cache
G4PVDataVector secDerivative
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)
virtual G4double GetLowEdgeEnergy(size_t binNumber) const
G4int operator==(const G4PhysicsVector &right) const
G4double GetLastValue() const
G4int operator!=(const G4PhysicsVector &right) const
G4PhysicsVector & operator=(const G4PhysicsVector &)
G4double GetLastEnergy() const
virtual G4bool Store(std::ofstream &fOut, G4bool ascii=false)
#define DBL_MAX
Definition: templates.hh:83