#include <ThetaImplicitEM.H>
|
| | ThetaImplicitEM ()=default |
| |
| | ~ThetaImplicitEM () override=default |
| |
| | ThetaImplicitEM (const ThetaImplicitEM &)=delete |
| |
| ThetaImplicitEM & | operator= (const ThetaImplicitEM &)=delete |
| |
| | ThetaImplicitEM (ThetaImplicitEM &&)=delete |
| |
| ThetaImplicitEM & | operator= (ThetaImplicitEM &&)=delete |
| |
| void | Define (WarpX *a_WarpX) override |
| | Read user-provided parameters that control the implicit solver. Allocate internal arrays for intermediate field values needed by the solver. More...
|
| |
| void | PrintParameters () const override |
| |
| void | OneStep (amrex::Real a_time, amrex::Real a_dt, int a_step) override |
| | Advance fields and particles by one time step using the specified implicit algorithm. More...
|
| |
| void | ComputeRHS (WarpXSolverVec &a_RHS, const WarpXSolverVec &a_E, amrex::Real a_time, amrex::Real a_dt, int a_nl_iter, bool a_from_jacobian) override |
| | Computes the RHS of the equation corresponding to the specified implicit algorithm. The discrete equations corresponding to numerical integration of ODEs are often written in the form U = b + RHS(U), where U is the variable to be solved for (e.g., the solution at the next time step), b is a constant (i.e., the solution from the previous time step), and RHS(U) is the right-hand-side of the equation. Iterative solvers, such as Picard and Newton, and higher-order Runge-Kutta methods, need to compute RHS(U) multiple times for different values of U. Thus, a routine that returns this value is needed. e.g., Ebar - E^n = cvac^2*0.5*dt*(curl(Bbar) - mu0*Jbar(Ebar,Bbar)) Here, U = Ebar, b = E^n, and the expression on the right is RHS(U). More...
|
| |
| amrex::Real | theta () const |
| |
| | ImplicitSolver ()=default |
| |
| virtual | ~ImplicitSolver ()=default |
| |
| | ImplicitSolver (const ImplicitSolver &)=delete |
| |
| ImplicitSolver & | operator= (const ImplicitSolver &)=delete |
| |
| | ImplicitSolver (ImplicitSolver &&)=delete |
| |
| ImplicitSolver & | operator= (ImplicitSolver &&)=delete |
| |
| bool | IsDefined () const |
| |
| void | GetParticleSolverParams (int &a_max_particle_iter, amrex::ParticleReal &a_particle_tol) const |
| |
|
| amrex::Real | m_theta = 0.5 |
| | Time-biasing parameter for fields used on RHS to advance system. More...
|
| |
| WarpXSolverVec | m_E |
| | Solver vectors to be used in the nonlinear solver to solve for the electric field E. The main logic for determining which variables should be WarpXSolverVec type is that it must have the same size and have the same centering of the data as the variable being solved for, which is E here. For example, if using a Yee grid then a container for curlB could be a WarpXSovlerVec, but magnetic field B should not be. More...
|
| |
| WarpXSolverVec | m_Eold |
| |
| amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > | m_Bold |
| | B is a derived variable from E. Need to save Bold to update B during the iterative nonlinear solve for E. Bold is owned here, but only used by WarpX. It is not used directly by the nonlinear solver, nor is it the same size as the solver vector (size E), and so it should not be WarpXSolverVec type. More...
|
| |
◆ ThetaImplicitEM() [1/3]
| ThetaImplicitEM::ThetaImplicitEM |
( |
| ) |
|
|
default |
◆ ~ThetaImplicitEM()
| ThetaImplicitEM::~ThetaImplicitEM |
( |
| ) |
|
|
overridedefault |
◆ ThetaImplicitEM() [2/3]
◆ ThetaImplicitEM() [3/3]
◆ ComputeRHS()
| void ThetaImplicitEM::ComputeRHS |
( |
WarpXSolverVec & |
a_RHS, |
|
|
const WarpXSolverVec & |
a_E, |
|
|
amrex::Real |
a_time, |
|
|
amrex::Real |
a_dt, |
|
|
int |
a_nl_iter, |
|
|
bool |
a_from_jacobian |
|
) |
| |
|
overridevirtual |
Computes the RHS of the equation corresponding to the specified implicit algorithm. The discrete equations corresponding to numerical integration of ODEs are often written in the form U = b + RHS(U), where U is the variable to be solved for (e.g., the solution at the next time step), b is a constant (i.e., the solution from the previous time step), and RHS(U) is the right-hand-side of the equation. Iterative solvers, such as Picard and Newton, and higher-order Runge-Kutta methods, need to compute RHS(U) multiple times for different values of U. Thus, a routine that returns this value is needed. e.g., Ebar - E^n = cvac^2*0.5*dt*(curl(Bbar) - mu0*Jbar(Ebar,Bbar)) Here, U = Ebar, b = E^n, and the expression on the right is RHS(U).
Implements ImplicitSolver.
◆ Define()
| void ThetaImplicitEM::Define |
( |
WarpX * |
a_WarpX | ) |
|
|
overridevirtual |
Read user-provided parameters that control the implicit solver. Allocate internal arrays for intermediate field values needed by the solver.
Implements ImplicitSolver.
◆ FinishFieldUpdate()
| void ThetaImplicitEM::FinishFieldUpdate |
( |
amrex::Real |
a_new_time | ) |
|
|
private |
Nonlinear solver is for the time-centered values of E. After the solver, need to use m_E and m_Eold to compute E^{n+1}.
◆ OneStep()
| void ThetaImplicitEM::OneStep |
( |
amrex::Real |
a_time, |
|
|
amrex::Real |
a_dt, |
|
|
int |
a_step |
|
) |
| |
|
overridevirtual |
Advance fields and particles by one time step using the specified implicit algorithm.
Implements ImplicitSolver.
◆ operator=() [1/2]
◆ operator=() [2/2]
◆ PrintParameters()
| void ThetaImplicitEM::PrintParameters |
( |
| ) |
const |
|
overridevirtual |
◆ theta()
| amrex::Real ThetaImplicitEM::theta |
( |
| ) |
const |
|
inline |
◆ UpdateWarpXFields()
| void ThetaImplicitEM::UpdateWarpXFields |
( |
const WarpXSolverVec & |
a_E, |
|
|
amrex::Real |
a_time, |
|
|
amrex::Real |
a_dt |
|
) |
| |
|
private |
Update the E and B fields owned by WarpX.
◆ m_Bold
B is a derived variable from E. Need to save Bold to update B during the iterative nonlinear solve for E. Bold is owned here, but only used by WarpX. It is not used directly by the nonlinear solver, nor is it the same size as the solver vector (size E), and so it should not be WarpXSolverVec type.
◆ m_E
Solver vectors to be used in the nonlinear solver to solve for the electric field E. The main logic for determining which variables should be WarpXSolverVec type is that it must have the same size and have the same centering of the data as the variable being solved for, which is E here. For example, if using a Yee grid then a container for curlB could be a WarpXSovlerVec, but magnetic field B should not be.
◆ m_Eold
◆ m_theta
| amrex::Real ThetaImplicitEM::m_theta = 0.5 |
|
private |
Time-biasing parameter for fields used on RHS to advance system.
The documentation for this class was generated from the following files:
- /home/docs/checkouts/readthedocs.org/user_builds/warpx/checkouts/24.08/Source/FieldSolver/ImplicitSolvers/ThetaImplicitEM.H
- /home/docs/checkouts/readthedocs.org/user_builds/warpx/checkouts/24.08/Source/FieldSolver/ImplicitSolvers/ThetaImplicitEM.cpp