MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Maxwell1_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_MAXWELL1_DECL_HPP
11#define MUELU_MAXWELL1_DECL_HPP
12
13#include "MueLu_ConfigDefs.hpp"
15
17
19#include "MueLu_Level_fwd.hpp"
25
26#include "Xpetra_Operator.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"
33
34namespace MueLu {
35
42template <class Scalar,
43 class LocalOrdinal,
44 class GlobalOrdinal,
45 class Node>
46class Maxwell1 : public VerboseObject, public Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
47#undef MUELU_MAXWELL1_SHORT
49
50 public:
51 using magnitudeType = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
52 using coordinateType = typename Teuchos::ScalarTraits<Scalar>::coordinateType;
53 using RealValuedMultiVector = typename Xpetra::MultiVector<coordinateType, LO, GO, NO>;
54
57 : Hierarchy11_(Teuchos::null)
58 , Hierarchy22_(Teuchos::null)
59 , HierarchyGmhd_(Teuchos::null)
61 }
62
64 Maxwell1(Teuchos::RCP<Hierarchy> H11, Teuchos::RCP<Hierarchy> H22)
65 : Hierarchy11_(H11)
66 , Hierarchy22_(H22)
67 , HierarchyGmhd_(Teuchos::null)
69 }
70
80 Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
81 const Teuchos::RCP<Matrix>& D0_Matrix,
82 const Teuchos::RCP<MultiVector>& Nullspace,
83 const Teuchos::RCP<RealValuedMultiVector>& Coords,
84 Teuchos::ParameterList& List,
85 bool ComputePrec = true)
87 RCP<Matrix> Kn_Matrix;
88 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, Teuchos::null, List);
89 resetMatrix(SM_Matrix, ComputePrec);
90 }
91
101 Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
102 const Teuchos::RCP<Matrix>& D0_Matrix,
103 const Teuchos::RCP<Matrix>& Kn_Matrix,
104 const Teuchos::RCP<MultiVector>& Nullspace,
105 const Teuchos::RCP<RealValuedMultiVector>& Coords,
106 Teuchos::ParameterList& List,
107 bool ComputePrec = true)
109 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, Teuchos::null, List);
110 resetMatrix(SM_Matrix, ComputePrec);
111 }
112
122 Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
123 const Teuchos::RCP<Matrix>& D0_Matrix,
124 const Teuchos::RCP<Matrix>& Kn_Matrix,
125 const Teuchos::RCP<MultiVector>& Nullspace,
126 const Teuchos::RCP<RealValuedMultiVector>& Coords,
127 const Teuchos::RCP<MultiVector>& Material,
128 Teuchos::ParameterList& List,
129 bool ComputePrec = true)
131 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, Material, List);
132 resetMatrix(SM_Matrix, ComputePrec);
133 }
134
145 Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
146 const Teuchos::RCP<Matrix>& D0_Matrix,
147 const Teuchos::RCP<Matrix>& Kn_Matrix,
148 const Teuchos::RCP<MultiVector>& Nullspace,
149 const Teuchos::RCP<RealValuedMultiVector>& Coords,
150 Teuchos::ParameterList& List, const Teuchos::RCP<Matrix>& GmhdA_Matrix,
151 bool ComputePrec = true)
153 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, Teuchos::null, List);
154 resetMatrix(SM_Matrix, ComputePrec);
155 GmhdA_Matrix_ = GmhdA_Matrix;
156 HierarchyGmhd_ = rcp(new Hierarchy("HierarchyGmhd"));
157 GMHDSetupHierarchy(List);
158 }
159
166 Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
167 Teuchos::ParameterList& List,
168 bool ComputePrec = true)
170 RCP<MultiVector> Nullspace = List.get<RCP<MultiVector> >("Nullspace", Teuchos::null);
171 RCP<RealValuedMultiVector> Coords = List.get<RCP<RealValuedMultiVector> >("Coordinates", Teuchos::null);
172 RCP<Matrix> D0_Matrix = List.get<RCP<Matrix> >("D0");
173 RCP<Matrix> Kn_Matrix;
174 if (List.isType<RCP<Matrix> >("Kn"))
175 Kn_Matrix = List.get<RCP<Matrix> >("Kn");
176
177 initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, Teuchos::null, List);
178
179 if (SM_Matrix != Teuchos::null)
180 resetMatrix(SM_Matrix, ComputePrec);
181 }
182
184 virtual ~Maxwell1() = default;
185
187 const Teuchos::RCP<const Map> getDomainMap() const;
188
190 const Teuchos::RCP<const Map> getRangeMap() const;
191
193 const Teuchos::RCP<Matrix>& getJacobian() const {
194 return SM_Matrix_;
195 }
196
198 void setParameters(Teuchos::ParameterList& list);
199
201 void compute(bool reuse = false);
202
204 void resetMatrix(Teuchos::RCP<Matrix> SM_Matrix_new, bool ComputePrec = true);
205
209 void apply(const MultiVector& X, MultiVector& Y,
210 Teuchos::ETransp mode = Teuchos::NO_TRANS,
211 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
212 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero()) const;
213
215 bool hasTransposeApply() const;
216
217 void describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel = Teuchos::VERB_HIGH) const;
218
220 void residual(const MultiVector& X,
221 const MultiVector& B,
222 MultiVector& R) const {
223 using STS = Teuchos::ScalarTraits<Scalar>;
224 R.update(STS::one(), B, STS::zero());
225 this->apply(X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
226 }
227
228 private:
230 Teuchos::RCP<Matrix> generate_kn() const;
231
233 void GMHDSetupHierarchy(Teuchos::ParameterList& List) const;
234
244 void initialize(const Teuchos::RCP<Matrix>& D0_Matrix,
245 const Teuchos::RCP<Matrix>& Kn_Matrix,
246 const Teuchos::RCP<MultiVector>& Nullspace,
247 const Teuchos::RCP<RealValuedMultiVector>& Coords,
248 const Teuchos::RCP<MultiVector>& Material,
249 Teuchos::ParameterList& List);
250
252 void applyInverseRefMaxwellAdditive(const MultiVector& RHS, MultiVector& X) const;
253
255 void applyInverseStandard(const MultiVector& RHS, MultiVector& X) const;
256
258 void allocateMemory(int numVectors) const;
259
261 void dump(const Matrix& A, std::string name) const;
262
264 void dump(const MultiVector& X, std::string name) const;
265
267 void dumpCoords(const RealValuedMultiVector& X, std::string name) const;
268
270 void dump(const Teuchos::ArrayRCP<bool>& v, std::string name) const;
271
273 void dump(const Kokkos::View<bool*, typename Node::device_type>& v, std::string name) const;
274
276 Teuchos::RCP<Teuchos::TimeMonitor> getTimer(std::string name, RCP<const Teuchos::Comm<int> > comm = Teuchos::null) const;
277
279 mutable Teuchos::ParameterList parameterList_, precList11_, precList22_;
280
282 Teuchos::RCP<Hierarchy> Hierarchy11_, Hierarchy22_, HierarchyGmhd_;
283
286
288 Kokkos::View<bool*, typename Node::device_type::memory_space> BCrowsKokkos_, BCcolsKokkos_, BCdomainKokkos_;
290 Teuchos::ArrayRCP<bool> BCrows_, BCcols_, BCdomain_;
292 Teuchos::RCP<MultiVector> Nullspace_;
294 Teuchos::RCP<RealValuedMultiVector> Coords_;
296 Teuchos::RCP<MultiVector> Material_;
300
302 typedef enum { MODE_STANDARD = 0,
307
309 RCP<Matrix> P11_;
310 mutable Teuchos::RCP<MultiVector> residualFine_, residual11c_, residual22_, update11c_, update22_;
311};
312
313} // namespace MueLu
314
315#define MUELU_MAXWELL1_SHORT
316#endif // MUELU_MAXWELL1_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form.
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
bool useKokkos_
Some options.
Teuchos::RCP< Matrix > generate_kn() const
Generates the Kn matrix.
Teuchos::RCP< Hierarchy > Hierarchy11_
Two hierarchies: one for the (1,1)-block, another for the (2,2)-block.
Teuchos::RCP< Hierarchy > HierarchyGmhd_
Teuchos::RCP< MultiVector > residualFine_
typename Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
Teuchos::ArrayRCP< bool > BCcols_
Teuchos::ArrayRCP< bool > BCrows_
Teuchos::RCP< Matrix > SM_Matrix_
Various matrices.
Teuchos::RCP< Hierarchy > Hierarchy22_
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)
Kokkos::View< bool *, typename Node::device_type::memory_space > BCrowsKokkos_
Vectors for BCs.
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)
void compute(bool reuse=false)
Setup the preconditioner.
void GMHDSetupHierarchy(Teuchos::ParameterList &List) const
Sets up hiearchy for GMHD matrices that include generalized Ohms law equations.
mode_type
Execution modes.
typename Teuchos::ScalarTraits< Scalar >::coordinateType coordinateType
void initialize(const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &Kn_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< RealValuedMultiVector > &Coords, const Teuchos::RCP< MultiVector > &Material, Teuchos::ParameterList &List)
Teuchos::ParameterList parameterList_
ParameterLists.
Teuchos::RCP< MultiVector > residual22_
Kokkos::View< bool *, typename Node::device_type::memory_space > BCdomainKokkos_
Teuchos::RCP< RealValuedMultiVector > Coords_
Coordinates.
Teuchos::RCP< Matrix > D0_Matrix_
virtual ~Maxwell1()=default
Destructor.
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
void applyInverseRefMaxwellAdditive(const MultiVector &RHS, MultiVector &X) const
apply RefMaxwell additive 2x2 style cycle
Teuchos::RCP< MultiVector > Nullspace_
Nullspace.
Teuchos::RCP< MultiVector > update11c_
Maxwell1(Teuchos::RCP< Hierarchy > H11, Teuchos::RCP< Hierarchy > H22)
Constructor with Hierarchies.
void applyInverseStandard(const MultiVector &RHS, MultiVector &X) const
apply standard Maxwell1 cycle
RCP< Matrix > P11_
Temporary memory (cached vectors for RefMaxwell-style)
Teuchos::RCP< MultiVector > Material_
Material.
void dumpCoords(const RealValuedMultiVector &X, std::string name) const
dump out real-valued multivector
Teuchos::RCP< MultiVector > update22_
Teuchos::ArrayRCP< bool > BCdomain_
Teuchos::ParameterList precList11_
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 residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
const Teuchos::RCP< Matrix > & getJacobian() const
Returns Jacobian matrix SM.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
Teuchos::RCP< Matrix > GmhdA_Matrix_
void dump(const Matrix &A, std::string name) const
dump out matrix
typename Xpetra::MultiVector< coordinateType, LO, GO, NO > RealValuedMultiVector
Maxwell1(const Teuchos::RCP< Matrix > &SM_Matrix, Teuchos::ParameterList &List, bool ComputePrec=true)
Teuchos::ParameterList precList22_
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
Teuchos::RCP< Matrix > Kn_Matrix_
Teuchos::RCP< MultiVector > residual11c_
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
get a (synced) timer
Kokkos::View< bool *, typename Node::device_type::memory_space > BCcolsKokkos_
void setParameters(Teuchos::ParameterList &list)
Set parameters.
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, const Teuchos::RCP< MultiVector > &Material, Teuchos::ParameterList &List, bool ComputePrec=true)
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 allocateMemory(int numVectors) const
allocate multivectors for solve
Verbose class for MueLu classes.
Namespace for MueLu classes and methods.