46 #ifndef MUELU_REFMAXWELL_DECL_HPP 
   47 #define MUELU_REFMAXWELL_DECL_HPP 
   61 #ifdef HAVE_MUELU_KOKKOS_REFACTOR 
   78 #include "MueLu_TrilinosSmoother.hpp" 
   79 #include "MueLu_Hierarchy.hpp" 
   88 #if defined(HAVE_MUELU_IFPACK2) && (!defined(HAVE_MUELU_EPETRA) || defined(HAVE_MUELU_INST_DOUBLE_INT_INT)) 
   90 #include "Ifpack2_Hiptmair.hpp" 
  125   class RefMaxwell : 
public VerboseObject, 
public Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> {
 
  127 #undef MUELU_REFMAXWELL_SHORT 
  173                bool ComputePrec = 
true)
 
  175       initialize(D0_Matrix,Ms_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
 
  197                bool ComputePrec = 
true)
 
  199       initialize(D0_Matrix,M1_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
 
  219       initialize(D0_Matrix,M1_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
 
  240       initialize(D0_Matrix,M1_Matrix,Teuchos::null,M1_Matrix,Nullspace,Coords,List);
 
  258       initialize(D0_Matrix,M1_Matrix,Teuchos::null,M1_Matrix,Nullspace,Coords,List);
 
  269                bool ComputePrec = 
true)
 
  276       if (List.
isType<Matrix>(
"Ms"))
 
  283       initialize(D0_Matrix,Ms_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
 
  285       if (SM_Matrix != Teuchos::null)
 
  307     void compute(
bool reuse=
false);
 
  329 #ifdef HAVE_MUELU_DEPRECATED_CODE 
  330     template <
class NewNode>
 
  380 #if defined(HAVE_MUELU_IFPACK2) && (!defined(HAVE_MUELU_EPETRA) || defined(HAVE_MUELU_INST_DOUBLE_INT_INT)) 
  388 #ifdef HAVE_MUELU_KOKKOS_REFACTOR 
  414 #define MUELU_REFMAXWELL_SHORT 
  415 #endif // MUELU_REFMAXWELL_DECL_HPP 
const Teuchos::RCP< Matrix > & getJacobian() const 
Returns Jacobian matrix SM. 
 
Teuchos::RCP< SmootherBase > PostSmoother_
 
Teuchos::RCP< const Import > Importer22_
 
RCP< Map< LocalOrdinal, GlobalOrdinal, Node2 > > clone(const Map< LocalOrdinal, GlobalOrdinal, Node1 > &map, const RCP< Node2 > &node2)
 
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::ArrayRCP< const bool > BCrows_
 
Teuchos::RCP< Teuchos::ParameterList > A22_RAP_reuse_data_
 
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. 
 
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MultiVector
 
Teuchos::RCP< MultiVector > D0res_
 
int BCrowcount_
Vectors for BCs. 
 
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::ArrayRCP< const bool > BCcols_
 
Teuchos::RCP< Matrix > A22_
 
Teuchos::RCP< MultiVector > D0resTmp_
 
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 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. 
 
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form...
 
Xpetra::MultiVector< magnitudeType, LO, GO, NO > RealValuedMultiVector
 
Teuchos::RCP< Matrix > M0inv_Matrix_
 
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_
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
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_
 
Teuchos::RCP< RealValuedMultiVector > CoordsH_
 
void applyInverse212(const MultiVector &RHS, MultiVector &X) const 
apply 2-1-2 algorithm for 2x2 solve 
 
Teuchos::RCP< SmootherBase > PreSmoother_
 
void applyInverse121(const MultiVector &RHS, MultiVector &X) const 
apply 1-2-1 algorithm for 2x2 solve 
 
void compute(bool reuse=false)
Setup the preconditioner. 
 
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)
 
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 
 
Teuchos::RCP< Matrix > AH_
 
void formCoarseMatrix()
Compute P11^{T}*A*P11 efficiently. 
 
Teuchos::RCP< MultiVector > P11xTmp_
 
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)