Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/global/management/src/G4PhysicsVector.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // --------------------------------------------------------------------
 34 
 35 #include "G4PhysicsVector.hh"
 36 #include <iomanip>
 37 
 38 // --------------------------------------------------------------
 39 G4PhysicsVector::G4PhysicsVector(G4bool val)
 40   : useSpline(val)
 41 {}
 42 
 43 // --------------------------------------------------------------------
 44 void G4PhysicsVector::Initialise()
 45 {
 46   if (1 < numberOfNodes)
 47   {
 48     idxmax = numberOfNodes - 2;
 49     edgeMin = binVector[0];
 50     edgeMax = binVector[idxmax + 1];
 51   }
 52 }
 53 
 54 // --------------------------------------------------------------
 55 G4bool G4PhysicsVector::Store(std::ofstream& fOut, G4bool ascii) const
 56 {
 57   // Ascii mode
 58   if (ascii)
 59   {
 60     fOut << *this;
 61     return true;
 62   }
 63   // Binary Mode
 64 
 65   // binning
 66   fOut.write((char*) (&edgeMin), sizeof edgeMin);
 67   fOut.write((char*) (&edgeMax), sizeof edgeMax);
 68   fOut.write((char*) (&numberOfNodes), sizeof numberOfNodes);
 69 
 70   // contents
 71   std::size_t size = dataVector.size();
 72   fOut.write((char*) (&size), sizeof size);
 73 
 74   auto value = new G4double[2 * size];
 75   for (std::size_t i = 0; i < size; ++i)
 76   {
 77     value[2 * i]     = binVector[i];
 78     value[2 * i + 1] = dataVector[i];
 79   }
 80   fOut.write((char*) (value), 2 * size * (sizeof(G4double)));
 81   delete[] value;
 82 
 83   return true;
 84 }
 85 
 86 // --------------------------------------------------------------
 87 G4bool G4PhysicsVector::Retrieve(std::ifstream& fIn, G4bool ascii)
 88 {
 89   // clear properties;
 90   dataVector.clear();
 91   binVector.clear();
 92   secDerivative.clear();
 93 
 94   // retrieve in ascii mode
 95   if (ascii)
 96   {
 97     // binning
 98     fIn >> edgeMin >> edgeMax >> numberOfNodes;
 99     if (fIn.fail() || numberOfNodes < 2)
100     {
101       return false;
102     }
103     // contents
104     G4int siz0 = 0;
105     fIn >> siz0;
106     if (siz0 < 2) { return false; }
107     auto siz = static_cast<std::size_t>(siz0);
108     if (fIn.fail() || siz != numberOfNodes)
109     {
110       return false;
111     }
112 
113     binVector.reserve(siz);
114     dataVector.reserve(siz);
115     G4double vBin, vData;
116 
117     for (std::size_t i = 0; i < siz; ++i)
118     {
119       vBin  = 0.;
120       vData = 0.;
121       fIn >> vBin >> vData;
122       if (fIn.fail())
123       {
124         return false;
125       }
126       binVector.push_back(vBin);
127       dataVector.push_back(vData);
128     }
129     Initialise();
130     return true;
131   }
132 
133   // retrieve in binary mode
134   // binning
135   fIn.read((char*) (&edgeMin), sizeof edgeMin);
136   fIn.read((char*) (&edgeMax), sizeof edgeMax);
137   fIn.read((char*) (&numberOfNodes), sizeof numberOfNodes);
138 
139   // contents
140   std::size_t size;
141   fIn.read((char*) (&size), sizeof size);
142 
143   auto value = new G4double[2 * size];
144   fIn.read((char*) (value), 2 * size * (sizeof(G4double)));
145   if (static_cast<G4int>(fIn.gcount()) != static_cast<G4int>(2 * size * (sizeof(G4double))))
146   {
147     delete[] value;
148     return false;
149   }
150 
151   binVector.reserve(size);
152   dataVector.reserve(size);
153   for (std::size_t i = 0; i < size; ++i)
154   {
155     binVector.push_back(value[2 * i]);
156     dataVector.push_back(value[2 * i + 1]);
157   }
158   delete[] value;
159 
160   Initialise();
161   return true;
162 }
163 
164 // --------------------------------------------------------------
165 void G4PhysicsVector::DumpValues(G4double unitE, G4double unitV) const
166 {
167   for (std::size_t i = 0; i < numberOfNodes; ++i)
168   {
169     G4cout << binVector[i] / unitE << "   " << dataVector[i] / unitV 
170            << G4endl;
171   }
172 }
173 
174 // --------------------------------------------------------------------
175 std::size_t G4PhysicsVector::FindBin(const G4double energy, 
176                                      std::size_t idx) const
177 {
178   if (idx + 1 < numberOfNodes && 
179       energy >= binVector[idx] && energy <= binVector[idx])
180   {
181     return idx;
182   } 
183   if (energy <= binVector[1])
184   {
185     return 0;
186   }
187   if (energy >= binVector[idxmax])
188   {
189     return idxmax;
190   }
191   return GetBin(energy); 
192 }
193 
194 // --------------------------------------------------------------------
195 void G4PhysicsVector::ScaleVector(const G4double factorE, 
196                                   const G4double factorV)
197 {
198   for (std::size_t i = 0; i < numberOfNodes; ++i)
199   {
200     binVector[i] *= factorE;
201     dataVector[i] *= factorV;
202   }
203   Initialise();
204 }
205 
206 // --------------------------------------------------------------------
207 void G4PhysicsVector::FillSecondDerivatives(const G4SplineType stype,
208               const G4double dir1,
209               const G4double dir2)
210 {
211   if (!useSpline) { return; }
212   // cannot compute derivatives for less than 5 points
213   const std::size_t nmin = (stype == G4SplineType::Base) ? 5 : 4;
214   if (nmin > numberOfNodes) 
215   {
216     if (0 < verboseLevel)
217     { 
218       G4cout << "### G4PhysicsVector: spline cannot be used for "
219        << numberOfNodes << " points - spline disabled" 
220        << G4endl;
221       DumpValues();
222     }
223     useSpline = false;
224     return;
225   }
226   // check energies of free vector
227   if (type == T_G4PhysicsFreeVector)
228   {
229     for (std::size_t i=0; i<=idxmax; ++i) 
230     {
231       if (binVector[i + 1] <= binVector[i])
232       {
233         if (0 < verboseLevel) 
234         {
235     G4cout << "### G4PhysicsVector: spline cannot be used, because "
236      << " E[" << i << "]=" << binVector[i]
237      << " >= E[" << i+1 << "]=" << binVector[i + 1]
238      << G4endl;
239     DumpValues();
240         }
241         useSpline = false;
242         return;
243       }
244     }
245   }
246 
247   // spline is possible
248   Initialise();
249   secDerivative.resize(numberOfNodes);
250 
251   if (1 < verboseLevel)
252   {
253     G4cout << "### G4PhysicsVector:: FillSecondDerivatives N=" 
254            << numberOfNodes << G4endl;
255     DumpValues();
256   }
257 
258   switch(stype) 
259   {
260     case G4SplineType::Base:
261       ComputeSecDerivative1();
262       break;
263 
264     case G4SplineType::FixedEdges:
265       ComputeSecDerivative2(dir1, dir2);
266       break;
267 
268     default:
269       ComputeSecDerivative0();
270   }
271 }
272 
273 // --------------------------------------------------------------
274 void G4PhysicsVector::ComputeSecDerivative0()
275 //  A simplified method of computation of second derivatives
276 {
277   std::size_t n = numberOfNodes - 1;
278 
279   for (std::size_t i = 1; i < n; ++i)
280   {
281     secDerivative[i] = 3.0 *
282       ((dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) -
283        (dataVector[i] - dataVector[i - 1]) /
284          (binVector[i] - binVector[i - 1])) /
285       (binVector[i + 1] - binVector[i - 1]);
286   }
287   secDerivative[n] = secDerivative[n - 1];
288   secDerivative[0] = secDerivative[1];
289 }
290 
291 // --------------------------------------------------------------
292 void G4PhysicsVector::ComputeSecDerivative1()
293 // Computation of second derivatives using "Not-a-knot" endpoint conditions
294 // B.I. Kvasov "Methods of shape-preserving spline approximation"
295 // World Scientific, 2000
296 {
297   std::size_t n = numberOfNodes - 1;
298   auto u = new G4double[n];
299   G4double p, sig;
300 
301   u[1] = ((dataVector[2] - dataVector[1]) / (binVector[2] - binVector[1]) -
302           (dataVector[1] - dataVector[0]) / (binVector[1] - binVector[0]));
303   u[1] = 6.0 * u[1] * (binVector[2] - binVector[1]) /
304          ((binVector[2] - binVector[0]) * (binVector[2] - binVector[0]));
305 
306   // Decomposition loop for tridiagonal algorithm. secDerivative[i]
307   // and u[i] are used for temporary storage of the decomposed factors.
308 
309   secDerivative[1] = (2.0 * binVector[1] - binVector[0] - binVector[2]) /
310                      (2.0 * binVector[2] - binVector[0] - binVector[1]);
311 
312   for(std::size_t i = 2; i < n - 1; ++i)
313   {
314     sig =
315       (binVector[i] - binVector[i - 1]) / (binVector[i + 1] - binVector[i - 1]);
316     p                = sig * secDerivative[i - 1] + 2.0;
317     secDerivative[i] = (sig - 1.0) / p;
318     u[i] =
319       (dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) -
320       (dataVector[i] - dataVector[i - 1]) / (binVector[i] - binVector[i - 1]);
321     u[i] =
322       (6.0 * u[i] / (binVector[i + 1] - binVector[i - 1])) - sig * u[i - 1] / p;
323   }
324 
325   sig =
326     (binVector[n - 1] - binVector[n - 2]) / (binVector[n] - binVector[n - 2]);
327   p = sig * secDerivative[n - 3] + 2.0;
328   u[n - 1] =
329     (dataVector[n] - dataVector[n - 1]) / (binVector[n] - binVector[n - 1]) -
330     (dataVector[n - 1] - dataVector[n - 2]) /
331       (binVector[n - 1] - binVector[n - 2]);
332   u[n - 1] = 6.0 * sig * u[n - 1] / (binVector[n] - binVector[n - 2]) -
333              (2.0 * sig - 1.0) * u[n - 2] / p;
334 
335   p = (1.0 + sig) + (2.0 * sig - 1.0) * secDerivative[n - 2];
336   secDerivative[n - 1] = u[n - 1] / p;
337 
338   // The back-substitution loop for the triagonal algorithm of solving
339   // a linear system of equations.
340 
341   for (std::size_t k = n - 2; k > 1; --k)
342   {
343     secDerivative[k] *=
344       (secDerivative[k + 1] - u[k] * (binVector[k + 1] - binVector[k - 1]) /
345                                 (binVector[k + 1] - binVector[k]));
346   }
347   secDerivative[n] =
348     (secDerivative[n - 1] - (1.0 - sig) * secDerivative[n - 2]) / sig;
349   sig = 1.0 - ((binVector[2] - binVector[1]) / (binVector[2] - binVector[0]));
350   secDerivative[1] *= (secDerivative[2] - u[1] / (1.0 - sig));
351   secDerivative[0] = (secDerivative[1] - sig * secDerivative[2]) / (1.0 - sig);
352 
353   delete[] u;
354 }
355 
356 // --------------------------------------------------------------
357 void G4PhysicsVector::ComputeSecDerivative2(G4double firstPointDerivative,
358                                             G4double endPointDerivative)
359 // A standard method of computation of second derivatives
360 // First derivatives at the first and the last point should be provided
361 // See for example W.H. Press et al. "Numerical recipes in C"
362 // Cambridge University Press, 1997.
363 {
364   std::size_t n = numberOfNodes - 1;
365   auto u = new G4double[n];
366   G4double p, sig, un;
367 
368   u[0] = (6.0 / (binVector[1] - binVector[0])) *
369          ((dataVector[1] - dataVector[0]) / (binVector[1] - binVector[0]) -
370           firstPointDerivative);
371 
372   secDerivative[0] = -0.5;
373 
374   // Decomposition loop for tridiagonal algorithm. secDerivative[i]
375   // and u[i] are used for temporary storage of the decomposed factors.
376 
377   for (std::size_t i = 1; i < n; ++i)
378   {
379     sig =
380       (binVector[i] - binVector[i - 1]) / (binVector[i + 1] - binVector[i - 1]);
381     p                = sig * (secDerivative[i - 1]) + 2.0;
382     secDerivative[i] = (sig - 1.0) / p;
383     u[i] =
384       (dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) -
385       (dataVector[i] - dataVector[i - 1]) / (binVector[i] - binVector[i - 1]);
386     u[i] =
387       6.0 * u[i] / (binVector[i + 1] - binVector[i - 1]) - sig * u[i - 1] / p;
388   }
389 
390   sig =
391     (binVector[n - 1] - binVector[n - 2]) / (binVector[n] - binVector[n - 2]);
392   p  = sig * secDerivative[n - 2] + 2.0;
393   un = (6.0 / (binVector[n] - binVector[n - 1])) *
394          (endPointDerivative - (dataVector[n] - dataVector[n - 1]) /
395                                  (binVector[n] - binVector[n - 1])) -
396        u[n - 1] / p;
397   secDerivative[n] = un / (secDerivative[n - 1] + 2.0);
398 
399   // The back-substitution loop for the triagonal algorithm of solving
400   // a linear system of equations.
401 
402   for (std::size_t k = n - 1; k > 0; --k)
403   {
404     secDerivative[k] *=
405       (secDerivative[k + 1] - u[k] * (binVector[k + 1] - binVector[k - 1]) /
406                                 (binVector[k + 1] - binVector[k]));
407   }
408   secDerivative[0] = 0.5 * (u[0] - secDerivative[1]);
409 
410   delete[] u;
411 }
412 
413 // --------------------------------------------------------------
414 std::ostream& operator<<(std::ostream& out, const G4PhysicsVector& pv)
415 {
416   // binning
417   G4long prec = out.precision();
418   out << std::setprecision(12) << pv.edgeMin << " " << pv.edgeMax << " "
419       << pv.numberOfNodes << G4endl;
420 
421   // contents
422   out << pv.dataVector.size() << G4endl;
423   for (std::size_t i = 0; i < pv.dataVector.size(); ++i)
424   {
425     out << pv.binVector[i] << "  " << pv.dataVector[i] << G4endl;
426   }
427   out.precision(prec);
428 
429   return out;
430 }
431 
432 //---------------------------------------------------------------
433 G4double G4PhysicsVector::GetEnergy(const G4double val) const
434 {
435   if (0 == numberOfNodes)
436   {
437     return 0.0;
438   }
439   if (1 == numberOfNodes || val <= dataVector[0])
440   {
441     return edgeMin;
442   }
443   if (val >= dataVector[numberOfNodes - 1])
444   {
445     return edgeMax;
446   }
447   std::size_t bin = std::lower_bound(dataVector.cbegin(), dataVector.cend(), val)
448                   - dataVector.cbegin() - 1;
449   if (bin > idxmax) { bin = idxmax; } 
450   G4double res = binVector[bin];
451   G4double del = dataVector[bin + 1] - dataVector[bin];
452   if (del > 0.0)
453   {
454     res += (val - dataVector[bin]) * (binVector[bin + 1] - res) / del;
455   }
456   return res;
457 }
458 
459 //---------------------------------------------------------------
460 void G4PhysicsVector::PrintPutValueError(std::size_t index, 
461                                          G4double val, 
462                                          const G4String& text)
463 {
464   G4ExceptionDescription ed;
465   ed << "Vector type: " << type << " length= " << numberOfNodes
466      << "; an attempt to put data at index= " << index
467      << " value= " << val << " in " << text;
468   G4Exception("G4PhysicsVector:", "gl0005", 
469               FatalException, ed, "Wrong operation");
470 }
471 
472 //---------------------------------------------------------------
473