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
38
39using namespace std;
40
41namespace Nektar
42{
43namespace LibUtilities
44{
45/**
46 * @class NekLinSysIterCGLoc
47 *
48 * Solves a linear system using iterative methods.
49 */
52 "ConjugateGradientLoc", NekLinSysIterCGLoc::create,
53 "NekLinSysIterCG solver in Local Space.");
54
57 const LibUtilities::CommSharedPtr &vRowComm, const int nDimen,
58 const NekSysKey &pKey)
59 : NekLinSysIter(pSession, vRowComm, nDimen, pKey)
60{
61 m_isLocal = true;
62}
63
65{
67}
68
70{
71}
72
73/**
74 *
75 */
77 const int nLocal, const Array<OneD, const NekDouble> &pInput,
78 Array<OneD, NekDouble> &pOutput, const int nDir, const NekDouble tol,
79 const NekDouble factor)
80{
81 boost::ignore_unused(tol, nDir);
82
83 m_tolerance = max(tol, 1.0E-16);
84 m_prec_factor = factor;
85
86 DoConjugateGradient(nLocal, pInput, pOutput);
87
88 return m_totalIterations;
89}
90
91/**  
92 * Solve a global linear system using the conjugate gradient method.  
93 * We solve only for the non-Dirichlet modes. The operator is evaluated  
94 * using an auxiliary function m_operator.DoNekSysLhsEval defined by the  
95 * specific solver. Distributed math routines are used to support  
96 * parallel execution of the solver.  
97 *  
98 * The implemented algorithm uses a reduced-communication reordering of  
99 * the standard PCG method (Demmel, Heath and Vorst, 1993)  
100 *  
101 * @param pInput Input residual of all DOFs.  
102 * @param pOutput Solution vector of all DOFs.  
103 */
105 const int nLocal, const Array<OneD, const NekDouble> &pInput,
106 Array<OneD, NekDouble> &pOutput)
107{
108 // Allocate array storage
109 Array<OneD, NekDouble> w_A(nLocal, 0.0);
110 Array<OneD, NekDouble> s_A(nLocal, 0.0);
111 Array<OneD, NekDouble> p_A(nLocal, 0.0);
112 Array<OneD, NekDouble> r_A(nLocal, 0.0);
113 Array<OneD, NekDouble> q_A(nLocal, 0.0);
114 Array<OneD, NekDouble> wk(nLocal, 0.0);
115
116 int k;
117 NekDouble alpha;
119 NekDouble rho;
120 NekDouble rho_new;
121 NekDouble mu;
122 NekDouble eps;
123 NekDouble min_resid;
124 Array<OneD, NekDouble> vExchange(3, 0.0);
125
126 // Copy initial residual from input
127 Vmath::Vcopy(nLocal, pInput, 1, r_A, 1);
128
129 // zero homogeneous out array ready for solution updates
130 // Should not be earlier in case input vector is same as
131 // output and above copy has been peformed
132 Vmath::Zero(nLocal, pOutput, 1);
133
134 // evaluate initial residual error for exit check
135 m_operator.DoAssembleLoc(r_A, wk, true);
136 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
137 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
138
139 eps = vExchange[2];
140
142 {
143 Set_Rhs_Magnitude(pInput);
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, w_A, true);
162 m_operator.DoNekSysLhsEval(w_A, s_A);
163
164 k = 0;
165
166 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
167 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
168 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
169
170 rho = vExchange[0];
171 mu = vExchange[1];
172 min_resid = m_rhs_magnitude;
173 beta = 0.0;
174 alpha = rho / mu;
176
177 // Continue until convergence
178 while (true)
179 {
180 if (k >= m_maxiter)
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(nLocal, beta, p_A, 1, w_A, 1, p_A, 1);
195 Vmath::Svtvp(nLocal, beta, q_A, 1, s_A, 1, q_A, 1);
196
197 // Update solution x_{k+1}
198 Vmath::Svtvp(nLocal, alpha, p_A, 1, pOutput, 1, pOutput, 1);
199
200 // Update residual vector r_{k+1}
201 Vmath::Svtvp(nLocal, -alpha, q_A, 1, r_A, 1, r_A, 1);
202
203 // Apply preconditioner
204 m_operator.DoNekSysPrecon(r_A, w_A, true);
205
206 // Perform the method-specific matrix-vector multiply operation.
207 m_operator.DoNekSysLhsEval(w_A, s_A);
208
209 // <r_{k+1}, w_{k+1}>
210 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
211 // <s_{k+1}, w_{k+1}>
212 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
213 // <r_{k+1}, r_{k+1}>
214 m_operator.DoAssembleLoc(r_A, wk, true);
215 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
216
217 // Perform inner-product exchanges
218 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
219
220 rho_new = vExchange[0];
221 mu = vExchange[1];
222 eps = vExchange[2];
223
225
226 // Test if norm is within tolerance
228 {
229 if (m_verbose && m_root)
230 {
231 cout << "CG iterations made = " << m_totalIterations
232 << " using tolerance of " << m_tolerance
233 << " (error = " << sqrt(eps / m_rhs_magnitude)
234 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
235 }
236 break;
237 }
238 min_resid = min(min_resid, eps);
239
240 // Compute search direction and solution coefficients
241 beta = rho_new / rho;
242 alpha = rho_new / (mu - rho_new * beta / alpha);
243 rho = rho_new;
244 k++;
245 }
246}
247} // namespace LibUtilities
248} // 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)
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
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.
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)
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 DoAssembleLoc(InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition: NekSys.h:163
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 Dot(int n, const T *w, const T *x)
dot (vector times vector): z = w*x
Definition: Vmath.cpp:1095
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