40 #include "petscversion.h"
44 namespace MultiRegions
48 "PETScMatMult",
"Sparse");
51 "PETScMatMult",
"Sparse",
54 "PETScMatMult",
"Shell",
65 const boost::weak_ptr<ExpList> &pExp,
66 const boost::shared_ptr<AssemblyMap> &pLocToGloMap)
78 PetscBool isInitialized;
79 PetscInitialized(&isInitialized);
82 PetscInitializeNoArguments();
86 MatCreate(PETSC_COMM_WORLD, &
m_matrix);
99 PetscBool isFinalized;
100 PetscFinalized(&isFinalized);
107 #ifdef NEKTAR_USE_MPI
109 MPI_Finalized(&mpiFinal);
110 isFinalized = isFinalized || mpiFinal ? PETSC_TRUE : PETSC_FALSE;
141 const int nHomDofs = pNumRows - pNumDir;
151 &pInput[pNumDir], INSERT_VALUES);
154 VecAssemblyBegin(
m_b);
155 VecAssemblyEnd (
m_b);
162 INSERT_VALUES, SCATTER_FORWARD);
164 INSERT_VALUES, SCATTER_FORWARD);
185 IS isGlobal, isLocal;
187 PETSC_COPY_VALUES, &isGlobal);
188 ISCreateStride (PETSC_COMM_SELF, nHomDofs, 0, 1, &isLocal);
191 VecCreate (PETSC_COMM_SELF, &
m_locVec);
192 VecSetSizes (
m_locVec, nHomDofs, PETSC_DECIDE);
199 ISDestroy(&isGlobal);
221 =
m_expList.lock()->GetSession()->GetComm();
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();
236 localCounts[rank] = nHomDofs;
239 for (n = 1; n < nProc; ++n)
241 localOffset[n] = localOffset[n-1] + localCounts[n-1];
244 int totHomDofs =
Vmath::Vsum(nProc, localCounts, 1);
245 vector<unsigned int> allUniIds(totHomDofs, 0);
248 for (n = 0; n < nHomDofs; ++n)
250 int gid = n + nDirDofs;
251 allUniIds[n + localOffset[rank]] = glo2uniMap[gid];
257 std::sort(allUniIds.begin(), allUniIds.end());
258 map<int,int> uniIdReorder;
261 for (cnt = n = 0; n < allUniIds.size(); ++n)
263 if (uniIdReorder.count(allUniIds[n]) > 0)
268 uniIdReorder[allUniIds[n]] = cnt++;
272 for (n = 0; n < nHomDofs; ++n)
274 int gid = n + nDirDofs;
275 int uniId = glo2uniMap[gid];
276 ASSERTL0(uniIdReorder.count(uniId) > 0,
"Error in ordering");
295 VecCreate (PETSC_COMM_WORLD, &
m_x);
297 VecSetFromOptions(
m_x);
308 ctx1->
nGlobal = ctx2->nGlobal = nGlobal;
309 ctx1->
nDir = ctx2->nDir = nDir;
310 ctx1->
linSys = ctx2->linSys =
this;
314 PETSC_DETERMINE, PETSC_DETERMINE,
316 MatShellSetOperation(
m_matrix, MATOP_MULT,
318 MatShellSetOperation(
m_matrix, MATOP_DESTROY,
322 PCCreate (PETSC_COMM_WORLD, &
m_pc);
323 #if PETSC_VERSION_GE(3,5,0)
328 PCSetType (
m_pc, PCSHELL);
331 PCShellSetContext(
m_pc, ctx2);
337 MatCreate (PETSC_COMM_WORLD, &
m_matrix);
340 PETSC_DETERMINE, PETSC_DETERMINE);
365 Vec &in, Vec &out,
ShellCtx *ctx,
bool precon)
367 const int nGlobal = ctx->
nGlobal;
368 const int nDir = ctx->
nDir;
369 const int nHomDofs = nGlobal - nDir;
376 INSERT_VALUES, SCATTER_FORWARD);
378 INSERT_VALUES, SCATTER_FORWARD);
384 PetscScalar *tmpLocIn;
385 VecGetArray (linSys->
m_locVec, &tmpLocIn);
387 VecRestoreArray(linSys->
m_locVec, &tmpLocIn);
392 linSys->
m_precon->DoPreconditioner(tmpIn, tmpOut);
401 &tmpOut[0], INSERT_VALUES);
402 VecAssemblyBegin(out);
403 VecAssemblyEnd (out);
419 Mat M, Vec in, Vec out)
423 MatShellGetContext(M, &ptr);
444 PC pc, Vec in, Vec out)
448 PCShellGetContext(pc, &ptr);
469 MatShellGetContext(M, &ptr);
487 PCShellGetContext(pc, &ptr);
503 KSPCreate(PETSC_COMM_WORLD, &
m_ksp);
505 m_ksp, tolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
506 KSPSetFromOptions(
m_ksp);
507 #if PETSC_VERSION_GE(3,5,0)
VecScatter m_ctx
PETSc scatter context that takes us between Nektar++ global ordering and PETSc vector ordering...
#define ASSERTL0(condition, msg)
GlobalLinSysPETSc * linSys
Pointer to the original calling object.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
PreconditionerSharedPtr m_precon
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
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.
static std::string matMult
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.
GlobalLinSysPETSc(const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExp, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
Mat m_matrix
PETSc matrix object.
int m_nLocal
Number of unique degrees of freedom on this process.
KSP m_ksp
KSP object that represents solver system.
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...
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 std::string matMultIds[]
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)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.