WarpX
Loading...
Searching...
No Matches
CurlCurlMLMGPC< T, Ops > Class Template Reference

Curl-curl Preconditioner. More...

#include <CurlCurlMLMGPC.H>

Inheritance diagram for CurlCurlMLMGPC< T, Ops >:
Preconditioner< T, Ops >

Public Types

using RT = typename T::value_type
 
- Public Types inherited from Preconditioner< T, Ops >
using RT = typename T::value_type
 

Public Member Functions

 CurlCurlMLMGPC ()=default
 Default constructor.
 
 ~CurlCurlMLMGPC () override=default
 Default destructor.
 
 CurlCurlMLMGPC (const CurlCurlMLMGPC &)=delete
 
CurlCurlMLMGPCoperator= (const CurlCurlMLMGPC &)=delete
 
 CurlCurlMLMGPC (CurlCurlMLMGPC &&) noexcept=delete
 
CurlCurlMLMGPCoperator= (CurlCurlMLMGPC &&) noexcept=delete
 
void Define (const T &, Ops *const) override
 Define the preconditioner.
 
void Update (const T &a_U) override
 Update the preconditioner.
 
void Apply (T &, const T &) override
 Apply (solve) the preconditioner given a RHS.
 
void printParameters () const override
 Print parameters.
 
bool IsDefined () const override
 Check if the nonlinear solver has been defined.
 
- Public Member Functions inherited from Preconditioner< T, Ops >
 Preconditioner ()=default
 Default constructor.
 
virtual ~Preconditioner ()=default
 Default destructor.
 
 Preconditioner (const Preconditioner &)=default
 
Preconditioneroperator= (const Preconditioner &)=default
 
 Preconditioner (Preconditioner &&) noexcept=default
 
Preconditioneroperator= (Preconditioner &&) noexcept=default
 
void CurTime (const RT a_time)
 Set the current time.
 
void CurTimeStep (const RT a_dt)
 Set the current time step size.
 

Protected Types

using MFArr = amrex::Array<amrex::MultiFab,3>
 

Protected Member Functions

void readParameters ()
 Read parameters.
 

Protected Attributes

bool m_is_defined = false
 
bool m_verbose = true
 
bool m_bottom_verbose = false
 
bool m_agglomeration = true
 
bool m_consolidation = true
 
bool m_use_gmres = false
 
bool m_use_gmres_pc = true
 
int m_max_iter = 10
 
int m_max_coarsening_level = 30
 
bool m_beta_scalar = true
 
RT m_atol = 1.0e-16
 
RT m_rtol = 1.0e-4
 
Ops * m_ops = nullptr
 
int m_num_amr_levels = 0
 
amrex::Vector< amrex::Geometrym_geom
 
amrex::Vector< amrex::BoxArraym_grids
 
amrex::Vector< amrex::DistributionMappingm_dmap
 
amrex::IntVect m_gv
 
const amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > * m_bcoefs = nullptr
 
amrex::Array< amrex::LinOpBCType, 3 > m_bc_lo
 
amrex::Array< amrex::LinOpBCType, 3 > m_bc_hi
 
std::unique_ptr< amrex::LPInfom_info
 
std::unique_ptr< amrex::MLCurlCurlm_curl_curl
 
std::unique_ptr< amrex::MLMGT< MFArr > > m_solver
 
std::unique_ptr< amrex::GMRESMLMGT< MFArr > > m_gmres_solver
 
- Protected Attributes inherited from Preconditioner< T, Ops >
RT m_time = 0.0
 
RT m_dt = 0.0
 

Detailed Description

template<class T, class Ops>
class CurlCurlMLMGPC< T, Ops >

Curl-curl Preconditioner.

Preconditioner that solves the curl-curl equation for the E-field, given a RHS. Uses AMReX's curl-curl linear operator and multigrid solver.

The equation solves for Eg in: curl ( alpha * curl ( Eg ) ) + beta * Eg = b where

  • alpha is a scalar
  • beta can either be a scalar that is constant in space or a MultiFab
  • Eg is the electric field.
  • b is a specified RHS with the same layout as Eg

This class is templated on a solution-type class T and an operator class Ops.

The Ops class must have the following function:

  • Return number of AMR levels
  • Return the amrex::Geometry object given an AMR level
  • Return hi and lo linear operator boundaries
  • Return the time step factor (theta) for the time integration scheme

The T class must have the following functions:

Member Typedef Documentation

◆ MFArr

template<class T, class Ops>
using CurlCurlMLMGPC< T, Ops >::MFArr = amrex::Array<amrex::MultiFab,3>
protected

◆ RT

template<class T, class Ops>
using CurlCurlMLMGPC< T, Ops >::RT = typename T::value_type

Constructor & Destructor Documentation

◆ CurlCurlMLMGPC() [1/3]

template<class T, class Ops>
CurlCurlMLMGPC< T, Ops >::CurlCurlMLMGPC ( )
default

Default constructor.

◆ ~CurlCurlMLMGPC()

template<class T, class Ops>
CurlCurlMLMGPC< T, Ops >::~CurlCurlMLMGPC ( )
overridedefault

Default destructor.

◆ CurlCurlMLMGPC() [2/3]

template<class T, class Ops>
CurlCurlMLMGPC< T, Ops >::CurlCurlMLMGPC ( const CurlCurlMLMGPC< T, Ops > & )
delete

◆ CurlCurlMLMGPC() [3/3]

template<class T, class Ops>
CurlCurlMLMGPC< T, Ops >::CurlCurlMLMGPC ( CurlCurlMLMGPC< T, Ops > && )
deletenoexcept

Member Function Documentation

◆ Apply()

template<class T, class Ops>
void CurlCurlMLMGPC< T, Ops >::Apply ( T & a_x,
const T & a_b )
overridevirtual

Apply (solve) the preconditioner given a RHS.

Given a right-hand-side b, solve: A x = b where A is the linear operator, in this case, the curl-curl operator: A x = curl (alpha * curl (x) ) + beta * x

Implements Preconditioner< T, Ops >.

◆ Define()

template<class T, class Ops>
void CurlCurlMLMGPC< T, Ops >::Define ( const T & a_U,
Ops * const a_ops )
overridevirtual

Define the preconditioner.

Implements Preconditioner< T, Ops >.

◆ IsDefined()

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::IsDefined ( ) const
inlinenodiscardoverridevirtual

Check if the nonlinear solver has been defined.

Implements Preconditioner< T, Ops >.

◆ operator=() [1/2]

template<class T, class Ops>
CurlCurlMLMGPC & CurlCurlMLMGPC< T, Ops >::operator= ( const CurlCurlMLMGPC< T, Ops > & )
delete

◆ operator=() [2/2]

template<class T, class Ops>
CurlCurlMLMGPC & CurlCurlMLMGPC< T, Ops >::operator= ( CurlCurlMLMGPC< T, Ops > && )
deletenoexcept

◆ printParameters()

template<class T, class Ops>
void CurlCurlMLMGPC< T, Ops >::printParameters ( ) const
overridevirtual

Print parameters.

Reimplemented from Preconditioner< T, Ops >.

◆ readParameters()

template<class T, class Ops>
void CurlCurlMLMGPC< T, Ops >::readParameters ( )
protected

Read parameters.

◆ Update()

template<class T, class Ops>
void CurlCurlMLMGPC< T, Ops >::Update ( const T & a_U)
overridevirtual

Update the preconditioner.

Implements Preconditioner< T, Ops >.

Member Data Documentation

◆ m_agglomeration

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::m_agglomeration = true
protected

◆ m_atol

template<class T, class Ops>
RT CurlCurlMLMGPC< T, Ops >::m_atol = 1.0e-16
protected

◆ m_bc_hi

template<class T, class Ops>
amrex::Array<amrex::LinOpBCType, 3> CurlCurlMLMGPC< T, Ops >::m_bc_hi
protected

◆ m_bc_lo

template<class T, class Ops>
amrex::Array<amrex::LinOpBCType, 3> CurlCurlMLMGPC< T, Ops >::m_bc_lo
protected

◆ m_bcoefs

template<class T, class Ops>
const amrex::Vector<amrex::Array<amrex::MultiFab*,3> >* CurlCurlMLMGPC< T, Ops >::m_bcoefs = nullptr
protected

◆ m_beta_scalar

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::m_beta_scalar = true
protected

◆ m_bottom_verbose

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::m_bottom_verbose = false
protected

◆ m_consolidation

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::m_consolidation = true
protected

◆ m_curl_curl

template<class T, class Ops>
std::unique_ptr<amrex::MLCurlCurl> CurlCurlMLMGPC< T, Ops >::m_curl_curl
protected

◆ m_dmap

template<class T, class Ops>
amrex::Vector<amrex::DistributionMapping> CurlCurlMLMGPC< T, Ops >::m_dmap
protected

◆ m_geom

template<class T, class Ops>
amrex::Vector<amrex::Geometry> CurlCurlMLMGPC< T, Ops >::m_geom
protected

◆ m_gmres_solver

template<class T, class Ops>
std::unique_ptr<amrex::GMRESMLMGT<MFArr> > CurlCurlMLMGPC< T, Ops >::m_gmres_solver
protected

◆ m_grids

template<class T, class Ops>
amrex::Vector<amrex::BoxArray> CurlCurlMLMGPC< T, Ops >::m_grids
protected

◆ m_gv

template<class T, class Ops>
amrex::IntVect CurlCurlMLMGPC< T, Ops >::m_gv
protected

◆ m_info

template<class T, class Ops>
std::unique_ptr<amrex::LPInfo> CurlCurlMLMGPC< T, Ops >::m_info
protected

◆ m_is_defined

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::m_is_defined = false
protected

◆ m_max_coarsening_level

template<class T, class Ops>
int CurlCurlMLMGPC< T, Ops >::m_max_coarsening_level = 30
protected

◆ m_max_iter

template<class T, class Ops>
int CurlCurlMLMGPC< T, Ops >::m_max_iter = 10
protected

◆ m_num_amr_levels

template<class T, class Ops>
int CurlCurlMLMGPC< T, Ops >::m_num_amr_levels = 0
protected

◆ m_ops

template<class T, class Ops>
Ops* CurlCurlMLMGPC< T, Ops >::m_ops = nullptr
protected

◆ m_rtol

template<class T, class Ops>
RT CurlCurlMLMGPC< T, Ops >::m_rtol = 1.0e-4
protected

◆ m_solver

template<class T, class Ops>
std::unique_ptr<amrex::MLMGT<MFArr> > CurlCurlMLMGPC< T, Ops >::m_solver
protected

◆ m_use_gmres

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::m_use_gmres = false
protected

◆ m_use_gmres_pc

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::m_use_gmres_pc = true
protected

◆ m_verbose

template<class T, class Ops>
bool CurlCurlMLMGPC< T, Ops >::m_verbose = true
protected

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