49 namespace MultiRegions
54 string PreconditionerBlock::className =
56 "Block", PreconditionerBlock::create,
"Block Diagonal Preconditioning");
64 PreconditionerBlock::PreconditionerBlock(
65 const std::shared_ptr<GlobalLinSys> &plinsys,
76 "Solver type not valid");
117 int i, j, k, n, cnt, gId;
118 int meshVertId, meshEdgeId, meshFaceId;
122 const int nExp = expList->GetExpSize();
123 const int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
126 PeriodicMap periodicVerts, periodicEdges, periodicFaces;
127 expList->GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
134 vector<map<int, vector<NekDouble>>> idMats(3);
138 vector<map<int, int>> gidMeshIds(3);
142 vector<map<int, int>> gidDofs(3);
154 for (cnt = n = 0; n < nExp; ++n)
156 exp = expList->GetExp(n);
160 m_linsys.lock()->GetStaticCondBlock(n)->GetBlock(0, 0);
164 for (i = 0; i < exp->GetNverts(); ++i)
166 meshVertId = exp->GetGeom()->GetVid(i);
167 int locId = exp->GetVertexMap(i);
170 gId = asmMap->GetLocalToGlobalMap(cnt + locId) - nDirBnd;
181 NekDouble vertVal = (*schurMat)(locId, locId);
185 auto gIt = idMats[0].find(gId);
187 if (gIt == idMats[0].end())
190 idMats[0][gId] = vector<NekDouble>(1, vertVal);
196 gIt->second[0] += vertVal;
203 auto pIt = periodicVerts.find(meshVertId);
204 if (pIt != periodicVerts.end())
206 for (j = 0; j < pIt->second.size(); ++j)
208 meshVertId = min(meshVertId, pIt->second[j].id);
214 gidMeshIds[0][gId] = meshVertId;
215 maxVertIds[0] = max(maxVertIds[0], meshVertId);
221 for (i = 0; i < exp->GetGeom()->GetNumEdges(); ++i)
223 meshEdgeId = exp->GetGeom()->GetEid(i);
231 auto pIt = periodicEdges.find(meshEdgeId);
232 if (pIt != periodicEdges.end())
234 pair<int, StdRegions::Orientation> idOrient =
237 meshEdgeId = idOrient.first;
238 edgeOrient = idOrient.second;
244 if (exp->GetGeom()->GetNumFaces())
247 ->GetEdgeInteriorToElementMap(i, bmap,
sign, edgeOrient);
249 ->GetEdgeInverseBoundaryMap(i);
253 exp->GetTraceInteriorToElementMap(i, bmap,
sign, edgeOrient);
255 ->GetTraceInverseBoundaryMap(i);
259 const int nEdgeCoeffs = bmap.size();
260 vector<NekDouble> tmpStore(nEdgeCoeffs * nEdgeCoeffs);
262 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
264 for (j = 0; j < nEdgeCoeffs; ++j)
272 asmMap->GetLocalToGlobalMap(cnt + bmap[j]) - nDirBnd);
284 for (k = 0; k < nEdgeCoeffs; ++k)
286 tmpStore[k + j * nEdgeCoeffs] =
287 sign1 *
sign[k] * (*schurMat)(bmap2[j], bmap2[k]);
296 gidDofs[1][gId] = nEdgeCoeffs;
299 auto gIt = idMats[1].find(gId);
300 if (gIt == idMats[1].end())
302 idMats[1][gId] = tmpStore;
306 ASSERTL1(tmpStore.size() == gIt->second.size(),
307 "Number of modes mismatch");
308 Vmath::Vadd(nEdgeCoeffs * nEdgeCoeffs, &gIt->second[0], 1,
309 &tmpStore[0], 1, &gIt->second[0], 1);
312 gidMeshIds[1][gId] = meshEdgeId;
313 maxVertIds[2] = max(maxVertIds[2], meshEdgeId);
314 maxVertIds[3] = max(maxVertIds[3], nEdgeCoeffs);
319 for (i = 0; i < exp->GetGeom()->GetNumFaces(); ++i)
321 meshFaceId = exp->GetGeom()->GetFid(i);
329 auto pIt = periodicFaces.find(meshFaceId);
330 if (pIt != periodicFaces.end())
332 meshFaceId = min(meshFaceId, pIt->second[0].id);
334 pIt->second[0].orient);
337 exp->GetTraceInteriorToElementMap(i, bmap,
sign, faceOrient);
339 ->GetTraceInverseBoundaryMap(i);
342 const int nFaceCoeffs = bmap.size();
343 vector<NekDouble> tmpStore(nFaceCoeffs * nFaceCoeffs);
345 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
347 for (j = 0; j < nFaceCoeffs; ++j)
350 asmMap->GetLocalToGlobalMap(cnt + bmap[j]) - nDirBnd);
362 for (k = 0; k < nFaceCoeffs; ++k)
364 tmpStore[k + j * nFaceCoeffs] =
365 sign1 *
sign[k] * (*schurMat)(bmap2[j], bmap2[k]);
374 gidDofs[2][gId] = nFaceCoeffs;
377 auto gIt = idMats[2].find(gId);
378 if (gIt == idMats[2].end())
380 idMats[2][gId] = tmpStore;
384 ASSERTL1(tmpStore.size() == gIt->second.size(),
385 "Number of modes mismatch");
386 Vmath::Vadd(nFaceCoeffs * nFaceCoeffs, &gIt->second[0], 1,
387 &tmpStore[0], 1, &gIt->second[0], 1);
390 gidMeshIds[2][gId] = meshFaceId;
391 maxVertIds[4] = max(maxVertIds[4], meshFaceId);
392 maxVertIds[5] = max(maxVertIds[5], nFaceCoeffs);
395 cnt += exp->GetNcoeffs();
400 m_comm = expList->GetSession()->GetComm()->GetRowComm();
405 vector<NekDouble> storageBuf;
406 vector<long> globalToUniversal;
408 for (i = 0, cnt = 1; i < 3; ++i)
410 const int maxDofs = maxVertIds[2 * i + 1];
415 for (
auto &gIt : idMats[i])
418 storageBuf.insert(storageBuf.end(), gIt.second.begin(),
422 ASSERTL1(gidMeshIds[i].count(gIt.first) > 0,
423 "Unable to find global ID " +
424 boost::lexical_cast<string>(gIt.first) +
426 meshVertId = gidMeshIds[i][gIt.first];
428 for (j = 0; j < gIt.second.size(); ++j)
430 globalToUniversal.push_back(cnt +
431 meshVertId * maxDofs * maxDofs + j);
438 cnt += (maxVertIds[2 * i] + 1) * maxDofs * maxDofs;
441 ASSERTL1(storageBuf.size() == globalToUniversal.size(),
442 "Storage buffer and global to universal map size does "
447 &globalToUniversal[0]);
452 expList->GetSession()->DefinesCmdLineArgument(
"verbose"));
459 n_blks[0] = idMats[0].size();
464 for (i = 1; i < 3; ++i)
466 for (
auto &gIt : idMats[i])
468 ASSERTL1(gidDofs[i].count(gIt.first) > 0,
469 "Unable to find number of degrees of freedom for "
471 boost::lexical_cast<string>(gIt.first));
473 n_blks[cnt++] = gidDofs[i][gIt.first];
488 for (
auto gIt = idMats[0].begin(); gIt != idMats[0].end(); ++gIt, ++cnt)
490 (*vertMat)(cnt, cnt) = 1.0 / storageData[cnt];
499 for (i = 1; i < 3; ++i)
501 for (
auto &gIt : idMats[i])
503 int nDofs = gidDofs[i][gIt.first];
508 for (j = 0; j < nDofs; ++j)
510 for (k = 0; k < nDofs; ++k)
512 (*tmp)(j, k) = storageData[k + j * nDofs + cnt];
516 cnt += nDofs * nDofs;
519 m_blkMat->SetBlock(cnt2, cnt2, tmp);
540 std::shared_ptr<MultiRegions::ExpList> expList =
541 ((
m_linsys.lock())->GetLocMat()).lock();
542 std::shared_ptr<MultiRegions::ExpList> trace = expList->GetTrace();
548 std::dynamic_pointer_cast<AssemblyMapDG>(
m_locToGloMap.lock());
550 int i, j, k, n, cnt, cnt2;
553 int nTrace = expList->GetTrace()->GetExpSize();
554 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
556 for (cnt = n = 0; n < nTrace; ++n)
563 cnt += trace->GetExp(n)->GetNcoeffs();
570 int maxTraceSize = -1;
571 map<int, int> arrayOffsets;
573 for (cnt = 0, n = nDir; n < nTrace; ++n)
575 int nCoeffs = trace->GetExp(n)->GetNcoeffs();
576 int nCoeffs2 = nCoeffs * nCoeffs;
578 arrayOffsets[n] = cnt;
581 if (maxTraceSize < nCoeffs2)
583 maxTraceSize = nCoeffs2;
588 m_comm = expList->GetSession()->GetComm()->GetRowComm();
595 for (cnt = n = 0; n < expList->GetExpSize(); ++n)
598 locExpansion = expList->GetExp(elmt);
601 asmMap->GetElmtToTrace()[elmt];
604 loc_mat = (
m_linsys.lock())->GetStaticCondBlock(n);
605 bnd_mat = loc_mat->GetBlock(0, 0);
607 int nFacets = locExpansion->GetNtraces();
609 for (cnt2 = i = 0; i < nFacets; ++i)
611 int nCoeffs = elmtToTraceMap[i]->GetNcoeffs();
612 int elmtId = elmtToTraceMap[i]->GetElmtId();
622 NekDouble *off = &tmpStore[arrayOffsets[elmtId]];
624 for (j = 0; j < nCoeffs; ++j)
626 NekDouble sign1 = asmMap->GetLocalToGlobalBndSign(cnt + j);
627 for (k = 0; k < nCoeffs; ++k)
629 NekDouble sign2 = asmMap->GetLocalToGlobalBndSign(cnt + k);
630 off[j * nCoeffs + k] +=
631 (*bnd_mat)(cnt2 + j, cnt2 + k) * sign1 * sign2;
642 for (cnt = 0, n = nDir; n < nTrace; ++n)
645 int nCoeffs = traceExp->GetNcoeffs();
646 int geomId = traceExp->GetGeom()->GetGlobalID();
648 for (i = 0; i < nCoeffs * nCoeffs; ++i)
650 uniIds[cnt++] = geomId * maxTraceSize + i + 1;
657 expList->GetSession()->DefinesCmdLineArgument(
"verbose"));
662 for (n = 0; n < nTrace - nDir; ++n)
664 n_blks[n] = trace->GetExp(n + nDir)->GetNcoeffs();
670 for (n = 0; n < nTrace - nDir; ++n)
672 int nCoeffs = trace->GetExp(n + nDir)->GetNcoeffs();
675 NekDouble *off = &tmpStore[arrayOffsets[n + nDir]];
677 for (i = 0; i < nCoeffs; ++i)
679 for (j = 0; j < nCoeffs; ++j)
681 (*tmp)(i, j) = off[i * nCoeffs + j];
699 int nNonDir = nGlobal - nDir;
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define sign(a, b)
return the sign(b)*a
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.
Describe a linear system.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
void BlockPreconditionerCG(void)
Construct a block preconditioner from for the continuous Galerkin system.
virtual void v_BuildPreconditioner()
virtual void v_InitObject()
void BlockPreconditionerHDG(void)
Construct a block preconditioner for the hybridized discontinuous Galerkin system.
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply preconditioner to pInput and store the result in pOutput.
DNekBlkMatSharedPtr m_blkMat
LibUtilities::CommSharedPtr m_comm
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
PreconFactory & GetPreconFactory()
pair< int, StdRegions::Orientation > DeterminePeriodicEdgeOrientId(int meshEdgeId, StdRegions::Orientation edgeOrient, const vector< PeriodicEntity > &periodicEdges)
Determine orientation of an edge to its periodic equivalents, as well as the ID of the representative...
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
StdRegions::Orientation DeterminePeriodicFaceOrient(StdRegions::Orientation faceOrient, StdRegions::Orientation perFaceOrient)
Determine relative orientation between two faces.
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.