MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_CGSolver_def.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_CGSOLVER_DEF_HPP
11#define MUELU_CGSOLVER_DEF_HPP
12
13#include <Xpetra_MatrixFactory.hpp>
14#include <Xpetra_MatrixFactory2.hpp>
15#include <Xpetra_MatrixMatrix.hpp>
16
17#include "MueLu_Utilities.hpp"
18#include "MueLu_Constraint.hpp"
19#include "MueLu_Monitor.hpp"
20
22
23namespace MueLu {
24
25template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28
29template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
30void CGSolver<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Iterate(const Matrix& Aref, const Constraint& C, const Matrix& P0, RCP<Matrix>& finalP) const {
31 // Note: this function matrix notations follow Saad's "Iterative methods", ed. 2, pg. 246
32 // So, X is the unknown prolongator, P's are conjugate directions, Z's are preconditioned P's
33 PrintMonitor m(*this, "CG iterations");
34
35 if (nIts_ == 0) {
36 finalP = MatrixFactory2::BuildCopy(rcpFromRef(P0));
37 return;
38 }
39
40 RCP<const Matrix> A = rcpFromRef(Aref);
41 ArrayRCP<const SC> D = Utilities::GetMatrixDiagonal_arcp(*A);
42 bool useTpetra = (A->getRowMap()->lib() == Xpetra::UseTpetra);
43
44 Teuchos::FancyOStream& mmfancy = this->GetOStream(Statistics2);
45
46 SC one = Teuchos::ScalarTraits<SC>::one();
47
48 RCP<Matrix> X, P, R, Z, AP;
49 RCP<Matrix> newX, tmpAP;
50#ifndef TWO_ARG_MATRIX_ADD
51 RCP<Matrix> newR, newP;
52#endif
53
54 SC oldRZ, newRZ, alpha, beta, app;
55
56 // T is used only for projecting onto
57 RCP<CrsMatrix> T_ = CrsMatrixFactory::Build(C.GetPattern());
58 T_->fillComplete(P0.getDomainMap(), P0.getRangeMap());
59 RCP<Matrix> T = rcp(new CrsMatrixWrap(T_));
60
61 // Initial P0 would only be used for multiplication
62 X = rcp_const_cast<Matrix>(rcpFromRef(P0));
63
64 tmpAP = MatrixMatrix::Multiply(*A, false, *X, false, mmfancy, true /*doFillComplete*/, true /*optimizeStorage*/);
65 C.Apply(*tmpAP, *T);
66
67 // R_0 = -A*X_0
68 R = Xpetra::MatrixFactory2<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(T);
69
70 R->scale(-one);
71
72 // Z_0 = M^{-1}R_0
73 Z = Xpetra::MatrixFactory2<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(R);
74 Utilities::MyOldScaleMatrix(*Z, D, true, true, false);
75
76 // P_0 = Z_0
77 P = Xpetra::MatrixFactory2<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(Z);
78
79 oldRZ = Utilities::Frobenius(*R, *Z);
80
81 for (size_t i = 0; i < nIts_; i++) {
82 // AP = constrain(A*P)
83 if (i == 0 || useTpetra) {
84 // Construct the MxM pattern from scratch
85 // This is done by default for Tpetra as the three argument version requires tmpAP
86 // to *not* be locally indexed which defeats the purpose
87 // TODO: need a three argument Tpetra version which allows reuse of already fill-completed matrix
88 tmpAP = MatrixMatrix::Multiply(*A, false, *P, false, mmfancy, true /*doFillComplete*/, true /*optimizeStorage*/);
89 } else {
90 // Reuse the MxM pattern
91 tmpAP = MatrixMatrix::Multiply(*A, false, *P, false, tmpAP, mmfancy, true /*doFillComplete*/, true /*optimizeStorage*/);
92 }
93 C.Apply(*tmpAP, *T);
94 AP = T;
95
96 app = Utilities::Frobenius(*AP, *P);
97 if (Teuchos::ScalarTraits<SC>::magnitude(app) < Teuchos::ScalarTraits<SC>::sfmin()) {
98 // It happens, for instance, if P = 0
99 // For example, if we use TentativePFactory for both nonzero pattern and initial guess
100 // I think it might also happen because of numerical breakdown, but we don't test for that yet
101 if (i == 0)
102 X = MatrixFactory2::BuildCopy(rcpFromRef(P0));
103 break;
104 }
105
106 // alpha = (R_i, Z_i)/(A*P_i, P_i)
107 alpha = oldRZ / app;
108 this->GetOStream(Runtime1, 1) << "alpha = " << alpha << std::endl;
109
110 // X_{i+1} = X_i + alpha*P_i
111#ifndef TWO_ARG_MATRIX_ADD
112 newX = Teuchos::null;
113 MatrixMatrix::TwoMatrixAdd(*P, false, alpha, *X, false, one, newX, mmfancy);
114 newX->fillComplete(P0.getDomainMap(), P0.getRangeMap());
115 X.swap(newX);
116#else
117 MatrixMatrix::TwoMatrixAdd(*P, false, alpha, *X, one);
118#endif
119
120 if (i == nIts_ - 1)
121 break;
122
123 // R_{i+1} = R_i - alpha*A*P_i
124#ifndef TWO_ARG_MATRIX_ADD
125 newR = Teuchos::null;
126 MatrixMatrix::TwoMatrixAdd(*AP, false, -alpha, *R, false, one, newR, mmfancy);
127 newR->fillComplete(P0.getDomainMap(), P0.getRangeMap());
128 R.swap(newR);
129#else
130 MatrixMatrix::TwoMatrixAdd(*AP, false, -alpha, *R, one);
131#endif
132
133 // Z_{i+1} = M^{-1} R_{i+1}
134 Z = MatrixFactory2::BuildCopy(R);
135 Utilities::MyOldScaleMatrix(*Z, D, true, true, false);
136
137 // beta = (R_{i+1}, Z_{i+1})/(R_i, Z_i)
138 newRZ = Utilities::Frobenius(*R, *Z);
139 beta = newRZ / oldRZ;
140
141 // P_{i+1} = Z_{i+1} + beta*P_i
142#ifndef TWO_ARG_MATRIX_ADD
143 newP = Teuchos::null;
144 MatrixMatrix::TwoMatrixAdd(*P, false, beta, *Z, false, one, newP, mmfancy);
145 newP->fillComplete(P0.getDomainMap(), P0.getRangeMap());
146 P.swap(newP);
147#else
148 MatrixMatrix::TwoMatrixAdd(*Z, false, one, *P, beta);
149#endif
150
151 oldRZ = newRZ;
152 }
153
154 finalP = X;
155}
156
157} // namespace MueLu
158
159#endif // ifndef MUELU_CGSOLVER_DECL_HPP
void Iterate(const Matrix &A, const Constraint &C, const Matrix &P0, RCP< Matrix > &P) const
Iterate.
size_t nIts_
Number of performed iterations.
Constraint space information for the potential prolongator.
RCP< const CrsGraph > GetPattern() const
void Apply(const Matrix &P, Matrix &Projected) const
Apply constraint.
static Scalar Frobenius(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &B)
static void MyOldScaleMatrix(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< const Scalar > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Runtime1
Description of what is happening (more verbose).