49 namespace MultiRegions
54 string PreconditionerBlock::className
57 PreconditionerBlock::create,
58 "Block Diagonal Preconditioning");
66 PreconditionerBlock::PreconditionerBlock(
67 const std::shared_ptr<GlobalLinSys> &plinsys,
79 "Solver type not valid");
120 int i, j, k, n, cnt, gId;
121 int meshVertId, meshEdgeId, meshFaceId;
125 const int nExp = expList->GetExpSize();
126 const int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
129 PeriodicMap periodicVerts, periodicEdges, periodicFaces;
130 expList->GetPeriodicEntities(
131 periodicVerts, periodicEdges, periodicFaces);
138 vector<map<int, vector<NekDouble> > > idMats(3);
142 vector<map<int, int> > gidMeshIds(3);
146 vector<map<int, int> > gidDofs(3);
158 for (cnt = n = 0; n < nExp; ++n)
160 exp = expList->GetExp(n);
164 m_linsys.lock()->GetStaticCondBlock(n)->GetBlock(0,0);
168 for (i = 0; i < exp->GetNverts(); ++i)
170 meshVertId = exp->GetGeom()->GetVid(i);
171 int locId = exp->GetVertexMap(i);
174 gId = asmMap->GetLocalToGlobalMap(
175 cnt + locId) - nDirBnd;
186 NekDouble vertVal = (*schurMat)(locId,locId);
190 auto gIt = idMats[0].find(gId);
192 if (gIt == idMats[0].end())
195 idMats[0][gId] = vector<NekDouble>(1, vertVal);
201 gIt->second[0] += vertVal;
208 auto pIt = periodicVerts.find(meshVertId);
209 if (pIt != periodicVerts.end())
211 for (j = 0; j < pIt->second.size(); ++j)
213 meshVertId = min(meshVertId, pIt->second[j].id);
219 gidMeshIds[0][gId] = meshVertId;
220 maxVertIds[0] = max(maxVertIds[0], meshVertId);
226 for (i = 0; i < exp->GetNedges(); ++i)
228 meshEdgeId = exp->GetGeom()->GetEid(i);
233 exp->GetGeom()->GetEorient(i);
237 auto pIt = periodicEdges.find(meshEdgeId);
238 if (pIt != periodicEdges.end())
240 pair<int, StdRegions::Orientation> idOrient =
242 meshEdgeId, edgeOrient, pIt->second);
243 meshEdgeId = idOrient.first;
244 edgeOrient = idOrient.second;
250 exp->GetEdgeInteriorMap(i, edgeOrient, bmap, sign);
251 bmap2 = exp->GetEdgeInverseBoundaryMap(i);
254 const int nEdgeCoeffs = bmap.num_elements();
255 vector<NekDouble> tmpStore(nEdgeCoeffs*nEdgeCoeffs);
257 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
259 for (j = 0; j < nEdgeCoeffs; ++j)
267 asmMap->GetLocalToGlobalMap(
281 for (k = 0; k < nEdgeCoeffs; ++k)
283 tmpStore[k+j*nEdgeCoeffs] =
284 sign1*sign[k]*(*schurMat)(bmap2[j], bmap2[k]);
293 gidDofs[1][gId] = nEdgeCoeffs;
296 auto gIt = idMats[1].find(gId);
297 if (gIt == idMats[1].end())
299 idMats[1][gId] = tmpStore;
303 ASSERTL1(tmpStore.size() == gIt->second.size(),
304 "Number of modes mismatch");
305 Vmath::Vadd(nEdgeCoeffs*nEdgeCoeffs, &gIt->second[0], 1,
306 &tmpStore[0], 1, &gIt->second[0], 1);
309 gidMeshIds[1][gId] = meshEdgeId;
310 maxVertIds[2] = max(maxVertIds[2], meshEdgeId);
311 maxVertIds[3] = max(maxVertIds[3], nEdgeCoeffs);
316 for (i = 0; i < exp->GetNfaces(); ++i)
318 meshFaceId = exp->GetGeom()->GetFid(i);
323 exp->GetGeom()->GetForient(i);
327 auto pIt = periodicFaces.find(meshFaceId);
328 if (pIt != periodicFaces.end())
330 meshFaceId = min(meshFaceId, pIt->second[0].id);
332 faceOrient, pIt->second[0].orient);
335 exp->GetFaceInteriorMap(i, faceOrient, bmap, sign);
336 bmap2 = exp->GetFaceInverseBoundaryMap(i);
339 const int nFaceCoeffs = bmap.num_elements();
340 vector<NekDouble> tmpStore(nFaceCoeffs*nFaceCoeffs);
342 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
344 for (j = 0; j < nFaceCoeffs; ++j)
347 asmMap->GetLocalToGlobalMap(cnt + bmap[j])
360 for (k = 0; k < nFaceCoeffs; ++k)
362 tmpStore[k+j*nFaceCoeffs] =
363 sign1*sign[k]*(*schurMat)(bmap2[j], bmap2[k]);
372 gidDofs[2][gId] = nFaceCoeffs;
375 auto gIt = idMats[2].find(gId);
376 if (gIt == idMats[2].end())
378 idMats[2][gId] = tmpStore;
382 ASSERTL1(tmpStore.size() == gIt->second.size(),
383 "Number of modes mismatch");
384 Vmath::Vadd(nFaceCoeffs*nFaceCoeffs, &gIt->second[0], 1,
385 &tmpStore[0], 1, &gIt->second[0], 1);
388 gidMeshIds[2][gId] = meshFaceId;
389 maxVertIds[4] = max(maxVertIds[4], meshFaceId);
390 maxVertIds[5] = max(maxVertIds[5], nFaceCoeffs);
393 cnt += exp->GetNcoeffs();
398 m_comm = expList->GetSession()->GetComm()->GetRowComm();
403 vector<NekDouble> storageBuf;
404 vector<long> globalToUniversal;
406 for (i = 0, cnt = 1; i < 3; ++i)
408 const int maxDofs = maxVertIds[2*i+1];
413 for (
auto &gIt : idMats[i])
416 storageBuf.insert(storageBuf.end(),
417 gIt.second.begin(), gIt.second.end());
420 ASSERTL1(gidMeshIds[i].count(gIt.first) > 0,
421 "Unable to find global ID " +
422 boost::lexical_cast<string>(gIt.first) +
424 meshVertId = gidMeshIds[i][gIt.first];
426 for (j = 0; j < gIt.second.size(); ++j)
428 globalToUniversal.push_back(
429 cnt + meshVertId*maxDofs*maxDofs + j);
436 cnt += (maxVertIds[2*i]+1)*maxDofs*maxDofs;
439 ASSERTL1(storageBuf.size() == globalToUniversal.size(),
440 "Storage buffer and global to universal map size does " 444 storageBuf.size(), &storageBuf[0]);
446 globalToUniversal.size(), &globalToUniversal[0]);
450 globalToUniversalMap,
m_comm,
451 expList->GetSession()->DefinesCmdLineArgument(
"verbose"));
456 1 + idMats[1].size() + idMats[2].size());
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 " 470 "global ID " + boost::lexical_cast<string>(
473 n_blks[cnt++] = gidDofs[i][gIt.first];
488 for (
auto gIt = idMats[0].begin(); gIt != idMats[0].end();
491 (*vertMat)(cnt, cnt) = 1.0/storageData[cnt];
500 for (i = 1; i < 3; ++i)
502 for (
auto &gIt : idMats[i])
504 int nDofs = gidDofs[i][gIt.first];
509 for (j = 0; j < nDofs; ++j)
511 for (k = 0; k < nDofs; ++k)
513 (*tmp)(j,k) = storageData[k+j*nDofs + cnt];
520 m_blkMat->SetBlock(cnt2, cnt2, tmp);
541 std::shared_ptr<MultiRegions::ExpList>
542 expList = ((
m_linsys.lock())->GetLocMat()).lock();
543 std::shared_ptr<MultiRegions::ExpList> trace = expList->GetTrace();
551 int i, j, k, n, cnt, cnt2;
554 int nTrace = expList->GetTrace()->GetExpSize();
557 for (cnt = n = 0; n < nTrace; ++n)
564 cnt += trace->GetExp(n)->GetNcoeffs();
571 int maxTraceSize = -1;
572 map<int, int> arrayOffsets;
574 for (cnt = 0, n = nDir; n < nTrace; ++n)
576 int nCoeffs = trace->GetExp(n)->GetNcoeffs();
577 int nCoeffs2 = nCoeffs * nCoeffs;
579 arrayOffsets[n] = cnt;
582 if (maxTraceSize < nCoeffs2)
584 maxTraceSize = nCoeffs2;
589 m_comm = expList->GetSession()->GetComm()->GetRowComm();
596 for (cnt = n = 0; n < expList->GetExpSize(); ++n)
599 locExpansion = expList->GetExp(elmt);
602 asmMap->GetElmtToTrace()[elmt];
605 loc_mat = (
m_linsys.lock())->GetStaticCondBlock(n);
606 bnd_mat = loc_mat->GetBlock(0,0);
608 int nFacets = locExpansion->GetNumBases() == 2 ?
609 locExpansion->GetNedges() : locExpansion->GetNfaces();
611 for (cnt2 = i = 0; i < nFacets; ++i)
613 int nCoeffs = elmtToTraceMap[i]->GetNcoeffs();
614 int elmtId = elmtToTraceMap[i]->GetElmtId ();
624 NekDouble *off = &tmpStore[arrayOffsets[elmtId]];
626 for (j = 0; j < nCoeffs; ++j)
628 NekDouble sign1 = asmMap->GetLocalToGlobalBndSign(
630 for (k = 0; k < nCoeffs; ++k)
632 NekDouble sign2 = asmMap->GetLocalToGlobalBndSign(
634 off[j*nCoeffs + k] +=
635 (*bnd_mat)(cnt2+j, cnt2+k) * sign1 * sign2;
646 for (cnt = 0, n = nDir; n < nTrace; ++n)
649 int nCoeffs = traceExp->GetNcoeffs();
650 int geomId = traceExp->GetGeom()->GetGlobalID();
652 for (i = 0; i < nCoeffs*nCoeffs; ++i)
654 uniIds[cnt++] = geomId * maxTraceSize + i + 1;
661 expList->GetSession()->DefinesCmdLineArgument(
"verbose"));
666 for (n = 0; n < nTrace - nDir; ++n)
668 n_blks[n] = trace->GetExp(n + nDir)->GetNcoeffs();
674 for (n = 0; n < nTrace - nDir; ++n)
676 int nCoeffs = trace->GetExp(n + nDir)->GetNcoeffs();
679 NekDouble *off = &tmpStore[arrayOffsets[n+nDir]];
681 for (i = 0; i < nCoeffs; ++i)
683 for (j = 0; j < nCoeffs; ++j)
685 (*tmp)(i,j) = off[i*nCoeffs + j];
704 int nNonDir = nGlobal-nDir;
virtual void v_BuildPreconditioner()
LibUtilities::CommSharedPtr m_comm
#define ASSERTL0(condition, msg)
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply preconditioner to pInput and store the result in pOutput.
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.
#define sign(a, b)
return the sign(b)*a
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
const std::weak_ptr< GlobalLinSys > m_linsys
DNekBlkMatSharedPtr m_blkMat
PreconFactory & GetPreconFactory()
std::shared_ptr< DNekMat > DNekMatSharedPtr
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
StdRegions::Orientation DeterminePeriodicFaceOrient(StdRegions::Orientation faceOrient, StdRegions::Orientation perFaceOrient)
Determine relative orientation between two faces.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void BlockPreconditionerHDG(void)
Construct a block preconditioner for the hybridized discontinuous Galerkin system.
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
int GetNumGlobalDirBndCoeffs() const
Returns the number of global Dirichlet boundary coefficients.
Describe a linear system.
virtual void v_InitObject()
std::shared_ptr< Expansion > ExpansionSharedPtr
void BlockPreconditionerCG(void)
Construct a block preconditioner from for the continuous Galerkin system.
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
std::weak_ptr< AssemblyMap > m_locToGloMap
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.