source: trunk/source/global/management/src/G4PhysicsVector.cc @ 1250

Last change on this file since 1250 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

File size: 13.8 KB
Line 
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,v 1.41 2009/12/07 09:21:27 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-03 $
29//
30//
31// --------------------------------------------------------------
32//      GEANT 4 class implementation file
33//
34//  G4PhysicsVector.cc
35//
36//  History:
37//    02 Dec. 1995, G.Cosmo : Structure created based on object model
38//    03 Mar. 1996, K.Amako : Implemented the 1st version
39//    01 Jul. 1996, K.Amako : Hidden bin from the user introduced
40//    12 Nov. 1998, K.Amako : A bug in GetVectorLength() fixed
41//    11 Nov. 2000, H.Kurashige : use STL vector for dataVector and binVector
42//    18 Jan. 2001, H.Kurashige : removed ptrNextTable
43//    09 Mar. 2001, H.Kurashige : added G4PhysicsVector type
44//    05 Sep. 2008, V.Ivanchenko : added protections for zero-length vector
45//    11 May  2009, A.Bagulya : added new implementation of methods
46//            ComputeSecondDerivatives - first derivatives at edge points
47//                                       should be provided by a user
48//            FillSecondDerivatives - default computation base on "not-a-knot"
49//                                    algorithm
50//    19 Jun. 2009, V.Ivanchenko : removed hidden bin
51// --------------------------------------------------------------
52
53#include "G4PhysicsVector.hh"
54#include <iomanip>
55
56// --------------------------------------------------------------
57
58G4PhysicsVector::G4PhysicsVector(G4bool spline)
59 : type(T_G4PhysicsVector),
60   edgeMin(0.), edgeMax(0.), numberOfNodes(0),
61   lastEnergy(-DBL_MAX), lastValue(0.), lastBin(0), useSpline(spline)
62{}
63
64// --------------------------------------------------------------
65
66G4PhysicsVector::~G4PhysicsVector() 
67{}
68
69// --------------------------------------------------------------
70
71G4PhysicsVector::G4PhysicsVector(const G4PhysicsVector& right)
72{
73  CopyData(right);
74}
75
76// --------------------------------------------------------------
77
78G4PhysicsVector& G4PhysicsVector::operator=(const G4PhysicsVector& right)
79{
80  if (&right==this)  { return *this; }
81  if (type != right.type)  { return *this; }
82
83  //DeleteData();
84  CopyData(right);
85
86  return *this;
87}
88
89// --------------------------------------------------------------
90
91G4int G4PhysicsVector::operator==(const G4PhysicsVector &right) const
92{
93  return (this == &right);
94}
95
96// --------------------------------------------------------------
97
98G4int G4PhysicsVector::operator!=(const G4PhysicsVector &right) const
99{
100  return (this != &right);
101}
102
103// --------------------------------------------------------------
104
105void G4PhysicsVector::DeleteData()
106{
107  secDerivative.clear();
108}
109
110// --------------------------------------------------------------
111
112void G4PhysicsVector::CopyData(const G4PhysicsVector& vec)
113{
114  type = vec.type;
115  edgeMin = vec.edgeMin;
116  edgeMax = vec.edgeMax;
117  numberOfNodes = vec.numberOfNodes;
118  lastEnergy = vec.lastEnergy;
119  lastValue = vec.lastValue;
120  lastBin = vec.lastBin;
121  dataVector = vec.dataVector;
122  binVector = vec.binVector;
123  useSpline = vec.useSpline;
124  comment = vec.comment;
125  secDerivative = vec.secDerivative;
126}
127
128// --------------------------------------------------------------
129
130G4double G4PhysicsVector::GetLowEdgeEnergy(size_t binNumber) const
131{
132  return binVector[binNumber];
133}
134
135// --------------------------------------------------------------
136
137G4bool G4PhysicsVector::Store(std::ofstream& fOut, G4bool ascii)
138{
139  // Ascii mode
140  if (ascii)
141  {
142    fOut << *this;
143    return true;
144  } 
145  // Binary Mode
146
147  // binning
148  fOut.write((char*)(&edgeMin), sizeof edgeMin);
149  fOut.write((char*)(&edgeMax), sizeof edgeMax);
150  fOut.write((char*)(&numberOfNodes), sizeof numberOfNodes);
151
152  // contents
153  size_t size = dataVector.size(); 
154  fOut.write((char*)(&size), sizeof size);
155
156  G4double* value = new G4double[2*size];
157  for(size_t i = 0; i < size; ++i)
158  {
159    value[2*i]  =  binVector[i];
160    value[2*i+1]=  dataVector[i];
161  }
162  fOut.write((char*)(value), 2*size*(sizeof (G4double)));
163  delete [] value;
164
165  return true;
166}
167
168// --------------------------------------------------------------
169
170G4bool G4PhysicsVector::Retrieve(std::ifstream& fIn, G4bool ascii)
171{
172  // clear properties;
173  lastEnergy=-DBL_MAX;
174  lastValue =0.;
175  lastBin   =0;
176  dataVector.clear();
177  binVector.clear();
178  secDerivative.clear();
179  comment = "";
180
181  // retrieve in ascii mode
182  if (ascii)
183  {
184    // binning
185    fIn >> edgeMin >> edgeMax >> numberOfNodes; 
186    if (fIn.fail())  { return false; }
187    // contents
188    size_t size=0;
189    fIn >> size;
190    if (fIn.fail())  { return false; }
191
192    binVector.reserve(size);
193    dataVector.reserve(size);
194    G4double vBin, vData;
195
196    for(size_t i = 0; i < size ; i++)
197    {
198      vBin = 0.;
199      vData= 0.;
200      fIn >> vBin >> vData;
201      if (fIn.fail())  { return false; }
202      binVector.push_back(vBin);
203      dataVector.push_back(vData);
204    }
205    return true ;
206  }
207
208  // retrieve in binary mode
209  // binning
210  fIn.read((char*)(&edgeMin), sizeof edgeMin);
211  fIn.read((char*)(&edgeMax), sizeof edgeMax);
212  fIn.read((char*)(&numberOfNodes), sizeof numberOfNodes ); 
213 
214  // contents
215  size_t size;
216  fIn.read((char*)(&size), sizeof size); 
217 
218  G4double* value = new G4double[2*size];
219  fIn.read((char*)(value),  2*size*(sizeof(G4double)) );
220  if (G4int(fIn.gcount()) != G4int(2*size*(sizeof(G4double))) )
221  {
222    delete [] value;
223    return false;
224  }
225
226  binVector.reserve(size);
227  dataVector.reserve(size);
228  for(size_t i = 0; i < size; ++i)
229  {
230    binVector.push_back(value[2*i]);
231    dataVector.push_back(value[2*i+1]);
232  }
233  delete [] value;
234  return true;
235}
236
237// --------------------------------------------------------------
238
239void 
240G4PhysicsVector::ScaleVector(G4double factorE, G4double factorV)
241{
242  size_t n = dataVector.size();
243  size_t i;
244  if(n > 0) { 
245    for(i=0; i<n; ++i) {
246      binVector[i]  *= factorE;
247      dataVector[i] *= factorV;
248    } 
249  }
250  n = secDerivative.size();
251  if(n > 0) { for(i=0; i<n; ++i) { secDerivative[i] *= factorV; } }
252
253  edgeMin *= factorE;
254  edgeMax *= factorE;
255  lastEnergy *= factorE;
256  lastValue  *= factorV;
257}
258
259// --------------------------------------------------------------
260
261void 
262G4PhysicsVector::ComputeSecondDerivatives(G4double firstPointDerivative, 
263                                          G4double endPointDerivative)
264  //  A standard method of computation of second derivatives
265  //  First derivatives at the first and the last point should be provided
266  //  See for example W.H. Press et al. "Numerical reciptes and C"
267  //  Cambridge University Press, 1997.
268{
269  if(4 > numberOfNodes)   // cannot compute derivatives for less than 4 bins
270  {
271    ComputeSecDerivatives();
272    return;
273  }
274
275  if(!SplinePossible()) { return; }
276
277  G4int n = numberOfNodes-1;
278
279  G4double* u = new G4double [n];
280 
281  G4double p, sig, un;
282
283  u[0] = (6.0/(binVector[1]-binVector[0]))
284    * ((dataVector[1]-dataVector[0])/(binVector[1]-binVector[0])
285       - firstPointDerivative);
286 
287  secDerivative[0] = - 0.5;
288
289  // Decomposition loop for tridiagonal algorithm. secDerivative[i]
290  // and u[i] are used for temporary storage of the decomposed factors.
291
292  for(G4int i=1; i<n; ++i)
293  {
294    sig = (binVector[i]-binVector[i-1]) / (binVector[i+1]-binVector[i-1]);
295    p = sig*secDerivative[i-1] + 2.0;
296    secDerivative[i] = (sig - 1.0)/p;
297    u[i] = (dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i])
298         - (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]);
299    u[i] = 6.0*u[i]/(binVector[i+1]-binVector[i-1]) - sig*u[i-1]/p;
300  }
301
302  sig = (binVector[n-1]-binVector[n-2]) / (binVector[n]-binVector[n-2]);
303  p = sig*secDerivative[n-2] + 2.0;
304  un = (6.0/(binVector[n]-binVector[n-1]))
305    *(endPointDerivative - 
306      (dataVector[n]-dataVector[n-1])/(binVector[n]-binVector[n-1])) - u[n-1]/p;
307  secDerivative[n] = un/(secDerivative[n-1] + 2.0);
308
309  // The back-substitution loop for the triagonal algorithm of solving
310  // a linear system of equations.
311   
312  for(G4int k=n-1; k>0; --k)
313  {
314    secDerivative[k] *= 
315      (secDerivative[k+1] - 
316       u[k]*(binVector[k+1]-binVector[k-1])/(binVector[k+1]-binVector[k]));
317  }
318  secDerivative[0] = 0.5*(u[0] - secDerivative[1]);
319
320  delete [] u;
321}
322
323// --------------------------------------------------------------
324
325void G4PhysicsVector::FillSecondDerivatives()
326  // Computation of second derivatives using "Not-a-knot" endpoint conditions
327  // B.I. Kvasov "Methods of shape-preserving spline approximation"
328  // World Scientific, 2000
329{ 
330  if(5 > numberOfNodes)  // cannot compute derivatives for less than 4 points
331  {
332    ComputeSecDerivatives();
333    return;
334  }
335
336  if(!SplinePossible()) { return; }
337 
338  G4int n = numberOfNodes-1;
339
340  //G4cout << "G4PhysicsVector::FillSecondDerivatives() n= " << n << G4endl;
341  // G4cout << *this << G4endl;
342
343  G4double* u = new G4double [n];
344 
345  G4double p, sig;
346
347  u[1] = ((dataVector[2]-dataVector[1])/(binVector[2]-binVector[1]) -
348          (dataVector[1]-dataVector[0])/(binVector[1]-binVector[0]));
349  u[1] = 6.0*u[1]*(binVector[2]-binVector[1])
350    / ((binVector[2]-binVector[0])*(binVector[2]-binVector[0]));
351 
352  // Decomposition loop for tridiagonal algorithm. secDerivative[i]
353  // and u[i] are used for temporary storage of the decomposed factors.
354
355  secDerivative[1] = (2.0*binVector[1]-binVector[0]-binVector[2])
356    / (2.0*binVector[2]-binVector[0]-binVector[1]);
357
358  for(G4int i=2; i<n-1; ++i)
359  {
360    sig = (binVector[i]-binVector[i-1]) / (binVector[i+1]-binVector[i-1]);
361    p = sig*secDerivative[i-1] + 2.0;
362    secDerivative[i] = (sig - 1.0)/p;
363    u[i] = (dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i])
364         - (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]);
365    u[i] = (6.0*u[i]/(binVector[i+1]-binVector[i-1])) - sig*u[i-1]/p;
366  }
367
368  sig = (binVector[n-1]-binVector[n-2]) / (binVector[n]-binVector[n-2]);
369  p = sig*secDerivative[n-3] + 2.0;
370  u[n-1] = (dataVector[n]-dataVector[n-1])/(binVector[n]-binVector[n-1])
371    - (dataVector[n-1]-dataVector[n-2])/(binVector[n-1]-binVector[n-2]);
372  u[n-1] = 6.0*sig*u[n-1]/(binVector[n]-binVector[n-2])
373    - (2.0*sig - 1.0)*u[n-2]/p; 
374
375  p = (1.0+sig) + (2.0*sig-1.0)*secDerivative[n-2];
376  secDerivative[n-1] = u[n-1]/p;
377
378  // The back-substitution loop for the triagonal algorithm of solving
379  // a linear system of equations.
380   
381  for(G4int k=n-2; k>1; --k)
382  {
383    secDerivative[k] *= 
384      (secDerivative[k+1] - 
385       u[k]*(binVector[k+1]-binVector[k-1])/(binVector[k+1]-binVector[k]));
386  }
387  secDerivative[n] = (secDerivative[n-1] - (1.0-sig)*secDerivative[n-2])/sig;
388  sig = 1.0 - ((binVector[2]-binVector[1])/(binVector[2]-binVector[0]));
389  secDerivative[1] *= (secDerivative[2] - u[1]/(1.0-sig));
390  secDerivative[0]  = (secDerivative[1] - sig*secDerivative[2])/(1.0-sig);
391
392  delete [] u;
393}
394
395// --------------------------------------------------------------
396
397void 
398G4PhysicsVector::ComputeSecDerivatives()
399  //  A simplified method of computation of second derivatives
400{
401  if(!SplinePossible())  { return; }
402
403  if(3 > numberOfNodes)  // cannot compute derivatives for less than 4 bins
404  {
405    useSpline = false;
406    return;
407  }
408
409  size_t n = numberOfNodes-1;
410
411  for(size_t i=1; i<n; ++i)
412  {
413    secDerivative[i] =
414      3.0*((dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i]) -
415           (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]))
416      /(binVector[i+1]-binVector[i-1]);
417  }
418  secDerivative[n] = secDerivative[n-1];
419  secDerivative[0] = secDerivative[1];
420}
421
422// --------------------------------------------------------------
423
424G4bool G4PhysicsVector::SplinePossible()
425  // Initialise second derivative array. If neighbor energy coincide
426  // or not ordered than spline cannot be applied
427{
428  if(!useSpline) return useSpline;
429  secDerivative.clear();
430  secDerivative.reserve(numberOfNodes);
431  for(size_t j=0; j<numberOfNodes; ++j)
432  {
433    secDerivative.push_back(0.0);
434    if(j > 0)
435    {
436      if(binVector[j]-binVector[j-1] <= 0.)  { useSpline = false; }
437    }
438  } 
439  return useSpline;
440}
441   
442// --------------------------------------------------------------
443
444std::ostream& operator<<(std::ostream& out, const G4PhysicsVector& pv)
445{
446  // binning
447  out << std::setprecision(12) << pv.edgeMin;
448  out <<" " << pv.edgeMax <<" "  << pv.numberOfNodes << G4endl; 
449
450  // contents
451  out << pv.dataVector.size() << G4endl; 
452  for(size_t i = 0; i < pv.dataVector.size(); i++)
453  {
454    out << std::setprecision(12) << pv.binVector[i] << "  "
455        << pv.dataVector[i] << G4endl;
456  }
457  return out;
458}
Note: See TracBrowser for help on using the repository browser.