|
WarpX
|
#include <ImplicitSolver.H>
Public Member Functions | |
| ImplicitSolver ()=default | |
| virtual | ~ImplicitSolver ()=default |
| ImplicitSolver (const ImplicitSolver &)=delete | |
| ImplicitSolver & | operator= (const ImplicitSolver &)=delete |
| ImplicitSolver (ImplicitSolver &&)=delete | |
| ImplicitSolver & | operator= (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::Geometry & | GetGeometry (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 | |
| WarpX * | m_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. | |
|
default |
|
virtualdefault |
|
delete |
|
delete |
| void ImplicitSolver::ComputeJfromMassMatrices | ( | const bool | a_J_from_MM_only | ) |
|
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.
|
nodiscardprotected |
Convert from WarpX FieldBoundaryType to amrex::LinOpBCType.
| void ImplicitSolver::CreateParticleAttributes | ( | ) | const |
|
protected |
Add J from particles included in mass matrices at start of each Newton step.
|
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.
|
inline |
|
nodiscard |
|
nodiscard |
|
nodiscard |
|
nodiscard |
|
nodiscard |
|
inlinenodiscard |
Return pointer to MultiFab array for mass matrix.
|
nodiscardpure virtual |
Implemented in SemiImplicitEM, StrangImplicitSpectralEM, and ThetaImplicitEM.
|
protected |
Initialize the Mass Matrices used for plasma response in nonlinear Newton solver.
|
inlinenodiscard |
|
inlinenodiscard |
|
pure virtual |
Advance fields and particles by one time step using the specified implicit algorithm.
Implemented in SemiImplicitEM, StrangImplicitSpectralEM, and ThetaImplicitEM.
|
delete |
|
delete |
|
protected |
parse nonlinear solver parameters (if one is used)
|
protected |
Perform operations needed before computing the Right Hand Side.
|
protected |
|
pure virtual |
Implemented in SemiImplicitEM, StrangImplicitSpectralEM, and ThetaImplicitEM.
|
protected |
Save E at start of each Newton step.
|
protected |
Scale mass matrices used for PC by c^2*mu0*theta*dt and add 1 to diagonal terms.
|
protected |
Communicate Mass Matrices used for PC and apply boundary conditions.
|
mutableprotected |
Time step.
|
protected |
|
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
|
protected |
Array of multifab pointers to mass matrix preconditioner.
|
protected |
Direction-dependent component number of mass matrices elements.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Nonlinear solver type and object.
|
protected |
Number of AMR levels.
|
protected |
whether to use suborbits for particles that fail to converge in m_max_particle_iterations iterations
|
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
|
protected |
whether the details of unconverged particles are printed out during the particle evolve loops
|
protected |
whether to skip the full Picard update of particles on the initial newton step
|
protected |
Time-biasing parameter for fields used on RHS to advance system.
|
protected |
Whether to use mass matrices for the implicit solver.
|
protected |
Whether to use mass matrices to compute J during linear stage of JFNK.
|
protected |
Whether to use mass matrices in the preconditioner.