|
| RefMaxwell () |
| Constructor. More...
|
|
| RefMaxwell (Teuchos::RCP< Hierarchy > HH, Teuchos::RCP< Hierarchy > H22) |
| Constructor with Hierarchies. More...
|
|
| RefMaxwell (const Teuchos::RCP< Matrix > &SM_Matrix, const Teuchos::RCP< Matrix > &Dk_1, const Teuchos::RCP< Matrix > &Dk_2, const Teuchos::RCP< Matrix > &D0, const Teuchos::RCP< Matrix > &M1_beta, const Teuchos::RCP< Matrix > &M1_alpha, const Teuchos::RCP< Matrix > &Mk_one, const Teuchos::RCP< Matrix > &Mk_1_one, const Teuchos::RCP< Matrix > &invMk_1_invBeta, const Teuchos::RCP< Matrix > &invMk_2_invAlpha, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< MultiVector > &Nullspace22, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, Teuchos::ParameterList &List, bool ComputePrec=true) |
|
| RefMaxwell (const Teuchos::RCP< Matrix > &SM_Matrix, const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Ms_Matrix, const Teuchos::RCP< Matrix > &M0inv_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, Teuchos::ParameterList &List, bool ComputePrec=true) |
|
| RefMaxwell (const Teuchos::RCP< Matrix > &SM_Matrix, const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &M0inv_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, Teuchos::ParameterList &List, bool ComputePrec=true) |
|
| RefMaxwell (const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &M0inv_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, Teuchos::ParameterList &List) |
|
| RefMaxwell (const Teuchos::RCP< Matrix > &SM_Matrix, const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, Teuchos::ParameterList &List, bool ComputePrec) |
|
| RefMaxwell (const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, Teuchos::ParameterList &List) |
|
| RefMaxwell (const Teuchos::RCP< Matrix > &SM_Matrix, Teuchos::ParameterList &List, bool ComputePrec=true) |
|
virtual | ~RefMaxwell () |
| Destructor. More...
|
|
const Teuchos::RCP< const Map > | getDomainMap () const |
| Returns the Xpetra::Map object associated with the domain of this operator. More...
|
|
const Teuchos::RCP< const Map > | getRangeMap () const |
| Returns the Xpetra::Map object associated with the range of this operator. More...
|
|
const Teuchos::RCP< Matrix > & | getJacobian () const |
| Returns Jacobian matrix SM. More...
|
|
void | setParameters (Teuchos::ParameterList &list) |
| Set parameters. More...
|
|
void | compute (bool reuse=false) |
| Setup the preconditioner. More...
|
|
void | resetMatrix (Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true) |
| Reset system matrix. More...
|
|
void | apply (const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const |
|
bool | hasTransposeApply () const |
| Indicates whether this operator supports applying the adjoint operator. More...
|
|
void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const |
|
void | residual (const MultiVector &X, const MultiVector &B, MultiVector &R) const |
| Compute a residual R = B - (*this) * X. More...
|
|
RCP< MultiVector > | buildNullspace (const int spaceNumber, const Kokkos::View< bool *, typename Node::device_type > &bcs, const bool applyBCs) |
| Builds a nullspace. More...
|
|
Teuchos::RCP< Matrix > | buildProjection (const int spaceNumber, const RCP< MultiVector > &EdgeNullspace) const |
| Builds a projection from a vector values space into a vector valued nodal space. More...
|
|
void | buildNodalProlongator (const Teuchos::RCP< Matrix > &A_nodal, Teuchos::RCP< Matrix > &P_nodal, Teuchos::RCP< MultiVector > &Nullspace_nodal, Teuchos::RCP< RealValuedMultiVector > &Coords_nodal) const |
|
RCP< Matrix > | buildVectorNodalProlongator (const Teuchos::RCP< Matrix > &P_nodal) const |
|
void | buildProlongator (const int spaceNumber, const Teuchos::RCP< Matrix > &A_nodal_Matrix, const RCP< MultiVector > &EdgeNullspace, Teuchos::RCP< Matrix > &edgeProlongator, Teuchos::RCP< MultiVector > &coarseEdgeNullspace, Teuchos::RCP< RealValuedMultiVector > &coarseNodalCoords) const |
|
RCP< Matrix > | buildAddon (const int spaceNumber) |
|
VerbLevel | GetVerbLevel () const |
| Get the verbosity level. More...
|
|
void | SetVerbLevel (const VerbLevel verbLevel) |
| Set the verbosity level of this object. More...
|
|
int | GetProcRankVerbose () const |
| Get proc rank used for printing. Do not use this information for any other purpose. More...
|
|
int | SetProcRankVerbose (int procRank) const |
| Set proc rank used for printing. More...
|
|
bool | IsPrint (MsgType type, int thisProcRankOnly=-1) const |
| Find out whether we need to print out information for a specific message type. More...
|
|
Teuchos::FancyOStream & | GetOStream (MsgType type, int thisProcRankOnly=0) const |
| Get an output stream for outputting the input message type. More...
|
|
Teuchos::FancyOStream & | GetBlackHole () const |
|
| VerboseObject () |
|
virtual | ~VerboseObject () |
| Destructor. More...
|
|
| VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) |
|
virtual const VerboseObject & | setVerbLevel (const EVerbosityLevel verbLevel) const |
|
virtual const VerboseObject & | setOverridingVerbLevel (const EVerbosityLevel verbLevel) const |
|
virtual EVerbosityLevel | getVerbLevel () const |
|
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
RCP< const ParameterList > | getValidVerboseObjectSublist () |
|
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
void | setupVerboseObjectSublist (ParameterList *paramList) |
|
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
void | readVerboseObjectSublist (ParameterList *paramList, RCP< FancyOStream > *oStream, EVerbosityLevel *verbLevel) |
|
void | readVerboseObjectSublist (ParameterList *paramList, VerboseObject< ObjectType > *verboseObject) |
|
virtual | ~VerboseObjectBase () |
|
| VerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) |
|
virtual const VerboseObjectBase & | setOStream (const RCP< FancyOStream > &oStream) const |
|
virtual const VerboseObjectBase & | setOverridingOStream (const RCP< FancyOStream > &oStream) const |
|
virtual VerboseObjectBase & | setLinePrefix (const std::string &linePrefix) |
|
virtual RCP< FancyOStream > | getOStream () const |
|
virtual RCP< FancyOStream > | getOverridingOStream () const |
|
virtual std::string | getLinePrefix () const |
|
virtual OSTab | getOSTab (const int tabs=1, const std::string &linePrefix="") const |
|
virtual void | removeEmptyProcessesInPlace (const RCP< const Map > &) |
|
virtual std::string | description () const |
|
void | describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
|
virtual | ~Describable () |
|
| LabeledObject () |
|
virtual | ~LabeledObject () |
|
virtual void | setObjectLabel (const std::string &objectLabel) |
|
virtual std::string | getObjectLabel () const |
|
DescribableStreamManipulatorState | describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default) |
|
std::ostream & | operator<< (std::ostream &os, const DescribableStreamManipulatorState &d) |
|
|
Teuchos::RCP
< Teuchos::ParameterList > | getValidParamterList () |
|
void | initialize (const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Ms_Matrix, const Teuchos::RCP< Matrix > &M0inv_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, Teuchos::ParameterList &List) |
|
void | initialize (const int k, const Teuchos::RCP< Matrix > &Dk_1, const Teuchos::RCP< Matrix > &Dk_2, const Teuchos::RCP< Matrix > &D0, const Teuchos::RCP< Matrix > &M1_beta, const Teuchos::RCP< Matrix > &M1_alpha, const Teuchos::RCP< Matrix > &Mk_one, const Teuchos::RCP< Matrix > &Mk_1_one, const Teuchos::RCP< Matrix > &invMk_1_invBeta, const Teuchos::RCP< Matrix > &invMk_2_invAlpha, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< MultiVector > &Nullspace22, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, Teuchos::ParameterList &List) |
|
void | determineSubHierarchyCommSizes (bool &doRebalancing, int &rebalanceStriding, int &numProcsCoarseA11, int &numProcsA22) |
| Determine how large the sub-communicators for the two hierarchies should be. More...
|
|
void | buildCoarse11Matrix () |
| Compute coarseA11 = P11^{T}*SM*P11 + addon efficiently. More...
|
|
void | rebalanceCoarse11Matrix (const int rebalanceStriding, const int numProcsCoarseA11) |
| rebalance the coarse A11 matrix, as well as P11, CoordsCoarse11 and Addon11 More...
|
|
void | build22Matrix (const bool reuse, const bool doRebalancing, const int rebalanceStriding, const int numProcsA22) |
| Setup A22 = D0^T SM D0 and rebalance it, as well as D0 and Coords_. More...
|
|
void | setupSubSolve (Teuchos::RCP< Hierarchy > &hierarchy, Teuchos::RCP< Operator > &thyraPrecOp, const Teuchos::RCP< Matrix > &A, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList ¶ms, std::string &label, const bool reuse, const bool isSingular=false) |
| Setup a subsolve. More...
|
|
void | setFineLevelSmoother11 () |
| Set the fine level smoother. More...
|
|
void | applyInverseAdditive (const MultiVector &RHS, MultiVector &X) const |
| apply additive algorithm for 2x2 solve More...
|
|
void | solveH (const MultiVector &RHS, MultiVector &X) const |
| apply solve to 1-1 block only More...
|
|
void | solve22 (const MultiVector &RHS, MultiVector &X) const |
| apply solve to 2-2 block only More...
|
|
void | allocateMemory (int numVectors) const |
| allocate multivectors for solve More...
|
|
void | dump (const RCP< Matrix > &A, std::string name) const |
| dump out matrix More...
|
|
void | dump (const RCP< MultiVector > &X, std::string name) const |
| dump out multivector More...
|
|
void | dumpCoords (const RCP< RealValuedMultiVector > &X, std::string name) const |
| dump out real-valued multivector More...
|
|
void | dump (const Teuchos::ArrayRCP< bool > &v, std::string name) const |
| dump out boolean ArrayView More...
|
|
void | dump (const Kokkos::View< bool *, typename Node::device_type > &v, std::string name) const |
| dump out boolean Kokkos::View More...
|
|
Teuchos::RCP
< Teuchos::TimeMonitor > | getTimer (std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const |
| get a (synced) timer More...
|
|
|
Teuchos::RCP< Hierarchy > | HierarchyCoarse11_ |
| Two hierarchies: one for the coarse (1,1)-block, another for the (2,2)-block. More...
|
|
Teuchos::RCP< Hierarchy > | Hierarchy22_ |
|
Teuchos::RCP< SmootherBase > | PreSmoother11_ |
|
Teuchos::RCP< SmootherBase > | PostSmoother11_ |
|
Teuchos::RCP< SmootherPrototype > | PreSmootherData11_ |
|
Teuchos::RCP< SmootherPrototype > | PostSmootherData11_ |
|
RCP< Operator > | thyraPrecOpH_ |
|
RCP< Operator > | thyraPrecOp22_ |
|
int | spaceNumber_ |
| The number of the space in the deRham complex. More...
|
|
std::string | solverName_ |
| The name of the solver. More...
|
|
size_t | dim_ |
| The spatial dimension. More...
|
|
Teuchos::RCP< Matrix > | SM_Matrix_ |
| The system that is getting preconditioned. More...
|
|
Teuchos::RCP< Matrix > | Dk_1_ |
| D_{k-1} matrix and its transpose. More...
|
|
Teuchos::RCP< Matrix > | Dk_1_T_ |
|
Teuchos::RCP< Matrix > | Dk_2_ |
| D_{k-2} matrix. More...
|
|
Teuchos::RCP< Matrix > | D0_ |
| D_0 matrix. More...
|
|
Teuchos::RCP< Matrix > | invMk_1_invBeta_ |
| inverse of mass matrices on (k-1)-th and (k-2)-th space with weights 1/beta and 1/alpha respectively More...
|
|
Teuchos::RCP< Matrix > | invMk_2_invAlpha_ |
|
Teuchos::RCP< Matrix > | Mk_one_ |
| mass matrices with unit weight on k-th and (k-1)-th spaces More...
|
|
Teuchos::RCP< Matrix > | Mk_1_one_ |
|
Teuchos::RCP< Matrix > | M1_beta_ |
| mass matrices on first space with weights beta and alpha respectively More...
|
|
Teuchos::RCP< Matrix > | M1_alpha_ |
|
Teuchos::RCP< Matrix > | P11_ |
| special prolongator for 11 block and its transpose More...
|
|
Teuchos::RCP< Matrix > | R11_ |
|
Teuchos::RCP< Matrix > | P22_ |
| special prolongator for 22 block and its transpose More...
|
|
Teuchos::RCP< Matrix > | R22_ |
|
Teuchos::RCP< Matrix > | coarseA11_ |
| coarse 11, 22 and coarse 22 blocks More...
|
|
Teuchos::RCP< Matrix > | A22_ |
|
Teuchos::RCP< Matrix > | coarseA22_ |
|
Teuchos::RCP< Matrix > | Addon11_ |
| the addon for the 11 block More...
|
|
Teuchos::RCP< Matrix > | Addon22_ |
| the addon for the 22 block More...
|
|
Teuchos::RCP< const Map > | DorigDomainMap_ |
|
Teuchos::RCP< const Import > | DorigImporter_ |
|
Kokkos::View< bool *, typename
Node::device_type::memory_space > | BCrows11_ |
| Vectors for BCs. More...
|
|
Kokkos::View< bool *, typename
Node::device_type::memory_space > | BCcols22_ |
|
Kokkos::View< bool *, typename
Node::device_type::memory_space > | BCdomain22_ |
|
int | globalNumberBoundaryUnknowns11_ |
|
int | globalNumberBoundaryUnknowns22_ |
|
Teuchos::RCP< MultiVector > | Nullspace11_ |
| Nullspace for (1.1) block. More...
|
|
Teuchos::RCP< MultiVector > | Nullspace22_ |
|
Teuchos::RCP
< RealValuedMultiVector > | NodalCoords_ |
| Coordinates. More...
|
|
Teuchos::RCP
< RealValuedMultiVector > | CoordsCoarse11_ |
|
Teuchos::RCP
< RealValuedMultiVector > | Coords22_ |
|
Teuchos::RCP< MultiVector > | NullspaceCoarse11_ |
| Nullspace for coarse (1,1) problem. More...
|
|
Teuchos::RCP< MultiVector > | CoarseNullspace22_ |
| Nullspace for coarse (2,2) problem. More...
|
|
Teuchos::RCP< const Import > | ImporterCoarse11_ |
| Importer to coarse (1,1) hierarchy. More...
|
|
Teuchos::RCP< const Import > | Importer22_ |
|
bool | Dk_1_T_R11_colMapsMatch_ |
|
bool | onlyBoundary11_ |
|
bool | onlyBoundary22_ |
|
Teuchos::ParameterList | parameterList_ |
| Parameter lists. More...
|
|
Teuchos::ParameterList | precList11_ |
|
Teuchos::ParameterList | precList22_ |
|
Teuchos::RCP
< Teuchos::ParameterList > | coarseA11_AP_reuse_data_ |
|
Teuchos::RCP
< Teuchos::ParameterList > | coarseA11_RAP_reuse_data_ |
|
Teuchos::RCP
< Teuchos::ParameterList > | A22_AP_reuse_data_ |
|
Teuchos::RCP
< Teuchos::ParameterList > | A22_RAP_reuse_data_ |
|
bool | disable_addon_ |
| Some options. More...
|
|
bool | disable_addon_22_ |
|
bool | dump_matrices_ |
|
bool | useKokkos_ |
|
bool | use_as_preconditioner_ |
|
bool | implicitTranspose_ |
|
bool | fuseProlongationAndUpdate_ |
|
bool | syncTimers_ |
|
bool | enable_reuse_ |
|
bool | skipFirst11Level_ |
|
bool | skipFirst22Level_ |
|
bool | applyBCsToAnodal_ |
|
bool | applyBCsToCoarse11_ |
|
bool | applyBCsTo22_ |
|
int | numItersCoarse11_ |
|
int | numIters22_ |
|
std::string | mode_ |
|
Teuchos::RCP< MultiVector > | P11res_ |
| Temporary memory. More...
|
|
Teuchos::RCP< MultiVector > | P11x_ |
|
Teuchos::RCP< MultiVector > | P11resSubComm_ |
|
Teuchos::RCP< MultiVector > | P11xSubComm_ |
|
Teuchos::RCP< MultiVector > | DresIntermediate_ |
|
Teuchos::RCP< MultiVector > | Dres_ |
|
Teuchos::RCP< MultiVector > | DxIntermediate_ |
|
Teuchos::RCP< MultiVector > | Dx_ |
|
Teuchos::RCP< MultiVector > | DresSubComm_ |
|
Teuchos::RCP< MultiVector > | DxSubComm_ |
|
Teuchos::RCP< MultiVector > | residual_ |
|
Teuchos::RCP< MultiVector > | P11resTmp_ |
|
Teuchos::RCP< MultiVector > | DresTmp_ |
|
Teuchos::RCP< MultiVector > | DTR11Tmp_ |
|
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
class MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form.
Let
M_k(gamma)
be the mass matrix with a weight coefficient gamma on the k-th space in the deRham complex
k=0 —> k=1 -—> k=2 -—> k=3
H(grad) —> H(curl) -—> H(div) -—> L^2
D_0 D_1 D_2
grad curl div
and let
D_k
be the discretized derivative from the k-th to the (k+1)-th space.
For example, in 3D, D_0 = grad, D_1 = curl and D_2 = div.
We want to solve the system
S e_k = b
where
S = M_k(alpha) + D_k^T * M_{k+1}(beta) * D_k.
In 3D and for k=1, S corresponds to a discretization of
alpha * identity + beta * curl curl.
In 3D and for k=2, S corresponds to
alpha * identity + beta * grad div.
We precondition S by constructing a block diagonal preconditioner for the equivalent 2x2 block system
e_k = a_k + D_{k-1}*p_{k-1}
( A11 M_k(alpha) * D_{k-1} ) ( a_k ) ( b ) ( } ( ) = ( ) ( D_{k-1}^T * M_k(alpha) A22 ) ( p_{k-1} ) ( D_{k-1}^T * b )
where
A11 = S + addon11 A22 = D_{k-1} ^T * S * D_{k-1} + addon22 addon11 = M_k(1) * D_{k-1} * M_{k-1}(1/beta)^-1 * D_{k-1}^T * M_k(1) addon22 = M_{k-1}(1) * D_{k-2} * M_{k-2}(1/alpha)^-1 * D_{k-2}^T * M_{k-1}(1)
We note that due to D_k * D_{k-1} = 0 we have that
A22 = D_{k-1} ^T * M_k(alpha) * D_{k-1} + addon22
The addon terms mimic the term that would be needed to augment to a Hodge Laplacian. In practice it seems that the addon terms are rarely strictly required. For k=1 (Maxwell) the addon22 term (in H(grad)) is always zero.
We cannot directly apply AMG to A11 and A22 in case they are not expressed in terms of nodal bases. Let Pi_k be the projection from the first order vectorial nodal finite element space into the k-th space. We will apply AMG to
Pi_k^T * A11 * Pi_k Pi_{k-1}^T * A22 * Pi_{k-1}
respectively. For k=1 A22 is already given in a nodal discretization and this step is omitted.
It can be beneficial to directly coarsen the projected diagonal blocks and skip any smoothing on them.
To obtain prolongators for the vectorial nodal problems, we construct the auxiliary operators
A11_nodal = D_{0}^T * M_1(beta) * D_{0} A22_nodal = D_{0}^T * M_1(alpha) * D_{0}
then construct typical nodal scalar HGrad prolongators.
We then replicate them into prolongators for vectorial H(grad) problem. Finally, we multiply from the left with the projection Pi_k onto the k-th FE space.
When k=1 this only needs to be done for the A11 block, as the A22 block is already given wrt a scalar nodal discretization.
The following input matrices need to be supplied by the user:
If addon11 is used we need:
- M_{k-1}(1/beta)^-1 - typically, this is the inverse of the lumped M_{k-1}(1/beta), not the actual inverse
- M_k(1)
If addon22 is used we need:
- M_{k-2}(1/alpha)^-1 - typically, this is the inverse of the lumped M_{k-2}(1/alpha), not the actual inverse
- M_{k-1}(1)
To construct the special prolongators we need
- D_{0}
- M_1(beta)
- M_1(alpha) If these mass matrices are not given, but M_1(1) is, then that matrix can be used instead. Alternatively, A11_nodal and A22_nodal can be passed directly.
We are using the following variable names:
| variable | matrix | note |-------------------------------—|--------------------—|-------------------— | SM
| S | | Dk_1
| D_{k-1} | same as D0 for k=1 | Dk_2
| D_{k-2} | | D0
| D_0 | same as Dk_1 for k=1 | Mk_one
| M_k(1) | | Mk_1_one
| M_{k-1}(1) | | M1_beta
| M_1(beta) | | M1_alpha
| M_1(alpha) | | invMk_1_invBeta
| M_{k-1}(1/beta)^{-1} | | invMk_2_invAlpha
| M_{k-2}(1/alpha)^{-1} |
For backwards compatibility the interfaces also allow
variable | matrix | note |
Ms | M_1(beta) | alias for M1_beta |
M1 | M_1(1) | alias for Mk_one when k=1 |
M0inv | M_0(1/beta) | alias for Mk_1_invBeta when k=1 |
Reference: P. Bochev, J. Hu, C. Siefert, and R. Tuminaro. "An algebraic multigrid approach based on
a compatible gauge reformulation of Maxwell's equations." SIAM Journal on Scientific Computing, 31(1), 557-583.
Parameter list options:
refmaxwell: mode
- a string
specifying the order of solve of the block system. Allowed values are: "additive" (default), "121", "212", "1", "2"
refmaxwell: disable addon
- bool
specifying whether the addon should be built for stabilization. Default: "true"
refmaxwell: use as preconditioner
- bool
specifying whether RefMaxwell is used as a preconditioner or as a solver.
refmaxwell: dump matrices
- bool
specifying whether the matrices should be dumped. Default: "false"
refmaxwell: 11list
and refmaxwell: 22list
- parameter list for the multigrid hierarchies on 11 and 22 blocks
refmaxwell: subsolves on subcommunicators
- bool
redistribute the two subsolves to disjoint sub-communicators (so that the additive solve can occur in parallel) Default: "false"
Definition at line 220 of file MueLu_RefMaxwell_decl.hpp.