Geant4-11
G4VelocityTable.hh
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// G4VelocityTable
27//
28// Class description:
29//
30// This class keeps a table of velocity as a function of the ratio
31// kinetic erngy and mass. G4VelocityTable is used by
32// G4Track::CalculateVelocity().
33
34// Author: Hisaya Kurashige, 17 August 2011
35// --------------------------------------------------------------------
36#ifndef G4VelocityTable_hh
37#define G4VelocityTable_hh 1
38
39#include <vector>
40#include <iostream>
41
42#include "globals.hh"
43#include "G4ios.hh"
45
47{
49
50 using G4VTDataVector = std::vector<G4double>;
51
52 public:
53
54 G4double Value(G4double theEnergy);
55 // Get the cross-section/energy-loss value corresponding to the
56 // given energy. An appropriate interpolation is used to calculate
57 // the value
58
60
61 static void SetVelocityTableProperties(G4double t_max, G4double t_min,
62 G4int nbin);
66
67 private:
68
71
73
74 std::size_t FindBinLocation(G4double theEnergy) const;
75 // Find the bin# in which theEnergy belongs - pure virtual function
76
77 inline G4double Interpolation() const;
78
79 // --------------------------------------------------------------------
80
81 G4double edgeMin = 0.0; // Energy of first point
82 G4double edgeMax = 0.0; // Energy of the last point
83
84 std::size_t numberOfNodes = 0;
85
86 G4VTDataVector dataVector; // Vector to keep the crossection/energyloss
87 G4VTDataVector binVector; // Vector to keep energy
88 G4VTDataVector secDerivative; // Vector to keep second derivatives
89
90 G4double dBin = 0.0; // Bin width - useful only for fixed binning
91 G4double baseBin = 0.0; // Set this in constructor for performance
92
93 G4double lastEnergy = -DBL_MAX; // Cache the last input value
94 G4double lastValue = 0.0; // Cache the last output value
95 std::size_t lastBin = 0; // Cache the last bin location
96
98
99 G4double maxT = 1000.0;
100 G4double minT = 0.0001;
102};
103
104// ----------------------
105// Inline methods
106// ----------------------
107
109{
110 // Linear interpolation is used to get the value. If the give energy
111 // is in the highest bin, no interpolation will be Done. Because
112 // there is an extra bin hidden from a user at locBin=numberOfBin,
113 // the following interpolation is valid even the current locBin=
114 // numberOfBin-1.
115
116 G4double intplFactor =
118 (binVector[lastBin + 1] - binVector[lastBin]); // Interpol. factor
119
120 return dataVector[lastBin] +
121 (dataVector[lastBin + 1] - dataVector[lastBin]) * intplFactor;
122}
123
124#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
static G4VelocityTable * GetVelocityTable()
static void SetVelocityTableProperties(G4double t_max, G4double t_min, G4int nbin)
std::size_t numberOfNodes
std::vector< G4double > G4VTDataVector
static G4ThreadLocal G4VelocityTable * theInstance
std::size_t lastBin
G4VTDataVector dataVector
G4double Value(G4double theEnergy)
static G4double GetMaxTOfVelocityTable()
G4VTDataVector secDerivative
G4VTDataVector binVector
static G4double GetMinTOfVelocityTable()
G4double Interpolation() const
void PrepareVelocityTable()
std::size_t FindBinLocation(G4double theEnergy) const
static G4int GetNbinOfVelocityTable()
#define DBL_MAX
Definition: templates.hh:62
#define G4ThreadLocal
Definition: tls.hh:77