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
38
39using namespace std;
40
41namespace Nektar
42{
43namespace LibUtilities
44{
45/**
46 * @class NekLinSysIterCG
47 *
48 * Solves a linear system using iterative methods.
49 */
52 "ConjugateGradient", NekLinSysIterCG::create,
53 "NekLinSysIterCG solver.");
54
57 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
58 const NekSysKey &pKey)
59 : NekLinSysIter(pSession, vRowComm, nDimen, pKey)
60{
61}
62
64{
66}
67
69{
70}
71
72/**
73 *
74 */
75int NekLinSysIterCG::v_SolveSystem(const int nGlobal,
76 const Array<OneD, const NekDouble> &pInput,
78 const int nDir, const NekDouble tol,
79 const NekDouble factor)
80{
81 m_tolerance = max(tol, 1.0E-16);
82 m_prec_factor = factor;
83
84 DoConjugateGradient(nGlobal, pInput, pOutput, nDir);
85
86 return m_totalIterations;
87}
88
89/**  
90 * Solve a global linear system using the conjugate gradient method.  
91 * We solve only for the non-Dirichlet modes. The operator is evaluated  
92 * using an auxiliary function m_operator.DoNekSysLhsEval defined by the  
93 * specific solver. Distributed math routines are used to support  
94 * parallel execution of the solver.  
95 *  
96 * The implemented algorithm uses a reduced-communication reordering of  
97 * the standard PCG method (Demmel, Heath and Vorst, 1993)  
98 *  
99 * @param pInput Input residual of all DOFs.  
100 * @param pOutput Solution vector of all DOFs.  
101 */
103 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
104 Array<OneD, NekDouble> &pOutput, const int nDir)
105{
106 // Get vector sizes
107 int nNonDir = nGlobal - nDir;
108
109 // Allocate array storage
110 Array<OneD, NekDouble> w_A(nGlobal, 0.0);
111 Array<OneD, NekDouble> s_A(nGlobal, 0.0);
112 Array<OneD, NekDouble> p_A(nNonDir, 0.0);
113 Array<OneD, NekDouble> r_A(nNonDir, 0.0);
114 Array<OneD, NekDouble> q_A(nNonDir, 0.0);
116
117 NekDouble alpha;
119 NekDouble rho;
120 NekDouble rho_new;
121 NekDouble mu;
122 NekDouble eps;
123 Array<OneD, NekDouble> vExchange(3, 0.0);
124
125 // Copy initial residual from input
126 Vmath::Vcopy(nNonDir, pInput + nDir, 1, r_A, 1);
127
128 // Zero homogeneous out array ready for solution updates
129 // Should not be earlier in case input vector is same as
130 // output and above copy has been peformed
131 Vmath::Zero(nNonDir, tmp = pOutput + nDir, 1);
132
133 // Evaluate initial residual error for exit check
134 vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
135
136 m_rowComm->AllReduce(vExchange[2], Nektar::LibUtilities::ReduceSum);
137
138 eps = vExchange[2];
139
141 {
142 NekVector<NekDouble> inGlob(nGlobal, pInput, eWrapper);
143 Set_Rhs_Magnitude(inGlob);
144 }
145
147
148 // If input residual is less than tolerance skip solve.
150 {
151 if (m_verbose && m_root)
152 {
153 cout << "CG iterations made = " << m_totalIterations
154 << " using tolerance of " << m_tolerance
155 << " (error = " << sqrt(eps / m_rhs_magnitude)
156 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
157 }
158 return;
159 }
160
161 m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
162
163 m_operator.DoNekSysLhsEval(w_A, s_A);
164
165 vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
166
167 vExchange[1] = Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
168
169 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
170
171 rho = vExchange[0];
172 mu = vExchange[1];
173 beta = 0.0;
174 alpha = rho / mu;
176
177 // Continue until convergence
178 while (true)
179 {
181 {
182 if (m_root)
183 {
184 cout << "CG iterations made = " << m_totalIterations
185 << " using tolerance of " << m_tolerance
186 << " (error = " << sqrt(eps / m_rhs_magnitude)
187 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
188 }
190 "Exceeded maximum number of iterations");
191 }
192
193 // Compute new search direction p_k, q_k
194 Vmath::Svtvp(nNonDir, beta, &p_A[0], 1, &w_A[nDir], 1, &p_A[0], 1);
195 Vmath::Svtvp(nNonDir, beta, &q_A[0], 1, &s_A[nDir], 1, &q_A[0], 1);
196
197 // Update solution x_{k+1}
198 Vmath::Svtvp(nNonDir, alpha, &p_A[0], 1, &pOutput[nDir], 1,
199 &pOutput[nDir], 1);
200
201 // Update residual vector r_{k+1}
202 Vmath::Svtvp(nNonDir, -alpha, &q_A[0], 1, &r_A[0], 1, &r_A[0], 1);
203
204 // Apply preconditioner
205 m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
206
207 // Perform the method-specific matrix-vector multiply operation.
208 m_operator.DoNekSysLhsEval(w_A, s_A);
209
210 // <r_{k+1}, w_{k+1}>
211 vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
212
213 // <s_{k+1}, w_{k+1}>
214 vExchange[1] =
215 Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
216
217 // <r_{k+1}, r_{k+1}>
218 vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
219
220 // Perform inner-product exchanges
221 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
222
223 rho_new = vExchange[0];
224 mu = vExchange[1];
225 eps = vExchange[2];
226
228
229 // Test if norm is within tolerance
231 {
232 if (m_verbose && m_root)
233 {
234 cout << "CG iterations made = " << m_totalIterations
235 << " using tolerance of " << m_tolerance
236 << " (error = " << sqrt(eps / m_rhs_magnitude)
237 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
238 }
239 break;
240 }
241
242 // Compute search direction and solution coefficients
243 beta = rho_new / rho;
244 alpha = rho_new / (mu - rho_new * beta / alpha);
245 rho = rho_new;
246 }
247}
248} // namespace LibUtilities
249} // namespace Nektar
#define ROOTONLY_NEKERROR(type, msg)
Definition: ErrorUtil.hpp:212
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
void DoConjugateGradient(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir)
Actual iterative solve.
virtual 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)
virtual void v_InitObject() override
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:306
LibUtilities::CommSharedPtr m_rowComm
Communicate.
Definition: NekSys.h:302
bool m_verbose
Verbose.
Definition: NekSys.h:308
NekDouble m_tolerance
Tolerance of iterative solver.
Definition: NekSys.h:300
NekSysOperators m_operator
Operators.
Definition: NekSys.h:311
int m_maxiter
Maximum iterations.
Definition: NekSys.h:298
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:150
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:143
std::shared_ptr< SessionReader > SessionReaderSharedPtr
NekLinSysIterFactory & GetNekLinSysIterFactory()
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:61
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
static const NekDouble kNekUnsetDouble
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.cpp:617
T Dot2(int n, const T *w, const T *x, const int *y)
dot2 (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:1142
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294