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

Preconditioner (wrapped as a Xpetra::Operator) for solving MultiPhysics PDEs. More...

#include <MueLu_MultiPhys_decl.hpp>

Inheritance diagram for MueLu::MultiPhys< 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

 MultiPhys ()
 Constructor. More...
 
 MultiPhys (const Teuchos::RCP< Matrix > &AmatMultiPhysics, const Teuchos::ArrayRCP< RCP< Matrix >> arrayOfAuxMatrices, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector >> arrayOfNullspaces, const Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector >> arrayOfCoords, const int nBlks, Teuchos::ParameterList &List, bool ComputePrec=true)
 
virtual ~MultiPhys ()
 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...
 
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 Member Functions

void initialize (const Teuchos::RCP< Matrix > &AmatMultiPhysics, const Teuchos::ArrayRCP< RCP< Matrix >> arrayOfAuxMatrices, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector >> arrayOfNullspaces, const Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector >> arrayOfCoords, const int nBlks, Teuchos::ParameterList &List)
 
void applyInverse (const MultiVector &RHS, MultiVector &X) const
 apply standard MultiPhys cycle More...
 
void allocateMemory (int numVectors) const
 allocate multivectors for solve 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::RCP< Matrix > AmatMultiphysics_
 Hierarchies: used to define P for (0,0)-block, .... (nBlks_-1,nBlks_-1) block. More...
 
Teuchos::RCP
< Teuchos::ParameterList
paramListMultiphysics_
 
Teuchos::RCP< HierarchyhierarchyMultiphysics_
 
Teuchos::ArrayRCP
< Teuchos::RCP
< Teuchos::ParameterList > > 
arrayOfParamLists_
 
Teuchos::ArrayRCP
< Teuchos::RCP< Hierarchy > > 
arrayOfHierarchies_
 
Teuchos::ArrayRCP
< Teuchos::RCP< Matrix > > 
arrayOfAuxMatrices_
 
Teuchos::ArrayRCP
< Teuchos::RCP< MultiVector > > 
arrayOfNullspaces_
 
Teuchos::ArrayRCP
< Teuchos::RCP
< RealValuedMultiVector > > 
arrayOfCoords_
 
int nBlks_
 
bool useKokkos_
 
bool enable_reuse_
 
bool syncTimers_
 

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::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >

Preconditioner (wrapped as a Xpetra::Operator) for solving MultiPhysics PDEs.

Definition at line 86 of file MueLu_MultiPhys_decl.hpp.

Member Typedef Documentation

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

Definition at line 91 of file MueLu_MultiPhys_decl.hpp.

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

Definition at line 92 of file MueLu_MultiPhys_decl.hpp.

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

Definition at line 93 of file MueLu_MultiPhys_decl.hpp.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 96 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::MultiPhys ( const Teuchos::RCP< Matrix > &  AmatMultiPhysics,
const Teuchos::ArrayRCP< RCP< Matrix >>  arrayOfAuxMatrices,
const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector >>  arrayOfNullspaces,
const Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector >>  arrayOfCoords,
const int  nBlks,
Teuchos::ParameterList List,
bool  ComputePrec = true 
)
inline

Constructor

Parameters
[in]AmatMultiPhysicsMultiphysics discretization matrix
[in]arrayOfAuxMatricesAuxiliary matrices used to generate subblock prolongators for multiphysics system
[in]arrayOfNullspacesNullspace multivectors used to generate subblock prolongators for multiphysics system
[in]arrayOfCoordsCoordinate multivectors used to generate subblock prolongators for multiphysics system
[in]nBlksnBlks x nBlks gives the block dimensions of the multiphysics operator
[in]ListParameter list
[in]ComputePrecIf true, compute the preconditioner immediately

Definition at line 112 of file MueLu_MultiPhys_decl.hpp.

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

Destructor.

Definition at line 129 of file MueLu_MultiPhys_decl.hpp.

Member Function Documentation

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

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

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

Set parameters.

Definition at line 91 of file MueLu_MultiPhys_def.hpp.

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

Setup the preconditioner.

Definition at line 122 of file MueLu_MultiPhys_def.hpp.

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

Reset system matrix.

Definition at line 230 of file MueLu_MultiPhys_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::MultiPhys< 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 242 of file MueLu_MultiPhys_def.hpp.

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

Indicates whether this operator supports applying the adjoint operator.

Definition at line 251 of file MueLu_MultiPhys_def.hpp.

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

Reimplemented from Teuchos::Describable.

Definition at line 278 of file MueLu_MultiPhys_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::MultiPhys< 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 160 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::initialize ( const Teuchos::RCP< Matrix > &  AmatMultiPhysics,
const Teuchos::ArrayRCP< RCP< Matrix >>  arrayOfAuxMatrices,
const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector >>  arrayOfNullspaces,
const Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector >>  arrayOfCoords,
const int  nBlks,
Teuchos::ParameterList List 
)
private

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

Parameters
[in]AmatMultiPhysicsMultiphysics discretization matrix
[in]arrayOfAuxMatricesArray of auxiliary matrices used to generate subblock prolongators for multiphysics system
[in]arrayOfNullspacesArray of nullspace multivectors used to generate subblock prolongators for multiphysics system
[in]arrayOfCoordsArray of coordinate multivectors used to generate subblock prolongators for multiphysics system
[in]nBlksnBlks x nBlks gives the block dimensions of the multiphysics operator
[in]ListParameter list

Definition at line 257 of file MueLu_MultiPhys_def.hpp.

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

apply standard MultiPhys cycle

Definition at line 237 of file MueLu_MultiPhys_def.hpp.

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

allocate multivectors for solve

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

Member Data Documentation

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

ParameterLists.

Definition at line 195 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Matrix> MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::AmatMultiphysics_
private

Hierarchies: used to define P for (0,0)-block, .... (nBlks_-1,nBlks_-1) block.

Definition at line 199 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Teuchos::ParameterList> MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::paramListMultiphysics_
private

Definition at line 200 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP<Hierarchy> MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hierarchyMultiphysics_
private

Definition at line 201 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP<Teuchos::RCP<Teuchos::ParameterList> > MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::arrayOfParamLists_
private

Definition at line 203 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP<Teuchos::RCP<Hierarchy> > MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::arrayOfHierarchies_
private

Definition at line 204 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP<Teuchos::RCP<Matrix> > MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::arrayOfAuxMatrices_
private

Definition at line 205 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP<Teuchos::RCP<MultiVector> > MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::arrayOfNullspaces_
private

Definition at line 206 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP<Teuchos::RCP<RealValuedMultiVector> > MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::arrayOfCoords_
private

Definition at line 207 of file MueLu_MultiPhys_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int MueLu::MultiPhys< Scalar, LocalOrdinal, GlobalOrdinal, Node >::nBlks_
private

Definition at line 209 of file MueLu_MultiPhys_decl.hpp.

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

Definition at line 210 of file MueLu_MultiPhys_decl.hpp.

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

Definition at line 210 of file MueLu_MultiPhys_decl.hpp.

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

Definition at line 210 of file MueLu_MultiPhys_decl.hpp.


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