MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_MultiPhys_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_MULTIPHYS_DECL_HPP
11#define MUELU_MULTIPHYS_DECL_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14#include "MueLu_BaseClass.hpp"
15
17
19#include "MueLu_TrilinosSmoother.hpp"
21#include "MueLu_Level_fwd.hpp"
26
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 MultiPhys : public VerboseObject, public Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
47#undef MUELU_MULTIPHYS_SHORT
49
50 public:
51 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitudeType;
52 typedef typename Teuchos::ScalarTraits<Scalar>::coordinateType coordinateType;
53 typedef typename Xpetra::MultiVector<coordinateType, LO, GO, NO> RealValuedMultiVector;
54
59 , nBlks_(0) {
60 }
61
74 MultiPhys(const Teuchos::RCP<Matrix>& AmatMultiPhysics,
75 const Teuchos::ArrayRCP<RCP<Matrix>> arrayOfAuxMatrices,
76 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfNullspaces,
77 const Teuchos::ArrayRCP<Teuchos::RCP<RealValuedMultiVector>> arrayOfCoords,
78 const int nBlks,
79 Teuchos::ParameterList& List,
80 bool ComputePrec = true,
81 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfMaterials = Teuchos::null,
82 bool OmitSubblockSmoother = true)
83 : AmatMultiphysics_(AmatMultiPhysics)
84 , arrayOfAuxMatrices_(arrayOfAuxMatrices)
85 , arrayOfNullspaces_(arrayOfNullspaces)
86 , arrayOfCoords_(arrayOfCoords)
87 , arrayOfMaterials_(arrayOfMaterials)
88 , OmitSubblockSmoother_(OmitSubblockSmoother)
89 , nBlks_(nBlks) {
90 initialize(AmatMultiPhysics, arrayOfAuxMatrices, arrayOfNullspaces, arrayOfCoords, nBlks, List, arrayOfMaterials);
91 compute(false);
92 }
93
95 virtual ~MultiPhys() {}
96
98 const Teuchos::RCP<const Map> getDomainMap() const;
99
101 const Teuchos::RCP<const Map> getRangeMap() const;
102
104 void setParameters(Teuchos::ParameterList& list);
105
107 void compute(bool reuse = false);
108
110 void resetMatrix(Teuchos::RCP<Matrix> SM_Matrix_new, bool ComputePrec = true);
111
115 void apply(const MultiVector& X, MultiVector& Y,
116 Teuchos::ETransp mode = Teuchos::NO_TRANS,
117 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
118 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero()) const;
119
121 bool hasTransposeApply() const;
122
123 void describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel = Teuchos::VERB_HIGH) const;
124
126 void residual(const MultiVector& X,
127 const MultiVector& B,
128 MultiVector& R) const {
129 using STS = Teuchos::ScalarTraits<Scalar>;
130 R.update(STS::one(), B, STS::zero());
131 this->apply(X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
132 }
133
135 Teuchos::RCP<Hierarchy> multiphysicsHierarchy() {
137 }
138
140 Teuchos::ArrayRCP<Teuchos::RCP<Hierarchy>> subblockHierarchies() {
141 return arrayOfHierarchies_;
142 }
143
144 private:
155 void initialize(const Teuchos::RCP<Matrix>& AmatMultiPhysics,
156 const Teuchos::ArrayRCP<RCP<Matrix>> arrayOfAuxMatrices,
157 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfNullspaces,
158 const Teuchos::ArrayRCP<Teuchos::RCP<RealValuedMultiVector>> arrayOfCoords,
159 const int nBlks,
160 Teuchos::ParameterList& List,
161 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfMaterials);
162
164 void applyInverse(const MultiVector& RHS, MultiVector& X) const;
165
167 void allocateMemory(int numVectors) const;
168
170 Teuchos::RCP<Teuchos::TimeMonitor> getTimer(std::string name, RCP<const Teuchos::Comm<int>> comm = Teuchos::null) const;
171
173 mutable Teuchos::ParameterList parameterList_;
174
176
177 Teuchos::RCP<Matrix> AmatMultiphysics_; // multiphysics discretization matrix
178 Teuchos::RCP<Teuchos::ParameterList> paramListMultiphysics_; // array of parameter lists directing MueLu's construct of subblock P operators
179 Teuchos::RCP<Hierarchy> hierarchyMultiphysics_; // multiphysics discretization matrix
180
181 Teuchos::ArrayRCP<Teuchos::RCP<Teuchos::ParameterList>> arrayOfParamLists_; // array of parameter lists directing MueLu's construct of subblock P operators
182 Teuchos::ArrayRCP<Teuchos::RCP<Hierarchy>> arrayOfHierarchies_;
183 Teuchos::ArrayRCP<Teuchos::RCP<Matrix>> arrayOfAuxMatrices_; // array of discretization/auxiliary matrices used to generate subblock prolongators
184 Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfNullspaces_; // array of nullspaces for smoothed aggregation.
185 Teuchos::ArrayRCP<Teuchos::RCP<RealValuedMultiVector>> arrayOfCoords_; // array of coordinates for smoothed aggregation/rebalancing.
186 Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfMaterials_; // array of materials for smoothed aggregation.
187
189
190 int nBlks_; // number of PDE sub-systems within multiphysics system
192};
193
194} // namespace MueLu
195
196#define MUELU_MULTIPHYS_SHORT
197#endif // MUELU_MULTIPHYS_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Preconditioner (wrapped as a Xpetra::Operator) for solving MultiPhysics PDEs.
Teuchos::RCP< Hierarchy > hierarchyMultiphysics_
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfNullspaces_
MultiPhys(const Teuchos::RCP< Matrix > &AmatMultiPhysics, const Teuchos::ArrayRCP< RCP< Matrix > > arrayOfAuxMatrices, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfNullspaces, const Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector > > arrayOfCoords, const int nBlks, Teuchos::ParameterList &List, bool ComputePrec=true, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfMaterials=Teuchos::null, bool OmitSubblockSmoother=true)
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
get a (synced) timer
Teuchos::ArrayRCP< Teuchos::RCP< Matrix > > arrayOfAuxMatrices_
void compute(bool reuse=false)
Setup the preconditioner.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
Teuchos::ParameterList parameterList_
ParameterLists.
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
void allocateMemory(int numVectors) const
allocate multivectors for solve
Teuchos::ScalarTraits< Scalar >::coordinateType coordinateType
void applyInverse(const MultiVector &RHS, MultiVector &X) const
apply standard MultiPhys cycle
void initialize(const Teuchos::RCP< Matrix > &AmatMultiPhysics, const Teuchos::ArrayRCP< RCP< Matrix > > arrayOfAuxMatrices, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfNullspaces, const Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector > > arrayOfCoords, const int nBlks, Teuchos::ParameterList &List, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfMaterials)
Xpetra::MultiVector< coordinateType, LO, GO, NO > RealValuedMultiVector
void setParameters(Teuchos::ParameterList &list)
Set parameters.
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
Teuchos::RCP< Matrix > AmatMultiphysics_
Hierarchies: used to define P for (0,0)-block, .... (nBlks_-1,nBlks_-1) block.
Teuchos::ArrayRCP< Teuchos::RCP< Hierarchy > > subblockHierarchies()
Return an array of hierarchies corresponding to each diagonal subblock.
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
virtual ~MultiPhys()
Destructor.
Teuchos::RCP< Hierarchy > multiphysicsHierarchy()
Return the multiphysics hierarchy.
Teuchos::ArrayRCP< Teuchos::RCP< Hierarchy > > arrayOfHierarchies_
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Teuchos::RCP< Teuchos::ParameterList > paramListMultiphysics_
Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector > > arrayOfCoords_
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfMaterials_
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::ArrayRCP< Teuchos::RCP< Teuchos::ParameterList > > arrayOfParamLists_
Verbose class for MueLu classes.
Namespace for MueLu classes and methods.