39#include "petscversion.h"
61 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
72 PetscBool isInitialized;
73 PetscInitialized(&isInitialized);
77 std::string commType =
78 m_expList.lock()->GetSession()->GetComm()->GetType();
79 if (commType.find(
"MPI") != std::string::npos)
82 std::static_pointer_cast<LibUtilities::CommMpi>(
83 m_expList.lock()->GetSession()->GetComm());
84 PETSC_COMM_WORLD = comm->GetComm();
87 PetscInitializeNoArguments();
91 MatCreate(PETSC_COMM_WORLD, &
m_matrix);
104 PetscBool isFinalized;
105 PetscFinalized(&isFinalized);
114 MPI_Finalized(&mpiFinal);
115 isFinalized = isFinalized || mpiFinal ? PETSC_TRUE : PETSC_FALSE;
144 const int nHomDofs = pNumRows - pNumDir;
153 locToGloMap->Assemble(pInput, Glo);
160 VecAssemblyBegin(
m_b);
166 KSPConvergedReason reason;
167 KSPGetConvergedReason(
m_ksp, &reason);
168 ASSERTL0(reason > 0,
"PETSc solver diverged, reason is: " +
169 std::string(KSPConvergedReasons[reason]));
180 locToGloMap->GlobalToLocal(Glo, pOutput);
196 IS isGlobal, isLocal;
198 PETSC_COPY_VALUES, &isGlobal);
199 ISCreateStride(PETSC_COMM_SELF, nHomDofs, 0, 1, &isLocal);
202 VecCreate(PETSC_COMM_SELF, &
m_locVec);
203 VecSetSizes(
m_locVec, nHomDofs, PETSC_DECIDE);
210 ISDestroy(&isGlobal);
232 m_expList.lock()->GetSession()->GetComm()->GetSpaceComm();
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();
247 localCounts[rank] = nHomDofs;
250 for (n = 1; n < nProc; ++n)
252 localOffset[n] = localOffset[n - 1] + localCounts[n - 1];
255 int totHomDofs =
Vmath::Vsum(nProc, localCounts, 1);
256 vector<unsigned int> allUniIds(totHomDofs, 0);
259 for (n = 0; n < nHomDofs; ++n)
261 int gid = n + nDirDofs;
262 allUniIds[n + localOffset[rank]] = glo2uniMap[gid];
268 std::sort(allUniIds.begin(), allUniIds.end());
269 map<int, int> uniIdReorder;
272 for (cnt = n = 0; n < allUniIds.size(); ++n)
274 if (uniIdReorder.count(allUniIds[n]) > 0)
279 uniIdReorder[allUniIds[n]] = cnt++;
283 for (n = 0; n < nHomDofs; ++n)
285 int gid = n + nDirDofs;
286 int uniId = glo2uniMap[gid];
287 ASSERTL0(uniIdReorder.count(uniId) > 0,
"Error in ordering");
306 VecCreate(PETSC_COMM_WORLD, &
m_x);
308 VecSetFromOptions(
m_x);
319 ctx1->
nGlobal = ctx2->nGlobal = nGlobal;
320 ctx1->
nDir = ctx2->nDir = nDir;
321 ctx1->
linSys = ctx2->linSys =
this;
325 PETSC_DETERMINE, (
void *)ctx1, &
m_matrix);
326 MatShellSetOperation(
m_matrix, MATOP_MULT,
328 MatShellSetOperation(
m_matrix, MATOP_DESTROY,
332 PCCreate(PETSC_COMM_WORLD, &
m_pc);
333#if PETSC_VERSION_GE(3, 5, 0)
338 PCSetType(
m_pc, PCSHELL);
341 PCShellSetContext(
m_pc, ctx2);
347 MatCreate(PETSC_COMM_WORLD, &
m_matrix);
377 const int nGlobal = ctx->
nGlobal;
378 const int nDir = ctx->
nDir;
379 const int nHomDofs = nGlobal - nDir;
385 VecScatterBegin(linSys->
m_ctx, in, linSys->
m_locVec, INSERT_VALUES,
387 VecScatterEnd(linSys->
m_ctx, in, linSys->
m_locVec, INSERT_VALUES,
394 PetscScalar *tmpLocIn;
395 VecGetArray(linSys->
m_locVec, &tmpLocIn);
397 VecRestoreArray(linSys->
m_locVec, &tmpLocIn);
402 linSys->
m_precon->DoPreconditioner(tmpIn, tmpOut);
410 VecSetValues(out, nHomDofs, &linSys->
m_reorderedMap[0], &tmpOut[0],
412 VecAssemblyBegin(out);
432 MatShellGetContext(M, &ptr);
456 PCShellGetContext(pc, &ptr);
477 MatShellGetContext(M, &ptr);
495 PCShellGetContext(pc, &ptr);
511 KSPCreate(PETSC_COMM_WORLD, &
m_ksp);
512 KSPSetTolerances(
m_ksp, tolerance, PETSC_DEFAULT, PETSC_DEFAULT,
514 KSPSetFromOptions(
m_ksp);
515#if PETSC_VERSION_GE(3, 5, 0)
#define ASSERTL0(condition, msg)
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.
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map.
Describe a linear system.
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
Mat m_matrix
PETSc matrix object.
PreconditionerSharedPtr m_precon
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 std::string matMult
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.
static std::string matMultIds[]
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.
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
void Zero(int n, T *x, const int incx)
Zero vector.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.