MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form. More...

#include <MueLu_RefMaxwell_decl.hpp>

Inheritance diagram for MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
MueLu::VerboseObject

Public Types

using magnitudeType = typename Teuchos::ScalarTraits< Scalar >::magnitudeType
 
using coordinateType = typename Teuchos::ScalarTraits< Scalar >::coordinateType
 
using RealValuedMultiVector = typename Xpetra::MultiVector< coordinateType, LO, GO, NO >
 

Public Member Functions

 RefMaxwell ()
 Constructor.
 
 RefMaxwell (Teuchos::RCP< Hierarchy > HH, Teuchos::RCP< Hierarchy > H22)
 Constructor with Hierarchies.
 
 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 > &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, const Teuchos::RCP< MultiVector > &Material_beta, const Teuchos::RCP< MultiVector > &Material_alpha, 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 > &Ms_Matrix, const Teuchos::RCP< Matrix > &M0inv_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace11, const Teuchos::RCP< RealValuedMultiVector > &NodalCoords, const Teuchos::RCP< MultiVector > &Material, 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 ()=default
 Destructor.
 
const Teuchos::RCP< const Map > getDomainMap () const
 Returns the Xpetra::Map object associated with the domain of this operator.
 
const Teuchos::RCP< const Map > getRangeMap () const
 Returns the Xpetra::Map object associated with the range of this operator.
 
const Teuchos::RCP< Matrix > & getJacobian () const
 Returns Jacobian matrix SM.
 
void setParameters (Teuchos::ParameterList &list)
 Set parameters.
 
void compute (bool reuse=false)
 Setup the preconditioner.
 
void resetMatrix (Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
 Reset system matrix.
 
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.
 
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.
 
RCP< MultiVector > buildNullspace (const int spaceNumber, const Kokkos::View< bool *, typename Node::device_type > &bcs, const bool applyBCs)
 Builds a nullspace.
 
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.
 
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)
 
- Public Member Functions inherited from MueLu::VerboseObject
 VerboseObject ()
 
virtual ~VerboseObject ()
 Destructor.
 
VerbLevel GetVerbLevel () const
 Get the verbosity level.
 
void SetVerbLevel (const VerbLevel verbLevel)
 Set the verbosity level of this object.
 
int GetProcRankVerbose () const
 Get proc rank used for printing. Do not use this information for any other purpose.
 
int SetProcRankVerbose (int procRank) const
 Set proc rank used for printing.
 
bool IsPrint (MsgType type, int thisProcRankOnly=-1) const
 Find out whether we need to print out information for a specific message type.
 
Teuchos::FancyOStream & GetOStream (MsgType type, int thisProcRankOnly=0) const
 Get an output stream for outputting the input message type.
 
Teuchos::FancyOStream & GetBlackHole () const
 

Private Member Functions

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, const Teuchos::RCP< MultiVector > &Material, 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, const Teuchos::RCP< MultiVector > &Material_beta, const Teuchos::RCP< MultiVector > &Material_alpha, 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.
 
void buildCoarse11Matrix ()
 Compute coarseA11 = P11^{T}*SM*P11 + addon efficiently.
 
void rebalanceCoarse11Matrix (const int rebalanceStriding, const int numProcsCoarseA11)
 rebalance the coarse A11 matrix, as well as P11, CoordsCoarse11 and Addon11
 
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_.
 
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, const Teuchos::RCP< MultiVector > &Material, Teuchos::ParameterList &params, std::string &label, const bool reuse, const bool isSingular=false)
 Setup a subsolve.
 
void setFineLevelSmoother11 ()
 Set the fine level smoother.
 
void applyInverseAdditive (const MultiVector &RHS, MultiVector &X) const
 apply additive algorithm for 2x2 solve
 
void solveH (const MultiVector &RHS, MultiVector &X) const
 apply solve to 1-1 block only
 
void solve22 (const MultiVector &RHS, MultiVector &X) const
 apply solve to 2-2 block only
 
void allocateMemory (int numVectors) const
 allocate multivectors for solve
 
void dump (const RCP< Matrix > &A, std::string name) const
 dump out matrix
 
void dump (const RCP< MultiVector > &X, std::string name) const
 dump out multivector
 
void dumpCoords (const RCP< RealValuedMultiVector > &X, std::string name) const
 dump out real-valued multivector
 
void dump (const Teuchos::ArrayRCP< bool > &v, std::string name) const
 dump out boolean ArrayView
 
void dump (const Kokkos::View< bool *, typename Node::device_type > &v, std::string name) const
 dump out boolean Kokkos::View
 
Teuchos::RCP< Teuchos::TimeMonitor > getTimer (std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
 get a (synced) timer
 

Private Attributes

Teuchos::RCP< HierarchyHierarchyCoarse11_
 Two hierarchies: one for the coarse (1,1)-block, another for the (2,2)-block.
 
Teuchos::RCP< HierarchyHierarchy22_
 
Teuchos::RCP< SmootherBasePreSmoother11_
 
Teuchos::RCP< SmootherBasePostSmoother11_
 
Teuchos::RCP< SmootherPrototypePreSmootherData11_
 
Teuchos::RCP< SmootherPrototypePostSmootherData11_
 
RCP< Operator > thyraPrecOpH_
 
RCP< Operator > thyraPrecOp22_
 
int spaceNumber_
 The number of the space in the deRham complex.
 
std::string solverName_
 The name of the solver.
 
size_t dim_
 The spatial dimension.
 
Teuchos::RCP< Matrix > SM_Matrix_
 The system that is getting preconditioned.
 
Teuchos::RCP< Matrix > Dk_1_
 D_{k-1} matrix and its transpose.
 
Teuchos::RCP< Matrix > Dk_1_T_
 
Teuchos::RCP< Matrix > Dk_2_
 D_{k-2} matrix.
 
Teuchos::RCP< Matrix > D0_
 D_0 matrix.
 
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
 
Teuchos::RCP< Matrix > invMk_2_invAlpha_
 
Teuchos::RCP< Matrix > Mk_one_
 mass matrices with unit weight on k-th and (k-1)-th spaces
 
Teuchos::RCP< Matrix > Mk_1_one_
 
Teuchos::RCP< Matrix > M1_beta_
 mass matrices on first space with weights beta and alpha respectively
 
Teuchos::RCP< Matrix > M1_alpha_
 
Teuchos::RCP< MultiVector > Material_beta_
 material for first space
 
Teuchos::RCP< MultiVector > Material_alpha_
 
Teuchos::RCP< Matrix > P11_
 special prolongator for 11 block and its transpose
 
Teuchos::RCP< Matrix > R11_
 
Teuchos::RCP< Matrix > P22_
 special prolongator for 22 block and its transpose
 
Teuchos::RCP< Matrix > R22_
 
Teuchos::RCP< Matrix > coarseA11_
 coarse 11, 22 and coarse 22 blocks
 
Teuchos::RCP< Matrix > A22_
 
Teuchos::RCP< Matrix > coarseA22_
 
Teuchos::RCP< Matrix > Addon11_
 the addon for the 11 block
 
Teuchos::RCP< Matrix > Addon22_
 the addon for the 22 block
 
Teuchos::RCP< const Map > DorigDomainMap_
 
Teuchos::RCP< const Import > DorigImporter_
 
Kokkos::View< bool *, typename Node::device_type::memory_space > BCrows11_
 Vectors for BCs.
 
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.
 
Teuchos::RCP< MultiVector > Nullspace22_
 
Teuchos::RCP< RealValuedMultiVectorNodalCoords_
 Coordinates.
 
Teuchos::RCP< RealValuedMultiVectorCoordsCoarse11_
 
Teuchos::RCP< RealValuedMultiVectorCoords22_
 
Teuchos::RCP< MultiVector > NullspaceCoarse11_
 Nullspace for coarse (1,1) problem.
 
Teuchos::RCP< MultiVector > CoarseNullspace22_
 Nullspace for coarse (2,2) problem.
 
Teuchos::RCP< const Import > ImporterCoarse11_
 Importer to coarse (1,1) hierarchy.
 
Teuchos::RCP< const Import > Importer22_
 
bool Dk_1_T_R11_colMapsMatch_
 
bool asyncTransfers_
 
bool onlyBoundary11_
 
bool onlyBoundary22_
 
Teuchos::ParameterList parameterList_
 Parameter lists.
 
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.
 
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.
 
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_
 
Teuchos::RCP< MultiVector > P11x_colmap_
 
Teuchos::RCP< MultiVector > Dx_colmap_
 

Additional Inherited Members

- Static Public Member Functions inherited from MueLu::VerboseObject
static void SetDefaultVerbLevel (const VerbLevel defaultVerbLevel)
 Set the default (global) verbosity level.
 
static VerbLevel GetDefaultVerbLevel ()
 Get the default (global) verbosity level.
 
static void SetMueLuOStream (const Teuchos::RCP< Teuchos::FancyOStream > &mueluOStream)
 
static void SetMueLuOFileStream (const std::string &filename)
 
static Teuchos::RCP< Teuchos::FancyOStream > GetMueLuOStream ()
 

Detailed Description

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:

  • S
  • D_{k-1}

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.

Member Typedef Documentation

◆ magnitudeType

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
using MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::magnitudeType = typename Teuchos::ScalarTraits<Scalar>::magnitudeType

Definition at line 225 of file MueLu_RefMaxwell_decl.hpp.

◆ coordinateType

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
using MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::coordinateType = typename Teuchos::ScalarTraits<Scalar>::coordinateType

Definition at line 226 of file MueLu_RefMaxwell_decl.hpp.

◆ RealValuedMultiVector

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
using MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::RealValuedMultiVector = typename Xpetra::MultiVector<coordinateType, LO, GO, NO>

Definition at line 227 of file MueLu_RefMaxwell_decl.hpp.

Constructor & Destructor Documentation

◆ RefMaxwell() [1/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::RefMaxwell ( )
inline

Constructor.

Definition at line 230 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [2/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::RefMaxwell ( Teuchos::RCP< Hierarchy HH,
Teuchos::RCP< Hierarchy H22 
)
inline

Constructor with Hierarchies.

Definition at line 238 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [3/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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 
)
inline

Definition at line 245 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [4/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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,
const Teuchos::RCP< MultiVector > &  Material_beta,
const Teuchos::RCP< MultiVector > &  Material_alpha,
Teuchos::ParameterList &  List,
bool  ComputePrec = true 
)
inline

Definition at line 272 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [5/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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 
)
inline

Constructor with Jacobian (with add on)

Parameters
[in]SM_MatrixJacobian
[in]D0_MatrixDiscrete Gradient
[in]Ms_MatrixEdge mass matrix for the nodal aggregates
[in]M0inv_MatrixInverse of lumped nodal mass matrix (add on only)
[in]M1_MatrixEdge mass matrix for the add on
[in]Nullspace11Null space (needed for periodic)
[in]NodalCoordsNodal coordinates
[in]ListParameter list
[in]ComputePrecIf true, compute the preconditioner immediately

Definition at line 313 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [6/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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,
const Teuchos::RCP< MultiVector > &  Material,
Teuchos::ParameterList &  List,
bool  ComputePrec = true 
)
inline

Definition at line 326 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [7/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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 
)
inline

Constructor with Jacobian (with add on)

Parameters
[in]SM_MatrixJacobian
[in]D0_MatrixDiscrete Gradient
[in]M0inv_MatrixInverse of lumped nodal mass matrix (add on only)
[in]M1_MatrixEdge mass matrix for the
[in]NullspaceNull space (needed for periodic)
[in]CoordsNodal coordinates
[in]ListParameter list
[in]ComputePrecIf true, compute the preconditioner immediately

Definition at line 351 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [8/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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 
)
inline

Constructor without Jacobian (with add on)

Parameters
[in]D0_MatrixDiscrete Gradient
[in]M0inv_MatrixInverse of lumped nodal mass matrix (add on only)
[in]M1_MatrixEdge mass matrix for the
[in]NullspaceNull space (needed for periodic)
[in]CoordsNodal coordinates
[in]ListParameter list

Definition at line 372 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [9/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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 
)
inline

Constructor with Jacobian (no add on)

Parameters
[in]SM_MatrixJacobian
[in]D0_MatrixDiscrete Gradient
[in]M1_MatrixEdge mass matrix for the
[in]NullspaceNull space (needed for periodic)
[in]CoordsNodal coordinates
[in]ListParameter list
[in]ComputePrecIf true, compute the preconditioner immediately

Definition at line 392 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [10/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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 
)
inline

Constructor without Jacobian (no add on)

Parameters
[in]D0_MatrixDiscrete Gradient
[in]M1_MatrixEdge mass matrix for the
[in]NullspaceNull space (needed for periodic)
[in]CoordsNodal coordinates
[in]ListParameter list

Definition at line 411 of file MueLu_RefMaxwell_decl.hpp.

◆ RefMaxwell() [11/11]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::RefMaxwell ( const Teuchos::RCP< Matrix > &  SM_Matrix,
Teuchos::ParameterList &  List,
bool  ComputePrec = true 
)

Constructor with parameter list

Parameters
[in]SM_MatrixJacobian
[in]ListParameter list
[in]ComputePrecIf true, compute the preconditioner immediately

Definition at line 2520 of file MueLu_RefMaxwell_def.hpp.

◆ ~RefMaxwell()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
virtual MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::~RefMaxwell ( )
virtualdefault

Destructor.

Member Function Documentation

◆ getDomainMap()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDomainMap ( ) const

Returns the Xpetra::Map object associated with the domain of this operator.

Definition at line 86 of file MueLu_RefMaxwell_def.hpp.

◆ getRangeMap()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRangeMap ( ) const

Returns the Xpetra::Map object associated with the range of this operator.

Definition at line 91 of file MueLu_RefMaxwell_def.hpp.

◆ getJacobian()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
const Teuchos::RCP< Matrix > & MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getJacobian ( ) const
inline

Returns Jacobian matrix SM.

Definition at line 440 of file MueLu_RefMaxwell_decl.hpp.

◆ setParameters()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setParameters ( Teuchos::ParameterList &  list)

Set parameters.

Definition at line 199 of file MueLu_RefMaxwell_def.hpp.

◆ compute()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::compute ( bool  reuse = false)

Setup the preconditioner.

Definition at line 293 of file MueLu_RefMaxwell_def.hpp.

◆ resetMatrix()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::resetMatrix ( Teuchos::RCP< Matrix >  SM_Matrix_new,
bool  ComputePrec = true 
)

Reset system matrix.

Definition at line 2106 of file MueLu_RefMaxwell_def.hpp.

◆ apply()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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

Returns in Y the result of a Xpetra::Operator applied to a Xpetra::MultiVector X.

Parameters
[in]X- MultiVector of dimension NumVectors to multiply with matrix.
[out]Y- MultiVector of dimension NumVectors containing result.

Definition at line 2453 of file MueLu_RefMaxwell_def.hpp.

◆ hasTransposeApply()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hasTransposeApply ( ) const

Indicates whether this operator supports applying the adjoint operator.

Definition at line 2515 of file MueLu_RefMaxwell_def.hpp.

◆ describe()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::VERB_HIGH 
) const

Definition at line 2836 of file MueLu_RefMaxwell_def.hpp.

◆ residual()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::residual ( const MultiVector &  X,
const MultiVector &  B,
MultiVector &  R 
) const
inline

Compute a residual R = B - (*this) * X.

Definition at line 467 of file MueLu_RefMaxwell_decl.hpp.

◆ getValidParamterList()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Teuchos::ParameterList > MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getValidParamterList ( )
private

Definition at line 97 of file MueLu_RefMaxwell_def.hpp.

◆ initialize() [1/2]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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,
const Teuchos::RCP< MultiVector > &  Material,
Teuchos::ParameterList &  List 
)
private

Initialize with matrices except the Jacobian (don't compute the preconditioner)

Note: This uses old notation that only makes sense for curl-curl problems.

Parameters
[in]D0_MatrixDiscrete Gradient
[in]Ms_MatrixEdge mass matrix for nodal aggregates
[in]M0inv_MatrixInverse of lumped nodal mass matrix (add on only)
[in]M1_MatrixEdge mass matrix for add on
[in]NullspaceNull space (needed for periodic)
[in]CoordsNodal coordinates
[in]ListParameter list

Definition at line 2605 of file MueLu_RefMaxwell_def.hpp.

◆ initialize() [2/2]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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,
const Teuchos::RCP< MultiVector > &  Material_beta,
const Teuchos::RCP< MultiVector > &  Material_alpha,
Teuchos::ParameterList &  List 
)
private

Initialize with matrices except the Jacobian (don't compute the preconditioner)

Parameters
[in]knumber of the space in the deRham sequence of the problem to be solved
[in]Dk_1Discrete derivative from (k-1)-th to k-th space
[in]Dk_2Discrete derivative from (k-2)-th to (k-1)-th space
[in]D0Discrete Gradient
[in]M1_betaMass matrix on 1-st space with weight beta for nodal aggregates
[in]M1_alphaMass matrix on 1-st space with weight alpha for nodal aggregates
[in]Mk_oneMass matrix on k-th space with unit weight for addon11
[in]Mk_1_oneMass matrix on (k-1)-th space with unit weight for addon22
[in]invMk_1_invBetaApproximate inverse of mass matrix on (k-1)-th space with weight 1/beta (addon11 only)
[in]invMk_2_invAlphaApproximate inverse of mass matrix on (k-2)-th space with weight 1/alpha (addon22 only)
[in]Nullspace11Null space (needed for periodic)
[in]Nullspace22Null space (needed for periodic)
[in]NodalCoordsNodal coordinates
[in]ListParameter list

Definition at line 2626 of file MueLu_RefMaxwell_def.hpp.

◆ determineSubHierarchyCommSizes()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::determineSubHierarchyCommSizes ( bool &  doRebalancing,
int &  rebalanceStriding,
int &  numProcsCoarseA11,
int &  numProcsA22 
)
private

Determine how large the sub-communicators for the two hierarchies should be.

Definition at line 611 of file MueLu_RefMaxwell_def.hpp.

◆ buildCoarse11Matrix()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buildCoarse11Matrix ( )
private

Compute coarseA11 = P11^{T}*SM*P11 + addon efficiently.

Definition at line 793 of file MueLu_RefMaxwell_def.hpp.

◆ rebalanceCoarse11Matrix()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::rebalanceCoarse11Matrix ( const int  rebalanceStriding,
const int  numProcsCoarseA11 
)
private

rebalance the coarse A11 matrix, as well as P11, CoordsCoarse11 and Addon11

Definition at line 879 of file MueLu_RefMaxwell_def.hpp.

◆ build22Matrix()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::build22Matrix ( const bool  reuse,
const bool  doRebalancing,
const int  rebalanceStriding,
const int  numProcsA22 
)
private

Setup A22 = D0^T SM D0 and rebalance it, as well as D0 and Coords_.

Definition at line 999 of file MueLu_RefMaxwell_def.hpp.

◆ buildNullspace()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buildNullspace ( const int  spaceNumber,
const Kokkos::View< bool *, typename Node::device_type > &  bcs,
const bool  applyBCs 
)

Builds a nullspace.

Definition at line 1412 of file MueLu_RefMaxwell_def.hpp.

◆ buildProjection()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buildProjection ( const int  spaceNumber,
const RCP< MultiVector > &  EdgeNullspace 
) const

Builds a projection from a vector values space into a vector valued nodal space.

Definition at line 1577 of file MueLu_RefMaxwell_def.hpp.

◆ buildNodalProlongator()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buildNodalProlongator ( const Teuchos::RCP< Matrix > &  A_nodal,
Teuchos::RCP< Matrix > &  P_nodal,
Teuchos::RCP< MultiVector > &  Nullspace_nodal,
Teuchos::RCP< RealValuedMultiVector > &  Coords_nodal 
) const

Setup an auxiliary nodal prolongator

Parameters
[in]A_nodal
[out]P_nodal
[out]Nullspace_nodal

Definition at line 1710 of file MueLu_RefMaxwell_def.hpp.

◆ buildVectorNodalProlongator()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buildVectorNodalProlongator ( const Teuchos::RCP< Matrix > &  P_nodal) const

Setup a vectorial nodal prolongator

Parameters
[in]P_nodal_scalar
[out]P_nodal_vector

Definition at line 1825 of file MueLu_RefMaxwell_def.hpp.

◆ buildProlongator()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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

Setup a special prolongator from vectorial nodal to edge or face discretization

Parameters
[in]spaceNumberthe type of target discretization
[in]A_nodal_Matrixused for the nodal aggregation
[in]EdgeNullspaceedge nullspace
[out]edgeProlongatoredge prolongator
[out]coarseEdgeNullspacecoarse edge nullspace
[out]coarseNodalCoordscoarse nodal coordinates

Definition at line 1890 of file MueLu_RefMaxwell_def.hpp.

◆ buildAddon()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buildAddon ( const int  spaceNumber)

Construct an addon matrix

Parameters
[in]spaceNumberthe type of target discretization

Definition at line 706 of file MueLu_RefMaxwell_def.hpp.

◆ setupSubSolve()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setupSubSolve ( Teuchos::RCP< Hierarchy > &  hierarchy,
Teuchos::RCP< Operator > &  thyraPrecOp,
const Teuchos::RCP< Matrix > &  A,
const Teuchos::RCP< MultiVector > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
const Teuchos::RCP< MultiVector > &  Material,
Teuchos::ParameterList &  params,
std::string &  label,
const bool  reuse,
const bool  isSingular = false 
)
private

Setup a subsolve.

Definition at line 2029 of file MueLu_RefMaxwell_def.hpp.

◆ setFineLevelSmoother11()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setFineLevelSmoother11 ( )
private

Set the fine level smoother.

Definition at line 1191 of file MueLu_RefMaxwell_def.hpp.

◆ applyInverseAdditive()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::applyInverseAdditive ( const MultiVector &  RHS,
MultiVector &  X 
) const
private

apply additive algorithm for 2x2 solve

Definition at line 2115 of file MueLu_RefMaxwell_def.hpp.

◆ solveH()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::solveH ( const MultiVector &  RHS,
MultiVector &  X 
) const
private

apply solve to 1-1 block only

Definition at line 2388 of file MueLu_RefMaxwell_def.hpp.

◆ solve22()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::solve22 ( const MultiVector &  RHS,
MultiVector &  X 
) const
private

apply solve to 2-2 block only

Definition at line 2419 of file MueLu_RefMaxwell_def.hpp.

◆ allocateMemory()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::allocateMemory ( int  numVectors) const
private

allocate multivectors for solve

Definition at line 1264 of file MueLu_RefMaxwell_def.hpp.

◆ dump() [1/4]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dump ( const RCP< Matrix > &  A,
std::string  name 
) const
private

dump out matrix

Definition at line 1340 of file MueLu_RefMaxwell_def.hpp.

◆ dump() [2/4]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dump ( const RCP< MultiVector > &  X,
std::string  name 
) const
private

dump out multivector

Definition at line 1348 of file MueLu_RefMaxwell_def.hpp.

◆ dumpCoords()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dumpCoords ( const RCP< RealValuedMultiVector > &  X,
std::string  name 
) const
private

dump out real-valued multivector

Definition at line 1356 of file MueLu_RefMaxwell_def.hpp.

◆ dump() [3/4]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dump ( const Teuchos::ArrayRCP< bool > &  v,
std::string  name 
) const
private

dump out boolean ArrayView

Definition at line 1364 of file MueLu_RefMaxwell_def.hpp.

◆ dump() [4/4]

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dump ( const Kokkos::View< bool *, typename Node::device_type > &  v,
std::string  name 
) const
private

dump out boolean Kokkos::View

Definition at line 1374 of file MueLu_RefMaxwell_def.hpp.

◆ getTimer()

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Teuchos::TimeMonitor > MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getTimer ( std::string  name,
RCP< const Teuchos::Comm< int > >  comm = Teuchos::null 
) const
private

get a (synced) timer

Definition at line 1388 of file MueLu_RefMaxwell_def.hpp.

Member Data Documentation

◆ HierarchyCoarse11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Hierarchy> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::HierarchyCoarse11_
private

Two hierarchies: one for the coarse (1,1)-block, another for the (2,2)-block.

Definition at line 639 of file MueLu_RefMaxwell_decl.hpp.

◆ Hierarchy22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Hierarchy> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Hierarchy22_
private

Definition at line 639 of file MueLu_RefMaxwell_decl.hpp.

◆ PreSmoother11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<SmootherBase> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::PreSmoother11_
private

Definition at line 640 of file MueLu_RefMaxwell_decl.hpp.

◆ PostSmoother11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<SmootherBase> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::PostSmoother11_
private

Definition at line 640 of file MueLu_RefMaxwell_decl.hpp.

◆ PreSmootherData11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<SmootherPrototype> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::PreSmootherData11_
private

Definition at line 641 of file MueLu_RefMaxwell_decl.hpp.

◆ PostSmootherData11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<SmootherPrototype> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::PostSmootherData11_
private

Definition at line 641 of file MueLu_RefMaxwell_decl.hpp.

◆ thyraPrecOpH_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP<Operator> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::thyraPrecOpH_
private

Definition at line 642 of file MueLu_RefMaxwell_decl.hpp.

◆ thyraPrecOp22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP<Operator> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::thyraPrecOp22_
private

Definition at line 642 of file MueLu_RefMaxwell_decl.hpp.

◆ spaceNumber_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::spaceNumber_
private

The number of the space in the deRham complex.

Definition at line 644 of file MueLu_RefMaxwell_decl.hpp.

◆ solverName_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
std::string MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::solverName_
private

The name of the solver.

Definition at line 646 of file MueLu_RefMaxwell_decl.hpp.

◆ dim_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dim_
private

The spatial dimension.

Definition at line 648 of file MueLu_RefMaxwell_decl.hpp.

◆ SM_Matrix_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SM_Matrix_
private

The system that is getting preconditioned.

Definition at line 650 of file MueLu_RefMaxwell_decl.hpp.

◆ Dk_1_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Dk_1_
private

D_{k-1} matrix and its transpose.

Definition at line 652 of file MueLu_RefMaxwell_decl.hpp.

◆ Dk_1_T_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Dk_1_T_
private

Definition at line 652 of file MueLu_RefMaxwell_decl.hpp.

◆ Dk_2_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Dk_2_
private

D_{k-2} matrix.

Definition at line 654 of file MueLu_RefMaxwell_decl.hpp.

◆ D0_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::D0_
private

D_0 matrix.

Definition at line 656 of file MueLu_RefMaxwell_decl.hpp.

◆ invMk_1_invBeta_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::invMk_1_invBeta_
private

inverse of mass matrices on (k-1)-th and (k-2)-th space with weights 1/beta and 1/alpha respectively

Definition at line 658 of file MueLu_RefMaxwell_decl.hpp.

◆ invMk_2_invAlpha_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::invMk_2_invAlpha_
private

Definition at line 658 of file MueLu_RefMaxwell_decl.hpp.

◆ Mk_one_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Mk_one_
private

mass matrices with unit weight on k-th and (k-1)-th spaces

Definition at line 660 of file MueLu_RefMaxwell_decl.hpp.

◆ Mk_1_one_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Mk_1_one_
private

Definition at line 660 of file MueLu_RefMaxwell_decl.hpp.

◆ M1_beta_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::M1_beta_
private

mass matrices on first space with weights beta and alpha respectively

Definition at line 662 of file MueLu_RefMaxwell_decl.hpp.

◆ M1_alpha_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::M1_alpha_
private

Definition at line 662 of file MueLu_RefMaxwell_decl.hpp.

◆ Material_beta_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Material_beta_
private

material for first space

Definition at line 664 of file MueLu_RefMaxwell_decl.hpp.

◆ Material_alpha_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Material_alpha_
private

Definition at line 664 of file MueLu_RefMaxwell_decl.hpp.

◆ P11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::P11_
private

special prolongator for 11 block and its transpose

Definition at line 666 of file MueLu_RefMaxwell_decl.hpp.

◆ R11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::R11_
private

Definition at line 666 of file MueLu_RefMaxwell_decl.hpp.

◆ P22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::P22_
private

special prolongator for 22 block and its transpose

Definition at line 668 of file MueLu_RefMaxwell_decl.hpp.

◆ R22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::R22_
private

Definition at line 668 of file MueLu_RefMaxwell_decl.hpp.

◆ coarseA11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::coarseA11_
private

coarse 11, 22 and coarse 22 blocks

Definition at line 670 of file MueLu_RefMaxwell_decl.hpp.

◆ A22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::A22_
private

Definition at line 670 of file MueLu_RefMaxwell_decl.hpp.

◆ coarseA22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::coarseA22_
private

Definition at line 670 of file MueLu_RefMaxwell_decl.hpp.

◆ Addon11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Addon11_
private

the addon for the 11 block

Definition at line 672 of file MueLu_RefMaxwell_decl.hpp.

◆ Addon22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Addon22_
private

the addon for the 22 block

Definition at line 674 of file MueLu_RefMaxwell_decl.hpp.

◆ DorigDomainMap_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<const Map> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DorigDomainMap_
private

Definition at line 675 of file MueLu_RefMaxwell_decl.hpp.

◆ DorigImporter_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<const Import> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DorigImporter_
private

Definition at line 676 of file MueLu_RefMaxwell_decl.hpp.

◆ BCrows11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Kokkos::View<bool *, typename Node::device_type::memory_space> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BCrows11_
private

Vectors for BCs.

Definition at line 678 of file MueLu_RefMaxwell_decl.hpp.

◆ BCcols22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Kokkos::View<bool *, typename Node::device_type::memory_space> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BCcols22_
private

Definition at line 678 of file MueLu_RefMaxwell_decl.hpp.

◆ BCdomain22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Kokkos::View<bool *, typename Node::device_type::memory_space> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BCdomain22_
private

Definition at line 678 of file MueLu_RefMaxwell_decl.hpp.

◆ globalNumberBoundaryUnknowns11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::globalNumberBoundaryUnknowns11_
private

Definition at line 679 of file MueLu_RefMaxwell_decl.hpp.

◆ globalNumberBoundaryUnknowns22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::globalNumberBoundaryUnknowns22_
private

Definition at line 679 of file MueLu_RefMaxwell_decl.hpp.

◆ Nullspace11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Nullspace11_
private

Nullspace for (1.1) block.

Definition at line 681 of file MueLu_RefMaxwell_decl.hpp.

◆ Nullspace22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Nullspace22_
private

Definition at line 681 of file MueLu_RefMaxwell_decl.hpp.

◆ NodalCoords_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<RealValuedMultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::NodalCoords_
private

Coordinates.

Definition at line 683 of file MueLu_RefMaxwell_decl.hpp.

◆ CoordsCoarse11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<RealValuedMultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CoordsCoarse11_
private

Definition at line 683 of file MueLu_RefMaxwell_decl.hpp.

◆ Coords22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<RealValuedMultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Coords22_
private

Definition at line 683 of file MueLu_RefMaxwell_decl.hpp.

◆ NullspaceCoarse11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::NullspaceCoarse11_
private

Nullspace for coarse (1,1) problem.

Definition at line 685 of file MueLu_RefMaxwell_decl.hpp.

◆ CoarseNullspace22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CoarseNullspace22_
private

Nullspace for coarse (2,2) problem.

Definition at line 687 of file MueLu_RefMaxwell_decl.hpp.

◆ ImporterCoarse11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<const Import> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ImporterCoarse11_
private

Importer to coarse (1,1) hierarchy.

Definition at line 689 of file MueLu_RefMaxwell_decl.hpp.

◆ Importer22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<const Import> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Importer22_
private

Definition at line 689 of file MueLu_RefMaxwell_decl.hpp.

◆ Dk_1_T_R11_colMapsMatch_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Dk_1_T_R11_colMapsMatch_
private

Definition at line 690 of file MueLu_RefMaxwell_decl.hpp.

◆ asyncTransfers_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::asyncTransfers_
private

Definition at line 690 of file MueLu_RefMaxwell_decl.hpp.

◆ onlyBoundary11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::onlyBoundary11_
private

Definition at line 691 of file MueLu_RefMaxwell_decl.hpp.

◆ onlyBoundary22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::onlyBoundary22_
private

Definition at line 691 of file MueLu_RefMaxwell_decl.hpp.

◆ parameterList_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ParameterList MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::parameterList_
private

Parameter lists.

Definition at line 693 of file MueLu_RefMaxwell_decl.hpp.

◆ precList11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ParameterList MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::precList11_
private

Definition at line 693 of file MueLu_RefMaxwell_decl.hpp.

◆ precList22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ParameterList MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::precList22_
private

Definition at line 693 of file MueLu_RefMaxwell_decl.hpp.

◆ coarseA11_AP_reuse_data_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Teuchos::ParameterList> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::coarseA11_AP_reuse_data_
private

Definition at line 694 of file MueLu_RefMaxwell_decl.hpp.

◆ coarseA11_RAP_reuse_data_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Teuchos::ParameterList> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::coarseA11_RAP_reuse_data_
private

Definition at line 694 of file MueLu_RefMaxwell_decl.hpp.

◆ A22_AP_reuse_data_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Teuchos::ParameterList> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::A22_AP_reuse_data_
private

Definition at line 695 of file MueLu_RefMaxwell_decl.hpp.

◆ A22_RAP_reuse_data_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Teuchos::ParameterList> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::A22_RAP_reuse_data_
private

Definition at line 695 of file MueLu_RefMaxwell_decl.hpp.

◆ disable_addon_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::disable_addon_
private

Some options.

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ disable_addon_22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::disable_addon_22_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ dump_matrices_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dump_matrices_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ useKokkos_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::useKokkos_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ use_as_preconditioner_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::use_as_preconditioner_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ implicitTranspose_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::implicitTranspose_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ fuseProlongationAndUpdate_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fuseProlongationAndUpdate_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ syncTimers_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::syncTimers_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ enable_reuse_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::enable_reuse_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ skipFirst11Level_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::skipFirst11Level_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ skipFirst22Level_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::skipFirst22Level_
private

Definition at line 697 of file MueLu_RefMaxwell_decl.hpp.

◆ applyBCsToAnodal_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::applyBCsToAnodal_
private

Definition at line 698 of file MueLu_RefMaxwell_decl.hpp.

◆ applyBCsToCoarse11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::applyBCsToCoarse11_
private

Definition at line 698 of file MueLu_RefMaxwell_decl.hpp.

◆ applyBCsTo22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::applyBCsTo22_
private

Definition at line 698 of file MueLu_RefMaxwell_decl.hpp.

◆ numItersCoarse11_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::numItersCoarse11_
private

Definition at line 699 of file MueLu_RefMaxwell_decl.hpp.

◆ numIters22_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::numIters22_
private

Definition at line 699 of file MueLu_RefMaxwell_decl.hpp.

◆ mode_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
std::string MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mode_
private

Definition at line 700 of file MueLu_RefMaxwell_decl.hpp.

◆ P11res_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::P11res_
mutableprivate

Temporary memory.

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ P11x_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::P11x_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ P11resSubComm_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::P11resSubComm_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ P11xSubComm_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::P11xSubComm_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ DresIntermediate_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DresIntermediate_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ Dres_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Dres_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ DxIntermediate_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DxIntermediate_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ Dx_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Dx_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ DresSubComm_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DresSubComm_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ DxSubComm_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DxSubComm_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ residual_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::residual_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ P11resTmp_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::P11resTmp_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ DresTmp_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DresTmp_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ DTR11Tmp_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DTR11Tmp_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ P11x_colmap_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::P11x_colmap_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.

◆ Dx_colmap_

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Dx_colmap_
private

Definition at line 703 of file MueLu_RefMaxwell_decl.hpp.


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