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{
73 DoConjugateGradient(nLocal, pInput, pOutput);
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 nLocal, const Array<OneD, const NekDouble> &pInput,
94{
95 // Allocate array storage
96 Array<OneD, NekDouble> w_A(nLocal, 0.0);
97 Array<OneD, NekDouble> s_A(nLocal, 0.0);
98 Array<OneD, NekDouble> p_A(nLocal, 0.0);
99 Array<OneD, NekDouble> r_A(nLocal, 0.0);
100 Array<OneD, NekDouble> q_A(nLocal, 0.0);
101 Array<OneD, NekDouble> wk(nLocal, 0.0);
102
103 NekDouble alpha;
105 NekDouble rho;
106 NekDouble rho_new;
107 NekDouble mu;
108 NekDouble eps;
109 Array<OneD, NekDouble> vExchange(3, 0.0);
110
111 // Copy initial residual from input
112 Vmath::Vcopy(nLocal, pInput, 1, r_A, 1);
113
114 // zero homogeneous out array ready for solution updates
115 // Should not be earlier in case input vector is same as
116 // output and above copy has been peformed
117 Vmath::Zero(nLocal, pOutput, 1);
118
119 // evaluate initial residual error for exit check
120 m_operator.DoAssembleLoc(r_A, wk, true);
121 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
122
123 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
124
125 eps = vExchange[2];
126
128 {
129 Set_Rhs_Magnitude(pInput);
130 }
131
133
134 // If input residual is less than tolerance skip solve.
136 {
137 if (m_verbose && m_root)
138 {
139 cout << "CG iterations made = " << m_totalIterations
140 << " using tolerance of " << m_NekLinSysTolerance
141 << " (error = " << sqrt(eps / m_rhs_magnitude)
142 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
143 }
144 return;
145 }
146
147 m_operator.DoNekSysPrecon(r_A, w_A, true);
148 m_operator.DoNekSysLhsEval(w_A, s_A);
149
150 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
151 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
152
153 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
154
155 rho = vExchange[0];
156 mu = vExchange[1];
157 beta = 0.0;
158 alpha = rho / mu;
160
161 // Continue until convergence
162 while (true)
163 {
165 {
166 if (m_root)
167 {
168 cout << "CG iterations made = " << m_totalIterations
169 << " using tolerance of " << m_NekLinSysTolerance
170 << " (error = " << sqrt(eps / m_rhs_magnitude)
171 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
172 }
174 "Exceeded maximum number of iterations");
175 }
176
177 // Compute new search direction p_k, q_k
178 Vmath::Svtvp(nLocal, beta, p_A, 1, w_A, 1, p_A, 1);
179 Vmath::Svtvp(nLocal, beta, q_A, 1, s_A, 1, q_A, 1);
180
181 // Update solution x_{k+1}
182 Vmath::Svtvp(nLocal, alpha, p_A, 1, pOutput, 1, pOutput, 1);
183
184 // Update residual vector r_{k+1}
185 Vmath::Svtvp(nLocal, -alpha, q_A, 1, r_A, 1, r_A, 1);
186
187 // Apply preconditioner
188 m_operator.DoNekSysPrecon(r_A, w_A, true);
189
190 // Perform the method-specific matrix-vector multiply operation.
191 m_operator.DoNekSysLhsEval(w_A, s_A);
192
193 // <r_{k+1}, w_{k+1}>
194 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
195
196 // <s_{k+1}, w_{k+1}>
197 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
198
199 // <r_{k+1}, r_{k+1}>
200 m_operator.DoAssembleLoc(r_A, wk, true);
201 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
202
203 // Perform inner-product exchanges
204 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
205
206 rho_new = vExchange[0];
207 mu = vExchange[1];
208 eps = vExchange[2];
209
211
212 // Test if norm is within tolerance
214 {
215 if (m_verbose && m_root)
216 {
217 cout << "CG iterations made = " << m_totalIterations
218 << " using tolerance of " << m_NekLinSysTolerance
219 << " (error = " << sqrt(eps / m_rhs_magnitude)
220 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
221 }
222 break;
223 }
224
225 // Compute search direction and solution coefficients
226 beta = rho_new / rho;
227 alpha = rho_new / (mu - rho_new * beta / alpha);
228 rho = rho_new;
229 }
230}
231} // 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)
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 Set_Rhs_Magnitude(const Array< OneD, NekDouble > &pIn)
LibUtilities::CommSharedPtr m_rowComm
Definition: NekSys.h:291
NekSysOperators m_operator
Definition: NekSys.h:298
void DoAssembleLoc(InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition: NekSys.h:162
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 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:294