MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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_fwd.hpp>

Inheritance diagram for MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
MueLu::VerboseObject Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > Teuchos::VerboseObject< VerboseObject > Teuchos::Describable Teuchos::VerboseObjectBase Teuchos::LabeledObject

Public Types

typedef Teuchos::ScalarTraits
< Scalar >::magnitudeType 
magnitudeType
 
typedef Xpetra::MultiVector
< magnitudeType, LO, GO, NO
RealValuedMultiVector
 

Public Member Functions

 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 > &D0_Matrix, const Teuchos::RCP< Matrix > &M0inv_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, 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 > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, 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 > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList &List, bool ComputePrec)
 
 RefMaxwell (const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList &List)
 
 RefMaxwell (const Teuchos::RCP< Matrix > &SM_Matrix, Teuchos::ParameterList &List, bool ComputePrec=true)
 
virtual ~RefMaxwell ()
 Destructor. More...
 
Teuchos::RCP< const MapgetDomainMap () const
 Returns the Xpetra::Map object associated with the domain of this operator. More...
 
Teuchos::RCP< const MapgetRangeMap () 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 ()
 Setup the preconditioner. More...
 
void buildProlongator ()
 Setup the prolongator for the (1,1)-block. More...
 
void formCoarseMatrix ()
 Compute P11^{T}*A*P11 efficiently. More...
 
void resetMatrix (Teuchos::RCP< Matrix > SM_Matrix_new)
 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
 
- Public Member Functions inherited from MueLu::VerboseObject
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::FancyOStreamGetOStream (MsgType type, int thisProcRankOnly=0) const
 Get an output stream for outputting the input message type. More...
 
Teuchos::FancyOStreamGetBlackHole () const
 
 VerboseObject ()
 
virtual ~VerboseObject ()
 Destructor. More...
 
- Public Member Functions inherited from Teuchos::VerboseObject< VerboseObject >
 VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null)
 
virtual const VerboseObjectsetVerbLevel (const EVerbosityLevel verbLevel) const
 
virtual const VerboseObjectsetOverridingVerbLevel (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)
 
- Public Member Functions inherited from Teuchos::VerboseObjectBase
virtual ~VerboseObjectBase ()
 
 VerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null)
 
virtual const VerboseObjectBasesetOStream (const RCP< FancyOStream > &oStream) const
 
virtual const VerboseObjectBasesetOverridingOStream (const RCP< FancyOStream > &oStream) const
 
virtual VerboseObjectBasesetLinePrefix (const std::string &linePrefix)
 
virtual RCP< FancyOStreamgetOStream () const
 
virtual RCP< FancyOStreamgetOverridingOStream () const
 
virtual std::string getLinePrefix () const
 
virtual OSTab getOSTab (const int tabs=1, const std::string &linePrefix="") const
 
- Public Member Functions inherited from Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >
virtual void removeEmptyProcessesInPlace (const RCP< const Map > &)
 
- Public Member Functions inherited from Teuchos::Describable
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)
 

Private Member Functions

void initialize (const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &M0inv_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList &List)
 
void solveH () const
 solve coarse (1,1) block More...
 
void solve22 () const
 solve (2,2) block More...
 
void applyInverseAdditive (const MultiVector &RHS, MultiVector &X) const
 apply additive algorithm for 2x2 solve More...
 
void applyInverse121 (const MultiVector &RHS, MultiVector &X) const
 apply 1-2-1 algorithm for 2x2 solve More...
 
void applyInverse212 (const MultiVector &RHS, MultiVector &X) const
 apply 2-1-2 algorithm for 2x2 solve More...
 
void applyInverse11only (const MultiVector &RHS, MultiVector &X) const
 apply solve to 1-1 block only More...
 

Private Attributes

Teuchos::RCP< HierarchyHierarchyH_
 Two hierarchies: one for the coarse (1,1)-block, another for the (2,2)-block. More...
 
Teuchos::RCP< HierarchyHierarchy22_
 
Teuchos::RCP< SmootherBasePreSmoother_
 
Teuchos::RCP< SmootherBasePostSmoother_
 
bool useHiptmairSmoothing_
 
Teuchos::RCP< Matrix > SM_Matrix_
 Various matrices. More...
 
Teuchos::RCP< Matrix > D0_Matrix_
 
Teuchos::RCP< Matrix > D0_T_Matrix_
 
Teuchos::RCP< Matrix > M0inv_Matrix_
 
Teuchos::RCP< Matrix > M1_Matrix_
 
Teuchos::RCP< Matrix > Ms_Matrix_
 
Teuchos::RCP< Matrix > A_nodal_Matrix_
 
Teuchos::RCP< Matrix > P11_
 
Teuchos::RCP< Matrix > R11_
 
Teuchos::RCP< Matrix > AH_
 
Teuchos::RCP< Matrix > A22_
 
Teuchos::ArrayRCP< const bool > BCrows_
 Vectors for BCs. More...
 
Teuchos::ArrayRCP< const bool > BCcols_
 
Teuchos::RCP< MultiVectorNullspace_
 Nullspace. More...
 
Teuchos::RCP
< RealValuedMultiVector
Coords_
 Coordinates. More...
 
Teuchos::RCP
< RealValuedMultiVector
CoordsH_
 
Teuchos::RCP< const Import > ImporterH_
 Importer to coarse (1,1) hierarchy. More...
 
Teuchos::RCP< const Import > Importer22_
 
Teuchos::ParameterList parameterList_
 Parameter lists. More...
 
Teuchos::ParameterList precList11_
 
Teuchos::ParameterList precList22_
 
Teuchos::ParameterList smootherList_
 
bool disable_addon_
 Some options. More...
 
bool dump_matrices_
 
bool useKokkos_
 
bool use_as_preconditioner_
 
std::string mode_
 
Teuchos::RCP< MultiVectorP11res_
 Temporary memory. More...
 
Teuchos::RCP< MultiVectorP11x_
 
Teuchos::RCP< MultiVectorD0res_
 
Teuchos::RCP< MultiVectorD0x_
 
Teuchos::RCP< MultiVectorresidual_
 
Teuchos::RCP< MultiVectorP11resTmp_
 
Teuchos::RCP< MultiVectorP11xTmp_
 
Teuchos::RCP< MultiVectorD0resTmp_
 
Teuchos::RCP< MultiVectorD0xTmp_
 

Additional Inherited Members

- Static Public Member Functions inherited from MueLu::VerboseObject
static void SetMueLuOStream (const Teuchos::RCP< Teuchos::FancyOStream > &mueluOStream)
 
static Teuchos::RCP
< Teuchos::FancyOStream
GetMueLuOStream ()
 
static void SetDefaultVerbLevel (const VerbLevel defaultVerbLevel)
 Set the default (global) verbosity level. More...
 
static VerbLevel GetDefaultVerbLevel ()
 Get the default (global) verbosity level. More...
 
- Static Public Member Functions inherited from Teuchos::VerboseObject< VerboseObject >
static void setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel)
 
static EVerbosityLevel getDefaultVerbLevel ()
 
- Static Public Member Functions inherited from Teuchos::VerboseObjectBase
static void setDefaultOStream (const RCP< FancyOStream > &defaultOStream)
 
static RCP< FancyOStreamgetDefaultOStream ()
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 
- Protected Member Functions inherited from Teuchos::VerboseObject< VerboseObject >
void initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null)
 
- Protected Member Functions inherited from Teuchos::VerboseObjectBase
void initializeVerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null)
 
virtual void informUpdatedVerbosityState () const
 

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.

This uses a 2x2 block reformulation.

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"
  • refmaxwell: disable addon - bool specifing whether the addon should be built for stabilization. Default: "true"
  • refmaxwell: use as preconditioner - bool specifing whether RefMaxwell is used as a preconditioner or as a solver.
  • refmaxwell: dump matrices - bool specifing whether the matrices should be dumped. Default: "false"
  • refmaxwell: prolongator compute algorithm - a string specifying the algorithm to build the prolongator. Allowed values are: "mat-mat" and "gustavson"
  • 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 54 of file MueLu_RefMaxwell_fwd.hpp.

Member Typedef Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::magnitudeType

Definition at line 131 of file MueLu_RefMaxwell_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Xpetra::MultiVector<magnitudeType,LO,GO,NO> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::RealValuedMultiVector

Definition at line 132 of file MueLu_RefMaxwell_decl.hpp.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 135 of file MueLu_RefMaxwell_decl.hpp.

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 144 of file MueLu_RefMaxwell_decl.hpp.

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 > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
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 163 of file MueLu_RefMaxwell_decl.hpp.

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 > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
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 190 of file MueLu_RefMaxwell_decl.hpp.

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 > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
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 210 of file MueLu_RefMaxwell_decl.hpp.

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 > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
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 235 of file MueLu_RefMaxwell_decl.hpp.

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

Constructor with parameter list

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

Definition at line 250 of file MueLu_RefMaxwell_decl.hpp.

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

Destructor.

Definition at line 272 of file MueLu_RefMaxwell_decl.hpp.

Member Function Documentation

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

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

Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 104 of file MueLu_RefMaxwell_def.hpp.

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

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

Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 110 of file MueLu_RefMaxwell_def.hpp.

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 281 of file MueLu_RefMaxwell_decl.hpp.

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 116 of file MueLu_RefMaxwell_def.hpp.

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

Setup the preconditioner.

Definition at line 146 of file MueLu_RefMaxwell_def.hpp.

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

Setup the prolongator for the (1,1)-block.

Definition at line 857 of file MueLu_RefMaxwell_def.hpp.

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

Compute P11^{T}*A*P11 efficiently.

Definition at line 1347 of file MueLu_RefMaxwell_def.hpp.

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

Reset system matrix.

Definition at line 1437 of file MueLu_RefMaxwell_def.hpp.

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
virtual

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.

Implements Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 1690 of file MueLu_RefMaxwell_def.hpp.

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 1766 of file MueLu_RefMaxwell_def.hpp.

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
virtual

Reimplemented from Teuchos::Describable.

Definition at line 1803 of file MueLu_RefMaxwell_def.hpp.

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 > &  M0inv_Matrix,
const Teuchos::RCP< Matrix > &  M1_Matrix,
const Teuchos::RCP< MultiVector > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
Teuchos::ParameterList List 
)
private

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

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

Definition at line 1773 of file MueLu_RefMaxwell_def.hpp.

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

solve coarse (1,1) block

Definition at line 1443 of file MueLu_RefMaxwell_def.hpp.

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

solve (2,2) block

Definition at line 1473 of file MueLu_RefMaxwell_def.hpp.

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 1502 of file MueLu_RefMaxwell_def.hpp.

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

apply 1-2-1 algorithm for 2x2 solve

Definition at line 1589 of file MueLu_RefMaxwell_def.hpp.

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

apply 2-1-2 algorithm for 2x2 solve

Definition at line 1629 of file MueLu_RefMaxwell_def.hpp.

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

apply solve to 1-1 block only

Definition at line 1668 of file MueLu_RefMaxwell_def.hpp.

Member Data Documentation

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

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

Definition at line 361 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 361 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 362 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 362 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 366 of file MueLu_RefMaxwell_decl.hpp.

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

Various matrices.

Definition at line 368 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 368 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 368 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 368 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 368 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 368 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 369 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 369 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 369 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 369 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 369 of file MueLu_RefMaxwell_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP<const bool> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BCrows_
private

Vectors for BCs.

Definition at line 375 of file MueLu_RefMaxwell_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP<const bool> MueLu::RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BCcols_
private

Definition at line 376 of file MueLu_RefMaxwell_decl.hpp.

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

Nullspace.

Definition at line 378 of file MueLu_RefMaxwell_decl.hpp.

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

Coordinates.

Definition at line 380 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 380 of file MueLu_RefMaxwell_decl.hpp.

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

Importer to coarse (1,1) hierarchy.

Definition at line 382 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 382 of file MueLu_RefMaxwell_decl.hpp.

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

Parameter lists.

Definition at line 384 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 384 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 384 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 384 of file MueLu_RefMaxwell_decl.hpp.

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

Some options.

Definition at line 386 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 386 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 386 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 386 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 387 of file MueLu_RefMaxwell_decl.hpp.

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 389 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 389 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 389 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 389 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 389 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 389 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 389 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 389 of file MueLu_RefMaxwell_decl.hpp.

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

Definition at line 389 of file MueLu_RefMaxwell_decl.hpp.


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