Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TExplicitEuler.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// G4TExplicitEuler
27//
28// Class description:
29//
30// Templated version of G4ExplicitEuler.
31// Adapted from G4G4TExplicitEuler class.
32//
33// Information from G4Explicit Euler:
34// ----------------------------------
35// Explicit Euler stepper for magnetic field: x_1 = x_0 + h * dx_0.
36// The simplistic approach to solving linear differential equations.
37// Take the current derivative and add it to the current position.
38
39// Author: Josh Xie (CERN, Google Summer of Code 2014), June 2014
40// Supervisors: Sandro Wenzel, John Apostolakis (CERN)
41// --------------------------------------------------------------------
42#ifndef G4TExplicitEuler_HH
43#define G4TExplicitEuler_HH
44
45#include "G4TMagErrorStepper.hh"
46#include "G4ThreeVector.hh"
47
48/**
49 * @brief G4TExplicitEuler is a templated version of G4ExplicitEuler stepper.
50 */
51
52template <class T_Equation, int N>
54 : public G4TMagErrorStepper<G4TExplicitEuler<T_Equation, N>, T_Equation, N>
55{
56 public:
57
58 static constexpr double IntegratorCorrection = 1.;
59
60 G4TExplicitEuler(T_Equation* EqRhs, G4int numberOfVariables = N)
61 : G4TMagErrorStepper<G4TExplicitEuler<T_Equation, N>, T_Equation, N>
62 (EqRhs, numberOfVariables), fEquation_Rhs(EqRhs)
63 {
64 if( numberOfVariables != N )
65 {
67 msg << "Equation has an incompatible number of variables." ;
68 msg << " template N = " << N << " equation-Nvar= "
69 << numberOfVariables;
70 G4Exception("G4TExplicitEuler: constructor", "GeomField0003",
72 }
73 }
74
75 ~G4TExplicitEuler() = default;
76
77 inline void DumbStepper(const G4double yIn[],
78 const G4double dydx[],
79 G4double h, G4double yOut[]) // override final
80 {
81 // Initialise time to t0, needed when it is not updated by the integration.
82 // yOut[7] = yIn[7]; // Better to set it to NaN; // TODO
83
84 for(G4int i = 0; i < N; ++i)
85 {
86 yOut[i] = yIn[i] + h * dydx[i]; // 1st and only Step
87 }
88
89 return;
90 }
91
92 G4int IntegratorOrder() const { return 1; }
93
94 private:
95
96 T_Equation* fEquation_Rhs;
97};
98
99#endif
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
~G4TExplicitEuler()=default
static constexpr double IntegratorCorrection
void DumbStepper(const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[])
G4int IntegratorOrder() const
G4TExplicitEuler(T_Equation *EqRhs, G4int numberOfVariables=N)
G4TMagErrorStepper(T_Equation *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
#define N
Definition crc32.c:57