Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GlobalLinSysPETSc.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File GlobalLinSys.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: GlobalLinSys definition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
38 
39 #include "petscis.h"
40 #include "petscversion.h"
41 
42 namespace Nektar
43 {
44  namespace MultiRegions
45  {
46  std::string GlobalLinSysPETSc::matMult =
48  "PETScMatMult", "Sparse");
49  std::string GlobalLinSysPETSc::matMultIds[] = {
51  "PETScMatMult", "Sparse",
54  "PETScMatMult", "Shell",
56  };
57 
58  /**
59  * @class GlobalLinSysPETSc
60  *
61  * Solves a linear system using PETSc.
62  */
64  const GlobalLinSysKey &pKey,
65  const boost::weak_ptr<ExpList> &pExp,
66  const boost::shared_ptr<AssemblyMap> &pLocToGloMap)
67  : GlobalLinSys(pKey, pExp, pLocToGloMap)
68  {
69  // Determine whether to use standard sparse matrix approach or
70  // shell.
71  m_matMult = pExp.lock()->GetSession()
72  ->GetSolverInfoAsEnum<PETScMatMult>(
73  "PETScMatMult");
74 
75  // Check PETSc is initialized. For some reason, this is needed on
76  // OS X as logging is not initialized properly in the call within
77  // CommMpi.
78  PetscBool isInitialized;
79  PetscInitialized(&isInitialized);
80  if (!isInitialized)
81  {
82  PetscInitializeNoArguments();
83  }
84 
85  // Create matrix
86  MatCreate(PETSC_COMM_WORLD, &m_matrix);
87  }
88 
89  /**
90  * @brief Clean up PETSc objects.
91  *
92  * Note that if SessionReader::Finalize is called before the end of the
93  * program, PETSc may have been finalized already, at which point we
94  * cannot deallocate our objects. If that's the case we do nothing and
95  * let the kernel clear up after us.
96  */
98  {
99  PetscBool isFinalized;
100  PetscFinalized(&isFinalized);
101 
102  // Sometimes, PetscFinalized returns false when (in fact) CommMpi's
103  // Finalise routine has been called. We therefore also need to check
104  // whether MPI has been finalised. This might arise from the
105  // additional call to PetscInitializeNoArguments in the constructor
106  // above.
107 #ifdef NEKTAR_USE_MPI
108  int mpiFinal = 0;
109  MPI_Finalized(&mpiFinal);
110  isFinalized = isFinalized || mpiFinal ? PETSC_TRUE : PETSC_FALSE;
111 #endif
112 
113  if (!isFinalized)
114  {
115  KSPDestroy(&m_ksp);
116  PCDestroy (&m_pc);
117  MatDestroy(&m_matrix);
118  VecDestroy(&m_x);
119  VecDestroy(&m_b);
120  VecDestroy(&m_locVec);
121  }
122  }
123 
124  /**
125  * @brief Solve linear system using PETSc.
126  *
127  * The general strategy being a PETSc solve is to:
128  *
129  * - Copy values into the PETSc vector #m_b
130  * - Solve the system #m_ksp and place result into #m_x.
131  * - Scatter results back into #m_locVec using #m_ctx scatter object.
132  * - Copy from #m_locVec to output array #pOutput.
133  */
135  const int pNumRows,
136  const Array<OneD,const NekDouble> &pInput,
137  Array<OneD, NekDouble> &pOutput,
138  const AssemblyMapSharedPtr &locToGloMap,
139  const int pNumDir)
140  {
141  const int nHomDofs = pNumRows - pNumDir;
142 
144  {
145  m_precon = CreatePrecon(locToGloMap);
146  m_precon->BuildPreconditioner();
147  }
148 
149  // Populate RHS vector from input
150  VecSetValues(m_b, nHomDofs, &m_reorderedMap[0],
151  &pInput[pNumDir], INSERT_VALUES);
152 
153  // Assemble RHS vector
154  VecAssemblyBegin(m_b);
155  VecAssemblyEnd (m_b);
156 
157  // Do system solve
158  KSPSolve(m_ksp, m_b, m_x);
159 
160  // Scatter results to local vector
161  VecScatterBegin(m_ctx, m_x, m_locVec,
162  INSERT_VALUES, SCATTER_FORWARD);
163  VecScatterEnd (m_ctx, m_x, m_locVec,
164  INSERT_VALUES, SCATTER_FORWARD);
165 
166  // Copy results into output vector
167  PetscScalar *tmp;
168  VecGetArray (m_locVec, &tmp);
169  Vmath::Vcopy (nHomDofs, tmp, 1, &pOutput[pNumDir], 1);
170  VecRestoreArray(m_locVec, &tmp);
171  }
172 
173  /**
174  * @brief Set up PETSc local (equivalent to Nektar++ global) and global
175  * (equivalent to universal) scatter maps.
176  *
177  * These maps are used in GlobalLinSysPETSc::v_SolveLinearSystem to
178  * scatter the solution vector back to each process.
179  */
181  {
182  const int nHomDofs = m_reorderedMap.size();
183 
184  // Create local and global numbering systems for vector
185  IS isGlobal, isLocal;
186  ISCreateGeneral (PETSC_COMM_SELF, nHomDofs, &m_reorderedMap[0],
187  PETSC_COPY_VALUES, &isGlobal);
188  ISCreateStride (PETSC_COMM_SELF, nHomDofs, 0, 1, &isLocal);
189 
190  // Create local vector for output
191  VecCreate (PETSC_COMM_SELF, &m_locVec);
192  VecSetSizes (m_locVec, nHomDofs, PETSC_DECIDE);
193  VecSetFromOptions(m_locVec);
194 
195  // Create scatter context
196  VecScatterCreate (m_x, isGlobal, m_locVec, isLocal, &m_ctx);
197 
198  // Clean up
199  ISDestroy(&isGlobal);
200  ISDestroy(&isLocal);
201  }
202 
203  /**
204  * @brief Calculate a reordering of universal IDs for PETSc.
205  *
206  * PETSc requires a unique, contiguous index of all global and universal
207  * degrees of freedom which represents its position inside the
208  * matrix. Presently Gs does not guarantee this, so this routine
209  * constructs a new universal mapping.
210  *
211  * @param glo2uniMap Global to universal map
212  * @param glo2unique Global to unique map
213  * @param pLocToGloMap Assembly map for this system
214  */
216  const Array<OneD, const int> &glo2uniMap,
217  const Array<OneD, const int> &glo2unique,
218  const AssemblyMapSharedPtr &pLocToGloMap)
219  {
221  = m_expList.lock()->GetSession()->GetComm();
222 
223  const int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
224  const int nHomDofs = glo2uniMap.num_elements() - nDirDofs;
225  const int nProc = vComm->GetSize();
226  const int rank = vComm->GetRank();
227 
228  int n, cnt;
229 
230  // Count number of unique degrees of freedom on each process.
231  m_nLocal = Vmath::Vsum(nHomDofs, glo2unique + nDirDofs, 1);
232  m_reorderedMap.resize(nHomDofs);
233 
234  // Reduce coefficient counts across all processors.
235  Array<OneD, int> localCounts(nProc, 0), localOffset(nProc, 0);
236  localCounts[rank] = nHomDofs;
237  vComm->AllReduce(localCounts, LibUtilities::ReduceSum);
238 
239  for (n = 1; n < nProc; ++n)
240  {
241  localOffset[n] = localOffset[n-1] + localCounts[n-1];
242  }
243 
244  int totHomDofs = Vmath::Vsum(nProc, localCounts, 1);
245  vector<unsigned int> allUniIds(totHomDofs, 0);
246 
247  // Assemble list of universal IDs
248  for (n = 0; n < nHomDofs; ++n)
249  {
250  int gid = n + nDirDofs;
251  allUniIds[n + localOffset[rank]] = glo2uniMap[gid];
252  }
253 
254  // Reduce this across processors so that each process has a list of
255  // all universal IDs.
256  vComm->AllReduce(allUniIds, LibUtilities::ReduceSum);
257  std::sort(allUniIds.begin(), allUniIds.end());
258  map<int,int> uniIdReorder;
259 
260  // Renumber starting from 0.
261  for (cnt = n = 0; n < allUniIds.size(); ++n)
262  {
263  if (uniIdReorder.count(allUniIds[n]) > 0)
264  {
265  continue;
266  }
267 
268  uniIdReorder[allUniIds[n]] = cnt++;
269  }
270 
271  // Populate reordering map.
272  for (n = 0; n < nHomDofs; ++n)
273  {
274  int gid = n + nDirDofs;
275  int uniId = glo2uniMap[gid];
276  ASSERTL0(uniIdReorder.count(uniId) > 0, "Error in ordering");
277  m_reorderedMap[n] = uniIdReorder[uniId];
278  }
279  }
280 
281  /**
282  * @brief Construct PETSc matrix and vector handles.
283  *
284  * @todo Preallocation should be done at this point, since presently
285  * matrix allocation takes a significant amount of time.
286  *
287  * @param nGlobal Number of global degrees of freedom in the system (on
288  * this processor)
289  * @param nDir Number of Dirichlet degrees of freedom (on this
290  * processor).
291  */
292  void GlobalLinSysPETSc::SetUpMatVec(int nGlobal, int nDir)
293  {
294  // CREATE VECTORS
295  VecCreate (PETSC_COMM_WORLD, &m_x);
296  VecSetSizes (m_x, m_nLocal, PETSC_DECIDE);
297  VecSetFromOptions(m_x);
298  VecDuplicate (m_x, &m_b);
299 
300  // CREATE MATRICES
302  {
303  // Create ShellCtx context object which will store the matrix
304  // size and a pointer to the linear system. We do this so that
305  // we can call a member function to the matrix-vector and
306  // preconditioning multiplication in a subclass.
307  ShellCtx *ctx1 = new ShellCtx(), *ctx2 = new ShellCtx();
308  ctx1->nGlobal = ctx2->nGlobal = nGlobal;
309  ctx1->nDir = ctx2->nDir = nDir;
310  ctx1->linSys = ctx2->linSys = this;
311 
312  // Set up MatShell object.
313  MatCreateShell (PETSC_COMM_WORLD, m_nLocal, m_nLocal,
314  PETSC_DETERMINE, PETSC_DETERMINE,
315  (void *)ctx1, &m_matrix);
316  MatShellSetOperation(m_matrix, MATOP_MULT,
317  (void(*)(void))DoMatrixMultiply);
318  MatShellSetOperation(m_matrix, MATOP_DESTROY,
319  (void(*)(void))DoDestroyMatCtx);
320 
321  // Create a PCShell to go alongside the MatShell.
322  PCCreate (PETSC_COMM_WORLD, &m_pc);
323 #if PETSC_VERSION_GE(3,5,0)
324  PCSetOperators (m_pc, m_matrix, m_matrix);
325 #else
326  PCSetOperators (m_pc, m_matrix, m_matrix, SAME_NONZERO_PATTERN);
327 #endif
328  PCSetType (m_pc, PCSHELL);
329  PCShellSetApply (m_pc, DoPreconditioner);
330  PCShellSetDestroy(m_pc, DoDestroyPCCtx);
331  PCShellSetContext(m_pc, ctx2);
332  }
333  else
334  {
335  // Otherwise we create a PETSc matrix and use MatSetFromOptions
336  // so that we can set various options on the command line.
337  MatCreate (PETSC_COMM_WORLD, &m_matrix);
338  MatSetType (m_matrix, MATAIJ);
339  MatSetSizes (m_matrix, m_nLocal, m_nLocal,
340  PETSC_DETERMINE, PETSC_DETERMINE);
341  MatSetFromOptions(m_matrix);
342  MatSetUp (m_matrix);
343  }
344  }
345 
346  /**
347  * @brief Perform either matrix multiplication or preconditioning using
348  * Nektar++ routines.
349  *
350  * This static function uses Nektar++ routines to calculate the
351  * matrix-vector product of @p M with @p in, storing the output in @p
352  * out.
353  *
354  * @todo There's a lot of scatters and copies that might possibly be
355  * eliminated to make this more efficient.
356  *
357  * @param in Input vector.
358  * @param out Output vector.
359  * @param ctx ShellCtx object that points to our instance of
360  * GlobalLinSysPETSc.
361  * @param precon If true, we apply a preconditioner, if false, we
362  * perform a matrix multiplication.
363  */
365  Vec &in, Vec &out, ShellCtx *ctx, bool precon)
366  {
367  const int nGlobal = ctx->nGlobal;
368  const int nDir = ctx->nDir;
369  const int nHomDofs = nGlobal - nDir;
370  GlobalLinSysPETSc *linSys = ctx->linSys;
371 
372  // Scatter from PETSc ordering to our local ordering. It's actually
373  // unclear whether this step might also do some communication in
374  // parallel, which is probably not ideal.
375  VecScatterBegin(linSys->m_ctx, in, linSys->m_locVec,
376  INSERT_VALUES, SCATTER_FORWARD);
377  VecScatterEnd (linSys->m_ctx, in, linSys->m_locVec,
378  INSERT_VALUES, SCATTER_FORWARD);
379 
380  // Temporary storage to pass to Nektar++
381  Array<OneD, NekDouble> tmpIn(nHomDofs), tmpOut(nHomDofs);
382 
383  // Get values from input vector and copy to tmpIn.
384  PetscScalar *tmpLocIn;
385  VecGetArray (linSys->m_locVec, &tmpLocIn);
386  Vmath::Vcopy (nHomDofs, tmpLocIn, 1, &tmpIn[0], 1);
387  VecRestoreArray(linSys->m_locVec, &tmpLocIn);
388 
389  // Do matrix multiply in Nektar++, store in tmpOut.
390  if (precon)
391  {
392  linSys->m_precon->DoPreconditioner(tmpIn, tmpOut);
393  }
394  else
395  {
396  linSys->v_DoMatrixMultiply(tmpIn, tmpOut);
397  }
398 
399  // Scatter back to PETSc ordering and put in out.
400  VecSetValues(out, nHomDofs, &linSys->m_reorderedMap[0],
401  &tmpOut[0], INSERT_VALUES);
402  VecAssemblyBegin(out);
403  VecAssemblyEnd (out);
404  }
405 
406  /**
407  * @brief Perform matrix multiplication using Nektar++ routines.
408  *
409  * This static function uses Nektar++ routines to calculate the
410  * matrix-vector product of @p M with @p in, storing the output in @p
411  * out.
412  *
413  * @param M Original MatShell matrix, which stores the ShellCtx
414  * object.
415  * @param in Input vector.
416  * @param out Output vector.
417  */
419  Mat M, Vec in, Vec out)
420  {
421  // Grab our shell context from M.
422  void *ptr;
423  MatShellGetContext(M, &ptr);
424  ShellCtx *ctx = (ShellCtx *)ptr;
425 
426  DoNekppOperation(in, out, ctx, false);
427 
428  // Must return 0, otherwise PETSc complains.
429  return 0;
430  }
431 
432  /**
433  * @brief Apply preconditioning using Nektar++ routines.
434  *
435  * This static function uses Nektar++ routines to apply the
436  * preconditioner stored in GlobalLinSysPETSc::m_precon from the context
437  * of @p pc to the vector @p in, storing the output in @p out.
438  *
439  * @param pc Preconditioner object that stores the ShellCtx.
440  * @param in Input vector.
441  * @param out Output vector.
442  */
444  PC pc, Vec in, Vec out)
445  {
446  // Grab our PCShell context from pc.
447  void *ptr;
448  PCShellGetContext(pc, &ptr);
449  ShellCtx *ctx = (ShellCtx *)ptr;
450 
451  DoNekppOperation(in, out, ctx, true);
452 
453  // Must return 0, otherwise PETSc complains.
454  return 0;
455  }
456 
457  /**
458  * @brief Destroy matrix shell context object.
459  *
460  * Note the matrix shell and preconditioner share a common context so
461  * this might have already been deallocated below, in which case we do
462  * nothing.
463  *
464  * @param M Matrix shell object
465  */
467  {
468  void *ptr;
469  MatShellGetContext(M, &ptr);
470  ShellCtx *ctx = (ShellCtx *)ptr;
471  delete ctx;
472  return 0;
473  }
474 
475  /**
476  * @brief Destroy preconditioner context object.
477  *
478  * Note the matrix shell and preconditioner share a common context so
479  * this might have already been deallocated above, in which case we do
480  * nothing.
481  *
482  * @param pc Preconditioner object
483  */
484  PetscErrorCode GlobalLinSysPETSc::DoDestroyPCCtx(PC pc)
485  {
486  void *ptr;
487  PCShellGetContext(pc, &ptr);
488  ShellCtx *ctx = (ShellCtx *)ptr;
489  delete ctx;
490  return 0;
491  }
492 
493  /**
494  * @brief Set up KSP solver object.
495  *
496  * This is reasonably generic setup -- most solver types can be changed
497  * using the .petscrc file.
498  *
499  * @param tolerance Residual tolerance to converge to.
500  */
502  {
503  KSPCreate(PETSC_COMM_WORLD, &m_ksp);
504  KSPSetTolerances(
505  m_ksp, tolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
506  KSPSetFromOptions(m_ksp);
507 #if PETSC_VERSION_GE(3,5,0)
508  KSPSetOperators(m_ksp, m_matrix, m_matrix);
509 #else
510  KSPSetOperators(m_ksp, m_matrix, m_matrix, SAME_NONZERO_PATTERN);
511 #endif
512 
514  {
515  KSPSetPC(m_ksp, m_pc);
516  }
517  }
518  }
519 }
VecScatter m_ctx
PETSc scatter context that takes us between Nektar++ global ordering and PETSc vector ordering...
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
GlobalLinSysPETSc * linSys
Pointer to the original calling object.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
void CalculateReordering(const Array< OneD, const int > &glo2uniMap, const Array< OneD, const int > &glo2unique, const AssemblyMapSharedPtr &pLocToGloMap)
Calculate a reordering of universal IDs for PETSc.
PC m_pc
PCShell for preconditioner.
void SetUpSolver(NekDouble tolerance)
Set up KSP solver object.
Internal struct for MatShell and PCShell calls to store current context for callback.
vector< int > m_reorderedMap
Reordering that takes universal IDs to a unique row in the PETSc matrix.
Vec m_x
PETSc vector objects used for local storage.
static PetscErrorCode DoPreconditioner(PC pc, Vec in, Vec out)
Apply preconditioning using Nektar++ routines.
int nGlobal
Number of global degrees of freedom.
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
GlobalLinSysPETSc(const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExp, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
int m_nLocal
Number of unique degrees of freedom on this process.
KSP m_ksp
KSP object that represents solver system.
double NekDouble
virtual void v_DoMatrixMultiply(const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)=0
virtual ~GlobalLinSysPETSc()
Clean up PETSc objects.
static PetscErrorCode DoMatrixMultiply(Mat M, Vec in, Vec out)
Perform matrix multiplication using Nektar++ routines.
Describe a linear system.
A PETSc global linear system.
void SetUpScatter()
Set up PETSc local (equivalent to Nektar++ global) and global (equivalent to universal) scatter maps...
A global linear system.
Definition: GlobalLinSys.h:74
int nDir
Number of Dirichlet degrees of freedom.
virtual void v_SolveLinearSystem(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
Solve linear system using PETSc.
static PetscErrorCode DoDestroyMatCtx(Mat M)
Destroy matrix shell context object.
void SetUpMatVec(int nGlobal, int nDir)
Construct PETSc matrix and vector handles.
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
static void DoNekppOperation(Vec &in, Vec &out, ShellCtx *ctx, bool precon)
Perform either matrix multiplication or preconditioning using Nektar++ routines.
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:723
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
static PetscErrorCode DoDestroyPCCtx(PC pc)
Destroy preconditioner context object.
PETScMatMult m_matMult
Enumerator to select matrix multiplication type.
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