WarpX
Loading...
Searching...
No Matches
ImplicitSolver Class Referenceabstract

#include <ImplicitSolver.H>

Inheritance diagram for ImplicitSolver:
SemiImplicitEM StrangImplicitSpectralEM ThetaImplicitEM

Public Member Functions

 ImplicitSolver ()=default
 
virtual ~ImplicitSolver ()=default
 
 ImplicitSolver (const ImplicitSolver &)=delete
 
ImplicitSolveroperator= (const ImplicitSolver &)=delete
 
 ImplicitSolver (ImplicitSolver &&)=delete
 
ImplicitSolveroperator= (ImplicitSolver &&)=delete
 
virtual void Define (WarpX *a_WarpX)=0
 Read user-provided parameters that control the implicit solver. Allocate internal arrays for intermediate field values needed by the solver.
 
bool IsDefined () const
 
virtual void PrintParameters () const =0
 
void CreateParticleAttributes () const
 
bool DoParticleSuborbits ()
 
virtual void OneStep (amrex::Real a_time, amrex::Real a_dt, int a_step)=0
 Advance fields and particles by one time step using the specified implicit algorithm.
 
virtual void ComputeRHS (WarpXSolverVec &a_RHS, const WarpXSolverVec &a_E, amrex::Real a_time, int a_nl_iter, bool a_from_jacobian)=0
 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).
 
int numAMRLevels () const
 
const amrex::GeometryGetGeometry (int) const
 
const amrex::Array< FieldBoundaryType, 3 > & GetFieldBoundaryLo () const
 
const amrex::Array< FieldBoundaryType, 3 > & GetFieldBoundaryHi () const
 
amrex::Array< amrex::LinOpBCType, 3 > GetLinOpBCLo () const
 
amrex::Array< amrex::LinOpBCType, 3 > GetLinOpBCHi () const
 
const amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > * GetMassMatricesCoeff () const
 Return pointer to MultiFab array for mass matrix.
 
virtual amrex::Real GetThetaForPC () const =0
 
void ComputeJfromMassMatrices (const bool a_J_from_MM_only)
 

Protected Member Functions

void parseNonlinearSolverParams (const amrex::ParmParse &pp)
 parse nonlinear solver parameters (if one is used)
 
void InitializeMassMatrices ()
 Initialize the Mass Matrices used for plasma response in nonlinear Newton solver.
 
void PrintBaseImplicitSolverParameters () const
 
void PreRHSOp (const amrex::Real a_cur_time, const int a_nl_iter, const bool a_from_jacobian)
 Perform operations needed before computing the Right Hand Side.
 
void SyncMassMatricesPCAndApplyBCs ()
 Communicate Mass Matrices used for PC and apply boundary conditions.
 
void SetMassMatricesForPC (const amrex::Real a_theta_dt)
 Scale mass matrices used for PC by c^2*mu0*theta*dt and add 1 to diagonal terms.
 
void CumulateJ ()
 Add J from particles included in mass matrices at start of each Newton step.
 
void SaveE ()
 Save E at start of each Newton step.
 
amrex::Array< amrex::LinOpBCType, 3 > convertFieldBCToLinOpBC (const amrex::Array< FieldBoundaryType, 3 > &) const
 Convert from WarpX FieldBoundaryType to amrex::LinOpBCType.
 

Protected Attributes

WarpXm_WarpX
 Pointer back to main WarpX class.
 
bool m_is_defined = false
 
int m_num_amr_levels = 1
 Number of AMR levels.
 
amrex::Real m_dt = 0.0
 Time step.
 
amrex::Real m_theta = 0.5
 Time-biasing parameter for fields used on RHS to advance system.
 
NonlinearSolverType m_nlsolver_type
 Nonlinear solver type and object.
 
std::unique_ptr< NonlinearSolver< WarpXSolverVec, ImplicitSolver > > m_nlsolver
 
amrex::ParticleReal m_particle_tolerance = amrex::ParticleReal(1.0e-10)
 tolerance used by the iterative method used to obtain a self-consistent update of the particle positions and velocities for given E and B on the grid
 
int m_max_particle_iterations = 21
 maximum iterations for the iterative method used to obtain a self-consistent update of the particle positions and velocities for given E and B on the grid
 
bool m_particle_suborbits = false
 whether to use suborbits for particles that fail to converge in m_max_particle_iterations iterations
 
bool m_print_unconverged_particle_details = false
 whether the details of unconverged particles are printed out during the particle evolve loops
 
bool m_skip_particle_picard_init = false
 whether to skip the full Picard update of particles on the initial newton step
 
bool m_use_mass_matrices = false
 Whether to use mass matrices for the implicit solver.
 
bool m_use_mass_matrices_jacobian = false
 Whether to use mass matrices to compute J during linear stage of JFNK.
 
bool m_use_mass_matrices_pc = false
 Whether to use mass matrices in the preconditioner.
 
amrex::IntVect m_ncomp_xx = amrex::IntVect{0}
 Direction-dependent component number of mass matrices elements.
 
amrex::IntVect m_ncomp_xy = amrex::IntVect{0}
 
amrex::IntVect m_ncomp_xz = amrex::IntVect{0}
 
amrex::IntVect m_ncomp_yx = amrex::IntVect{0}
 
amrex::IntVect m_ncomp_yy = amrex::IntVect{0}
 
amrex::IntVect m_ncomp_yz = amrex::IntVect{0}
 
amrex::IntVect m_ncomp_zx = amrex::IntVect{0}
 
amrex::IntVect m_ncomp_zy = amrex::IntVect{0}
 
amrex::IntVect m_ncomp_zz = amrex::IntVect{0}
 
amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > m_mmpc_mfarrvec
 Array of multifab pointers to mass matrix preconditioner.
 

Constructor & Destructor Documentation

◆ ImplicitSolver() [1/3]

ImplicitSolver::ImplicitSolver ( )
default

◆ ~ImplicitSolver()

virtual ImplicitSolver::~ImplicitSolver ( )
virtualdefault

◆ ImplicitSolver() [2/3]

ImplicitSolver::ImplicitSolver ( const ImplicitSolver & )
delete

◆ ImplicitSolver() [3/3]

ImplicitSolver::ImplicitSolver ( ImplicitSolver && )
delete

Member Function Documentation

◆ ComputeJfromMassMatrices()

void ImplicitSolver::ComputeJfromMassMatrices ( const bool a_J_from_MM_only)

◆ ComputeRHS()

virtual void ImplicitSolver::ComputeRHS ( WarpXSolverVec & a_RHS,
const WarpXSolverVec & a_E,
amrex::Real a_time,
int a_nl_iter,
bool a_from_jacobian )
pure virtual

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).

Implemented in SemiImplicitEM, StrangImplicitSpectralEM, and ThetaImplicitEM.

◆ convertFieldBCToLinOpBC()

Array< LinOpBCType, 3 > ImplicitSolver::convertFieldBCToLinOpBC ( const amrex::Array< FieldBoundaryType, 3 > & a_fbc) const
nodiscardprotected

Convert from WarpX FieldBoundaryType to amrex::LinOpBCType.

◆ CreateParticleAttributes()

void ImplicitSolver::CreateParticleAttributes ( ) const

◆ CumulateJ()

void ImplicitSolver::CumulateJ ( )
protected

Add J from particles included in mass matrices at start of each Newton step.

◆ Define()

virtual void ImplicitSolver::Define ( WarpX * a_WarpX)
pure virtual

Read user-provided parameters that control the implicit solver. Allocate internal arrays for intermediate field values needed by the solver.

Implemented in SemiImplicitEM, StrangImplicitSpectralEM, and ThetaImplicitEM.

◆ DoParticleSuborbits()

bool ImplicitSolver::DoParticleSuborbits ( )
inline

◆ GetFieldBoundaryHi()

const Array< FieldBoundaryType, 3 > & ImplicitSolver::GetFieldBoundaryHi ( ) const
nodiscard

◆ GetFieldBoundaryLo()

const Array< FieldBoundaryType, 3 > & ImplicitSolver::GetFieldBoundaryLo ( ) const
nodiscard

◆ GetGeometry()

const Geometry & ImplicitSolver::GetGeometry ( int a_lvl) const
nodiscard

◆ GetLinOpBCHi()

Array< LinOpBCType, 3 > ImplicitSolver::GetLinOpBCHi ( ) const
nodiscard

◆ GetLinOpBCLo()

Array< LinOpBCType, 3 > ImplicitSolver::GetLinOpBCLo ( ) const
nodiscard

◆ GetMassMatricesCoeff()

const amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > * ImplicitSolver::GetMassMatricesCoeff ( ) const
inlinenodiscard

Return pointer to MultiFab array for mass matrix.

◆ GetThetaForPC()

virtual amrex::Real ImplicitSolver::GetThetaForPC ( ) const
nodiscardpure virtual

◆ InitializeMassMatrices()

void ImplicitSolver::InitializeMassMatrices ( )
protected

Initialize the Mass Matrices used for plasma response in nonlinear Newton solver.

◆ IsDefined()

bool ImplicitSolver::IsDefined ( ) const
inlinenodiscard

◆ numAMRLevels()

int ImplicitSolver::numAMRLevels ( ) const
inlinenodiscard

◆ OneStep()

virtual void ImplicitSolver::OneStep ( amrex::Real a_time,
amrex::Real a_dt,
int a_step )
pure virtual

Advance fields and particles by one time step using the specified implicit algorithm.

Implemented in SemiImplicitEM, StrangImplicitSpectralEM, and ThetaImplicitEM.

◆ operator=() [1/2]

ImplicitSolver & ImplicitSolver::operator= ( const ImplicitSolver & )
delete

◆ operator=() [2/2]

ImplicitSolver & ImplicitSolver::operator= ( ImplicitSolver && )
delete

◆ parseNonlinearSolverParams()

void ImplicitSolver::parseNonlinearSolverParams ( const amrex::ParmParse & pp)
protected

parse nonlinear solver parameters (if one is used)

◆ PreRHSOp()

void ImplicitSolver::PreRHSOp ( const amrex::Real a_cur_time,
const int a_nl_iter,
const bool a_from_jacobian )
protected

Perform operations needed before computing the Right Hand Side.

◆ PrintBaseImplicitSolverParameters()

void ImplicitSolver::PrintBaseImplicitSolverParameters ( ) const
protected

◆ PrintParameters()

virtual void ImplicitSolver::PrintParameters ( ) const
pure virtual

◆ SaveE()

void ImplicitSolver::SaveE ( )
protected

Save E at start of each Newton step.

◆ SetMassMatricesForPC()

void ImplicitSolver::SetMassMatricesForPC ( const amrex::Real a_theta_dt)
protected

Scale mass matrices used for PC by c^2*mu0*theta*dt and add 1 to diagonal terms.

◆ SyncMassMatricesPCAndApplyBCs()

void ImplicitSolver::SyncMassMatricesPCAndApplyBCs ( )
protected

Communicate Mass Matrices used for PC and apply boundary conditions.

Member Data Documentation

◆ m_dt

amrex::Real ImplicitSolver::m_dt = 0.0
mutableprotected

Time step.

◆ m_is_defined

bool ImplicitSolver::m_is_defined = false
protected

◆ m_max_particle_iterations

int ImplicitSolver::m_max_particle_iterations = 21
protected

maximum iterations for the iterative method used to obtain a self-consistent update of the particle positions and velocities for given E and B on the grid

◆ m_mmpc_mfarrvec

amrex::Vector<amrex::Array<amrex::MultiFab*,3> > ImplicitSolver::m_mmpc_mfarrvec
protected

Array of multifab pointers to mass matrix preconditioner.

◆ m_ncomp_xx

amrex::IntVect ImplicitSolver::m_ncomp_xx = amrex::IntVect{0}
protected

Direction-dependent component number of mass matrices elements.

◆ m_ncomp_xy

amrex::IntVect ImplicitSolver::m_ncomp_xy = amrex::IntVect{0}
protected

◆ m_ncomp_xz

amrex::IntVect ImplicitSolver::m_ncomp_xz = amrex::IntVect{0}
protected

◆ m_ncomp_yx

amrex::IntVect ImplicitSolver::m_ncomp_yx = amrex::IntVect{0}
protected

◆ m_ncomp_yy

amrex::IntVect ImplicitSolver::m_ncomp_yy = amrex::IntVect{0}
protected

◆ m_ncomp_yz

amrex::IntVect ImplicitSolver::m_ncomp_yz = amrex::IntVect{0}
protected

◆ m_ncomp_zx

amrex::IntVect ImplicitSolver::m_ncomp_zx = amrex::IntVect{0}
protected

◆ m_ncomp_zy

amrex::IntVect ImplicitSolver::m_ncomp_zy = amrex::IntVect{0}
protected

◆ m_ncomp_zz

amrex::IntVect ImplicitSolver::m_ncomp_zz = amrex::IntVect{0}
protected

◆ m_nlsolver

std::unique_ptr<NonlinearSolver<WarpXSolverVec,ImplicitSolver> > ImplicitSolver::m_nlsolver
protected

◆ m_nlsolver_type

NonlinearSolverType ImplicitSolver::m_nlsolver_type
protected

Nonlinear solver type and object.

◆ m_num_amr_levels

int ImplicitSolver::m_num_amr_levels = 1
protected

Number of AMR levels.

◆ m_particle_suborbits

bool ImplicitSolver::m_particle_suborbits = false
protected

whether to use suborbits for particles that fail to converge in m_max_particle_iterations iterations

◆ m_particle_tolerance

amrex::ParticleReal ImplicitSolver::m_particle_tolerance = amrex::ParticleReal(1.0e-10)
protected

tolerance used by the iterative method used to obtain a self-consistent update of the particle positions and velocities for given E and B on the grid

◆ m_print_unconverged_particle_details

bool ImplicitSolver::m_print_unconverged_particle_details = false
protected

whether the details of unconverged particles are printed out during the particle evolve loops

◆ m_skip_particle_picard_init

bool ImplicitSolver::m_skip_particle_picard_init = false
protected

whether to skip the full Picard update of particles on the initial newton step

◆ m_theta

amrex::Real ImplicitSolver::m_theta = 0.5
protected

Time-biasing parameter for fields used on RHS to advance system.

◆ m_use_mass_matrices

bool ImplicitSolver::m_use_mass_matrices = false
protected

Whether to use mass matrices for the implicit solver.

◆ m_use_mass_matrices_jacobian

bool ImplicitSolver::m_use_mass_matrices_jacobian = false
protected

Whether to use mass matrices to compute J during linear stage of JFNK.

◆ m_use_mass_matrices_pc

bool ImplicitSolver::m_use_mass_matrices_pc = false
protected

Whether to use mass matrices in the preconditioner.

◆ m_WarpX

WarpX* ImplicitSolver::m_WarpX
protected

Pointer back to main WarpX class.


The documentation for this class was generated from the following files: