Nektar++
NekLinSysIterCG.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NekLinSysIterCG.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: NekLinSysIterCG definition
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38using namespace std;
39
41{
42/**
43 * @class NekLinSysIterCG
44 *
45 * Solves a linear system using iterative methods.
46 */
49 "ConjugateGradient", NekLinSysIterCG::create,
50 "NekLinSysIterCG solver.");
51
54 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
55 const NekSysKey &pKey)
56 : NekLinSysIter(pSession, vRowComm, nDimen, pKey)
57{
58}
59
61{
63}
64
65/**
66 *
67 */
68int NekLinSysIterCG::v_SolveSystem(const int nGlobal,
69 const Array<OneD, const NekDouble> &pInput,
71 const int nDir, const NekDouble tol,
72 [[maybe_unused]] const NekDouble factor)
73{
74 m_tolerance = max(tol, 1.0E-16);
75
76 DoConjugateGradient(nGlobal, pInput, pOutput, nDir);
77
78 return m_totalIterations;
79}
80
81/**  
82 * Solve a global linear system using the conjugate gradient method.  
83 * We solve only for the non-Dirichlet modes. The operator is evaluated  
84 * using an auxiliary function m_operator.DoNekSysLhsEval defined by the  
85 * specific solver. Distributed math routines are used to support  
86 * parallel execution of the solver.  
87 *  
88 * The implemented algorithm uses a reduced-communication reordering of  
89 * the standard PCG method (Demmel, Heath and Vorst, 1993)  
90 *  
91 * @param pInput Input residual of all DOFs.  
92 * @param pOutput Solution vector of all DOFs.  
93 */
95 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
96 Array<OneD, NekDouble> &pOutput, const int nDir)
97{
98 // Get vector sizes
99 int nNonDir = nGlobal - nDir;
100
101 // Allocate array storage
102 Array<OneD, NekDouble> w_A(nGlobal, 0.0);
103 Array<OneD, NekDouble> s_A(nGlobal, 0.0);
104 Array<OneD, NekDouble> p_A(nNonDir, 0.0);
105 Array<OneD, NekDouble> r_A(nNonDir, 0.0);
106 Array<OneD, NekDouble> q_A(nNonDir, 0.0);
108
109 NekDouble alpha;
111 NekDouble rho;
112 NekDouble rho_new;
113 NekDouble mu;
114 NekDouble eps;
115 Array<OneD, NekDouble> vExchange(3, 0.0);
116
117 // Copy initial residual from input
118 Vmath::Vcopy(nNonDir, pInput + nDir, 1, r_A, 1);
119
120 // Zero homogeneous out array ready for solution updates
121 // Should not be earlier in case input vector is same as
122 // output and above copy has been peformed
123 Vmath::Zero(nNonDir, tmp = pOutput + nDir, 1);
124
125 // Evaluate initial residual error for exit check
126 vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
127
128 m_rowComm->AllReduce(vExchange[2], Nektar::LibUtilities::ReduceSum);
129
130 eps = vExchange[2];
131
133 {
134 NekVector<NekDouble> inGlob(nGlobal, pInput, eWrapper);
135 Set_Rhs_Magnitude(inGlob);
136 }
137
139
140 // If input residual is less than tolerance skip solve.
142 {
143 if (m_verbose && m_root)
144 {
145 cout << "CG iterations made = " << m_totalIterations
146 << " using tolerance of " << m_tolerance
147 << " (error = " << sqrt(eps / m_rhs_magnitude)
148 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
149 }
150 return;
151 }
152
153 m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
154 m_operator.DoNekSysLhsEval(w_A, s_A);
155
156 vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
157 vExchange[1] = Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
158
159 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
160
161 rho = vExchange[0];
162 mu = vExchange[1];
163 beta = 0.0;
164 alpha = rho / mu;
166
167 // Continue until convergence
168 while (true)
169 {
171 {
172 if (m_root)
173 {
174 cout << "CG iterations made = " << m_totalIterations
175 << " using tolerance of " << m_tolerance
176 << " (error = " << sqrt(eps / m_rhs_magnitude)
177 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
178 }
180 "Exceeded maximum number of iterations");
181 }
182
183 // Compute new search direction p_k, q_k
184 Vmath::Svtvp(nNonDir, beta, &p_A[0], 1, &w_A[nDir], 1, &p_A[0], 1);
185 Vmath::Svtvp(nNonDir, beta, &q_A[0], 1, &s_A[nDir], 1, &q_A[0], 1);
186
187 // Update solution x_{k+1}
188 Vmath::Svtvp(nNonDir, alpha, &p_A[0], 1, &pOutput[nDir], 1,
189 &pOutput[nDir], 1);
190
191 // Update residual vector r_{k+1}
192 Vmath::Svtvp(nNonDir, -alpha, &q_A[0], 1, &r_A[0], 1, &r_A[0], 1);
193
194 // Apply preconditioner
195 m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
196
197 // Perform the method-specific matrix-vector multiply operation.
198 m_operator.DoNekSysLhsEval(w_A, s_A);
199
200 // <r_{k+1}, w_{k+1}>
201 vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
202
203 // <s_{k+1}, w_{k+1}>
204 vExchange[1] =
205 Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
206
207 // <r_{k+1}, r_{k+1}>
208 vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
209
210 // Perform inner-product exchanges
211 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
212
213 rho_new = vExchange[0];
214 mu = vExchange[1];
215 eps = vExchange[2];
216
218
219 // Test if norm is within tolerance
221 {
222 if (m_verbose && m_root)
223 {
224 cout << "CG iterations made = " << m_totalIterations
225 << " using tolerance of " << m_tolerance
226 << " (error = " << sqrt(eps / m_rhs_magnitude)
227 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
228 }
229 break;
230 }
231
232 // Compute search direction and solution coefficients
233 beta = rho_new / rho;
234 alpha = rho_new / (mu - rho_new * beta / alpha);
235 rho = rho_new;
236 }
237}
238} // namespace Nektar::LibUtilities
#define ROOTONLY_NEKERROR(type, msg)
Definition: ErrorUtil.hpp:205
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
void DoConjugateGradient(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir)
Actual iterative solve.
int v_SolveSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir, const NekDouble tol, const NekDouble factor) override
NekLinSysIterCG(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
Constructor for full direct matrix solve.
static NekLinSysIterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
NekDouble m_rhs_magnitude
Dot product of rhs to normalise stopping criterion.
void Set_Rhs_Magnitude(const NekVector< NekDouble > &pIn)
Array< OneD, int > m_map
Global to universal unique map.
bool m_root
Root if parallel.
Definition: NekSys.h:297
LibUtilities::CommSharedPtr m_rowComm
Communicate.
Definition: NekSys.h:293
bool m_verbose
Verbose.
Definition: NekSys.h:299
NekDouble m_tolerance
Tolerance of iterative solver.
Definition: NekSys.h:291
NekSysOperators m_operator
Operators.
Definition: NekSys.h:302
int m_maxiter
Maximum iterations.
Definition: NekSys.h:289
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:148
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:141
std::shared_ptr< SessionReader > SessionReaderSharedPtr
NekLinSysIterFactory & GetNekLinSysIterFactory()
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
static const NekDouble kNekUnsetDouble
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
T Dot2(int n, const T *w, const T *x, const int *y)
dot product
Definition: Vmath.hpp:790
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294