46 #ifndef MUELU_REFMAXWELL_DECL_HPP 47 #define MUELU_REFMAXWELL_DECL_HPP 58 #include "MueLu_TrilinosSmoother.hpp" 59 #include "MueLu_Hierarchy.hpp" 61 #include "Xpetra_Map_fwd.hpp" 62 #include "Xpetra_Matrix_fwd.hpp" 63 #include "Xpetra_MatrixFactory_fwd.hpp" 64 #include "Xpetra_MultiVectorFactory_fwd.hpp" 65 #include "Xpetra_CrsMatrixWrap_fwd.hpp" 66 #include "Xpetra_BlockedCrsMatrix_fwd.hpp" 67 #include "Xpetra_ExportFactory_fwd.hpp" 83 template <
class Scalar,
87 class RefMaxwell :
public Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> {
89 #undef MUELU_REFMAXWELL_SHORT 104 RefMaxwell(Teuchos::RCP<Hierarchy> H11, Teuchos::RCP<Hierarchy> H22) :
124 const Teuchos::RCP<Matrix> & D0_Matrix,
125 const Teuchos::RCP<Matrix> & M0inv_Matrix,
126 const Teuchos::RCP<Matrix> & M1_Matrix,
127 const Teuchos::RCP<MultiVector> & Nullspace,
128 const Teuchos::RCP<MultiVector> & Coords,
129 Teuchos::ParameterList& List,
130 bool ComputePrec =
true)
132 initialize(D0_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
151 const Teuchos::RCP<Matrix> & M0inv_Matrix,
152 const Teuchos::RCP<Matrix> & M1_Matrix,
153 const Teuchos::RCP<MultiVector> & Nullspace,
154 const Teuchos::RCP<MultiVector> & Coords,
157 initialize(D0_Matrix,M0inv_Matrix,M1_Matrix,Nullspace,Coords,List);
171 const Teuchos::RCP<Matrix> & D0_Matrix,
172 const Teuchos::RCP<Matrix> & M1_Matrix,
173 const Teuchos::RCP<MultiVector> & Nullspace,
174 const Teuchos::RCP<MultiVector> & Coords,
175 Teuchos::ParameterList& List,
176 bool ComputePrec =
true)
178 initialize(D0_Matrix,Teuchos::null,M1_Matrix,Nullspace,Coords,List);
196 const Teuchos::RCP<Matrix> & M1_Matrix,
197 const Teuchos::RCP<MultiVector> & Nullspace,
198 const Teuchos::RCP<MultiVector> & Coords,
201 initialize(D0_Matrix,Teuchos::null,M1_Matrix,Nullspace,Coords,List);
226 void resetMatrix(Teuchos::RCP<Matrix> SM_Matrix_new);
240 void apply (
const MultiVector& X, MultiVector& Y,
241 Teuchos::ETransp mode = Teuchos::NO_TRANS,
242 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
243 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero())
const;
248 template <
class NewNode>
249 Teuchos::RCP< RefMaxwell<Scalar, LocalOrdinal, GlobalOrdinal, NewNode> >
250 clone (
const RCP<NewNode>& new_node)
const {
259 std::vector<LocalOrdinal>& dirichletRows) {
260 dirichletRows.resize(0);
261 for(
size_t i=0; i<A->getNodeNumRows(); i++) {
262 Teuchos::ArrayView<const LocalOrdinal> indices;
263 Teuchos::ArrayView<const Scalar> values;
264 A->getLocalRowView(i,indices,values);
266 for (
int j=0; j<indices.size(); j++) {
272 if (Teuchos::ScalarTraits<Scalar>::magnitude(values[j]) > 1.0e-16) {
276 if (nnz == 1 || nnz == 2) {
277 dirichletRows.push_back(i);
283 std::vector<LocalOrdinal>& dirichletRows,
284 std::vector<LocalOrdinal>& dirichletCols) {
285 Teuchos::RCP<const Map> domMap = A->getDomainMap();
286 Teuchos::RCP<const Map> colMap = A->getColMap();
287 Teuchos::RCP<Export> exporter = ExportFactory::Build(colMap,domMap);
288 Teuchos::RCP<MultiVector> myColsToZero = MultiVectorFactory::Build(colMap,1);
289 Teuchos::RCP<MultiVector> globalColsToZero = MultiVectorFactory::Build(domMap,1);
290 myColsToZero->putScalar((Scalar)0.0);
291 globalColsToZero->putScalar((Scalar)0.0);
292 for(
size_t i=0; i<dirichletRows.size(); i++) {
293 Teuchos::ArrayView<const LocalOrdinal> indices;
294 Teuchos::ArrayView<const Scalar> values;
295 A->getLocalRowView(dirichletRows[i],indices,values);
296 for(
int j=0; j<indices.size(); j++)
297 myColsToZero->replaceLocalValue(indices[j],0,(Scalar)1.0);
299 globalColsToZero->doExport(*myColsToZero,*exporter,Xpetra::ADD);
300 myColsToZero->doImport(*globalColsToZero,*exporter,Xpetra::INSERT);
301 Teuchos::ArrayRCP<const Scalar> myCols = myColsToZero->getData(0);
302 dirichletCols.resize(colMap->getNodeNumElements());
303 for(
size_t i=0; i<colMap->getNodeNumElements(); i++) {
304 if(Teuchos::ScalarTraits<Scalar>::magnitude(myCols[i])>0.0)
312 std::vector<LocalOrdinal>& dirichletRows) {
313 for(
size_t i=0; i<dirichletRows.size(); i++) {
314 Teuchos::ArrayView<const LocalOrdinal> indices;
315 Teuchos::ArrayView<const Scalar> values;
316 A->getLocalRowView(dirichletRows[i],indices,values);
317 std::vector<Scalar> vec;
318 vec.resize(indices.size());
319 Teuchos::ArrayView<Scalar> zerovalues(vec);
320 for(
int j=0; j<indices.size(); j++)
321 zerovalues[j]=(Scalar)1.0e-32;
322 A->replaceLocalValues(dirichletRows[i],indices,zerovalues);
327 std::vector<LocalOrdinal>& dirichletCols) {
328 for(
size_t i=0; i<A->getNodeNumRows(); i++) {
329 Teuchos::ArrayView<const LocalOrdinal> indices;
330 Teuchos::ArrayView<const Scalar> values;
331 A->getLocalRowView(i,indices,values);
332 std::vector<Scalar> vec;
333 vec.resize(indices.size());
334 Teuchos::ArrayView<Scalar> zerovalues(vec);
335 for(
int j=0; j<indices.size(); j++) {
336 if(dirichletCols[indices[j]]==1)
337 zerovalues[j]=(Scalar)1.0e-32;
339 zerovalues[j]=values[j];
341 A->replaceLocalValues(i,indices,zerovalues);
346 Teuchos::RCP<const Map> rowMap = A->getRowMap();
347 RCP<Matrix> DiagMatrix = MatrixFactory::Build(rowMap,1);
348 RCP<Matrix> NewMatrix = MatrixFactory::Build(rowMap,1);
349 for(
size_t i=0; i<A->getNodeNumRows(); i++) {
350 Teuchos::ArrayView<const LocalOrdinal> indices;
351 Teuchos::ArrayView<const Scalar> values;
352 A->getLocalRowView(i,indices,values);
354 for (
int j=0; j<indices.size(); j++) {
355 if (Teuchos::ScalarTraits<Scalar>::magnitude(values[j]) > tol) {
359 Scalar one = (Scalar)1.0;
360 Scalar zero = (Scalar)0.0;
361 GlobalOrdinal row = rowMap->getGlobalElement(i);
363 DiagMatrix->insertGlobalValues(row,
364 Teuchos::ArrayView<GlobalOrdinal>(&row,1),
365 Teuchos::ArrayView<Scalar>(&one,1));
368 DiagMatrix->insertGlobalValues(row,
369 Teuchos::ArrayView<GlobalOrdinal>(&row,1),
370 Teuchos::ArrayView<Scalar>(&zero,1));
373 DiagMatrix->fillComplete();
376 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
377 Xpetra::MatrixMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TwoMatrixAdd(*DiagMatrix,
false,(Scalar)1.0,*A,
false,(Scalar)1.0,NewMatrix,*out);
378 NewMatrix->fillComplete();
392 void initialize(
const Teuchos::RCP<Matrix> & D0_Matrix,
393 const Teuchos::RCP<Matrix> & M0inv_Matrix,
394 const Teuchos::RCP<Matrix> & M1_Matrix,
395 const Teuchos::RCP<MultiVector> & Nullspace,
396 const Teuchos::RCP<MultiVector> & Coords,
397 Teuchos::ParameterList& List);
420 #define MUELU_REFMAXWELL_SHORT 421 #endif // MUELU_REFMAXWELL_DECL_HPP
RefMaxwell(Teuchos::RCP< Hierarchy > H11, Teuchos::RCP< Hierarchy > H22)
Constructor with Hierarchies.
Teuchos::RCP< Matrix > Ms_Matrix_
void initialize(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< MultiVector > &Coords, Teuchos::ParameterList &List)
void Apply_BCsToMatrixRows(Teuchos::RCP< Matrix > &A, std::vector< LocalOrdinal > &dirichletRows)
Teuchos::RCP< Matrix > D0_Matrix_
Teuchos::RCP< Matrix > SM_Matrix_
Various matrices.
void findDirichletCols(Teuchos::RCP< Matrix > A, std::vector< LocalOrdinal > &dirichletRows, std::vector< LocalOrdinal > &dirichletCols)
std::vector< LocalOrdinal > BCcols_
Teuchos::ParameterList smootherList_
std::vector< LocalOrdinal > BCrows_
Vectors for BCs.
Teuchos::RCP< MultiVector > Coords_
Teuchos::RCP< Matrix > A22_
Teuchos::RCP< Hierarchy > HierarchySmoother_
Namespace for MueLu classes and methods.
void findDirichletRows(Teuchos::RCP< Matrix > A, std::vector< LocalOrdinal > &dirichletRows)
void setParameters(Teuchos::ParameterList &list)
Set parameters.
void applyInverseAdditive(const MultiVector &RHS, MultiVector &X) const
apply additive algorithm for 2x2 solve
void buildProlongator()
Setup the prolongator for the (1,1)-block.
Teuchos::RCP< Matrix > TMT_Agg_Matrix_
Preconditioner (wrapped as a Tpetra::Operator) for Maxwell's equations in curl-curl form...
Teuchos::RCP< Matrix > M0inv_Matrix_
void Apply_BCsToMatrixCols(Teuchos::RCP< Matrix > &A, std::vector< LocalOrdinal > &dirichletCols)
void applyInverse121(const MultiVector &RHS, MultiVector &X) const
apply 1-2-1 algorithm for 2x2 solve
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< MultiVector > &Coords, Teuchos::ParameterList &List, bool ComputePrec=true)
Teuchos::ParameterList precList11_
void compute()
Setup the preconditioner.
Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
Teuchos::RCP< Matrix > P11_
Teuchos::RCP< Hierarchy > Hierarchy22_
Teuchos::RCP< Matrix > M1_Matrix_
Teuchos::RCP< Matrix > A11_
Teuchos::RCP< Hierarchy > Hierarchy11_
Two hierarchies: one for the (1,1)-block, another for the (2,2)-block.
bool disable_addon_
Some options.
Teuchos::RCP< Level > TopLevel_
Top Level.
virtual ~RefMaxwell()
Destructor.
Teuchos::ParameterList parameterList_
Parameter lists.
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
void applyInverse212(const MultiVector &RHS, MultiVector &X) const
apply 2-1-2 algorithm for 2x2 solve
Teuchos::RCP< RefMaxwell< Scalar, LocalOrdinal, GlobalOrdinal, NewNode > > clone(const RCP< NewNode > &new_node) const
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< MultiVector > &Coords, Teuchos::ParameterList &List)
RefMaxwell(const Teuchos::RCP< Matrix > &D0_Matrix, const Teuchos::RCP< Matrix > &M1_Matrix, const Teuchos::RCP< MultiVector > &Nullspace, const Teuchos::RCP< MultiVector > &Coords, Teuchos::ParameterList &List)
void formCoarseMatrix()
Compute P11^{T}*A*P11 efficiently.
void Remove_Zeroed_Rows(Teuchos::RCP< Matrix > &A, double tol=1.0e-14)
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::ParameterList precList22_
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< MultiVector > &Coords, Teuchos::ParameterList &List, bool ComputePrec=true)
Teuchos::RCP< MultiVector > Nullspace_
Nullspace.
Teuchos::RCP< Matrix > TMT_Matrix_
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new)
Reset system matrix.
Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.