49 namespace MultiRegions
54 string PreconditionerBlock::className =
56 "Block", PreconditionerBlock::create,
"Block Diagonal Preconditioning");
64 PreconditionerBlock::PreconditionerBlock(
65 const std::shared_ptr<GlobalLinSys> &plinsys,
79 "Solver type not valid");
133 int i, j, k, n, cnt, gId;
134 int meshVertId, meshEdgeId, meshFaceId;
138 const int nExp = expList->GetExpSize();
139 const int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
142 PeriodicMap periodicVerts, periodicEdges, periodicFaces;
143 expList->GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
153 vector<map<int, vector<NekDouble>>> idMats(nStorage);
157 vector<map<int, int>> gidMeshIds(3);
161 vector<map<int, int>> gidDofs(nStorage);
173 for (cnt = n = 0; n < nExp; ++n)
175 exp = expList->GetExp(n);
185 auto tmpMat =
m_linsys.lock()->GetBlock(n);
188 int nBndDofs = exp->NumBndryCoeffs();
189 int nIntDofs = exp->GetNcoeffs() - nBndDofs;
193 exp->GetBoundaryMap(bndMap);
194 exp->GetInteriorMap(intMap);
202 for (
int i = 0; i < nBndDofs; ++i)
204 for (
int j = 0; j < nBndDofs; ++j)
206 (*bndryMat)(i, j) = (*tmpMat)(bndMap[i], bndMap[j]);
214 vector<NekDouble> tmpStore(nIntDofs * nIntDofs);
215 for (
int i = 0; i < nIntDofs; ++i)
217 for (
int j = 0; j < nIntDofs; ++j)
219 tmpStore[j + i * nIntDofs] =
220 (*tmpMat)(intMap[i], intMap[j]);
225 idMats[3][n] = tmpStore;
226 gidDofs[3][n] = nIntDofs;
234 schurMat =
m_linsys.lock()->GetStaticCondBlock(n)->GetBlock(0, 0);
239 for (i = 0; i < exp->GetNverts(); ++i)
241 meshVertId = exp->GetGeom()->GetVid(i);
242 int locId = exp->GetVertexMap(i);
245 gId = asmMap->GetLocalToGlobalMap(cnt + locId) - nDirBnd;
256 NekDouble vertVal = (*schurMat)(locId, locId);
260 auto gIt = idMats[0].find(gId);
262 if (gIt == idMats[0].end())
265 idMats[0][gId] = vector<NekDouble>(1, vertVal);
271 gIt->second[0] += vertVal;
278 auto pIt = periodicVerts.find(meshVertId);
279 if (pIt != periodicVerts.end())
281 for (j = 0; j < pIt->second.size(); ++j)
283 meshVertId = min(meshVertId, pIt->second[j].id);
289 gidMeshIds[0][gId] = meshVertId;
290 maxVertIds[0] = max(maxVertIds[0], meshVertId);
296 for (i = 0; i < exp->GetGeom()->GetNumEdges(); ++i)
298 meshEdgeId = exp->GetGeom()->GetEid(i);
306 auto pIt = periodicEdges.find(meshEdgeId);
307 if (pIt != periodicEdges.end())
309 pair<int, StdRegions::Orientation> idOrient =
312 meshEdgeId = idOrient.first;
313 edgeOrient = idOrient.second;
319 if (exp->GetGeom()->GetNumFaces())
322 ->GetEdgeInteriorToElementMap(i, bmap,
sign, edgeOrient);
324 ->GetEdgeInverseBoundaryMap(i);
328 exp->GetTraceInteriorToElementMap(i, bmap,
sign, edgeOrient);
330 ->GetTraceInverseBoundaryMap(i);
334 const int nEdgeCoeffs = bmap.size();
335 vector<NekDouble> tmpStore(nEdgeCoeffs * nEdgeCoeffs);
337 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
339 for (j = 0; j < nEdgeCoeffs; ++j)
347 asmMap->GetLocalToGlobalMap(cnt + bmap[j]) - nDirBnd);
359 for (k = 0; k < nEdgeCoeffs; ++k)
361 tmpStore[k + j * nEdgeCoeffs] =
362 sign1 *
sign[k] * (*schurMat)(bmap2[j], bmap2[k]);
371 gidDofs[1][gId] = nEdgeCoeffs;
374 auto gIt = idMats[1].find(gId);
375 if (gIt == idMats[1].end())
377 idMats[1][gId] = tmpStore;
381 ASSERTL1(tmpStore.size() == gIt->second.size(),
382 "Number of modes mismatch");
383 Vmath::Vadd(nEdgeCoeffs * nEdgeCoeffs, &gIt->second[0], 1,
384 &tmpStore[0], 1, &gIt->second[0], 1);
387 gidMeshIds[1][gId] = meshEdgeId;
388 maxVertIds[2] = max(maxVertIds[2], meshEdgeId);
389 maxVertIds[3] = max(maxVertIds[3], nEdgeCoeffs);
394 for (i = 0; i < exp->GetGeom()->GetNumFaces(); ++i)
396 meshFaceId = exp->GetGeom()->GetFid(i);
404 auto pIt = periodicFaces.find(meshFaceId);
405 if (pIt != periodicFaces.end())
407 meshFaceId = min(meshFaceId, pIt->second[0].id);
409 pIt->second[0].orient);
412 exp->GetTraceInteriorToElementMap(i, bmap,
sign, faceOrient);
414 ->GetTraceInverseBoundaryMap(i);
417 const int nFaceCoeffs = bmap.size();
418 vector<NekDouble> tmpStore(nFaceCoeffs * nFaceCoeffs);
420 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
422 for (j = 0; j < nFaceCoeffs; ++j)
425 asmMap->GetLocalToGlobalMap(cnt + bmap[j]) - nDirBnd);
437 for (k = 0; k < nFaceCoeffs; ++k)
439 tmpStore[k + j * nFaceCoeffs] =
440 sign1 *
sign[k] * (*schurMat)(bmap2[j], bmap2[k]);
449 gidDofs[2][gId] = nFaceCoeffs;
452 auto gIt = idMats[2].find(gId);
453 if (gIt == idMats[2].end())
455 idMats[2][gId] = tmpStore;
459 ASSERTL1(tmpStore.size() == gIt->second.size(),
460 "Number of modes mismatch");
461 Vmath::Vadd(nFaceCoeffs * nFaceCoeffs, &gIt->second[0], 1,
462 &tmpStore[0], 1, &gIt->second[0], 1);
465 gidMeshIds[2][gId] = meshFaceId;
466 maxVertIds[4] = max(maxVertIds[4], meshFaceId);
467 maxVertIds[5] = max(maxVertIds[5], nFaceCoeffs);
470 cnt += exp->GetNcoeffs();
475 m_comm = expList->GetSession()->GetComm()->GetRowComm();
480 vector<NekDouble> storageBuf;
481 vector<long> globalToUniversal;
483 for (i = 0, cnt = 1; i < 3; ++i)
485 const int maxDofs = maxVertIds[2 * i + 1];
490 for (
auto &gIt : idMats[i])
493 storageBuf.insert(storageBuf.end(), gIt.second.begin(),
497 ASSERTL1(gidMeshIds[i].count(gIt.first) > 0,
498 "Unable to find global ID " +
499 boost::lexical_cast<string>(gIt.first) +
501 meshVertId = gidMeshIds[i][gIt.first];
503 for (j = 0; j < gIt.second.size(); ++j)
505 globalToUniversal.push_back(cnt +
506 meshVertId * maxDofs * maxDofs + j);
513 cnt += (maxVertIds[2 * i] + 1) * maxDofs * maxDofs;
516 ASSERTL1(storageBuf.size() == globalToUniversal.size(),
517 "Storage buffer and global to universal map size does "
522 &globalToUniversal[0]);
527 expList->GetSession()->DefinesCmdLineArgument(
"verbose"));
531 int nblksSize = 1 + idMats[1].size() + idMats[2].size();
534 nblksSize += idMats[3].size();
539 n_blks[0] = idMats[0].size();
544 for (i = 1; i < nStorage; ++i)
546 for (
auto &gIt : idMats[i])
548 ASSERTL1(gidDofs[i].count(gIt.first) > 0,
549 "Unable to find number of degrees of freedom for "
551 boost::lexical_cast<string>(gIt.first));
553 n_blks[cnt++] = gidDofs[i][gIt.first];
568 for (
auto gIt = idMats[0].begin(); gIt != idMats[0].end(); ++gIt, ++cnt)
570 (*vertMat)(cnt, cnt) = 1.0 / storageData[cnt];
579 for (i = 1; i < nStorage; ++i)
581 for (
auto &gIt : idMats[i])
583 int nDofs = gidDofs[i][gIt.first];
588 for (j = 0; j < nDofs; ++j)
590 for (k = 0; k < nDofs; ++k)
596 (*tmp)(j, k) = gIt.second[k + j * nDofs];
601 (*tmp)(j, k) = storageData[k + j * nDofs + cnt];
606 cnt += nDofs * nDofs;
609 m_blkMat->SetBlock(cnt2, cnt2, tmp);
630 std::shared_ptr<MultiRegions::ExpList> expList =
631 ((
m_linsys.lock())->GetLocMat()).lock();
632 std::shared_ptr<MultiRegions::ExpList> trace = expList->GetTrace();
638 std::dynamic_pointer_cast<AssemblyMapDG>(
m_locToGloMap.lock());
640 int i, j, k, n, cnt, cnt2;
643 int nTrace = expList->GetTrace()->GetExpSize();
644 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
646 for (cnt = n = 0; n < nTrace; ++n)
653 cnt += trace->GetExp(n)->GetNcoeffs();
660 int maxTraceSize = -1;
661 map<int, int> arrayOffsets;
663 for (cnt = 0, n = nDir; n < nTrace; ++n)
665 int nCoeffs = trace->GetExp(n)->GetNcoeffs();
666 int nCoeffs2 = nCoeffs * nCoeffs;
668 arrayOffsets[n] = cnt;
671 if (maxTraceSize < nCoeffs2)
673 maxTraceSize = nCoeffs2;
678 m_comm = expList->GetSession()->GetComm()->GetRowComm();
685 for (cnt = n = 0; n < expList->GetExpSize(); ++n)
688 locExpansion = expList->GetExp(elmt);
691 asmMap->GetElmtToTrace()[elmt];
694 loc_mat = (
m_linsys.lock())->GetStaticCondBlock(n);
695 bnd_mat = loc_mat->GetBlock(0, 0);
697 int nFacets = locExpansion->GetNtraces();
699 for (cnt2 = i = 0; i < nFacets; ++i)
701 int nCoeffs = elmtToTraceMap[i]->GetNcoeffs();
702 int elmtId = elmtToTraceMap[i]->GetElmtId();
712 NekDouble *off = &tmpStore[arrayOffsets[elmtId]];
714 for (j = 0; j < nCoeffs; ++j)
716 NekDouble sign1 = asmMap->GetLocalToGlobalBndSign(cnt + j);
717 for (k = 0; k < nCoeffs; ++k)
719 NekDouble sign2 = asmMap->GetLocalToGlobalBndSign(cnt + k);
720 off[j * nCoeffs + k] +=
721 (*bnd_mat)(cnt2 + j, cnt2 + k) * sign1 * sign2;
732 for (cnt = 0, n = nDir; n < nTrace; ++n)
735 int nCoeffs = traceExp->GetNcoeffs();
736 int geomId = traceExp->GetGeom()->GetGlobalID();
738 for (i = 0; i < nCoeffs * nCoeffs; ++i)
740 uniIds[cnt++] = geomId * maxTraceSize + i + 1;
747 expList->GetSession()->DefinesCmdLineArgument(
"verbose"));
752 for (n = 0; n < nTrace - nDir; ++n)
754 n_blks[n] = trace->GetExp(n + nDir)->GetNcoeffs();
760 for (n = 0; n < nTrace - nDir; ++n)
762 int nCoeffs = trace->GetExp(n + nDir)->GetNcoeffs();
765 NekDouble *off = &tmpStore[arrayOffsets[n + nDir]];
767 for (i = 0; i < nCoeffs; ++i)
769 for (j = 0; j < nCoeffs; ++j)
771 (*tmp)(i, j) = off[i * nCoeffs + j];
791 int nGlobal =
m_isFull ? asmMap->GetNumGlobalCoeffs()
792 : asmMap->GetNumGlobalBndCoeffs();
793 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
794 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.
void BlockPreconditionerHDG(void)
Construct a block preconditioner for the hybridized discontinuous Galerkin system.
virtual void v_InitObject() override
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
Apply preconditioner to pInput and store the result in pOutput.
DNekBlkMatSharedPtr m_blkMat
virtual void v_BuildPreconditioner() override
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.