MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

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

#include <MueLu_Maxwell1_decl.hpp>

Inheritance diagram for MueLu::Maxwell1< 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 Teuchos::ScalarTraits
< Scalar >::coordinateType 
coordinateType
 
typedef Xpetra::MultiVector
< coordinateType, LO, GO, NO
RealValuedMultiVector
 

Public Member Functions

 Maxwell1 ()
 Constructor. More...
 
 Maxwell1 (Teuchos::RCP< Hierarchy > H11, Teuchos::RCP< Hierarchy > H22)
 Constructor with Hierarchies. More...
 
 Maxwell1 (const Teuchos::RCP< Matrix > &SM_Matrix, const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList &List, bool ComputePrec=true)
 
 Maxwell1 (const Teuchos::RCP< Matrix > &SM_Matrix, const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Kn_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList &List, bool ComputePrec=true)
 
 Maxwell1 (const Teuchos::RCP< Matrix > &SM_Matrix, const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Kn_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList &List, const Teuchos::RCP< Matrix > &GmhdA_Matrix, bool ComputePrec=true)
 
 Maxwell1 (const Teuchos::RCP< Matrix > &SM_Matrix, Teuchos::ParameterList &List, bool ComputePrec=true)
 
virtual ~Maxwell1 ()
 Destructor. More...
 
const Teuchos::RCP< const MapgetDomainMap () const
 Returns the Xpetra::Map object associated with the domain of this operator. More...
 
const 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 (bool reuse=false)
 Setup the preconditioner. More...
 
void resetMatrix (Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
 Reset system matrix. More...
 
void apply (const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
 
bool hasTransposeApply () const
 Indicates whether this operator supports applying the adjoint operator. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
 
void residual (const MultiVector &X, const MultiVector &B, MultiVector &R) const
 Compute a residual R = B - (*this) * X. More...
 
- 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 Types

enum  mode_type { MODE_STANDARD = 0, MODE_REFMAXWELL, MODE_EDGE_ONLY, MODE_GMHD_STANDARD }
 Execution modes. More...
 

Private Member Functions

Teuchos::RCP< Matrix > generate_kn () const
 Generates the Kn matrix. More...
 
void GMHDSetupHierarchy (Teuchos::ParameterList &List) const
 Sets up hiearchy for GMHD matrices that include generalized Ohms law equations. More...
 
void initialize (const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Kn_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, Teuchos::ParameterList &List)
 
void applyInverseRefMaxwellAdditive (const MultiVector &RHS, MultiVector &X) const
 apply RefMaxwell additive 2x2 style cycle More...
 
void applyInverseStandard (const MultiVector &RHS, MultiVector &X) const
 apply standard Maxwell1 cycle More...
 
void allocateMemory (int numVectors) const
 allocate multivectors for solve More...
 
void dump (const Matrix &A, std::string name) const
 dump out matrix More...
 
void dump (const MultiVector &X, std::string name) const
 dump out multivector More...
 
void dumpCoords (const RealValuedMultiVector &X, std::string name) const
 dump out real-valued multivector More...
 
void dump (const Teuchos::ArrayRCP< bool > &v, std::string name) const
 dump out boolean ArrayView More...
 
void dump (const Kokkos::View< bool *, typename Node::device_type > &v, std::string name) const
 dump out boolean Kokkos::View More...
 
Teuchos::RCP
< Teuchos::TimeMonitor
getTimer (std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
 get a (synced) timer More...
 

Private Attributes

Teuchos::ParameterList parameterList_
 ParameterLists. More...
 
Teuchos::ParameterList precList11_
 
Teuchos::ParameterList precList22_
 
Teuchos::RCP< HierarchyHierarchy11_
 Two hierarchies: one for the (1,1)-block, another for the (2,2)-block. More...
 
Teuchos::RCP< HierarchyHierarchy22_
 
Teuchos::RCP< HierarchyHierarchyGmhd_
 
Teuchos::RCP< Matrix > SM_Matrix_
 Various matrices. More...
 
Teuchos::RCP< Matrix > D0_Matrix_
 
Teuchos::RCP< Matrix > Kn_Matrix_
 
Teuchos::RCP< Matrix > GmhdA_Matrix_
 
Kokkos::View< bool *, typename
Node::device_type::memory_space > 
BCrowsKokkos_
 Vectors for BCs. More...
 
Kokkos::View< bool *, typename
Node::device_type::memory_space > 
BCcolsKokkos_
 
Kokkos::View< bool *, typename
Node::device_type::memory_space > 
BCdomainKokkos_
 
int BCedges_
 
int BCnodes_
 
Teuchos::ArrayRCP< bool > BCrows_
 
Teuchos::ArrayRCP< bool > BCcols_
 
Teuchos::ArrayRCP< bool > BCdomain_
 
Teuchos::RCP< MultiVectorNullspace_
 Nullspace. More...
 
Teuchos::RCP
< RealValuedMultiVector
Coords_
 Coordinates. More...
 
bool useKokkos_
 Some options. More...
 
bool allEdgesBoundary_
 
bool allNodesBoundary_
 
bool dump_matrices_
 
bool enable_reuse_
 
bool syncTimers_
 
bool applyBCsTo22_
 
mode_type mode_
 
RCP< Matrix > P11_
 Temporary memory (cached vectors for RefMaxwell-style) More...
 
Teuchos::RCP< MultiVectorresidualFine_
 
Teuchos::RCP< MultiVectorresidual11c_
 
Teuchos::RCP< MultiVectorresidual22_
 
Teuchos::RCP< MultiVectorupdate11c_
 
Teuchos::RCP< MultiVectorupdate22_
 

Additional Inherited Members

- Static Public Member Functions inherited from MueLu::VerboseObject
static void SetMueLuOStream (const Teuchos::RCP< Teuchos::FancyOStream > &mueluOStream)
 
static void SetMueLuOFileStream (const std::string &filename)
 
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::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >

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

Definition at line 80 of file MueLu_Maxwell1_decl.hpp.

Member Typedef Documentation

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

Definition at line 85 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Teuchos::ScalarTraits<Scalar>::coordinateType MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::coordinateType

Definition at line 86 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 87 of file MueLu_Maxwell1_decl.hpp.

Member Enumeration Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
enum MueLu::Maxwell1::mode_type
private

Execution modes.

Enumerator
MODE_STANDARD 
MODE_REFMAXWELL 
MODE_EDGE_ONLY 
MODE_GMHD_STANDARD 

Definition at line 310 of file MueLu_Maxwell1_decl.hpp.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 90 of file MueLu_Maxwell1_decl.hpp.

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

Constructor with Hierarchies.

Definition at line 98 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Maxwell1 ( const Teuchos::RCP< Matrix > &  SM_Matrix,
const Teuchos::RCP< Matrix > &  D0_Matrix,
const Teuchos::RCP< MultiVector > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
Teuchos::ParameterList List,
bool  ComputePrec = true 
)
inline

Constructor with Jacobian

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

Definition at line 114 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Maxwell1 ( const Teuchos::RCP< Matrix > &  SM_Matrix,
const Teuchos::RCP< Matrix > &  D0_Matrix,
const Teuchos::RCP< Matrix > &  Kn_Matrix,
const Teuchos::RCP< MultiVector > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
Teuchos::ParameterList List,
bool  ComputePrec = true 
)
inline

Constructor with Jacobian and nodal matrix

Parameters
[in]SM_MatrixJacobian
[in]D0_MatrixDiscrete Gradient
[in]Kn_MatrixNodal Laplacian
[in]CoordsNodal coordinates
[in]ListParameter list
[in]ComputePrecIf true, compute the preconditioner immediately

Definition at line 135 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Maxwell1 ( const Teuchos::RCP< Matrix > &  SM_Matrix,
const Teuchos::RCP< Matrix > &  D0_Matrix,
const Teuchos::RCP< Matrix > &  Kn_Matrix,
const Teuchos::RCP< MultiVector > &  Nullspace,
const Teuchos::RCP< RealValuedMultiVector > &  Coords,
Teuchos::ParameterList List,
const Teuchos::RCP< Matrix > &  GmhdA_Matrix,
bool  ComputePrec = true 
)
inline

Gmhd GMHD Constructor with Jacobian and nodal matrix AND Gmhd matrix

Parameters
[in]SM_MatrixJacobian
[in]D0_MatrixDiscrete Gradient
[in]Kn_MatrixNodal Laplacian
[in]CoordsNodal coordinates
[in]ListParameter list
[in]GmhdA_MatrixGmhd matrix including generalized Ohms law equations
[in]ComputePrecIf true, compute the preconditioner immediately

Definition at line 157 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Maxwell1 ( 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 178 of file MueLu_Maxwell1_decl.hpp.

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

Destructor.

Definition at line 196 of file MueLu_Maxwell1_decl.hpp.

Member Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > MueLu::Maxwell1< 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 82 of file MueLu_Maxwell1_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > MueLu::Maxwell1< 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 87 of file MueLu_Maxwell1_def.hpp.

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

Returns Jacobian matrix SM.

Definition at line 205 of file MueLu_Maxwell1_decl.hpp.

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

Set parameters.

Definition at line 92 of file MueLu_Maxwell1_def.hpp.

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

Setup the preconditioner.

Definition at line 258 of file MueLu_Maxwell1_def.hpp.

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

Reset system matrix.

Definition at line 679 of file MueLu_Maxwell1_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::Maxwell1< 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 747 of file MueLu_Maxwell1_def.hpp.

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

Indicates whether this operator supports applying the adjoint operator.

Definition at line 763 of file MueLu_Maxwell1_def.hpp.

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

Reimplemented from Teuchos::Describable.

Definition at line 843 of file MueLu_Maxwell1_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::Maxwell1< 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 232 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::generate_kn ( ) const
private

Generates the Kn matrix.

Definition at line 586 of file MueLu_Maxwell1_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GMHDSetupHierarchy ( Teuchos::ParameterList List) const
private

Sets up hiearchy for GMHD matrices that include generalized Ohms law equations.

Definition at line 239 of file MueLu_Maxwell1_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::initialize ( const Teuchos::RCP< Matrix > &  D0_Matrix,
const Teuchos::RCP< Matrix > &  Kn_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]Kn_MatrixKn nodal matrix
[in]NullspaceNull space (needed for periodic)
[in]CoordsNodal coordinates
[in]ListParameter list

Definition at line 769 of file MueLu_Maxwell1_def.hpp.

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

apply RefMaxwell additive 2x2 style cycle

Definition at line 688 of file MueLu_Maxwell1_def.hpp.

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

apply standard Maxwell1 cycle

Definition at line 742 of file MueLu_Maxwell1_def.hpp.

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

allocate multivectors for solve

Definition at line 597 of file MueLu_Maxwell1_def.hpp.

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

dump out matrix

Definition at line 618 of file MueLu_Maxwell1_def.hpp.

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

dump out multivector

Definition at line 626 of file MueLu_Maxwell1_def.hpp.

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

dump out real-valued multivector

Definition at line 634 of file MueLu_Maxwell1_def.hpp.

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

dump out boolean ArrayView

Definition at line 642 of file MueLu_Maxwell1_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::Maxwell1< 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 652 of file MueLu_Maxwell1_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Teuchos::TimeMonitor > MueLu::Maxwell1< 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 664 of file MueLu_Maxwell1_def.hpp.

Member Data Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ParameterList MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::parameterList_
mutableprivate

ParameterLists.

Definition at line 289 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ParameterList MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::precList11_
mutableprivate

Definition at line 289 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ParameterList MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::precList22_
mutableprivate

Definition at line 289 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Hierarchy> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Hierarchy11_
private

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

Definition at line 292 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 292 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Hierarchy> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::HierarchyGmhd_
private

Definition at line 292 of file MueLu_Maxwell1_decl.hpp.

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

Various matrices.

Definition at line 295 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 295 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Kn_Matrix_
private

Definition at line 295 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GmhdA_Matrix_
private

Definition at line 295 of file MueLu_Maxwell1_decl.hpp.

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

Vectors for BCs.

Definition at line 298 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 298 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 298 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BCedges_
private

Definition at line 299 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BCnodes_
private

Definition at line 299 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 300 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 300 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP<bool> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BCdomain_
private

Definition at line 300 of file MueLu_Maxwell1_decl.hpp.

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

Nullspace.

Definition at line 302 of file MueLu_Maxwell1_decl.hpp.

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

Coordinates.

Definition at line 304 of file MueLu_Maxwell1_decl.hpp.

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

Some options.

Definition at line 306 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::allEdgesBoundary_
private

Definition at line 306 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::allNodesBoundary_
private

Definition at line 306 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 306 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 306 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 306 of file MueLu_Maxwell1_decl.hpp.

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

Definition at line 307 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
mode_type MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mode_
private

Definition at line 314 of file MueLu_Maxwell1_decl.hpp.

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

Temporary memory (cached vectors for RefMaxwell-style)

Definition at line 317 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::residualFine_
mutableprivate

Definition at line 318 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::residual11c_
mutableprivate

Definition at line 318 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::residual22_
mutableprivate

Definition at line 318 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::update11c_
mutableprivate

Definition at line 318 of file MueLu_Maxwell1_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<MultiVector> MueLu::Maxwell1< Scalar, LocalOrdinal, GlobalOrdinal, Node >::update22_
mutableprivate

Definition at line 318 of file MueLu_Maxwell1_decl.hpp.


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