46 #ifndef MUELU_REFMAXWELL_DECL_HPP
47 #define MUELU_REFMAXWELL_DECL_HPP
62 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
80 #include "MueLu_TrilinosSmoother.hpp"
81 #include "MueLu_Hierarchy.hpp"
83 #include "Xpetra_Map_fwd.hpp"
84 #include "Xpetra_Matrix_fwd.hpp"
85 #include "Xpetra_MatrixFactory_fwd.hpp"
86 #include "Xpetra_MultiVectorFactory_fwd.hpp"
87 #include "Xpetra_VectorFactory_fwd.hpp"
88 #include "Xpetra_CrsMatrixWrap_fwd.hpp"
90 #if defined(HAVE_MUELU_IFPACK2) && (!defined(HAVE_MUELU_EPETRA))
91 #define MUELU_REFMAXWELL_CAN_USE_HIPTMAIR
93 #include "Ifpack2_Hiptmair.hpp"
128 class RefMaxwell :
public VerboseObject,
public Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> {
130 #undef MUELU_REFMAXWELL_SHORT
177 bool ComputePrec =
true)
179 initialize(D0_Matrix,Ms_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
201 bool ComputePrec =
true)
203 initialize(D0_Matrix,M1_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
223 initialize(D0_Matrix,M1_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
244 initialize(D0_Matrix,M1_Matrix,Teuchos::null,M1_Matrix,Nullspace,Coords,List);
262 initialize(D0_Matrix,M1_Matrix,Teuchos::null,M1_Matrix,Nullspace,Coords,List);
273 bool ComputePrec =
true)
287 initialize(D0_Matrix,Ms_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
289 if (SM_Matrix != Teuchos::null)
311 void compute(
bool reuse=
false);
325 void apply (
const MultiVector& X, MultiVector& Y,
337 const MultiVector & B,
338 MultiVector & R)
const {
340 R.update(STS::one(),B,STS::zero());
375 void solveH(
const MultiVector& RHS, MultiVector& X)
const;
378 void solve22(
const MultiVector& RHS, MultiVector& X)
const;
384 void dump(
const Matrix& A, std::string name)
const;
387 void dump(
const MultiVector& X, std::string name)
const;
395 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
407 xpImporter->setDistributorParameters(matvecParams);
410 xpExporter->setDistributorParameters(matvecParams);
416 #if defined(MUELU_REFMAXWELL_CAN_USE_HIPTMAIR)
424 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
450 #define MUELU_REFMAXWELL_SHORT
451 #endif // MUELU_REFMAXWELL_DECL_HPP
const Teuchos::RCP< Matrix > & getJacobian() const
Returns Jacobian matrix SM.
Teuchos::RCP< SmootherBase > PostSmoother_
Teuchos::RCP< const Import > Importer22_
void setMatvecParams(Matrix &A, RCP< ParameterList > matvecParams)
set parameters
Teuchos::RCP< Matrix > D0_T_Matrix_
MueLu::DefaultLocalOrdinal LocalOrdinal
void solve22(const MultiVector &RHS, MultiVector &X) const
apply solve to 2-2 block only
Teuchos::RCP< Teuchos::ParameterList > A22_RAP_reuse_data_
int BCedges_
Vectors for BCs.
Teuchos::RCP< Matrix > Ms_Matrix_
Teuchos::RCP< MultiVector > D0xTmp_
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)
Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
Teuchos::RCP< Matrix > D0_Matrix_
T & get(const std::string &name, T def_value)
Teuchos::RCP< Matrix > SM_Matrix_
Various matrices.
Teuchos::RCP< MultiVector > D0res_
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)
Teuchos::ParameterList smootherList_
Teuchos::RCP< Matrix > A22_
Teuchos::RCP< MultiVector > D0resTmp_
Teuchos::RCP< MultiVector > D0TR11Tmp_
void solveH(const MultiVector &RHS, MultiVector &X) const
apply solve to 1-1 block only
Teuchos::RCP< Matrix > Addon_Matrix_
bool useHiptmairSmoothing_
void setParameters(Teuchos::ParameterList &list)
Set parameters.
void dump(const Matrix &A, std::string name) const
dump out matrix
void buildProlongator()
Setup the prolongator for the (1,1)-block.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
Teuchos::RCP< MultiVector > P11res_
Temporary memory.
Static class that holds the complete list of valid MueLu parameters.
Teuchos::RCP< Matrix > M0inv_Matrix_
void allocateMemory(int numVectors) const
allocate multivectors for solve
Teuchos::RCP< MultiVector > residual_
Teuchos::RCP< MultiVector > D0x_
void applyInverseAdditive(const MultiVector &RHS, MultiVector &X) const
apply additive algorithm for 2x2 solve
Teuchos::RCP< MultiVector > P11resTmp_
MueLu::DefaultScalar Scalar
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Teuchos::RCP< Teuchos::ParameterList > AH_AP_reuse_data_
Teuchos::ParameterList precList11_
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)
Teuchos::RCP< Matrix > P11_
Teuchos::RCP< Matrix > A_nodal_Matrix_
Teuchos::RCP< Hierarchy > Hierarchy22_
Teuchos::RCP< Teuchos::ParameterList > AH_RAP_reuse_data_
Teuchos::RCP< Matrix > M1_Matrix_
bool fuseProlongationAndUpdate_
Teuchos::RCP< RealValuedMultiVector > CoordsH_
Teuchos::ArrayRCP< bool > BCrows_
Teuchos::ArrayRCP< bool > BCcols_
Xpetra::MultiVector< coordinateType, LO, GO, NO > RealValuedMultiVector
void applyInverse212(const MultiVector &RHS, MultiVector &X) const
apply 2-1-2 algorithm for 2x2 solve
Teuchos::RCP< SmootherBase > PreSmoother_
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
get a (synced) timer
void applyInverse121(const MultiVector &RHS, MultiVector &X) const
apply 1-2-1 algorithm for 2x2 solve
void compute(bool reuse=false)
Setup the preconditioner.
Teuchos::ArrayRCP< bool > BCdomain_
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
bool disable_addon_
Some options.
void initialize(const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Ms_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)
virtual ~RefMaxwell()
Destructor.
Teuchos::ParameterList parameterList_
Parameter lists.
Teuchos::RCP< RealValuedMultiVector > Coords_
Coordinates.
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)
void dumpCoords(const RealValuedMultiVector &X, std::string name) const
dump out real-valued multivector
bool use_as_preconditioner_
bool isType(const std::string &name) const
Teuchos::RCP< MultiVector > P11x_
Teuchos::RCP< Matrix > R11_
Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
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
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
Teuchos::RCP< Matrix > AH_
void formCoarseMatrix()
Compute P11^{T}*A*P11 efficiently.
Teuchos::RCP< MultiVector > P11xTmp_
bool D0_T_R11_colMapsMatch_
Teuchos::ScalarTraits< Scalar >::coordinateType coordinateType
Teuchos::RCP< Hierarchy > HierarchyH_
Two hierarchies: one for the coarse (1,1)-block, another for the (2,2)-block.
RefMaxwell(Teuchos::RCP< Hierarchy > HH, Teuchos::RCP< Hierarchy > H22)
Constructor with Hierarchies.
Teuchos::RCP< Teuchos::ParameterList > A22_AP_reuse_data_
Teuchos::ParameterList precList22_
Teuchos::RCP< MultiVector > Nullspace_
Nullspace.
RefMaxwell(const Teuchos::RCP< Matrix > &SM_Matrix, Teuchos::ParameterList &List, bool ComputePrec=true)
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Teuchos::RCP< const Import > ImporterH_
Importer to coarse (1,1) hierarchy.
RefMaxwell(const Teuchos::RCP< Matrix > &SM_Matrix, const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Ms_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)