Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GlobalLinSysIterative.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: GlobalLinSysIterative.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: GlobalLinSysIterative definition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 namespace Nektar
39 {
40  namespace MultiRegions
41  {
42  /**
43  * @class GlobalLinSysIterative
44  *
45  * Solves a linear system using iterative methods.
46  */
47 
48  /// Constructor for full direct matrix solve.
50  const GlobalLinSysKey &pKey,
51  const boost::weak_ptr<ExpList> &pExpList,
52  const boost::shared_ptr<AssemblyMap>
53  &pLocToGloMap)
54  : GlobalLinSys(pKey, pExpList, pLocToGloMap),
55  m_rhs_magnitude(NekConstants::kNekUnsetDouble),
56  m_rhs_mag_sm(0.9),
58  m_totalIterations(0),
59  m_useProjection(false),
60  m_numPrevSols(0)
61  {
63  = pExpList.lock()->GetSession();
64 
65  m_tolerance = pLocToGloMap->GetIterativeTolerance();
66  m_maxiter = pLocToGloMap->GetMaxIterations();
67 
68  LibUtilities::CommSharedPtr vComm = m_expList.lock()->GetComm()->GetRowComm();
69  m_root = (vComm->GetRank())? false : true;
70  m_verbose = (vSession->DefinesCmdLineArgument("verbose"))? true :false;
71 
72  int successiveRHS;
73 
74  if((successiveRHS = pLocToGloMap->GetSuccessiveRHS()))
75  {
76  m_prevLinSol.set_capacity(successiveRHS);
77  m_useProjection = true;
78  }
79  else
80  {
81  m_useProjection = false;
82  }
83  }
84 
86  {
87  }
88 
89 
90  /**
91  *
92  */
94  const int nGlobal,
95  const Array<OneD,const NekDouble> &pInput,
96  Array<OneD, NekDouble> &pOutput,
97  const AssemblyMapSharedPtr &plocToGloMap,
98  const int nDir)
99  {
100  if (m_useProjection)
101  {
102  DoAconjugateProjection(nGlobal, pInput, pOutput, plocToGloMap, nDir);
103  }
104  else
105  {
106  // applying plain Conjugate Gradient
107  DoConjugateGradient(nGlobal, pInput, pOutput, plocToGloMap, nDir);
108  }
109  }
110 
111 
112  /**
113  * This method implements A-conjugate projection technique
114  * in order to speed up successive linear solves with
115  * right-hand sides arising from time-dependent discretisations.
116  * (P.F.Fischer, Comput. Methods Appl. Mech. Engrg. 163, 1998)
117  */
119  const int nGlobal,
120  const Array<OneD,const NekDouble> &pInput,
121  Array<OneD, NekDouble> &pOutput,
122  const AssemblyMapSharedPtr &plocToGloMap,
123  const int nDir)
124  {
125  // Get the communicator for performing data exchanges
127  = m_expList.lock()->GetComm()->GetRowComm();
128 
129  // Get vector sizes
130  int nNonDir = nGlobal - nDir;
132 
133  if (0 == m_numPrevSols)
134  {
135  // no previous solutions found, call CG
136 
137  DoConjugateGradient(nGlobal, pInput, pOutput, plocToGloMap, nDir);
138 
139  UpdateKnownSolutions(nGlobal, pOutput, nDir);
140  }
141  else
142  {
143  // Create NekVector wrappers for linear algebra operations
144  NekVector<NekDouble> b (nNonDir, pInput + nDir, eWrapper);
145  NekVector<NekDouble> x (nNonDir, tmp = pOutput + nDir, eWrapper);
146 
147  // check the input vector (rhs) is not zero
148 
149  NekDouble rhsNorm = Vmath::Dot2(nNonDir,
150  pInput + nDir,
151  pInput + nDir,
152  m_map + nDir);
153 
154  vComm->AllReduce(rhsNorm, Nektar::LibUtilities::ReduceSum);
155 
156  if (rhsNorm < NekConstants::kNekZeroTol)
157  {
158  Array<OneD, NekDouble> tmp = pOutput+nDir;
159  Vmath::Zero(nNonDir, tmp, 1);
160  return;
161  }
162 
163  // Allocate array storage
164  Array<OneD, NekDouble> px_s (nGlobal, 0.0);
165  Array<OneD, NekDouble> pb_s (nGlobal, 0.0);
166  Array<OneD, NekDouble> tmpAx_s (nGlobal, 0.0);
167  Array<OneD, NekDouble> tmpx_s (nGlobal, 0.0);
168 
169  NekVector<NekDouble> pb (nNonDir, tmp = pb_s + nDir, eWrapper);
170  NekVector<NekDouble> px (nNonDir, tmp = px_s + nDir, eWrapper);
171  NekVector<NekDouble> tmpAx (nNonDir, tmp = tmpAx_s + nDir, eWrapper);
172  NekVector<NekDouble> tmpx (nNonDir, tmp = tmpx_s + nDir, eWrapper);
173 
174 
175  // notation follows the paper cited:
176  // \alpha_i = \tilda{x_i}^T b^n
177  // projected x, px = \sum \alpha_i \tilda{x_i}
178 
179  Array<OneD, NekDouble> alpha (m_prevLinSol.size(), 0.0);
180  for (int i = 0; i < m_prevLinSol.size(); i++)
181  {
182  alpha[i] = Vmath::Dot2(nNonDir,
183  m_prevLinSol[i],
184  pInput + nDir,
185  m_map + nDir);
186  }
187  vComm->AllReduce(alpha, Nektar::LibUtilities::ReduceSum);
188 
189  for (int i = 0; i < m_prevLinSol.size(); i++)
190  {
191  if (alpha[i] < NekConstants::kNekZeroTol)
192  {
193  continue;
194  }
195 
196  NekVector<NekDouble> xi (nNonDir, m_prevLinSol[i], eWrapper);
197  px += alpha[i] * xi;
198  }
199 
200  // pb = b^n - A px
201  Vmath::Vcopy(nNonDir,
202  pInput.get() + nDir, 1,
203  pb_s.get() + nDir, 1);
204 
205  v_DoMatrixMultiply(px_s, tmpAx_s);
206 
207  pb -= tmpAx;
208 
209 
210  // solve the system with projected rhs
211  DoConjugateGradient(nGlobal, pb_s, tmpx_s, plocToGloMap, nDir);
212 
213 
214  // remainder solution + projection of previous solutions
215  x = tmpx + px;
216 
217  // save the auxiliary solution to prev. known solutions
218  UpdateKnownSolutions(nGlobal, tmpx_s, nDir);
219  }
220  }
221 
222 
223  /**
224  * Calculating A-norm of an input vector,
225  * A-norm(x) := sqrt( < x, Ax > )
226  */
228  const int nGlobal,
229  const Array<OneD,const NekDouble> &in,
230  const int nDir)
231  {
232  // Get the communicator for performing data exchanges
234  = m_expList.lock()->GetComm()->GetRowComm();
235 
236  // Get vector sizes
237  int nNonDir = nGlobal - nDir;
238 
239  // Allocate array storage
240  Array<OneD, NekDouble> tmpAx_s (nGlobal, 0.0);
241 
242  v_DoMatrixMultiply(in, tmpAx_s);
243 
244  NekDouble anorm_sq = Vmath::Dot2(nNonDir,
245  in + nDir,
246  tmpAx_s + nDir,
247  m_map + nDir);
248  vComm->AllReduce(anorm_sq, Nektar::LibUtilities::ReduceSum);
249  return std::sqrt(anorm_sq);
250  }
251 
252  /**
253  * Updates the storage of previously known solutions.
254  * Performs normalisation of input vector wrt A-norm.
255  */
257  const int nGlobal,
258  const Array<OneD,const NekDouble> &newX,
259  const int nDir)
260  {
261  // Get vector sizes
262  int nNonDir = nGlobal - nDir;
263 
264  // Get the communicator for performing data exchanges
266  = m_expList.lock()->GetComm()->GetRowComm();
267 
268  // Check the solution is non-zero
269  NekDouble solNorm = Vmath::Dot2(nNonDir,
270  newX + nDir,
271  newX + nDir,
272  m_map + nDir);
273  vComm->AllReduce(solNorm, Nektar::LibUtilities::ReduceSum);
274 
275  if (solNorm < NekConstants::kNekZeroTol)
276  {
277  return;
278  }
279 
280 
281  // Allocate array storage
282  Array<OneD, NekDouble> tmpAx_s (nGlobal, 0.0);
283  Array<OneD, NekDouble> px_s (nGlobal, 0.0);
284  Array<OneD, NekDouble> tmp1, tmp2;
285 
286  // Create NekVector wrappers for linear algebra operations
287  NekVector<NekDouble> px (nNonDir, tmp1 = px_s + nDir, eWrapper);
288  NekVector<NekDouble> tmpAx (nNonDir, tmp2 = tmpAx_s + nDir, eWrapper);
289 
290 
291  // calculating \tilda{x} - sum \alpha_i\tilda{x}_i
292 
293  Vmath::Vcopy(nNonDir,
294  tmp1 = newX + nDir, 1,
295  tmp2 = px_s + nDir, 1);
296 
297  if (m_prevLinSol.size() > 0)
298  {
299  v_DoMatrixMultiply(newX, tmpAx_s);
300  }
301 
302  Array<OneD, NekDouble> alpha (m_prevLinSol.size(), 0.0);
303  for (int i = 0; i < m_prevLinSol.size(); i++)
304  {
305  alpha[i] = Vmath::Dot2(nNonDir,
306  m_prevLinSol[i],
307  tmpAx_s + nDir,
308  m_map + nDir);
309  }
310  vComm->AllReduce(alpha, Nektar::LibUtilities::ReduceSum);
311 
312  for (int i = 0; i < m_prevLinSol.size(); i++)
313  {
314  if (alpha[i] < NekConstants::kNekZeroTol)
315  {
316  continue;
317  }
318 
319  NekVector<NekDouble> xi (nNonDir, m_prevLinSol[i], eWrapper);
320  px -= alpha[i] * xi;
321  }
322 
323 
324  // Some solutions generated by CG are identical zeros, see
325  // solutions generated for Test_Tet_equitri.xml (IncNavierStokesSolver).
326  // Not going to store identically zero solutions.
327 
328  NekDouble anorm = CalculateAnorm(nGlobal, px_s, nDir);
329  if (anorm < NekConstants::kNekZeroTol)
330  {
331  return;
332  }
333 
334  // normalisation of new solution
335  Vmath::Smul(nNonDir, 1.0/anorm, px_s.get() + nDir, 1, px_s.get() + nDir, 1);
336 
337  // updating storage with non-Dirichlet-dof part of new solution vector
338  m_prevLinSol.push_back(px_s + nDir);
339  m_numPrevSols++;
340  }
341 
342 
343 
344  /**  
345  * Solve a global linear system using the conjugate gradient method.  
346  * We solve only for the non-Dirichlet modes. The operator is evaluated  
347  * using an auxiliary function v_DoMatrixMultiply defined by the  
348  * specific solver. Distributed math routines are used to support  
349  * parallel execution of the solver.  
350  *  
351  * The implemented algorithm uses a reduced-communication reordering of  
352  * the standard PCG method (Demmel, Heath and Vorst, 1993)  
353  *  
354  * @param pInput Input residual of all DOFs.  
355  * @param pOutput Solution vector of all DOFs.  
356  */
358  const int nGlobal,
359  const Array<OneD,const NekDouble> &pInput,
360  Array<OneD, NekDouble> &pOutput,
361  const AssemblyMapSharedPtr &plocToGloMap,
362  const int nDir)
363  {
364  if (!m_precon)
365  {
366  v_UniqueMap();
367  m_precon = CreatePrecon(plocToGloMap);
368  m_precon->BuildPreconditioner();
369  }
370 
371  // Get the communicator for performing data exchanges
373  = m_expList.lock()->GetComm()->GetRowComm();
374 
375  // Get vector sizes
376  int nNonDir = nGlobal - nDir;
377 
378  // Allocate array storage
379  Array<OneD, NekDouble> w_A (nGlobal, 0.0);
380  Array<OneD, NekDouble> s_A (nGlobal, 0.0);
381  Array<OneD, NekDouble> p_A (nNonDir, 0.0);
382  Array<OneD, NekDouble> r_A (nNonDir, 0.0);
383  Array<OneD, NekDouble> q_A (nNonDir, 0.0);
385 
386  // Create NekVector wrappers for linear algebra operations
387  NekVector<NekDouble> in (nNonDir,pInput + nDir, eWrapper);
388  NekVector<NekDouble> out(nNonDir,tmp = pOutput + nDir,eWrapper);
389  NekVector<NekDouble> w (nNonDir,tmp = w_A + nDir, eWrapper);
390  NekVector<NekDouble> s (nNonDir,tmp = s_A + nDir, eWrapper);
391  NekVector<NekDouble> p (nNonDir,p_A, eWrapper);
392  NekVector<NekDouble> r (nNonDir,r_A, eWrapper);
393  NekVector<NekDouble> q (nNonDir,q_A, eWrapper);
394 
395  int k;
396  NekDouble alpha, beta, rho, rho_new, mu, eps, min_resid;
397  Array<OneD, NekDouble> vExchange(3,0.0);
398 
399  // Copy initial residual from input
400  r = in;
401  // zero homogeneous out array ready for solution updates
402  // Should not be earlier in case input vector is same as
403  // output and above copy has been peformed
404  Vmath::Zero(nNonDir,tmp = pOutput + nDir,1);
405 
406 
407  // evaluate initial residual error for exit check
408  vExchange[2] = Vmath::Dot2(nNonDir,
409  r_A,
410  r_A,
411  m_map + nDir);
412 
413  vComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
414 
415  eps = vExchange[2];
416 
418  {
419  NekVector<NekDouble> inGlob (nGlobal, pInput, eWrapper);
420  Set_Rhs_Magnitude(inGlob);
421  }
422 
423  // If input residual is less than tolerance skip solve.
425  {
426  if (m_verbose && m_root)
427  {
428  cout << "CG iterations made = " << m_totalIterations
429  << " using tolerance of " << m_tolerance
430  << " (error = " << sqrt(eps/m_rhs_magnitude)
431  << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")"
432  << endl;
433  }
434  return;
435  }
436 
437  m_totalIterations = 1;
438  m_precon->DoPreconditioner(r_A, tmp = w_A + nDir);
439 
440  v_DoMatrixMultiply(w_A, s_A);
441 
442  k = 0;
443 
444  vExchange[0] = Vmath::Dot2(nNonDir,
445  r_A,
446  w_A + nDir,
447  m_map + nDir);
448 
449  vExchange[1] = Vmath::Dot2(nNonDir,
450  s_A + nDir,
451  w_A + nDir,
452  m_map + nDir);
453 
454  vComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
455 
456  rho = vExchange[0];
457  mu = vExchange[1];
458  min_resid = m_rhs_magnitude;
459  beta = 0.0;
460  alpha = rho/mu;
461 
462  // Continue until convergence
463  while (true)
464  {
465  if(k >= m_maxiter)
466  {
467  if (m_root)
468  {
469  cout << "CG iterations made = " << m_totalIterations
470  << " using tolerance of " << m_tolerance
471  << " (error = " << sqrt(eps/m_rhs_magnitude)
472  << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")"
473  << endl;
474  ASSERTL0(false,
475  "Exceeded maximum number of iterations");
476  }
477  }
478 
479  // Compute new search direction p_k, q_k
480  Vmath::Svtvp(nNonDir, beta, &p_A[0], 1, &w_A[nDir], 1, &p_A[0], 1);
481  Vmath::Svtvp(nNonDir, beta, &q_A[0], 1, &s_A[nDir], 1, &q_A[0], 1);
482 
483  // Update solution x_{k+1}
484  Vmath::Svtvp(nNonDir, alpha, &p_A[0], 1, &pOutput[nDir], 1, &pOutput[nDir], 1);
485 
486  // Update residual vector r_{k+1}
487  Vmath::Svtvp(nNonDir, -alpha, &q_A[0], 1, &r_A[0], 1, &r_A[0], 1);
488 
489  // Apply preconditioner
490  m_precon->DoPreconditioner(r_A, tmp = w_A + nDir);
491 
492  // Perform the method-specific matrix-vector multiply operation.
493  v_DoMatrixMultiply(w_A, s_A);
494 
495  // <r_{k+1}, w_{k+1}>
496  vExchange[0] = Vmath::Dot2(nNonDir,
497  r_A,
498  w_A + nDir,
499  m_map + nDir);
500  // <s_{k+1}, w_{k+1}>
501  vExchange[1] = Vmath::Dot2(nNonDir,
502  s_A + nDir,
503  w_A + nDir,
504  m_map + nDir);
505 
506  // <r_{k+1}, r_{k+1}>
507  vExchange[2] = Vmath::Dot2(nNonDir,
508  r_A,
509  r_A,
510  m_map + nDir);
511 
512  // Perform inner-product exchanges
513  vComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
514 
515  rho_new = vExchange[0];
516  mu = vExchange[1];
517  eps = vExchange[2];
518 
520  // test if norm is within tolerance
521  if (eps < m_tolerance * m_tolerance * m_rhs_magnitude)
522  {
523  if (m_verbose && m_root)
524  {
525  cout << "CG iterations made = " << m_totalIterations
526  << " using tolerance of " << m_tolerance
527  << " (error = " << sqrt(eps/m_rhs_magnitude)
528  << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")"
529  << endl;
530  }
531  break;
532  }
533  min_resid = min(min_resid, eps);
534 
535  // Compute search direction and solution coefficients
536  beta = rho_new/rho;
537  alpha = rho_new/(mu - rho_new*beta/alpha);
538  rho = rho_new;
539  k++;
540  }
541  }
542 
544  const NekVector<NekDouble> &pIn)
545  {
546  Array<OneD, NekDouble> vExchange(1);
547  vExchange[0] = Vmath::Dot(pIn.GetDimension(),&pIn[0],&pIn[0]);
548 
549  m_expList.lock()->GetComm()->GetRowComm()->AllReduce(
551 
552  // To ensure that very different rhs values are not being
553  // used in subsequent solvers such as the velocit solve in
554  // INC NS. If this works we then need to work out a better
555  // way to control this.
556  NekDouble new_rhs_mag = (vExchange[0] > 1e-6)? vExchange[0] : 1.0;
557 
559  {
560  m_rhs_magnitude = new_rhs_mag;
561  }
562  else
563  {
565  (1.0-m_rhs_mag_sm)*new_rhs_mag);
566  }
567  }
568 
569  }
570 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
NekDouble m_rhs_mag_sm
cnt to how many times rhs_magnitude is called
void Set_Rhs_Magnitude(const NekVector< NekDouble > &pIn)
void UpdateKnownSolutions(const int pGlobalBndDofs, const Array< OneD, const NekDouble > &pSolution, const int pNumDirBndDofs)
bool m_root
Provide verbose output and root if parallel.
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:471
bool m_useProjection
Whether to apply projection technique.
boost::circular_buffer< Array< OneD, NekDouble > > m_prevLinSol
Storage for solutions to previous linear problems.
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)=0
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
NekDouble CalculateAnorm(const int nGlobal, const Array< OneD, const NekDouble > &in, const int nDir)
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
static const NekDouble kNekZeroTol
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
static PreconditionerSharedPtr NullPreconditionerSharedPtr
NekDouble m_tolerance
Tolerance of iterative solver.
int m_numPrevSols
Total counter of previous solutions.
Array< OneD, int > m_map
Global to universal unique map.
double NekDouble
GlobalLinSysIterative(const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
static const NekDouble kNekUnsetDouble
Describe a linear system.
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:900
unsigned int GetDimension() const
Returns the number of dimensions for the point.
Definition: NekVector.cpp:212
void DoAconjugateProjection(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
A-conjugate projection technique.
NekDouble m_rhs_magnitude
dot product of rhs to normalise stopping criterion
A global linear system.
Definition: GlobalLinSys.h:74
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:940
virtual void v_SolveLinearSystem(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
Solve the matrix system.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void DoConjugateGradient(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
Actual iterative solve.
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map...
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129