Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // 26 // G4PhysicsVector class implementation 26 // G4PhysicsVector class implementation 27 // 27 // 28 // Authors: 28 // Authors: 29 // - 02 Dec. 1995, G.Cosmo: Structure created 29 // - 02 Dec. 1995, G.Cosmo: Structure created based on object model 30 // - 03 Mar. 1996, K.Amako: Implemented the 1s 30 // - 03 Mar. 1996, K.Amako: Implemented the 1st version 31 // Revisions: 31 // Revisions: 32 // - 11 Nov. 2000, H.Kurashige: Use STL vector 32 // - 11 Nov. 2000, H.Kurashige: Use STL vector for dataVector and binVector >> 33 // - 02 Apr. 2008, A.Bagulya: Added SplineInterpolation() and SetSpline() >> 34 // - 19 Jun. 2009, V.Ivanchenko: Removed hidden bin >> 35 // - 15 Mar. 2019 M.Novak: added Value method with the known log-energy value >> 36 // that can avoid the log call in case of log-vectors 33 // ------------------------------------------- 37 // -------------------------------------------------------------------- 34 38 35 #include "G4PhysicsVector.hh" 39 #include "G4PhysicsVector.hh" 36 #include <iomanip> 40 #include <iomanip> 37 41 38 // ------------------------------------------- 42 // -------------------------------------------------------------- >> 43 39 G4PhysicsVector::G4PhysicsVector(G4bool val) 44 G4PhysicsVector::G4PhysicsVector(G4bool val) 40 : useSpline(val) 45 : useSpline(val) 41 {} 46 {} 42 47 43 // ------------------------------------------- << 48 // -------------------------------------------------------------- 44 void G4PhysicsVector::Initialise() << 49 >> 50 G4PhysicsVector::~G4PhysicsVector() {} >> 51 >> 52 // -------------------------------------------------------------- >> 53 >> 54 G4PhysicsVector::G4PhysicsVector(const G4PhysicsVector& right) >> 55 { >> 56 invdBin = right.invdBin; >> 57 baseBin = right.baseBin; >> 58 verboseLevel = right.verboseLevel; >> 59 >> 60 DeleteData(); >> 61 CopyData(right); >> 62 } >> 63 >> 64 // -------------------------------------------------------------- >> 65 >> 66 G4PhysicsVector& G4PhysicsVector::operator=(const G4PhysicsVector& right) >> 67 { >> 68 if(&right == this) >> 69 { >> 70 return *this; >> 71 } >> 72 invdBin = right.invdBin; >> 73 baseBin = right.baseBin; >> 74 verboseLevel = right.verboseLevel; >> 75 >> 76 DeleteData(); >> 77 CopyData(right); >> 78 return *this; >> 79 } >> 80 >> 81 // -------------------------------------------------------------- >> 82 >> 83 G4bool G4PhysicsVector::operator==(const G4PhysicsVector& right) const >> 84 { >> 85 return (this == &right); >> 86 } >> 87 >> 88 // -------------------------------------------------------------- >> 89 >> 90 G4bool G4PhysicsVector::operator!=(const G4PhysicsVector& right) const >> 91 { >> 92 return (this != &right); >> 93 } >> 94 >> 95 // -------------------------------------------------------------- >> 96 >> 97 void G4PhysicsVector::DeleteData() >> 98 { >> 99 useSpline = false; >> 100 secDerivative.clear(); >> 101 } >> 102 >> 103 // -------------------------------------------------------------- >> 104 >> 105 void G4PhysicsVector::CopyData(const G4PhysicsVector& vec) 45 { 106 { 46 if (1 < numberOfNodes) << 107 type = vec.type; >> 108 edgeMin = vec.edgeMin; >> 109 edgeMax = vec.edgeMax; >> 110 numberOfNodes = vec.numberOfNodes; >> 111 useSpline = vec.useSpline; >> 112 >> 113 std::size_t i; >> 114 dataVector.resize(numberOfNodes); >> 115 for(i = 0; i < numberOfNodes; ++i) 47 { 116 { 48 idxmax = numberOfNodes - 2; << 117 dataVector[i] = (vec.dataVector)[i]; 49 edgeMin = binVector[0]; << 118 } 50 edgeMax = binVector[idxmax + 1]; << 119 binVector.resize(numberOfNodes); >> 120 for(i = 0; i < numberOfNodes; ++i) >> 121 { >> 122 binVector[i] = (vec.binVector)[i]; >> 123 } >> 124 if(0 < (vec.secDerivative).size()) >> 125 { >> 126 secDerivative.resize(numberOfNodes); >> 127 for(i = 0; i < numberOfNodes; ++i) >> 128 { >> 129 secDerivative[i] = (vec.secDerivative)[i]; >> 130 } 51 } 131 } 52 } 132 } 53 133 54 // ------------------------------------------- 134 // -------------------------------------------------------------- >> 135 >> 136 G4double G4PhysicsVector::GetLowEdgeEnergy(std::size_t binNumber) const >> 137 { >> 138 return binVector[binNumber]; >> 139 } >> 140 >> 141 // -------------------------------------------------------------- >> 142 55 G4bool G4PhysicsVector::Store(std::ofstream& f 143 G4bool G4PhysicsVector::Store(std::ofstream& fOut, G4bool ascii) const 56 { 144 { 57 // Ascii mode 145 // Ascii mode 58 if (ascii) << 146 if(ascii) 59 { 147 { 60 fOut << *this; 148 fOut << *this; 61 return true; 149 return true; 62 } 150 } 63 // Binary Mode 151 // Binary Mode 64 152 65 // binning 153 // binning 66 fOut.write((char*) (&edgeMin), sizeof edgeMi 154 fOut.write((char*) (&edgeMin), sizeof edgeMin); 67 fOut.write((char*) (&edgeMax), sizeof edgeMa 155 fOut.write((char*) (&edgeMax), sizeof edgeMax); 68 fOut.write((char*) (&numberOfNodes), sizeof 156 fOut.write((char*) (&numberOfNodes), sizeof numberOfNodes); 69 157 70 // contents 158 // contents 71 std::size_t size = dataVector.size(); 159 std::size_t size = dataVector.size(); 72 fOut.write((char*) (&size), sizeof size); 160 fOut.write((char*) (&size), sizeof size); 73 161 74 auto value = new G4double[2 * size]; << 162 G4double* value = new G4double[2 * size]; 75 for (std::size_t i = 0; i < size; ++i) << 163 for(std::size_t i = 0; i < size; ++i) 76 { 164 { 77 value[2 * i] = binVector[i]; 165 value[2 * i] = binVector[i]; 78 value[2 * i + 1] = dataVector[i]; 166 value[2 * i + 1] = dataVector[i]; 79 } 167 } 80 fOut.write((char*) (value), 2 * size * (size 168 fOut.write((char*) (value), 2 * size * (sizeof(G4double))); 81 delete[] value; 169 delete[] value; 82 170 83 return true; 171 return true; 84 } 172 } 85 173 86 // ------------------------------------------- 174 // -------------------------------------------------------------- >> 175 87 G4bool G4PhysicsVector::Retrieve(std::ifstream 176 G4bool G4PhysicsVector::Retrieve(std::ifstream& fIn, G4bool ascii) 88 { 177 { 89 // clear properties; 178 // clear properties; 90 dataVector.clear(); 179 dataVector.clear(); 91 binVector.clear(); 180 binVector.clear(); 92 secDerivative.clear(); 181 secDerivative.clear(); 93 182 94 // retrieve in ascii mode 183 // retrieve in ascii mode 95 if (ascii) << 184 if(ascii) 96 { 185 { 97 // binning 186 // binning 98 fIn >> edgeMin >> edgeMax >> numberOfNodes 187 fIn >> edgeMin >> edgeMax >> numberOfNodes; 99 if (fIn.fail() || numberOfNodes < 2) << 188 if(fIn.fail() || numberOfNodes < 2) 100 { 189 { 101 return false; 190 return false; 102 } 191 } 103 // contents 192 // contents 104 G4int siz0 = 0; << 193 G4int siz = 0; 105 fIn >> siz0; << 194 fIn >> siz; 106 if (siz0 < 2) { return false; } << 195 if(fIn.fail() || siz != G4int(numberOfNodes)) 107 auto siz = static_cast<std::size_t>(siz0); << 108 if (fIn.fail() || siz != numberOfNodes) << 109 { 196 { 110 return false; 197 return false; 111 } 198 } 112 199 113 binVector.reserve(siz); 200 binVector.reserve(siz); 114 dataVector.reserve(siz); 201 dataVector.reserve(siz); 115 G4double vBin, vData; 202 G4double vBin, vData; 116 203 117 for (std::size_t i = 0; i < siz; ++i) << 204 for(G4int i = 0; i < siz; ++i) 118 { 205 { 119 vBin = 0.; 206 vBin = 0.; 120 vData = 0.; 207 vData = 0.; 121 fIn >> vBin >> vData; 208 fIn >> vBin >> vData; 122 if (fIn.fail()) << 209 if(fIn.fail()) 123 { 210 { 124 return false; 211 return false; 125 } 212 } 126 binVector.push_back(vBin); 213 binVector.push_back(vBin); 127 dataVector.push_back(vData); 214 dataVector.push_back(vData); 128 } 215 } 129 Initialise(); << 216 >> 217 // to remove any inconsistency >> 218 numberOfNodes = siz; >> 219 edgeMin = binVector[0]; >> 220 edgeMax = binVector[numberOfNodes - 1]; 130 return true; 221 return true; 131 } 222 } 132 223 133 // retrieve in binary mode 224 // retrieve in binary mode 134 // binning 225 // binning 135 fIn.read((char*) (&edgeMin), sizeof edgeMin) 226 fIn.read((char*) (&edgeMin), sizeof edgeMin); 136 fIn.read((char*) (&edgeMax), sizeof edgeMax) 227 fIn.read((char*) (&edgeMax), sizeof edgeMax); 137 fIn.read((char*) (&numberOfNodes), sizeof nu 228 fIn.read((char*) (&numberOfNodes), sizeof numberOfNodes); 138 229 139 // contents 230 // contents 140 std::size_t size; 231 std::size_t size; 141 fIn.read((char*) (&size), sizeof size); 232 fIn.read((char*) (&size), sizeof size); 142 233 143 auto value = new G4double[2 * size]; << 234 G4double* value = new G4double[2 * size]; 144 fIn.read((char*) (value), 2 * size * (sizeof 235 fIn.read((char*) (value), 2 * size * (sizeof(G4double))); 145 if (static_cast<G4int>(fIn.gcount()) != stat << 236 if(G4int(fIn.gcount()) != G4int(2 * size * (sizeof(G4double)))) 146 { 237 { 147 delete[] value; 238 delete[] value; 148 return false; 239 return false; 149 } 240 } 150 241 151 binVector.reserve(size); 242 binVector.reserve(size); 152 dataVector.reserve(size); 243 dataVector.reserve(size); 153 for (std::size_t i = 0; i < size; ++i) << 244 for(std::size_t i = 0; i < size; ++i) 154 { 245 { 155 binVector.push_back(value[2 * i]); 246 binVector.push_back(value[2 * i]); 156 dataVector.push_back(value[2 * i + 1]); 247 dataVector.push_back(value[2 * i + 1]); 157 } 248 } 158 delete[] value; 249 delete[] value; 159 250 160 Initialise(); << 251 // to remove any inconsistency >> 252 numberOfNodes = size; >> 253 edgeMin = binVector[0]; >> 254 edgeMax = binVector[numberOfNodes - 1]; >> 255 161 return true; 256 return true; 162 } 257 } 163 258 164 // ------------------------------------------- 259 // -------------------------------------------------------------- >> 260 165 void G4PhysicsVector::DumpValues(G4double unit 261 void G4PhysicsVector::DumpValues(G4double unitE, G4double unitV) const 166 { 262 { 167 for (std::size_t i = 0; i < numberOfNodes; + << 263 for(std::size_t i = 0; i < numberOfNodes; ++i) 168 { 264 { 169 G4cout << binVector[i] / unitE << " " << << 265 G4cout << binVector[i] / unitE << " " << dataVector[i] / unitV << G4endl; 170 << G4endl; << 171 } 266 } 172 } 267 } 173 268 174 // ------------------------------------------- 269 // -------------------------------------------------------------------- 175 std::size_t G4PhysicsVector::FindBin(const G4d << 176 std::size << 177 { << 178 if (idx + 1 < numberOfNodes && << 179 energy >= binVector[idx] && energy <= bi << 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 270 194 // ------------------------------------------- << 271 void G4PhysicsVector::ScaleVector(G4double factorE, G4double factorV) 195 void G4PhysicsVector::ScaleVector(const G4doub << 196 const G4doub << 197 { 272 { 198 for (std::size_t i = 0; i < numberOfNodes; + << 273 std::size_t n = dataVector.size(); >> 274 std::size_t i; >> 275 for(i = 0; i < n; ++i) 199 { 276 { 200 binVector[i] *= factorE; 277 binVector[i] *= factorE; 201 dataVector[i] *= factorV; 278 dataVector[i] *= factorV; 202 } 279 } 203 Initialise(); << 280 secDerivative.clear(); >> 281 >> 282 edgeMin = binVector[0]; >> 283 edgeMax = binVector[n - 1]; 204 } 284 } 205 285 206 // ------------------------------------------- << 286 // -------------------------------------------------------------- 207 void G4PhysicsVector::FillSecondDerivatives(co << 287 208 const G4double dir1, << 288 void G4PhysicsVector::ComputeSecondDerivatives(G4double firstPointDerivative, 209 const G4double dir2) << 289 G4double endPointDerivative) 210 { << 290 // A standard method of computation of second derivatives 211 if (!useSpline) { return; } << 291 // First derivatives at the first and the last point should be provided 212 // cannot compute derivatives for less than << 292 // See for example W.H. Press et al. "Numerical recipes in C" 213 const std::size_t nmin = (stype == G4SplineT << 293 // Cambridge University Press, 1997. 214 if (nmin > numberOfNodes) << 294 { 215 { << 295 if(4 > numberOfNodes) // cannot compute derivatives for less than 4 bins 216 if (0 < verboseLevel) << 296 { 217 { << 297 ComputeSecDerivatives(); 218 G4cout << "### G4PhysicsVector: spline c << 219 << numberOfNodes << " points - spline d << 220 << G4endl; << 221 DumpValues(); << 222 } << 223 useSpline = false; << 224 return; 298 return; 225 } 299 } 226 // check energies of free vector << 300 227 if (type == T_G4PhysicsFreeVector) << 301 if(!SplinePossible()) 228 { 302 { 229 for (std::size_t i=0; i<=idxmax; ++i) << 303 return; 230 { << 231 if (binVector[i + 1] <= binVector[i]) << 232 { << 233 if (0 < verboseLevel) << 234 { << 235 G4cout << "### G4PhysicsVector: spline can << 236 << " E[" << i << "]=" << binVector[i] << 237 << " >= E[" << i+1 << "]=" << binVector[i << 238 << G4endl; << 239 DumpValues(); << 240 } << 241 useSpline = false; << 242 return; << 243 } << 244 } << 245 } 304 } 246 305 247 // spline is possible << 306 useSpline = true; 248 Initialise(); << 249 secDerivative.resize(numberOfNodes); << 250 307 251 if (1 < verboseLevel) << 308 G4int n = numberOfNodes - 1; >> 309 >> 310 G4double* u = new G4double[n]; >> 311 >> 312 G4double p, sig, un; >> 313 >> 314 u[0] = (6.0 / (binVector[1] - binVector[0])) * >> 315 ((dataVector[1] - dataVector[0]) / (binVector[1] - binVector[0]) - >> 316 firstPointDerivative); >> 317 >> 318 secDerivative[0] = -0.5; >> 319 >> 320 // Decomposition loop for tridiagonal algorithm. secDerivative[i] >> 321 // and u[i] are used for temporary storage of the decomposed factors. >> 322 >> 323 for(G4int i = 1; i < n; ++i) 252 { 324 { 253 G4cout << "### G4PhysicsVector:: FillSecon << 325 sig = 254 << numberOfNodes << G4endl; << 326 (binVector[i] - binVector[i - 1]) / (binVector[i + 1] - binVector[i - 1]); 255 DumpValues(); << 327 p = sig * (secDerivative[i - 1]) + 2.0; >> 328 secDerivative[i] = (sig - 1.0) / p; >> 329 u[i] = >> 330 (dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) - >> 331 (dataVector[i] - dataVector[i - 1]) / (binVector[i] - binVector[i - 1]); >> 332 u[i] = >> 333 6.0 * u[i] / (binVector[i + 1] - binVector[i - 1]) - sig * u[i - 1] / p; 256 } 334 } 257 335 258 switch(stype) << 336 sig = 259 { << 337 (binVector[n - 1] - binVector[n - 2]) / (binVector[n] - binVector[n - 2]); 260 case G4SplineType::Base: << 338 p = sig * secDerivative[n - 2] + 2.0; 261 ComputeSecDerivative1(); << 339 un = (6.0 / (binVector[n] - binVector[n - 1])) * 262 break; << 340 (endPointDerivative - (dataVector[n] - dataVector[n - 1]) / >> 341 (binVector[n] - binVector[n - 1])) - >> 342 u[n - 1] / p; >> 343 secDerivative[n] = un / (secDerivative[n - 1] + 2.0); 263 344 264 case G4SplineType::FixedEdges: << 345 // The back-substitution loop for the triagonal algorithm of solving 265 ComputeSecDerivative2(dir1, dir2); << 346 // a linear system of equations. 266 break; << 267 347 268 default: << 348 for(G4int k = n - 1; k > 0; --k) 269 ComputeSecDerivative0(); << 349 { >> 350 secDerivative[k] *= >> 351 (secDerivative[k + 1] - u[k] * (binVector[k + 1] - binVector[k - 1]) / >> 352 (binVector[k + 1] - binVector[k])); 270 } 353 } >> 354 secDerivative[0] = 0.5 * (u[0] - secDerivative[1]); >> 355 >> 356 delete[] u; 271 } 357 } 272 358 273 // ------------------------------------------- 359 // -------------------------------------------------------------- 274 void G4PhysicsVector::ComputeSecDerivative0() << 360 275 // A simplified method of computation of seco << 361 void G4PhysicsVector::FillSecondDerivatives() 276 { 362 { 277 std::size_t n = numberOfNodes - 1; << 363 // Computation of second derivatives using "Not-a-knot" endpoint conditions >> 364 // B.I. Kvasov "Methods of shape-preserving spline approximation" >> 365 // World Scientific, 2000 278 366 279 for (std::size_t i = 1; i < n; ++i) << 367 if(5 > numberOfNodes) // cannot compute derivatives for less than 4 points 280 { 368 { 281 secDerivative[i] = 3.0 * << 369 ComputeSecDerivatives(); 282 ((dataVector[i + 1] - dataVector[i]) / ( << 370 return; 283 (dataVector[i] - dataVector[i - 1]) / << 284 (binVector[i] - binVector[i - 1])) / << 285 (binVector[i + 1] - binVector[i - 1]); << 286 } 371 } 287 secDerivative[n] = secDerivative[n - 1]; << 288 secDerivative[0] = secDerivative[1]; << 289 } << 290 372 291 // ------------------------------------------- << 373 if(!SplinePossible()) 292 void G4PhysicsVector::ComputeSecDerivative1() << 374 { 293 // Computation of second derivatives using "No << 375 return; 294 // B.I. Kvasov "Methods of shape-preserving sp << 376 } 295 // World Scientific, 2000 << 377 296 { << 378 useSpline = true; 297 std::size_t n = numberOfNodes - 1; << 379 298 auto u = new G4double[n]; << 380 G4int n = numberOfNodes - 1; >> 381 >> 382 G4double* u = new G4double[n]; >> 383 299 G4double p, sig; 384 G4double p, sig; 300 385 301 u[1] = ((dataVector[2] - dataVector[1]) / (b 386 u[1] = ((dataVector[2] - dataVector[1]) / (binVector[2] - binVector[1]) - 302 (dataVector[1] - dataVector[0]) / (b 387 (dataVector[1] - dataVector[0]) / (binVector[1] - binVector[0])); 303 u[1] = 6.0 * u[1] * (binVector[2] - binVecto 388 u[1] = 6.0 * u[1] * (binVector[2] - binVector[1]) / 304 ((binVector[2] - binVector[0]) * (bin 389 ((binVector[2] - binVector[0]) * (binVector[2] - binVector[0])); 305 390 306 // Decomposition loop for tridiagonal algori 391 // Decomposition loop for tridiagonal algorithm. secDerivative[i] 307 // and u[i] are used for temporary storage o 392 // and u[i] are used for temporary storage of the decomposed factors. 308 393 309 secDerivative[1] = (2.0 * binVector[1] - bin 394 secDerivative[1] = (2.0 * binVector[1] - binVector[0] - binVector[2]) / 310 (2.0 * binVector[2] - bin 395 (2.0 * binVector[2] - binVector[0] - binVector[1]); 311 396 312 for(std::size_t i = 2; i < n - 1; ++i) << 397 for(G4int i = 2; i < n - 1; ++i) 313 { 398 { 314 sig = 399 sig = 315 (binVector[i] - binVector[i - 1]) / (bin 400 (binVector[i] - binVector[i - 1]) / (binVector[i + 1] - binVector[i - 1]); 316 p = sig * secDerivative[i - 401 p = sig * secDerivative[i - 1] + 2.0; 317 secDerivative[i] = (sig - 1.0) / p; 402 secDerivative[i] = (sig - 1.0) / p; 318 u[i] = 403 u[i] = 319 (dataVector[i + 1] - dataVector[i]) / (b 404 (dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) - 320 (dataVector[i] - dataVector[i - 1]) / (b 405 (dataVector[i] - dataVector[i - 1]) / (binVector[i] - binVector[i - 1]); 321 u[i] = 406 u[i] = 322 (6.0 * u[i] / (binVector[i + 1] - binVec 407 (6.0 * u[i] / (binVector[i + 1] - binVector[i - 1])) - sig * u[i - 1] / p; 323 } 408 } 324 409 325 sig = 410 sig = 326 (binVector[n - 1] - binVector[n - 2]) / (b 411 (binVector[n - 1] - binVector[n - 2]) / (binVector[n] - binVector[n - 2]); 327 p = sig * secDerivative[n - 3] + 2.0; 412 p = sig * secDerivative[n - 3] + 2.0; 328 u[n - 1] = 413 u[n - 1] = 329 (dataVector[n] - dataVector[n - 1]) / (bin 414 (dataVector[n] - dataVector[n - 1]) / (binVector[n] - binVector[n - 1]) - 330 (dataVector[n - 1] - dataVector[n - 2]) / 415 (dataVector[n - 1] - dataVector[n - 2]) / 331 (binVector[n - 1] - binVector[n - 2]); 416 (binVector[n - 1] - binVector[n - 2]); 332 u[n - 1] = 6.0 * sig * u[n - 1] / (binVector 417 u[n - 1] = 6.0 * sig * u[n - 1] / (binVector[n] - binVector[n - 2]) - 333 (2.0 * sig - 1.0) * u[n - 2] / p; 418 (2.0 * sig - 1.0) * u[n - 2] / p; 334 419 335 p = (1.0 + sig) + (2.0 * sig - 1.0) * secDer << 420 p = (1.0 + sig) + (2.0 * sig - 1.0) * secDerivative[n - 2]; 336 secDerivative[n - 1] = u[n - 1] / p; 421 secDerivative[n - 1] = u[n - 1] / p; 337 422 338 // The back-substitution loop for the triago 423 // The back-substitution loop for the triagonal algorithm of solving 339 // a linear system of equations. 424 // a linear system of equations. 340 425 341 for (std::size_t k = n - 2; k > 1; --k) << 426 for(G4int k = n - 2; k > 1; --k) 342 { 427 { 343 secDerivative[k] *= 428 secDerivative[k] *= 344 (secDerivative[k + 1] - u[k] * (binVecto 429 (secDerivative[k + 1] - u[k] * (binVector[k + 1] - binVector[k - 1]) / 345 (binVector[k + 430 (binVector[k + 1] - binVector[k])); 346 } 431 } 347 secDerivative[n] = 432 secDerivative[n] = 348 (secDerivative[n - 1] - (1.0 - sig) * secD 433 (secDerivative[n - 1] - (1.0 - sig) * secDerivative[n - 2]) / sig; 349 sig = 1.0 - ((binVector[2] - binVector[1]) / 434 sig = 1.0 - ((binVector[2] - binVector[1]) / (binVector[2] - binVector[0])); 350 secDerivative[1] *= (secDerivative[2] - u[1] 435 secDerivative[1] *= (secDerivative[2] - u[1] / (1.0 - sig)); 351 secDerivative[0] = (secDerivative[1] - sig * 436 secDerivative[0] = (secDerivative[1] - sig * secDerivative[2]) / (1.0 - sig); 352 437 353 delete[] u; 438 delete[] u; 354 } 439 } 355 440 356 // ------------------------------------------- 441 // -------------------------------------------------------------- 357 void G4PhysicsVector::ComputeSecDerivative2(G4 << 442 358 G4 << 443 void G4PhysicsVector::ComputeSecDerivatives() 359 // A standard method of computation of second << 360 // First derivatives at the first and the last << 361 // See for example W.H. Press et al. "Numerica << 362 // Cambridge University Press, 1997. << 363 { 444 { 364 std::size_t n = numberOfNodes - 1; << 445 // A simplified method of computation of second derivatives 365 auto u = new G4double[n]; << 366 G4double p, sig, un; << 367 446 368 u[0] = (6.0 / (binVector[1] - binVector[0])) << 447 if(3 > numberOfNodes) // cannot compute derivatives for less than 4 bins 369 ((dataVector[1] - dataVector[0]) / (b << 448 { 370 firstPointDerivative); << 449 useSpline = false; >> 450 return; >> 451 } 371 452 372 secDerivative[0] = -0.5; << 453 if(!SplinePossible()) >> 454 { >> 455 return; >> 456 } 373 457 374 // Decomposition loop for tridiagonal algori << 458 useSpline = true; 375 // and u[i] are used for temporary storage o << 459 >> 460 std::size_t n = numberOfNodes - 1; 376 461 377 for (std::size_t i = 1; i < n; ++i) << 462 for(std::size_t i = 1; i < n; ++i) 378 { 463 { 379 sig = << 464 secDerivative[i] = 380 (binVector[i] - binVector[i - 1]) / (bin << 465 3.0 * 381 p = sig * (secDerivative[i << 466 ((dataVector[i + 1] - dataVector[i]) / (binVector[i + 1] - binVector[i]) - 382 secDerivative[i] = (sig - 1.0) / p; << 467 (dataVector[i] - dataVector[i - 1]) / 383 u[i] = << 468 (binVector[i] - binVector[i - 1])) / 384 (dataVector[i + 1] - dataVector[i]) / (b << 469 (binVector[i + 1] - binVector[i - 1]); 385 (dataVector[i] - dataVector[i - 1]) / (b << 386 u[i] = << 387 6.0 * u[i] / (binVector[i + 1] - binVect << 388 } 470 } >> 471 secDerivative[n] = secDerivative[n - 1]; >> 472 secDerivative[0] = secDerivative[1]; >> 473 } 389 474 390 sig = << 475 // -------------------------------------------------------------- 391 (binVector[n - 1] - binVector[n - 2]) / (b << 392 p = sig * secDerivative[n - 2] + 2.0; << 393 un = (6.0 / (binVector[n] - binVector[n - 1] << 394 (endPointDerivative - (dataVector[n] << 395 (binVector[n] << 396 u[n - 1] / p; << 397 secDerivative[n] = un / (secDerivative[n - 1 << 398 476 399 // The back-substitution loop for the triago << 477 G4bool G4PhysicsVector::SplinePossible() 400 // a linear system of equations. << 478 { >> 479 // Initialise second derivative array. If neighbor energy coincide >> 480 // or not ordered than spline cannot be applied 401 481 402 for (std::size_t k = n - 1; k > 0; --k) << 482 G4bool result = true; >> 483 for(std::size_t j = 1; j < numberOfNodes; ++j) 403 { 484 { 404 secDerivative[k] *= << 485 if(binVector[j] <= binVector[j - 1]) 405 (secDerivative[k + 1] - u[k] * (binVecto << 486 { 406 (binVector[k + << 487 result = false; >> 488 useSpline = false; >> 489 secDerivative.clear(); >> 490 break; >> 491 } 407 } 492 } 408 secDerivative[0] = 0.5 * (u[0] - secDerivati << 493 secDerivative.resize(numberOfNodes, 0.0); 409 << 494 return result; 410 delete[] u; << 411 } 495 } 412 496 413 // ------------------------------------------- 497 // -------------------------------------------------------------- >> 498 414 std::ostream& operator<<(std::ostream& out, co 499 std::ostream& operator<<(std::ostream& out, const G4PhysicsVector& pv) 415 { 500 { 416 // binning 501 // binning 417 G4long prec = out.precision(); << 502 G4int prec = out.precision(); 418 out << std::setprecision(12) << pv.edgeMin < 503 out << std::setprecision(12) << pv.edgeMin << " " << pv.edgeMax << " " 419 << pv.numberOfNodes << G4endl; 504 << pv.numberOfNodes << G4endl; 420 505 421 // contents 506 // contents 422 out << pv.dataVector.size() << G4endl; 507 out << pv.dataVector.size() << G4endl; 423 for (std::size_t i = 0; i < pv.dataVector.si << 508 for(std::size_t i = 0; i < pv.dataVector.size(); ++i) 424 { 509 { 425 out << pv.binVector[i] << " " << pv.dataV 510 out << pv.binVector[i] << " " << pv.dataVector[i] << G4endl; 426 } 511 } 427 out.precision(prec); << 512 out << std::setprecision(prec); 428 513 429 return out; 514 return out; 430 } 515 } 431 516 432 //-------------------------------------------- 517 //--------------------------------------------------------------- 433 G4double G4PhysicsVector::GetEnergy(const G4do << 518 >> 519 G4double G4PhysicsVector::Value(G4double theEnergy, std::size_t& lastIdx) const 434 { 520 { 435 if (0 == numberOfNodes) << 521 G4double y; >> 522 if(theEnergy <= edgeMin) 436 { 523 { 437 return 0.0; << 524 lastIdx = 0; >> 525 y = dataVector[0]; 438 } 526 } 439 if (1 == numberOfNodes || val <= dataVector[ << 527 else if(theEnergy >= edgeMax) 440 { 528 { 441 return edgeMin; << 529 lastIdx = numberOfNodes - 1; >> 530 y = dataVector[lastIdx]; 442 } 531 } 443 if (val >= dataVector[numberOfNodes - 1]) << 532 else 444 { 533 { 445 return edgeMax; << 534 lastIdx = FindBin(theEnergy, lastIdx); >> 535 y = Interpolation(lastIdx, theEnergy); 446 } 536 } 447 std::size_t bin = std::lower_bound(dataVecto << 537 return y; 448 - dataVector.cbegin() - 1; << 538 } 449 if (bin > idxmax) { bin = idxmax; } << 539 >> 540 //--------------------------------------------------------------- >> 541 >> 542 G4double G4PhysicsVector::FindLinearEnergy(G4double rand) const >> 543 { >> 544 if(1 >= numberOfNodes) >> 545 { >> 546 return 0.0; >> 547 } >> 548 G4double y = rand * dataVector[numberOfNodes - 1]; >> 549 std::size_t bin = std::lower_bound(dataVector.begin(), dataVector.end(), y) - >> 550 dataVector.begin() - 1; >> 551 bin = std::min(bin, numberOfNodes - 2); 450 G4double res = binVector[bin]; 552 G4double res = binVector[bin]; 451 G4double del = dataVector[bin + 1] - dataVec 553 G4double del = dataVector[bin + 1] - dataVector[bin]; 452 if (del > 0.0) << 554 if(del > 0.0) 453 { 555 { 454 res += (val - dataVector[bin]) * (binVecto << 556 res += (y - dataVector[bin]) * (binVector[bin + 1] - res) / del; 455 } 557 } 456 return res; 558 return res; 457 } 559 } 458 560 459 //-------------------------------------------- 561 //--------------------------------------------------------------- 460 void G4PhysicsVector::PrintPutValueError(std:: << 562 461 G4dou << 563 void G4PhysicsVector::PrintPutValueError(std::size_t index) 462 const << 463 { 564 { 464 G4ExceptionDescription ed; 565 G4ExceptionDescription ed; 465 ed << "Vector type: " << type << " length= " << 566 ed << "Vector type " << type << " length= " << numberOfNodes 466 << "; an attempt to put data at index= " << 567 << " an attempt to put data at index= " << index; 467 << " value= " << val << " in " << text; << 568 G4Exception("G4PhysicsVector::PutValue()", "gl0005", FatalException, ed, 468 G4Exception("G4PhysicsVector:", "gl0005", << 569 "Memory overwritten"); 469 FatalException, ed, "Wrong opera << 470 } 570 } 471 << 472 //-------------------------------------------- << 473 571