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)
72{
73 DoConjugateGradient(nGlobal, pInput, pOutput, nDir);
74
75 return m_totalIterations;
76}
77
78/**  
79 * Solve a global linear system using the conjugate gradient method.  
80 * We solve only for the non-Dirichlet modes. The operator is evaluated  
81 * using an auxiliary function m_operator.DoNekSysLhsEval defined by the  
82 * specific solver. Distributed math routines are used to support  
83 * parallel execution of the solver.  
84 *  
85 * The implemented algorithm uses a reduced-communication reordering of  
86 * the standard PCG method (Demmel, Heath and Vorst, 1993)  
87 *  
88 * @param pInput Input residual of all DOFs.  
89 * @param pOutput Solution vector of all DOFs.  
90 */
92 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
93 Array<OneD, NekDouble> &pOutput, const int nDir)
94{
95 // Get vector sizes
96 int nNonDir = nGlobal - nDir;
97
98 // Allocate array storage
99 Array<OneD, NekDouble> w_A(nGlobal, 0.0);
100 Array<OneD, NekDouble> s_A(nGlobal, 0.0);
101 Array<OneD, NekDouble> p_A(nNonDir, 0.0);
102 Array<OneD, NekDouble> r_A(nNonDir, 0.0);
103 Array<OneD, NekDouble> q_A(nNonDir, 0.0);
105
106 NekDouble alpha;
108 NekDouble rho;
109 NekDouble rho_new;
110 NekDouble mu;
111 NekDouble eps;
112 Array<OneD, NekDouble> vExchange(3, 0.0);
113
114 // Copy initial residual from input
115 Vmath::Vcopy(nNonDir, pInput + nDir, 1, r_A, 1);
116
117 // Zero homogeneous out array ready for solution updates
118 // Should not be earlier in case input vector is same as
119 // output and above copy has been peformed
120 Vmath::Zero(nNonDir, tmp = pOutput + nDir, 1);
121
122 // Evaluate initial residual error for exit check
123 vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
124
125 m_rowComm->AllReduce(vExchange[2], Nektar::LibUtilities::ReduceSum);
126
127 eps = vExchange[2];
128
130 {
131 Set_Rhs_Magnitude(pInput);
132 }
133
135
136 // If input residual is less than tolerance skip solve.
138 {
139 if (m_verbose && m_root)
140 {
141 cout << "CG iterations made = " << m_totalIterations
142 << " using tolerance of " << m_NekLinSysTolerance
143 << " (error = " << sqrt(eps / m_rhs_magnitude)
144 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
145 }
146 return;
147 }
148
149 m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
150 m_operator.DoNekSysLhsEval(w_A, s_A);
151
152 vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
153 vExchange[1] = Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
154
155 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
156
157 rho = vExchange[0];
158 mu = vExchange[1];
159 beta = 0.0;
160 alpha = rho / mu;
162
163 // Continue until convergence
164 while (true)
165 {
167 {
168 if (m_root)
169 {
170 cout << "CG iterations made = " << m_totalIterations
171 << " using tolerance of " << m_NekLinSysTolerance
172 << " (error = " << sqrt(eps / m_rhs_magnitude)
173 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
174 }
176 "Exceeded maximum number of iterations");
177 }
178
179 // Compute new search direction p_k, q_k
180 Vmath::Svtvp(nNonDir, beta, &p_A[0], 1, &w_A[nDir], 1, &p_A[0], 1);
181 Vmath::Svtvp(nNonDir, beta, &q_A[0], 1, &s_A[nDir], 1, &q_A[0], 1);
182
183 // Update solution x_{k+1}
184 Vmath::Svtvp(nNonDir, alpha, &p_A[0], 1, &pOutput[nDir], 1,
185 &pOutput[nDir], 1);
186
187 // Update residual vector r_{k+1}
188 Vmath::Svtvp(nNonDir, -alpha, &q_A[0], 1, &r_A[0], 1, &r_A[0], 1);
189
190 // Apply preconditioner
191 m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
192
193 // Perform the method-specific matrix-vector multiply operation.
194 m_operator.DoNekSysLhsEval(w_A, s_A);
195
196 // <r_{k+1}, w_{k+1}>
197 vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
198
199 // <s_{k+1}, w_{k+1}>
200 vExchange[1] =
201 Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
202
203 // <r_{k+1}, r_{k+1}>
204 vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
205
206 // Perform inner-product exchanges
207 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
208
209 rho_new = vExchange[0];
210 mu = vExchange[1];
211 eps = vExchange[2];
212
214
215 // Test if norm is within tolerance
217 {
218 if (m_verbose && m_root)
219 {
220 cout << "CG iterations made = " << m_totalIterations
221 << " using tolerance of " << m_NekLinSysTolerance
222 << " (error = " << sqrt(eps / m_rhs_magnitude)
223 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
224 }
225 break;
226 }
227
228 // Compute search direction and solution coefficients
229 beta = rho_new / rho;
230 alpha = rho_new / (mu - rho_new * beta / alpha);
231 rho = rho_new;
232 }
233}
234} // 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.
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) 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:291
NekSysOperators m_operator
Definition: NekSys.h:298
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:149
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:142
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
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294