Nektar++
NekLinSysIterCGLoc.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NekLinSysIterCGLoc.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: NekLinSysIterCGLoc definition
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38using namespace std;
39
41{
42/**
43 * @class NekLinSysIterCGLoc
44 *
45 * Solves a linear system using iterative methods.
46 */
49 "ConjugateGradientLoc", NekLinSysIterCGLoc::create,
50 "NekLinSysIterCG solver in Local Space.");
51
54 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
55 const NekSysKey &pKey)
56 : NekLinSysIter(pSession, vRowComm, nDimen, pKey)
57{
58 m_isLocal = true;
59}
60
62{
64}
65
66/**
67 *
68 */
70 const int nLocal, const Array<OneD, const NekDouble> &pInput,
71 Array<OneD, NekDouble> &pOutput, [[maybe_unused]] const int nDir,
72 const NekDouble tol, [[maybe_unused]] const NekDouble factor)
73{
74 m_tolerance = max(tol, 1.0E-16);
75
76 DoConjugateGradient(nLocal, pInput, pOutput);
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 nLocal, const Array<OneD, const NekDouble> &pInput,
97{
98 // Allocate array storage
99 Array<OneD, NekDouble> w_A(nLocal, 0.0);
100 Array<OneD, NekDouble> s_A(nLocal, 0.0);
101 Array<OneD, NekDouble> p_A(nLocal, 0.0);
102 Array<OneD, NekDouble> r_A(nLocal, 0.0);
103 Array<OneD, NekDouble> q_A(nLocal, 0.0);
104 Array<OneD, NekDouble> wk(nLocal, 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(nLocal, pInput, 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(nLocal, pOutput, 1);
121
122 // evaluate initial residual error for exit check
123 m_operator.DoAssembleLoc(r_A, wk, true);
124 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
125
126 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
127
128 eps = vExchange[2];
129
131 {
132 Set_Rhs_Magnitude(pInput);
133 }
134
136
137 // If input residual is less than tolerance skip solve.
139 {
140 if (m_verbose && m_root)
141 {
142 cout << "CG iterations made = " << m_totalIterations
143 << " using tolerance of " << m_tolerance
144 << " (error = " << sqrt(eps / m_rhs_magnitude)
145 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
146 }
147 return;
148 }
149
150 m_operator.DoNekSysPrecon(r_A, w_A, true);
151 m_operator.DoNekSysLhsEval(w_A, s_A);
152
153 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
154 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
155
156 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
157
158 rho = vExchange[0];
159 mu = vExchange[1];
160 beta = 0.0;
161 alpha = rho / mu;
163
164 // Continue until convergence
165 while (true)
166 {
168 {
169 if (m_root)
170 {
171 cout << "CG iterations made = " << m_totalIterations
172 << " using tolerance of " << m_tolerance
173 << " (error = " << sqrt(eps / m_rhs_magnitude)
174 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
175 }
177 "Exceeded maximum number of iterations");
178 }
179
180 // Compute new search direction p_k, q_k
181 Vmath::Svtvp(nLocal, beta, p_A, 1, w_A, 1, p_A, 1);
182 Vmath::Svtvp(nLocal, beta, q_A, 1, s_A, 1, q_A, 1);
183
184 // Update solution x_{k+1}
185 Vmath::Svtvp(nLocal, alpha, p_A, 1, pOutput, 1, pOutput, 1);
186
187 // Update residual vector r_{k+1}
188 Vmath::Svtvp(nLocal, -alpha, q_A, 1, r_A, 1, r_A, 1);
189
190 // Apply preconditioner
191 m_operator.DoNekSysPrecon(r_A, w_A, true);
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::Dot(nLocal, r_A, w_A);
198
199 // <s_{k+1}, w_{k+1}>
200 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
201
202 // <r_{k+1}, r_{k+1}>
203 m_operator.DoAssembleLoc(r_A, wk, true);
204 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
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_tolerance
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.
Definition: NekFactory.hpp:197
void DoConjugateGradient(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
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
static NekLinSysIterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
NekLinSysIterCGLoc(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
Constructor for full direct matrix solve.
NekDouble m_rhs_magnitude
Dot product of rhs to normalise stopping criterion.
void Set_Rhs_Magnitude(const NekVector< NekDouble > &pIn)
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 DoAssembleLoc(InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition: NekSys.h:161
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 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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294