MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_Maxwell1_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_MAXWELL1_DECL_HPP
47 #define MUELU_MAXWELL1_DECL_HPP
48 
49 #include "MueLu_ConfigDefs.hpp"
50 #include "MueLu_BaseClass.hpp"
51 
53 
54 #include "MueLu_Utilities_fwd.hpp"
55 #include "MueLu_Level_fwd.hpp"
56 #include "MueLu_Hierarchy_fwd.hpp"
57 #include "MueLu_RAPFactory_fwd.hpp"
58 #include "MueLu_PerfUtils_fwd.hpp"
60 
61 #include "Xpetra_Map_fwd.hpp"
62 #include "Xpetra_Matrix_fwd.hpp"
67 
68 namespace MueLu {
69 
76 template <class Scalar,
77  class LocalOrdinal,
78  class GlobalOrdinal,
79  class Node>
80 class Maxwell1 : public VerboseObject, public Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
81 #undef MUELU_MAXWELL1_SHORT
82 #include "MueLu_UseShortNames.hpp"
83 
84  public:
88 
91  : Hierarchy11_(Teuchos::null)
92  , Hierarchy22_(Teuchos::null)
93  , HierarchyGmhd_(Teuchos::null)
94  , mode_(MODE_STANDARD) {
95  }
96 
99  : Hierarchy11_(H11)
100  , Hierarchy22_(H22)
101  , HierarchyGmhd_(Teuchos::null)
102  , mode_(MODE_STANDARD) {
103  }
104 
114  Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
115  const Teuchos::RCP<Matrix>& D0_Matrix,
116  const Teuchos::RCP<MultiVector>& Nullspace,
119  bool ComputePrec = true)
120  : mode_(MODE_STANDARD) {
121  RCP<Matrix> Kn_Matrix;
122  initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, List);
123  resetMatrix(SM_Matrix, ComputePrec);
124  }
125 
135  Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
136  const Teuchos::RCP<Matrix>& D0_Matrix,
137  const Teuchos::RCP<Matrix>& Kn_Matrix,
138  const Teuchos::RCP<MultiVector>& Nullspace,
141  bool ComputePrec = true)
142  : mode_(MODE_STANDARD) {
143  initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, List);
144  resetMatrix(SM_Matrix, ComputePrec);
145  }
146 
157  Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
158  const Teuchos::RCP<Matrix>& D0_Matrix,
159  const Teuchos::RCP<Matrix>& Kn_Matrix,
160  const Teuchos::RCP<MultiVector>& Nullspace,
162  Teuchos::ParameterList& List, const Teuchos::RCP<Matrix>& GmhdA_Matrix,
163  bool ComputePrec = true)
165  initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, List);
166  resetMatrix(SM_Matrix, ComputePrec);
167  GmhdA_Matrix_ = GmhdA_Matrix;
168  HierarchyGmhd_ = rcp(new Hierarchy("HierarchyGmhd"));
169  GMHDSetupHierarchy(List);
170  }
171 
178  Maxwell1(const Teuchos::RCP<Matrix>& SM_Matrix,
180  bool ComputePrec = true)
181  : mode_(MODE_STANDARD) {
182  RCP<MultiVector> Nullspace = List.get<RCP<MultiVector> >("Nullspace", Teuchos::null);
183  RCP<RealValuedMultiVector> Coords = List.get<RCP<RealValuedMultiVector> >("Coordinates", Teuchos::null);
184  RCP<Matrix> D0_Matrix = List.get<RCP<Matrix> >("D0");
185  RCP<Matrix> Kn_Matrix;
186  if (List.isType<RCP<Matrix> >("Kn"))
187  Kn_Matrix = List.get<RCP<Matrix> >("Kn");
188 
189  initialize(D0_Matrix, Kn_Matrix, Nullspace, Coords, List);
190 
191  if (SM_Matrix != Teuchos::null)
192  resetMatrix(SM_Matrix, ComputePrec);
193  }
194 
196  virtual ~Maxwell1() {}
197 
200 
202  const Teuchos::RCP<const Map> getRangeMap() const;
203 
206  return SM_Matrix_;
207  }
208 
211 
213  void compute(bool reuse = false);
214 
216  void resetMatrix(Teuchos::RCP<Matrix> SM_Matrix_new, bool ComputePrec = true);
217 
221  void apply(const MultiVector& X, MultiVector& Y,
225 
227  bool hasTransposeApply() const;
228 
230 
232  void residual(const MultiVector& X,
233  const MultiVector& B,
234  MultiVector& R) const {
235  using STS = Teuchos::ScalarTraits<Scalar>;
236  R.update(STS::one(), B, STS::zero());
237  this->apply(X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
238  }
239 
240  private:
243 
245  void GMHDSetupHierarchy(Teuchos::ParameterList& List) const;
246 
255  void initialize(const Teuchos::RCP<Matrix>& D0_Matrix,
256  const Teuchos::RCP<Matrix>& Kn_Matrix,
257  const Teuchos::RCP<MultiVector>& Nullspace,
259  Teuchos::ParameterList& List);
260 
262  void applyInverseRefMaxwellAdditive(const MultiVector& RHS, MultiVector& X) const;
263 
265  void applyInverseStandard(const MultiVector& RHS, MultiVector& X) const;
266 
268  void allocateMemory(int numVectors) const;
269 
271  void dump(const Matrix& A, std::string name) const;
272 
274  void dump(const MultiVector& X, std::string name) const;
275 
277  void dumpCoords(const RealValuedMultiVector& X, std::string name) const;
278 
280  void dump(const Teuchos::ArrayRCP<bool>& v, std::string name) const;
281 
283  void dump(const Kokkos::View<bool*, typename Node::device_type>& v, std::string name) const;
284 
286  Teuchos::RCP<Teuchos::TimeMonitor> getTimer(std::string name, RCP<const Teuchos::Comm<int> > comm = Teuchos::null) const;
287 
290 
293 
296 
298  Kokkos::View<bool*, typename Node::device_type::memory_space> BCrowsKokkos_, BCcolsKokkos_, BCdomainKokkos_;
308 
310  typedef enum { MODE_STANDARD = 0,
315 
319 };
320 
321 } // namespace MueLu
322 
323 #define MUELU_MAXWELL1_SHORT
324 #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_
MueLu::DefaultNode Node
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&#39;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_
Maxwell1()
Constructor.
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.