MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Constraint_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_CONSTRAINT_DECL_HPP
11#define MUELU_CONSTRAINT_DECL_HPP
12
13#include "Teuchos_ScalarTraits.hpp"
14
15#include <Xpetra_MultiVector_fwd.hpp>
16#include <Xpetra_MultiVectorFactory.hpp>
17#include <Xpetra_Matrix.hpp>
18#include <Xpetra_CrsGraph_fwd.hpp>
19#include <Xpetra_MatrixFactory_fwd.hpp>
20
21#include "MueLu_ConfigDefs.hpp"
22#include "MueLu_BaseClass.hpp"
24
25#ifdef HAVE_MUELU_BELOS
26#include <BelosLinearProblem.hpp>
27#include <BelosSolverFactory.hpp>
28#include <BelosXpetraAdapter.hpp>
29#endif
30
31namespace MueLu {
32
34
74template <class Scalar = DefaultScalar,
77 class Node = DefaultNode>
79 : public Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
80 public BaseClass {
81#undef MUELU_CONSTRAINT_SHORT
83 public:
89 using MagnitudeType = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
90
91 Constraint() = default;
92
93 virtual MagnitudeType ResidualNorm(const RCP<const Matrix> P) const = 0;
94
96 virtual const RCP<const Map> getDomainMap() const;
97
99 virtual const RCP<const Map> getRangeMap() const;
100
102
103
105 virtual void apply(const MultiVector& P,
106 MultiVector& Projected,
107 Teuchos::ETransp mode = Teuchos::NO_TRANS,
108 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
109 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero()) const;
110
112 void residual(const MultiVector& X,
113 const MultiVector& B,
114 MultiVector& R) const;
115
117
118 RCP<const CrsGraph> GetPattern() const;
119
120 void SetPattern(RCP<const CrsGraph>& Ppattern);
121
122 void SetConstraintsMatrix(RCP<Matrix>& X);
123
124 RCP<Matrix> GetConstraintMatrix();
125
126 virtual typename CrsGraph::local_graph_type FindBlocks(RCP<const CrsGraph>& XXt);
127
128 void AssignMatrixEntriesToVector(const Matrix& P, const RCP<const CrsGraph>& pattern, MultiVector& vecP) const;
129
130 void AssignMatrixEntriesToVector(const Matrix& P, MultiVector& vecP) const;
131
132 RCP<Matrix> GetMatrixWithEntriesFromVector(MultiVector& vecP) const;
133
134 void LeastSquaresSolve(const MultiVector& B, MultiVector& C) const;
135
136 protected:
137 void PrepareLeastSquaresSolve(const std::string& solverType, bool detect_singular_blocks = false);
138
139 private:
141 RCP<Matrix> X_;
142
144 RCP<const CrsGraph> Ppattern_;
145
146 std::string solverType_;
147
148#ifdef HAVE_MUELU_BELOS
149 using MV = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
150 using OP = Belos::OperatorT<MV>;
151 RCP<Belos::LinearProblem<Scalar, MV, OP>> problem_;
152 RCP<Belos::SolverManager<Scalar, MV, OP>> solver_;
153#endif
154
156 void PrepareLeastSquaresSolveBelos(bool detect_singular_blocks);
157
159 void LeastSquaresSolveBelos(const MultiVector& B, MultiVector& C) const;
160
162 RCP<Matrix> invXXt_;
163
165 void PrepareLeastSquaresSolveDirect(bool detect_singular_blocks);
166
168 void LeastSquaresSolveDirect(const MultiVector& B, MultiVector& C) const;
169
170 // temporary memory
171 RCP<MultiVector> temp1_, temp2_, temp3_;
172};
173
174} // namespace MueLu
175
176#define MUELU_CONSTRAINT_SHORT
177#endif
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Base class for MueLu classes.
Constraint space information for the potential prolongator.
void LeastSquaresSolveDirect(const MultiVector &B, MultiVector &C) const
Direct solve of least-squares problem.
RCP< Matrix > GetConstraintMatrix()
virtual const RCP< const Map > getDomainMap() const
The Map associated with the domain of this operator, which must be compatible with X....
void PrepareLeastSquaresSolveDirect(bool detect_singular_blocks)
Prepare direct solution of least-squares problem.
Belos::OperatorT< MV > OP
RCP< Belos::SolverManager< Scalar, MV, OP > > solver_
RCP< const CrsGraph > Ppattern_
Nonzero sparsity pattern.
virtual void apply(const MultiVector &P, MultiVector &Projected, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Apply constraint.
RCP< Matrix > X_
The constraints matrix.
void PrepareLeastSquaresSolveBelos(bool detect_singular_blocks)
Prepare least-squares solve using Belos.
virtual CrsGraph::local_graph_type FindBlocks(RCP< const CrsGraph > &XXt)
void PrepareLeastSquaresSolve(const std::string &solverType, bool detect_singular_blocks=false)
typename Teuchos::ScalarTraits< Scalar >::magnitudeType MagnitudeType
RCP< Belos::LinearProblem< Scalar, MV, OP > > problem_
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
RCP< Matrix > invXXt_
Inverse of X*X^T.
virtual MagnitudeType ResidualNorm(const RCP< const Matrix > P) const =0
Constraint()=default
RCP< MultiVector > temp2_
void LeastSquaresSolve(const MultiVector &B, MultiVector &C) const
void SetConstraintsMatrix(RCP< Matrix > &X)
void LeastSquaresSolveBelos(const MultiVector &B, MultiVector &C) const
Perform least-squares solve using Belos.
RCP< Matrix > GetMatrixWithEntriesFromVector(MultiVector &vecP) const
RCP< const CrsGraph > GetPattern() const
RCP< MultiVector > temp3_
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MV
RCP< MultiVector > temp1_
void SetPattern(RCP< const CrsGraph > &Ppattern)
virtual const RCP< const Map > getRangeMap() const
The Map associated with the range of this operator, which must be compatible with Y....
void AssignMatrixEntriesToVector(const Matrix &P, const RCP< const CrsGraph > &pattern, MultiVector &vecP) const
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar