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