Nektar++
Loading...
Searching...
No Matches
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
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{
73 DoConjugateGradient(nLocal, pInput, pOutput);
74
75 return m_totalIterations;
76}
77
78void NekLinSysIterCGLoc::v_DoIterate(const int nGlobal,
79 const Array<OneD, NekDouble> &rhs,
81 [[maybe_unused]] const int nDir,
82 NekDouble &err, int &iter)
83{
84 DoConjugateGradient(nGlobal, rhs, x);
85 iter = m_totalIterations;
86 err = m_finalError;
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 nLocal, const Array<OneD, const NekDouble> &pInput,
104 Array<OneD, NekDouble> &pOutput)
105{
106 // Allocate array storage
107 Array<OneD, NekDouble> w_A(nLocal, 0.0);
108 Array<OneD, NekDouble> s_A(nLocal, 0.0);
109 Array<OneD, NekDouble> p_A(nLocal, 0.0);
110 Array<OneD, NekDouble> r_A(nLocal, 0.0);
111 Array<OneD, NekDouble> q_A(nLocal, 0.0);
112 Array<OneD, NekDouble> wk(nLocal, 0.0);
113
114 NekDouble alpha;
116 NekDouble rho;
117 NekDouble rho_new;
118 NekDouble mu;
119 NekDouble eps;
120 Array<OneD, NekDouble> vExchange(3, 0.0);
121
122 // Copy initial residual from input
123 Vmath::Vcopy(nLocal, pInput, 1, r_A, 1);
124
125 // zero homogeneous out array ready for solution updates
126 // Should not be earlier in case input vector is same as
127 // output and above copy has been peformed
128 Vmath::Zero(nLocal, pOutput, 1);
129
130 // evaluate initial residual error for exit check
131 m_operator.DoAssembleLoc(r_A, wk, true);
132 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
133
134 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
135
136 eps = vExchange[2];
137
139 {
140 Set_Rhs_Magnitude(pInput);
141 }
142
144
145 // If input residual is less than tolerance skip solve.
147 {
149 if (m_verbose && m_root)
150 {
151 cout << "CG iterations made = " << m_totalIterations
152 << " using tolerance of " << m_NekLinSysTolerance
153 << " (error = " << m_finalError
154 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
155 }
156 return;
157 }
158
159 m_operator.DoNekSysPrecon(r_A, w_A, true);
160 m_operator.DoNekSysLhsEval(w_A, s_A);
161
162 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
163 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
164
165 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
166
167 rho = vExchange[0];
168 mu = vExchange[1];
169 beta = 0.0;
170 alpha = rho / mu;
172
173 // Continue until convergence
174 while (true)
175 {
177 {
179 if (m_root)
180 {
181 cout << "CG iterations made = " << m_totalIterations
182 << " using tolerance of " << m_NekLinSysTolerance
183 << " (error = " << m_finalError
184 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")"
185 << " WARNING: Exceeded maxIt" << endl;
186 }
187 break;
188 }
189
190 // Compute new search direction p_k, q_k
191 Vmath::Svtvp(nLocal, beta, p_A, 1, w_A, 1, p_A, 1);
192 Vmath::Svtvp(nLocal, beta, q_A, 1, s_A, 1, q_A, 1);
193
194 // Update solution x_{k+1}
195 Vmath::Svtvp(nLocal, alpha, p_A, 1, pOutput, 1, pOutput, 1);
196
197 // Update residual vector r_{k+1}
198 Vmath::Svtvp(nLocal, -alpha, q_A, 1, r_A, 1, r_A, 1);
199
200 // Apply preconditioner
201 m_operator.DoNekSysPrecon(r_A, w_A, true);
202
203 // Perform the method-specific matrix-vector multiply operation.
204 m_operator.DoNekSysLhsEval(w_A, s_A);
205
206 // <r_{k+1}, w_{k+1}>
207 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
208
209 // <s_{k+1}, w_{k+1}>
210 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
211
213 {
214 // <r_{k+1}, r_{k+1}>
215 m_operator.DoAssembleLoc(r_A, wk, true);
216 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
217 }
218 // Perform inner-product exchanges
219 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
220
221 rho_new = vExchange[0];
222 mu = vExchange[1];
223 eps = vExchange[2];
224
226
227 // Test if norm is within tolerance
229 {
231 if (m_verbose && m_root)
232 {
233 cout << "CG iterations made = " << m_totalIterations
234 << " using tolerance of " << m_NekLinSysTolerance
235 << " (error = " << m_finalError
236 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
237 }
238 break;
239 }
240
241 // Compute search direction and solution coefficients
242 beta = rho_new / rho;
243 alpha = rho_new / (mu - rho_new * beta / alpha);
244 rho = rho_new;
245 }
246}
247} // 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)
Actual iterative solve.
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.
int v_SolveSystem(const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir) override
void v_DoIterate(const int nGlobal, const Array< OneD, NekDouble > &rhs, Array< OneD, NekDouble > &x, const int nDir, NekDouble &err, int &iter) override
void Set_Rhs_Magnitude(const Array< OneD, NekDouble > &pIn)
LibUtilities::CommSharedPtr m_rowComm
Definition NekSys.h:299
NekSysOperators m_operator
Definition NekSys.h:306
void DoAssembleLoc(InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition NekSys.h:168
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 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