10 #ifndef MUELU_MAXWELL1_DECL_HPP
11 #define MUELU_MAXWELL1_DECL_HPP
46 #undef MUELU_MAXWELL1_SHORT
84 bool ComputePrec =
true)
87 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, List);
106 bool ComputePrec =
true)
108 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, List);
128 bool ComputePrec =
true)
130 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, List);
145 bool ComputePrec =
true)
154 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, List);
156 if (SM_Matrix != Teuchos::null)
178 void compute(
bool reuse =
false);
201 R.
update(STS::one(), B, STS::zero());
236 void dump(
const Matrix& A, std::string name)
const;
248 void dump(
const Kokkos::View<bool*, typename Node::device_type>& v, std::string name)
const;
288 #define MUELU_MAXWELL1_SHORT
289 #endif // MUELU_MAXWELL1_DECL_HPP
const Teuchos::RCP< Matrix > & getJacobian() const
Returns Jacobian matrix SM.
Teuchos::RCP< MultiVector > Nullspace_
Nullspace.
Teuchos::RCP< Matrix > GmhdA_Matrix_
Xpetra::MultiVector< coordinateType, LO, GO, NO > RealValuedMultiVector
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
MueLu::DefaultLocalOrdinal LocalOrdinal
void dump(const Matrix &A, std::string name) const
dump out matrix
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)
Teuchos::RCP< Hierarchy > Hierarchy11_
Two hierarchies: one for the (1,1)-block, another for the (2,2)-block.
Kokkos::View< bool *, typename Node::device_type::memory_space > BCdomainKokkos_
Teuchos::RCP< Hierarchy > HierarchyGmhd_
T & get(const std::string &name, T def_value)
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
bool useKokkos_
Some options.
void GMHDSetupHierarchy(Teuchos::ParameterList &List) const
Sets up hiearchy for GMHD matrices that include generalized Ohms law equations.
Teuchos::ArrayRCP< bool > BCcols_
Teuchos::RCP< Matrix > D0_Matrix_
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)
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
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)
Teuchos::RCP< MultiVector > residual11c_
Teuchos::ParameterList precList22_
Kokkos::View< bool *, typename Node::device_type::memory_space > BCcolsKokkos_
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
Teuchos::RCP< MultiVector > update22_
Verbose class for MueLu classes.
virtual void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)=0
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form...
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
Kokkos::View< bool *, typename Node::device_type::memory_space > BCrowsKokkos_
Vectors for BCs.
void setParameters(Teuchos::ParameterList &list)
Set parameters.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::ArrayRCP< bool > BCdomain_
Teuchos::ParameterList precList11_
Teuchos::RCP< Hierarchy > Hierarchy22_
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Teuchos::ScalarTraits< Scalar >::coordinateType coordinateType
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
get a (synced) timer
mode_type
Execution modes.
void applyInverseRefMaxwellAdditive(const MultiVector &RHS, MultiVector &X) const
apply RefMaxwell additive 2x2 style cycle
void allocateMemory(int numVectors) const
allocate multivectors for solve
Teuchos::RCP< Matrix > Kn_Matrix_
RCP< Matrix > P11_
Temporary memory (cached vectors for RefMaxwell-style)
Teuchos::ArrayRCP< bool > BCrows_
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 applyInverseStandard(const MultiVector &RHS, MultiVector &X) const
apply standard Maxwell1 cycle
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
Maxwell1(Teuchos::RCP< Hierarchy > H11, Teuchos::RCP< Hierarchy > H22)
Constructor with Hierarchies.
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)
void compute(bool reuse=false)
Setup the preconditioner.
bool isType(const std::string &name) const
Teuchos::RCP< MultiVector > update11c_
Maxwell1(const Teuchos::RCP< Matrix > &SM_Matrix, Teuchos::ParameterList &List, bool ComputePrec=true)
virtual ~Maxwell1()
Destructor.
Teuchos::RCP< MultiVector > residual22_
Teuchos::RCP< Matrix > generate_kn() const
Generates the Kn matrix.
Teuchos::RCP< RealValuedMultiVector > Coords_
Coordinates.
void dumpCoords(const RealValuedMultiVector &X, std::string name) const
dump out real-valued multivector
Teuchos::RCP< Matrix > SM_Matrix_
Various matrices.
Teuchos::RCP< MultiVector > residualFine_
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Teuchos::ParameterList parameterList_
ParameterLists.