Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4TSimpleRunge.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// G4TSimpleRunge
27//
28// Class description:
29//
30// Templated version of G4SimpleRunge
31//
32//
33// Created: Josh Xie (supported by Google Summer of Code 2014 )
34// Supervisors: Sandro Wenzel, John Apostolakis (CERN)
35// Adapted from G4G4TSimpleRunge class
36// --------------------------------------------------------------------
37#ifndef G4TSimpleRunge_HH
38#define G4TSimpleRunge_HH
39
40#include <cassert>
41#include "G4TMagErrorStepper.hh"
42#include "G4ThreeVector.hh"
43
44template <class T_Equation, int N>
46 : public G4TMagErrorStepper<G4TSimpleRunge<T_Equation, N>, T_Equation, N>
47{
48 public: // with description
49 static constexpr double IntegratorCorrection = 1. / ((1 << 2) - 1);
50
51 G4TSimpleRunge(T_Equation* EqRhs, G4int numberOfVariables = 6)
52 : G4TMagErrorStepper<G4TSimpleRunge<T_Equation, N>, T_Equation, N>(
53 EqRhs, numberOfVariables)
54 , fNumberOfVariables(numberOfVariables)
55 , fEquation_Rhs(EqRhs)
56
57 {
58 // default GetNumberOfStateVariables() == 12
59 assert(this->GetNumberOfStateVariables() <= 12);
60 }
61
63
64 inline void RightHandSide(G4double y[],
65 G4double dydx[])
66 {
67 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
68 }
69
70 inline void DumbStepper(const G4double yIn[],
71 const G4double dydx[],
72 G4double h, G4double yOut[]) // override final
73 {
74 // Initialise time to t0, needed when it is not updated by the integration.
75 yTemp[7] = yOut[7] = yIn[7]; // Better to set it to NaN; // TODO
76
77 for(G4int i = 0; i < N; ++i)
78 {
79 yTemp[i] = yIn[i] + 0.5 * h * dydx[i];
80 }
81
82 this->RightHandSide(yTemp, dydxTemp);
83
84 for(G4int i = 0; i < N; ++i)
85 {
86 yOut[i] = yIn[i] + h * (dydxTemp[i]);
87 }
88 }
89
90 public: // without description
91 inline G4int IntegratorOrder() const { return 2; }
92
93 private:
94 G4int fNumberOfVariables;
95 G4double dydxTemp[N > 12 ? N : 12];
96 G4double yTemp[N > 12 ? N : 12];
97
98 T_Equation* fEquation_Rhs;
99 // scratch space
100};
101
102#endif /* G4TSimpleRunge_HH */
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4int GetNumberOfStateVariables() const
static constexpr double IntegratorCorrection
G4int IntegratorOrder() const
void DumbStepper(const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[])
G4TSimpleRunge(T_Equation *EqRhs, G4int numberOfVariables=6)
void RightHandSide(G4double y[], G4double dydx[])
#define N
Definition crc32.c:57