39 #include "petscversion.h"
43 namespace MultiRegions
54 const boost::weak_ptr<ExpList> &pExp,
55 const boost::shared_ptr<AssemblyMap> &pLocToGloMap)
61 PetscInitializeNoArguments();
64 MatCreate(PETSC_COMM_WORLD, &
m_matrix);
88 const int nHomDofs = pNumRows - pNumDir;
92 &pInput[pNumDir], INSERT_VALUES);
95 VecAssemblyBegin(
m_b);
103 INSERT_VALUES, SCATTER_FORWARD);
105 INSERT_VALUES, SCATTER_FORWARD);
126 IS isGlobal, isLocal;
128 PETSC_COPY_VALUES, &isGlobal);
129 ISCreateStride (PETSC_COMM_SELF, nHomDofs, 0, 1, &isLocal);
132 VecCreate (PETSC_COMM_SELF, &
m_locVec);
133 VecSetSizes (
m_locVec, nHomDofs, PETSC_DECIDE);
140 ISDestroy(&isGlobal);
158 =
m_expList.lock()->GetSession()->GetComm();
160 const int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
161 const int nHomDofs = glo2uniMap.num_elements() - nDirDofs;
162 const int nProc = vComm->GetSize();
163 const int rank = vComm->GetRank();
173 localCounts[rank] = nHomDofs;
176 for (n = 1; n < nProc; ++n)
178 localOffset[n] = localOffset[n-1] + localCounts[n-1];
181 int totHomDofs =
Vmath::Vsum(nProc, localCounts, 1);
182 vector<unsigned int> allUniIds(totHomDofs, 0);
185 for (n = 0; n < nHomDofs; ++n)
187 int gid = n + nDirDofs;
188 allUniIds[n + localOffset[rank]] = glo2uniMap[gid];
194 std::sort(allUniIds.begin(), allUniIds.end());
195 map<int,int> uniIdReorder;
198 for (cnt = n = 0; n < allUniIds.size(); ++n)
200 if (uniIdReorder.count(allUniIds[n]) > 0)
205 uniIdReorder[allUniIds[n]] = cnt++;
209 for (n = 0; n < nHomDofs; ++n)
211 int gid = n + nDirDofs;
212 int uniId = glo2uniMap[gid];
213 ASSERTL0(uniIdReorder.count(uniId) > 0,
"Error in ordering");
227 VecCreate (PETSC_COMM_WORLD, &
m_x);
229 VecSetFromOptions(
m_x);
233 MatCreate (PETSC_COMM_WORLD, &
m_matrix);
236 PETSC_DETERMINE, PETSC_DETERMINE);
249 KSPCreate(PETSC_COMM_WORLD, &
m_ksp);
251 m_ksp, tolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
252 KSPSetFromOptions(
m_ksp);
253 #if PETSC_VERSION_GE(3,5,0)
#define ASSERTL0(condition, msg)
void SetUpMatVec()
Construct PETSc matrix and vector handles.
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.
void SetUpSolver(NekDouble tolerance)
Set up KSP solver object.
vector< int > m_reorderedMap
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.
virtual ~GlobalLinSysPETSc()
Describe a linear system.
void SetUpScatter()
Set up PETSc local (equivalent to Nektar++ global) and global (equivalent to universal) scatter maps...
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.
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)
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.