35 #include <boost/core/ignore_unused.hpp>
47 namespace MultiRegions
59 string GlobalLinSysPETScStaticCond::className =
61 "PETScStaticCond", GlobalLinSysPETScStaticCond::create,
62 "PETSc static condensation.");
64 string GlobalLinSysPETScStaticCond::className2 =
66 "PETScMultiLevelStaticCond", GlobalLinSysPETScStaticCond::create,
67 "PETSc multi-level static condensation.");
88 GlobalLinSysPETScStaticCond::GlobalLinSysPETScStaticCond(
90 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
97 "This constructor is only valid when using static "
100 pLocToGloMap->GetGlobalSysSolnType(),
101 "The local to global map is not set up for the requested "
113 const std::shared_ptr<AssemblyMap> &pLocToGloMap,
143 int n, n_exp =
m_expList.lock()->GetNumElmts();
150 for (n = 0; n < n_exp; ++n)
156 m_precon->TransformedSchurCompl(n, cnt, mat);
158 cnt += mat->GetRows();
174 int i, j, n, cnt, gid1, gid2, loc_lda;
177 const int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
195 pLocToGloMap->GetGlobalToUniversalBndMapUnique(),
199 SetUpMatVec(pLocToGloMap->GetNumGlobalBndCoeffs(), nDirDofs);
205 SetUpSolver(pLocToGloMap->GetIterativeTolerance());
215 for (n = cnt = 0; n <
m_schurCompl->GetNumberOfBlockRows(); ++n)
218 loc_lda = loc_mat->GetRows();
220 for (i = 0; i < loc_lda; ++i)
222 gid1 = pLocToGloMap->GetLocalToGlobalBndMap(cnt + i) - nDirDofs;
223 sign1 = pLocToGloMap->GetLocalToGlobalBndSign(cnt + i);
227 for (j = 0; j < loc_lda; ++j)
229 gid2 = pLocToGloMap->GetLocalToGlobalBndMap(cnt + j) -
231 sign2 = pLocToGloMap->GetLocalToGlobalBndSign(cnt + j);
235 value = sign1 * sign2 * (*loc_mat)(i, j);
236 MatSetValue(
m_matrix, gid1ro, gid2ro, value,
246 MatAssemblyBegin(
m_matrix, MAT_FINAL_ASSEMBLY);
247 MatAssemblyEnd(
m_matrix, MAT_FINAL_ASSEMBLY);
255 unsigned int nbdry = localMat->GetRows();
256 unsigned int nblks = 1;
257 unsigned int esize[1] = {nbdry};
260 nblks, nblks, esize, esize);
261 schurComplBlock->SetBlock(0, 0, localMat);
263 return schurComplBlock;
269 boost::ignore_unused(F_bnd);
286 m_precon->DoTransformBasisToLowEnergy(pInOut);
292 m_precon->DoTransformCoeffsFromLowEnergy(pInOut);
298 m_precon->DoTransformCoeffsToLowEnergy(pInput, pOutput);
317 int nLocBndDofs = asmMap->GetNumLocalBndCoeffs();
318 int nDirDofs = asmMap->GetNumGlobalDirBndCoeffs();
321 asmMap->GlobalToLocalBnd(input, in.GetPtr(), nDirDofs);
322 out = (*m_schurCompl) * in;
323 asmMap->AssembleBnd(out.
GetPtr(), output, nDirDofs);
331 const std::shared_ptr<AssemblyMap> &l2gMap)
335 mkey, pExpList, pSchurCompl, pBinvD, pC, pInvD, l2gMap,
m_precon);
336 sys->Initialise(l2gMap);
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map.
void Initialise(const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Describe a linear system.
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
A PETSc global linear system.
PETScMatMult m_matMult
Enumerator to select matrix multiplication type.
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.
Mat m_matrix
PETSc matrix object.
PreconditionerSharedPtr m_precon
void SetUpSolver(NekDouble tolerance)
Set up KSP solver object.
void SetUpMatVec(int nGlobal, int nDir)
Construct PETSc matrix and vector handles.
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.
virtual void v_CoeffsFwdTransform(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
virtual void v_AssembleSchurComplement(std::shared_ptr< AssemblyMap > locToGloMap) override
Assemble the Schur complement matrix.
virtual GlobalLinSysStaticCondSharedPtr v_Recurse(const GlobalLinSysKey &mkey, const std::weak_ptr< ExpList > &pExpList, const DNekScalBlkMatSharedPtr pSchurCompl, const DNekScalBlkMatSharedPtr pBinvD, const DNekScalBlkMatSharedPtr pC, const DNekScalBlkMatSharedPtr pInvD, const std::shared_ptr< AssemblyMap > &locToGloMap) override
virtual ~GlobalLinSysPETScStaticCond()
virtual void v_CoeffsBwdTransform(Array< OneD, NekDouble > &pInOut) override
virtual void v_InitObject() override
virtual void v_DoMatrixMultiply(const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output) override
Apply matrix-vector multiplication using local approach and the assembly map.
virtual DNekScalBlkMatSharedPtr v_GetStaticCondBlock(unsigned int n) override
Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided ...
virtual void v_PreSolve(int scLevel, Array< OneD, NekDouble > &F_bBnd) override
GlobalLinSysPETScStaticCond(const GlobalLinSysKey &mkey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &locToGloMap)
Constructor for full direct matrix solve.
virtual void v_BasisFwdTransform(Array< OneD, NekDouble > &pInOut) override
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
std::weak_ptr< AssemblyMap > m_locToGloMap
Local to global map.
void SetupTopLevel(const std::shared_ptr< AssemblyMap > &locToGloMap)
Set up the storage for the Schur complement or the top level of the multi-level Schur complement.
DNekScalBlkMatSharedPtr m_BinvD
Block matrix.
DNekScalBlkMatSharedPtr m_C
Block matrix.
DNekScalBlkMatSharedPtr m_invD
Block matrix.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
Array< OneD, DataType > & GetPtr()
@ ePETScMultiLevelStaticCond
std::shared_ptr< GlobalLinSysStaticCond > GlobalLinSysStaticCondSharedPtr
GlobalLinSysFactory & GetGlobalLinSysFactory()
std::shared_ptr< Preconditioner > PreconditionerSharedPtr
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
std::shared_ptr< GlobalLinSysPETScStaticCond > GlobalLinSysPETScStaticCondSharedPtr
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr