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 
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 namespace LibUtilities
44 {
45 /**
46  * @class NekLinSysIterCG
47  *
48  * Solves a linear system using iterative methods.
49  */
50 string NekLinSysIterCG::className =
52  "ConjugateGradient", NekLinSysIterCG::create,
53  "NekLinSysIterCG solver.");
54 
55 NekLinSysIterCG::NekLinSysIterCG(
57  const LibUtilities::CommSharedPtr &vComm, const int nDimen,
58  const NekSysKey &pKey)
59  : NekLinSysIter(pSession, vComm, nDimen, pKey)
60 {
61 }
62 
64 {
66 }
67 
69 {
70 }
71 
72 /**
73  *
74  */
75 int NekLinSysIterCG::v_SolveSystem(const int nGlobal,
76  const Array<OneD, const NekDouble> &pInput,
77  Array<OneD, NekDouble> &pOutput,
78  const int nDir, const NekDouble tol,
79  const NekDouble factor)
80 {
81  boost::ignore_unused(tol);
82 
83  m_tolerance = max(tol, 1.0E-16);
84  m_prec_factor = factor;
85 
86  DoConjugateGradient(nGlobal, pInput, pOutput, nDir);
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 nGlobal, const Array<OneD, const NekDouble> &pInput,
106  Array<OneD, NekDouble> &pOutput, const int nDir)
107 {
108  // Get vector sizes
109  int nNonDir = nGlobal - nDir;
110 
111  // Allocate array storage
112  Array<OneD, NekDouble> w_A(nGlobal, 0.0);
113  Array<OneD, NekDouble> s_A(nGlobal, 0.0);
114  Array<OneD, NekDouble> p_A(nNonDir, 0.0);
115  Array<OneD, NekDouble> r_A(nNonDir, 0.0);
116  Array<OneD, NekDouble> q_A(nNonDir, 0.0);
118 
119  // Create NekVector wrappers for linear algebra operations
120  NekVector<NekDouble> in(nNonDir, pInput + nDir, eWrapper);
121  NekVector<NekDouble> out(nNonDir, tmp = pOutput + nDir, eWrapper);
122  NekVector<NekDouble> w(nNonDir, tmp = w_A + nDir, eWrapper);
123  NekVector<NekDouble> s(nNonDir, tmp = s_A + nDir, eWrapper);
124  NekVector<NekDouble> p(nNonDir, p_A, eWrapper);
125  NekVector<NekDouble> r(nNonDir, r_A, eWrapper);
126  NekVector<NekDouble> q(nNonDir, q_A, eWrapper);
127 
128  int k;
129  NekDouble alpha;
130  NekDouble beta;
131  NekDouble rho;
132  NekDouble rho_new;
133  NekDouble mu;
134  NekDouble eps;
135  NekDouble min_resid;
136  Array<OneD, NekDouble> vExchange(3, 0.0);
137 
138  // Copy initial residual from input
139  r = in;
140  // zero homogeneous out array ready for solution updates
141  // Should not be earlier in case input vector is same as
142  // output and above copy has been peformed
143  Vmath::Zero(nNonDir, tmp = pOutput + nDir, 1);
144 
145  // evaluate initial residual error for exit check
146  vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
147 
148  m_Comm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
149 
150  eps = vExchange[2];
151 
153  {
154  NekVector<NekDouble> inGlob(nGlobal, pInput, eWrapper);
155  Set_Rhs_Magnitude(inGlob);
156  }
157 
158  m_totalIterations = 0;
159 
160  // If input residual is less than tolerance skip solve.
162  {
163  if (m_verbose && m_root)
164  {
165  cout << "CG iterations made = " << m_totalIterations
166  << " using tolerance of " << m_tolerance
167  << " (error = " << sqrt(eps / m_rhs_magnitude)
168  << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
169  }
170  return;
171  }
172 
173  m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
174 
175  m_operator.DoNekSysLhsEval(w_A, s_A);
176 
177  k = 0;
178 
179  vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
180 
181  vExchange[1] = Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
182 
183  m_Comm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
184 
185  rho = vExchange[0];
186  mu = vExchange[1];
187  min_resid = m_rhs_magnitude;
188  beta = 0.0;
189  alpha = rho / mu;
190  m_totalIterations = 1;
191 
192  // Continue until convergence
193  while (true)
194  {
195  if (k >= m_maxiter)
196  {
197  if (m_root)
198  {
199  cout << "CG iterations made = " << m_totalIterations
200  << " using tolerance of " << m_tolerance
201  << " (error = " << sqrt(eps / m_rhs_magnitude)
202  << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
203  }
205  "Exceeded maximum number of iterations");
206  }
207 
208  // Compute new search direction p_k, q_k
209  Vmath::Svtvp(nNonDir, beta, &p_A[0], 1, &w_A[nDir], 1, &p_A[0], 1);
210  Vmath::Svtvp(nNonDir, beta, &q_A[0], 1, &s_A[nDir], 1, &q_A[0], 1);
211 
212  // Update solution x_{k+1}
213  Vmath::Svtvp(nNonDir, alpha, &p_A[0], 1, &pOutput[nDir], 1,
214  &pOutput[nDir], 1);
215 
216  // Update residual vector r_{k+1}
217  Vmath::Svtvp(nNonDir, -alpha, &q_A[0], 1, &r_A[0], 1, &r_A[0], 1);
218 
219  // Apply preconditioner
220  m_operator.DoNekSysPrecon(r_A, tmp = w_A + nDir);
221 
222  // Perform the method-specific matrix-vector multiply operation.
223  m_operator.DoNekSysLhsEval(w_A, s_A);
224 
225  // <r_{k+1}, w_{k+1}>
226  vExchange[0] = Vmath::Dot2(nNonDir, r_A, w_A + nDir, m_map + nDir);
227  // <s_{k+1}, w_{k+1}>
228  vExchange[1] =
229  Vmath::Dot2(nNonDir, s_A + nDir, w_A + nDir, m_map + nDir);
230 
231  // <r_{k+1}, r_{k+1}>
232  vExchange[2] = Vmath::Dot2(nNonDir, r_A, r_A, m_map + nDir);
233 
234  // Perform inner-product exchanges
235  m_Comm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
236 
237  rho_new = vExchange[0];
238  mu = vExchange[1];
239  eps = vExchange[2];
240 
242 
243  // Test if norm is within tolerance
245  {
246  if (m_verbose && m_root)
247  {
248  cout << "CG iterations made = " << m_totalIterations
249  << " using tolerance of " << m_tolerance
250  << " (error = " << sqrt(eps / m_rhs_magnitude)
251  << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
252  }
253  break;
254  }
255  min_resid = min(min_resid, eps);
256 
257  // Compute search direction and solution coefficients
258  beta = rho_new / rho;
259  alpha = rho_new / (mu - rho_new * beta / alpha);
260  rho = rho_new;
261  k++;
262  }
263 }
264 } // namespace LibUtilities
265 } // namespace Nektar
#define ROOTONLY_NEKERROR(type, msg)
Definition: ErrorUtil.hpp:213
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
void DoConjugateGradient(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir)
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)
NekDouble m_rhs_magnitude
Dot product of rhs to normalise stopping criterion.
void Set_Rhs_Magnitude(const NekVector< NekDouble > &pIn)
Array< OneD, int > m_map
Global to universal unique map.
bool m_root
Root if parallel.
Definition: NekSys.h:289
bool m_verbose
Verbose.
Definition: NekSys.h:291
NekDouble m_tolerance
Tolerance of iterative solver.
Definition: NekSys.h:283
NekSysOperators m_operator
Operators.
Definition: NekSys.h:294
LibUtilities::CommSharedPtr m_Comm
Communicate.
Definition: NekSys.h:285
int m_maxiter
Maximum iterations.
Definition: NekSys.h:281
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:142
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:135
std::shared_ptr< SessionReader > SessionReaderSharedPtr
NekLinSysIterFactory & GetNekLinSysIterFactory()
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
static const NekDouble kNekUnsetDouble
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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:565
T Dot2(int n, const T *w, const T *x, const int *y)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:1084
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267