Nektar++
Loading...
Searching...
No Matches
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
64
65/**
66 *
67 */
68int NekLinSysIterCG::v_SolveSystem(const int nGlobal,
69 const Array<OneD, const NekDouble> &pInput,
71 const int nDir)
72{
73 DoConjugateGradient(nGlobal, pInput, pOutput, nDir);
74
75 return m_totalIterations;
76}
77
78void NekLinSysIterCG::v_DoIterate(const int nGlobal,
79 const Array<OneD, NekDouble> &rhs,
80 Array<OneD, NekDouble> &x, const int nDir,
81 NekDouble &err, int &iter)
82{
83 DoConjugateGradient(nGlobal, rhs, x, nDir);
84 iter = m_totalIterations;
85 err = m_finalError;
86}
87
88/**  
89 * Solve a global linear system using the conjugate gradient method.  
90 * We solve only for the non-Dirichlet modes. The operator is evaluated  
91 * using an auxiliary function m_operator.DoNekSysLhsEval defined by the  
92 * specific solver. Distributed math routines are used to support  
93 * parallel execution of the solver.  
94 *  
95 * The implemented algorithm uses a reduced-communication reordering of  
96 * the standard PCG method (Demmel, Heath and Vorst, 1993)  
97 *  
98 * @param pInput Input residual of all DOFs.  
99 * @param pOutput Solution vector of all DOFs.  
100 */
102 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
103 Array<OneD, NekDouble> &pOutput, const int nDir)
104{
105 // Get vector sizes
106 int nNonDir = nGlobal - nDir;
107
108 // Allocate array storage
109 Array<OneD, NekDouble> w_A(nGlobal, 0.0);
110 Array<OneD, NekDouble> s_A(nGlobal, 0.0);
111 Array<OneD, NekDouble> p_A(nNonDir, 0.0);
112 Array<OneD, NekDouble> r_A(nNonDir, 0.0);
113 Array<OneD, NekDouble> q_A(nNonDir, 0.0);
115
116 NekDouble alpha;
118 NekDouble rho;
119 NekDouble rho_new;
120 NekDouble mu;
121 NekDouble eps;
122 Array<OneD, NekDouble> vExchange(3, 0.0);
123
124 // Copy initial residual from input
125 Vmath::Vcopy(nNonDir, pInput + nDir, 1, r_A, 1);
126
127 // Zero homogeneous out array ready for solution updates
128 // Should not be earlier in case input vector is same as
129 // output and above copy has been peformed
130 Vmath::Zero(nNonDir, tmp = pOutput + nDir, 1);
131
132 // Evaluate initial residual error for exit check
133 vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
134
135 m_rowComm->AllReduce(vExchange[2], Nektar::LibUtilities::ReduceSum);
136
137 eps = vExchange[2];
138
140 {
141 Set_Rhs_Magnitude(pInput);
142 }
143
145
146 // If input residual is less than tolerance skip solve.
148 {
150 if (m_verbose && m_root)
151 {
152 cout << "CG iterations made = " << m_totalIterations
153 << " using tolerance of " << m_NekLinSysTolerance
154 << " (error = " << m_finalError
155 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
156 }
157 return;
158 }
159
160 m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
161 m_operator.DoNekSysLhsEval(w_A, s_A);
162
163 vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
164 vExchange[1] = Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
165
166 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
167
168 rho = vExchange[0];
169 mu = vExchange[1];
170 beta = 0.0;
171 alpha = rho / mu;
173
174 // Continue until convergence
175 while (true)
176 {
178 {
180 if (m_root)
181 {
182 cout << "CG iterations made = " << m_totalIterations
183 << " using tolerance of " << m_NekLinSysTolerance
184 << " (error = " << m_finalError
185 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")"
186 << " WARNING: Exceeded maxIt" << endl;
187 }
188 break;
189 }
190
191 // Compute new search direction p_k, q_k
192 Vmath::Svtvp(nNonDir, beta, &p_A[0], 1, &w_A[nDir], 1, &p_A[0], 1);
193 Vmath::Svtvp(nNonDir, beta, &q_A[0], 1, &s_A[nDir], 1, &q_A[0], 1);
194
195 // Update solution x_{k+1}
196 Vmath::Svtvp(nNonDir, alpha, &p_A[0], 1, &pOutput[nDir], 1,
197 &pOutput[nDir], 1);
198
199 // Update residual vector r_{k+1}
200 Vmath::Svtvp(nNonDir, -alpha, &q_A[0], 1, &r_A[0], 1, &r_A[0], 1);
201
202 // Apply preconditioner
203 m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
204
205 // Perform the method-specific matrix-vector multiply operation.
206 m_operator.DoNekSysLhsEval(w_A, s_A);
207
208 if (m_mapIsOnes)
209 {
210 // <r_{k+1}, w_{k+1}>
211 vExchange[0] = Vmath::Dot(nNonDir, r_A, w_A + nDir);
212
213 // <s_{k+1}, w_{k+1}>
214 vExchange[1] = Vmath::Dot(nNonDir, s_A + nDir, w_A + nDir);
215
217 {
218 // <r_{k+1}, r_{k+1}>
219 vExchange[2] = Vmath::Dot(nNonDir, r_A, r_A);
220 }
221 }
222 else
223 {
224 // <r_{k+1}, w_{k+1}>
225 vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
226
227 // <s_{k+1}, w_{k+1}>
228 vExchange[1] =
229 Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
230
232 {
233 // <r_{k+1}, r_{k+1}>
234 vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
235 }
236 }
237 // Perform inner-product exchanges
238 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
239
240 rho_new = vExchange[0];
241 mu = vExchange[1];
242 eps = vExchange[2];
243
245
246 // Test if norm is within tolerance
248 {
250 if (m_verbose && m_root)
251 {
252 cout << "CG iterations made = " << m_totalIterations
253 << " using tolerance of " << m_NekLinSysTolerance
254 << " (error = " << m_finalError
255 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
256 }
257 break;
258 }
259
260 // Compute search direction and solution coefficients
261 beta = rho_new / rho;
262 alpha = rho_new / (mu - rho_new * beta / alpha);
263 rho = rho_new;
264 }
265}
266} // namespace Nektar::LibUtilities
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
void DoConjugateGradient(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir)
Actual iterative solve.
void v_DoIterate(const int nGlobal, const Array< OneD, NekDouble > &rhs, Array< OneD, NekDouble > &x, const int nDir, NekDouble &err, int &iter) override
int v_SolveSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir) 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)
void Set_Rhs_Magnitude(const Array< OneD, NekDouble > &pIn)
Array< OneD, int > m_map
Global to universal unique map.
LibUtilities::CommSharedPtr m_rowComm
Definition NekSys.h:299
NekSysOperators m_operator
Definition NekSys.h:306
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition NekSys.h:155
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition NekSys.h:148
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
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
T Dot(int n, const T *w, const T *x)
dot product
Definition Vmath.hpp:761
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
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290