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