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->GetGeom()->GetNumEdges(); ++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 if(exp->GetGeom()->GetNumFaces())
253 GetEdgeInteriorToElementMap(i, bmap,
sign, edgeOrient);
255 GetEdgeInverseBoundaryMap(i);
259 exp->GetTraceInteriorToElementMap(i, bmap,
sign, edgeOrient);
261 GetTraceInverseBoundaryMap(i);
265 const int nEdgeCoeffs = bmap.size();
266 vector<NekDouble> tmpStore(nEdgeCoeffs*nEdgeCoeffs);
268 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
270 for (j = 0; j < nEdgeCoeffs; ++j)
278 asmMap->GetLocalToGlobalMap(
292 for (k = 0; k < nEdgeCoeffs; ++k)
294 tmpStore[k+j*nEdgeCoeffs] =
295 sign1*
sign[k]*(*schurMat)(bmap2[j], bmap2[k]);
304 gidDofs[1][gId] = nEdgeCoeffs;
307 auto gIt = idMats[1].find(gId);
308 if (gIt == idMats[1].end())
310 idMats[1][gId] = tmpStore;
314 ASSERTL1(tmpStore.size() == gIt->second.size(),
315 "Number of modes mismatch");
316 Vmath::Vadd(nEdgeCoeffs*nEdgeCoeffs, &gIt->second[0], 1,
317 &tmpStore[0], 1, &gIt->second[0], 1);
320 gidMeshIds[1][gId] = meshEdgeId;
321 maxVertIds[2] = max(maxVertIds[2], meshEdgeId);
322 maxVertIds[3] = max(maxVertIds[3], nEdgeCoeffs);
327 for (i = 0; i < exp->GetGeom()->GetNumFaces(); ++i)
329 meshFaceId = exp->GetGeom()->GetFid(i);
334 exp->GetGeom()->GetForient(i);
338 auto pIt = periodicFaces.find(meshFaceId);
339 if (pIt != periodicFaces.end())
341 meshFaceId = min(meshFaceId, pIt->second[0].id);
343 faceOrient, pIt->second[0].orient);
346 exp->GetTraceInteriorToElementMap(i, bmap,
sign, faceOrient);
348 ->GetTraceInverseBoundaryMap(i);
351 const int nFaceCoeffs = bmap.size();
352 vector<NekDouble> tmpStore(nFaceCoeffs*nFaceCoeffs);
354 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
356 for (j = 0; j < nFaceCoeffs; ++j)
359 asmMap->GetLocalToGlobalMap(cnt + bmap[j])
372 for (k = 0; k < nFaceCoeffs; ++k)
374 tmpStore[k+j*nFaceCoeffs] =
375 sign1*
sign[k]*(*schurMat)(bmap2[j], bmap2[k]);
384 gidDofs[2][gId] = nFaceCoeffs;
387 auto gIt = idMats[2].find(gId);
388 if (gIt == idMats[2].end())
390 idMats[2][gId] = tmpStore;
394 ASSERTL1(tmpStore.size() == gIt->second.size(),
395 "Number of modes mismatch");
396 Vmath::Vadd(nFaceCoeffs*nFaceCoeffs, &gIt->second[0], 1,
397 &tmpStore[0], 1, &gIt->second[0], 1);
400 gidMeshIds[2][gId] = meshFaceId;
401 maxVertIds[4] = max(maxVertIds[4], meshFaceId);
402 maxVertIds[5] = max(maxVertIds[5], nFaceCoeffs);
405 cnt += exp->GetNcoeffs();
410 m_comm = expList->GetSession()->GetComm()->GetRowComm();
415 vector<NekDouble> storageBuf;
416 vector<long> globalToUniversal;
418 for (i = 0, cnt = 1; i < 3; ++i)
420 const int maxDofs = maxVertIds[2*i+1];
425 for (
auto &gIt : idMats[i])
428 storageBuf.insert(storageBuf.end(),
429 gIt.second.begin(), gIt.second.end());
432 ASSERTL1(gidMeshIds[i].count(gIt.first) > 0,
433 "Unable to find global ID " +
434 boost::lexical_cast<string>(gIt.first) +
436 meshVertId = gidMeshIds[i][gIt.first];
438 for (j = 0; j < gIt.second.size(); ++j)
440 globalToUniversal.push_back(
441 cnt + meshVertId*maxDofs*maxDofs + j);
448 cnt += (maxVertIds[2*i]+1)*maxDofs*maxDofs;
451 ASSERTL1(storageBuf.size() == globalToUniversal.size(),
452 "Storage buffer and global to universal map size does "
456 storageBuf.size(), &storageBuf[0]);
458 globalToUniversal.size(), &globalToUniversal[0]);
462 globalToUniversalMap,
m_comm,
463 expList->GetSession()->DefinesCmdLineArgument(
"verbose"));
468 1 + idMats[1].size() + idMats[2].size());
471 n_blks[0] = idMats[0].size();
476 for (i = 1; i < 3; ++i)
478 for (
auto &gIt : idMats[i])
480 ASSERTL1(gidDofs[i].count(gIt.first) > 0,
481 "Unable to find number of degrees of freedom for "
482 "global ID " + boost::lexical_cast<string>(
485 n_blks[cnt++] = gidDofs[i][gIt.first];
500 for (
auto gIt = idMats[0].begin(); gIt != idMats[0].end();
503 (*vertMat)(cnt, cnt) = 1.0/storageData[cnt];
512 for (i = 1; i < 3; ++i)
514 for (
auto &gIt : idMats[i])
516 int nDofs = gidDofs[i][gIt.first];
521 for (j = 0; j < nDofs; ++j)
523 for (k = 0; k < nDofs; ++k)
525 (*tmp)(j,k) = storageData[k+j*nDofs + cnt];
532 m_blkMat->SetBlock(cnt2, cnt2, tmp);
553 std::shared_ptr<MultiRegions::ExpList>
554 expList = ((
m_linsys.lock())->GetLocMat()).lock();
555 std::shared_ptr<MultiRegions::ExpList> trace = expList->GetTrace();
563 int i, j, k, n, cnt, cnt2;
566 int nTrace = expList->GetTrace()->GetExpSize();
569 for (cnt = n = 0; n < nTrace; ++n)
576 cnt += trace->GetExp(n)->GetNcoeffs();
583 int maxTraceSize = -1;
584 map<int, int> arrayOffsets;
586 for (cnt = 0, n = nDir; n < nTrace; ++n)
588 int nCoeffs = trace->GetExp(n)->GetNcoeffs();
589 int nCoeffs2 = nCoeffs * nCoeffs;
591 arrayOffsets[n] = cnt;
594 if (maxTraceSize < nCoeffs2)
596 maxTraceSize = nCoeffs2;
601 m_comm = expList->GetSession()->GetComm()->GetRowComm();
608 for (cnt = n = 0; n < expList->GetExpSize(); ++n)
611 locExpansion = expList->GetExp(elmt);
614 asmMap->GetElmtToTrace()[elmt];
617 loc_mat = (
m_linsys.lock())->GetStaticCondBlock(n);
618 bnd_mat = loc_mat->GetBlock(0,0);
620 int nFacets = locExpansion->GetNtraces();
622 for (cnt2 = i = 0; i < nFacets; ++i)
624 int nCoeffs = elmtToTraceMap[i]->GetNcoeffs();
625 int elmtId = elmtToTraceMap[i]->GetElmtId ();
635 NekDouble *off = &tmpStore[arrayOffsets[elmtId]];
637 for (j = 0; j < nCoeffs; ++j)
639 NekDouble sign1 = asmMap->GetLocalToGlobalBndSign(
641 for (k = 0; k < nCoeffs; ++k)
643 NekDouble sign2 = asmMap->GetLocalToGlobalBndSign(
645 off[j*nCoeffs + k] +=
646 (*bnd_mat)(cnt2+j, cnt2+k) * sign1 * sign2;
657 for (cnt = 0, n = nDir; n < nTrace; ++n)
660 int nCoeffs = traceExp->GetNcoeffs();
661 int geomId = traceExp->GetGeom()->GetGlobalID();
663 for (i = 0; i < nCoeffs*nCoeffs; ++i)
665 uniIds[cnt++] = geomId * maxTraceSize + i + 1;
672 expList->GetSession()->DefinesCmdLineArgument(
"verbose"));
677 for (n = 0; n < nTrace - nDir; ++n)
679 n_blks[n] = trace->GetExp(n + nDir)->GetNcoeffs();
685 for (n = 0; n < nTrace - nDir; ++n)
687 int nCoeffs = trace->GetExp(n + nDir)->GetNcoeffs();
690 NekDouble *off = &tmpStore[arrayOffsets[n+nDir]];
692 for (i = 0; i < nCoeffs; ++i)
694 for (j = 0; j < nCoeffs; ++j)
696 (*tmp)(i,j) = off[i*nCoeffs + j];
715 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.
int GetNumGlobalDirBndCoeffs() const
Returns the number of global Dirichlet boundary coefficients.
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.