Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TSimpleHeum.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// G4TSimpleHeum
27//
28// Class description:
29//
30// Templated version of G4SimpleHeum
31//
32// Created: Josh Xie (supported by Google Summer of Code 2014 )
33// Supervisors: Sandro Wenzel, John Apostolakis (CERN)
34// --------------------------------------------------------------------
35// Adapted from G4G4TSimpleHeum class
36// Original desription:
37//
38// Simple Heum stepper for magnetic field:
39// x_1 = x_0 +
40// h * 1/4 * dx(t0,x0) +
41// 3/4 * dx(t0+2/3*h, x0+2/3*h*(dx(t0+h/3,x0+h/3*dx(t0,x0))))
42//
43// Third order solver.
44// Created: W.Wander <[email protected]>, 12/09/1997
45// --------------------------------------------------------------------
46
47#ifndef TSIMPLEHEUM_HH
48#define TSIMPLEHEUM_HH
49
50#include <cassert>
51#include "G4TMagErrorStepper.hh"
52#include "G4ThreeVector.hh"
53
54template <class T_Equation, unsigned int N>
56 : public G4TMagErrorStepper<G4TSimpleHeum<T_Equation, N>, T_Equation, N>
57{
58 public: // with description
59 constexpr static unsigned int gIntegratorOrder = 3;
60 static constexpr double IntegratorCorrection= 1.0 /
61 ((1<<gIntegratorOrder) - 1);
62
63 G4TSimpleHeum(T_Equation* EqRhs, unsigned int numberOfVariables = 6);
64
66 // Constructor and destructor.
67
68 inline void RightHandSide(G4double y[],
69 G4double dydx[])
70 {
71 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
72 }
73
74 inline void DumbStepper(const G4double yIn[],
75 const G4double dydx[],
76 G4double h, G4double yOut[]); // override final
77
78 public: // without description
80
81 private:
82 G4int fNumberOfVariables;
83
84 G4double dydxTemp[N];
85 G4double dydxTemp2[N];
86 G4double yTemp[N];
87 G4double yTemp2[N];
88 // scratch space
89
90 T_Equation* fEquation_Rhs;
91};
92
93template <class T_Equation, unsigned int N >
95 unsigned int numberOfVariables )
96 : G4TMagErrorStepper<G4TSimpleHeum<T_Equation, N>, T_Equation, N>(
97 EqRhs, numberOfVariables)
98 , fNumberOfVariables(numberOfVariables)
99 , fEquation_Rhs(EqRhs)
100{
101 assert(fNumberOfVariables == N);
102 if( dynamic_cast<G4EquationOfMotion*>(EqRhs) == nullptr )
103 {
104 G4Exception("G4TSimpleHeum: constructor", "GeomField0001",
105 FatalException, "Equation is not an G4EquationOfMotion.");
106 }
107}
108
109template <class T_Equation, unsigned int N >
110inline void
112 const G4double dydx[],
113 G4double h, G4double yOut[])
114{
115 for(unsigned int i = 0; i < N; ++i)
116 {
117 yTemp[i] = yIn[i] + (1.0 / 3.0) * h * dydx[i];
118 }
119
120 this->RightHandSide(yTemp, dydxTemp);
121
122 for(unsigned int i = 0; i < N; ++i)
123 {
124 yTemp2[i] = yIn[i] + (2.0 / 3.0) * h * dydxTemp[i];
125 }
126
127 this->RightHandSide(yTemp2, dydxTemp2);
128
129 for(unsigned int i = 0; i < N; ++i)
130 {
131 yOut[i] = yIn[i] + h * (0.25 * dydx[i] + 0.75 * dydxTemp2[i]);
132 }
133
134 if(fNumberOfVariables == 12)
135 {
136 this->NormalisePolarizationVector(yOut);
137 }
138}
139
140#endif /* TSIMPLEHEUM_HH */
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
static constexpr double IntegratorCorrection
G4TSimpleHeum(T_Equation *EqRhs, unsigned int numberOfVariables=6)
G4int IntegratorOrder() const
void DumbStepper(const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[])
static constexpr unsigned int gIntegratorOrder
void RightHandSide(G4double y[], G4double dydx[])