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