10 #ifndef MUELU_MULTIPHYS_DECL_HPP
11 #define MUELU_MULTIPHYS_DECL_HPP
19 #include "MueLu_TrilinosSmoother.hpp"
47 #undef MUELU_MULTIPHYS_SHORT
78 bool ComputePrec =
true)
84 initialize(AmatMultiPhysics, arrayOfAuxMatrices, arrayOfNullspaces, arrayOfCoords, nBlks, List);
101 void compute(
bool reuse =
false);
124 R.
update(STS::one(), B, STS::zero());
175 #define MUELU_MULTIPHYS_SHORT
176 #endif // MUELU_MULTIPHYS_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
Preconditioner (wrapped as a Xpetra::Operator) for solving MultiPhysics PDEs.
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
Xpetra::MultiVector< coordinateType, LO, GO, NO > RealValuedMultiVector
Teuchos::ArrayRCP< Teuchos::RCP< Matrix > > arrayOfAuxMatrices_
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
virtual ~MultiPhys()
Destructor.
Verbose class for MueLu classes.
void setParameters(Teuchos::ParameterList &list)
Set parameters.
virtual void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)=0
Teuchos::ScalarTraits< Scalar >::coordinateType coordinateType
Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector > > arrayOfCoords_
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)
void applyInverse(const MultiVector &RHS, MultiVector &X) const
apply standard MultiPhys cycle
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
Teuchos::RCP< Matrix > AmatMultiphysics_
Hierarchies: used to define P for (0,0)-block, .... (nBlks_-1,nBlks_-1) block.
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Teuchos::ArrayRCP< Teuchos::RCP< Teuchos::ParameterList > > arrayOfParamLists_
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int >> comm=Teuchos::null) const
get a (synced) timer
Teuchos::ParameterList parameterList_
ParameterLists.
Teuchos::RCP< Hierarchy > hierarchyMultiphysics_
Teuchos::ArrayRCP< Teuchos::RCP< Hierarchy > > arrayOfHierarchies_
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
Teuchos::RCP< Teuchos::ParameterList > paramListMultiphysics_
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
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 resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
void compute(bool reuse=false)
Setup the preconditioner.
Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfNullspaces_
void allocateMemory(int numVectors) const
allocate multivectors for solve