10#ifndef MUELU_MULTIPHYS_DECL_HPP
11#define MUELU_MULTIPHYS_DECL_HPP
19#include "MueLu_TrilinosSmoother.hpp"
27#include "Xpetra_Map_fwd.hpp"
28#include "Xpetra_Matrix_fwd.hpp"
29#include "Xpetra_MatrixFactory_fwd.hpp"
30#include "Xpetra_MultiVectorFactory_fwd.hpp"
31#include "Xpetra_VectorFactory_fwd.hpp"
32#include "Xpetra_CrsMatrixWrap_fwd.hpp"
47#undef MUELU_MULTIPHYS_SHORT
51 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
magnitudeType;
52 typedef typename Teuchos::ScalarTraits<Scalar>::coordinateType
coordinateType;
74 MultiPhys(
const Teuchos::RCP<Matrix>& AmatMultiPhysics,
75 const Teuchos::ArrayRCP<RCP<Matrix>> arrayOfAuxMatrices,
76 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfNullspaces,
77 const Teuchos::ArrayRCP<Teuchos::RCP<RealValuedMultiVector>> arrayOfCoords,
79 Teuchos::ParameterList& List,
80 bool ComputePrec =
true,
81 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfMaterials = Teuchos::null,
82 bool OmitSubblockSmoother =
true)
90 initialize(AmatMultiPhysics, arrayOfAuxMatrices, arrayOfNullspaces, arrayOfCoords, nBlks, List, arrayOfMaterials);
107 void compute(
bool reuse =
false);
110 void resetMatrix(Teuchos::RCP<Matrix> SM_Matrix_new,
bool ComputePrec =
true);
115 void apply(
const MultiVector& X, MultiVector& Y,
116 Teuchos::ETransp mode = Teuchos::NO_TRANS,
117 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
118 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero())
const;
123 void describe(Teuchos::FancyOStream& out,
const Teuchos::EVerbosityLevel verbLevel = Teuchos::VERB_HIGH)
const;
127 const MultiVector& B,
128 MultiVector& R)
const {
129 using STS = Teuchos::ScalarTraits<Scalar>;
130 R.update(STS::one(), B, STS::zero());
131 this->
apply(X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
155 void initialize(
const Teuchos::RCP<Matrix>& AmatMultiPhysics,
156 const Teuchos::ArrayRCP<RCP<Matrix>> arrayOfAuxMatrices,
157 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfNullspaces,
158 const Teuchos::ArrayRCP<Teuchos::RCP<RealValuedMultiVector>> arrayOfCoords,
160 Teuchos::ParameterList& List,
161 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfMaterials);
164 void applyInverse(
const MultiVector& RHS, MultiVector& X)
const;
196#define MUELU_MULTIPHYS_SHORT
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Preconditioner (wrapped as a Xpetra::Operator) for solving MultiPhysics PDEs.
Teuchos::RCP< Hierarchy > hierarchyMultiphysics_
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfNullspaces_
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, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfMaterials=Teuchos::null, bool OmitSubblockSmoother=true)
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
get a (synced) timer
Teuchos::ArrayRCP< Teuchos::RCP< Matrix > > arrayOfAuxMatrices_
void compute(bool reuse=false)
Setup the preconditioner.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
Teuchos::ParameterList parameterList_
ParameterLists.
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
void allocateMemory(int numVectors) const
allocate multivectors for solve
Teuchos::ScalarTraits< Scalar >::coordinateType coordinateType
void applyInverse(const MultiVector &RHS, MultiVector &X) const
apply standard MultiPhys cycle
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, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfMaterials)
Xpetra::MultiVector< coordinateType, LO, GO, NO > RealValuedMultiVector
void setParameters(Teuchos::ParameterList &list)
Set parameters.
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
Teuchos::RCP< Matrix > AmatMultiphysics_
Hierarchies: used to define P for (0,0)-block, .... (nBlks_-1,nBlks_-1) block.
Teuchos::ArrayRCP< Teuchos::RCP< Hierarchy > > subblockHierarchies()
Return an array of hierarchies corresponding to each diagonal subblock.
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
virtual ~MultiPhys()
Destructor.
Teuchos::RCP< Hierarchy > multiphysicsHierarchy()
Return the multiphysics hierarchy.
Teuchos::ArrayRCP< Teuchos::RCP< Hierarchy > > arrayOfHierarchies_
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Teuchos::RCP< Teuchos::ParameterList > paramListMultiphysics_
Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector > > arrayOfCoords_
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfMaterials_
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::ArrayRCP< Teuchos::RCP< Teuchos::ParameterList > > arrayOfParamLists_
bool OmitSubblockSmoother_
Verbose class for MueLu classes.
Namespace for MueLu classes and methods.