Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysicsFreeVector.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// G4PhysicsFreeVector class implementation
27//
28// Authors:
29// - 02 Dec. 1995, G.Cosmo: Structure created based on object model
30// - 06 Jun. 1996, K.Amako: Implemented the 1st version
31// Revisions:
32// - 11 Nov. 2000, H.Kurashige: Use STL vector for dataVector and binVector
33// - 25 Aug. 2021, V.Ivanchenko updated for Geant4 11.0
34// --------------------------------------------------------------------
35
37#include "G4Exp.hh"
38
39// --------------------------------------------------------------------
43
44// --------------------------------------------------------------------
46 : G4PhysicsFreeVector(static_cast<std::size_t>(length), false)
47{}
48
49// --------------------------------------------------------------------
51 : G4PhysicsVector(spline)
52{
53 numberOfNodes = length;
54
55 if (0 < length) {
56 binVector.resize(numberOfNodes, 0.0);
57 dataVector.resize(numberOfNodes, 0.0);
58 }
59 Initialise();
60}
61
62// --------------------------------------------------------------------
64 G4double, G4bool spline)
65 : G4PhysicsFreeVector(length, spline)
66{}
67
68// --------------------------------------------------------------------
69G4PhysicsFreeVector::G4PhysicsFreeVector(const std::vector<G4double>& energies,
70 const std::vector<G4double>& values,
71 G4bool spline)
72 : G4PhysicsVector(spline)
73{
74 numberOfNodes = energies.size();
75
76 if (numberOfNodes != values.size())
77 {
78 G4ExceptionDescription ed;
79 ed << "The size of energy vector " << numberOfNodes << " != " << values.size();
80 G4Exception("G4PhysicsFreeVector constructor: ","glob04", FatalException, ed);
81 }
82
83 binVector = energies;
84 dataVector = values;
85 Initialise();
86}
87
88// --------------------------------------------------------------------
90 const G4double* values,
91 std::size_t length,
92 G4bool spline)
93 : G4PhysicsVector(spline)
94{
95 numberOfNodes = length;
96
97 if (0 < numberOfNodes)
98 {
101
102 for(std::size_t i = 0; i < numberOfNodes; ++i)
103 {
104 binVector[i] = energies[i];
105 dataVector[i] = values[i];
106 }
107 }
108 Initialise();
109}
110
111// --------------------------------------------------------------------
112void G4PhysicsFreeVector::PutValues(const std::size_t index,
113 const G4double e,
114 const G4double value)
115{
116 if (index >= numberOfNodes)
117 {
118 PrintPutValueError(index, value, "G4PhysicsFreeVector::PutValues ");
119 return;
120 }
121 binVector[index] = e;
122 dataVector[index] = value;
123 if(index == 0)
124 {
125 edgeMin = e;
126 }
127 else if(numberOfNodes == index + 1)
128 {
129 edgeMax = e;
130 }
131}
132
133// --------------------------------------------------------------------
135 const G4double value)
136{
137 auto binLoc = std::lower_bound(binVector.cbegin(), binVector.cend(), energy);
138 auto dataLoc = dataVector.cbegin();
139 dataLoc += binLoc - binVector.cbegin();
140
141 binVector.insert(binLoc, energy);
142 dataVector.insert(dataLoc, value);
143
145 Initialise();
146}
147
148// --------------------------------------------------------------------
150{
151 // check if log search is applicable
152 if (0 >= n || edgeMin <= 0.0 || edgeMin == edgeMax || numberOfNodes < 3)
153 {
154 return;
155 }
156 nLogNodes = static_cast<std::size_t>(static_cast<G4int>(numberOfNodes)/n);
157 if (nLogNodes < 3) { nLogNodes = 3; }
158 scale.resize(nLogNodes, 0);
159 imax1 = nLogNodes - 2;
160 iBin1 = (imax1 + 1) / G4Log(edgeMax/edgeMin);
162 scale[0] = 0;
163 scale[imax1 + 1] = idxmax;
164 std::size_t j = 0;
165 for (std::size_t i = 1; i <= imax1; ++i)
166 {
168 for (; j <= idxmax; ++j)
169 {
170 if (binVector[j] <= e && e < binVector[j+1])
171 {
172 scale[i] = j;
173 break;
174 }
175 }
176 }
177}
178
179// --------------------------------------------------------------------
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void InsertValues(const G4double energy, const G4double value)
void PutValues(const std::size_t index, const G4double energy, const G4double value)
void EnableLogBinSearch(const G4int n=1)
G4PhysicsFreeVector(G4bool spline=false)
std::size_t nLogNodes
void PrintPutValueError(std::size_t index, G4double value, const G4String &text)
std::size_t numberOfNodes
std::vector< G4double > dataVector
std::vector< G4double > binVector
virtual void Initialise()
std::vector< std::size_t > scale