46using namespace LibUtilities;
56 "LowEnergy Preconditioning");
66 const std::shared_ptr<GlobalLinSys> &plinsys,
78 "Solver type not valid");
80 std::shared_ptr<MultiRegions::ExpList> expList =
81 ((
m_linsys.lock())->GetLocMat()).lock();
83 m_comm = expList->GetComm();
87 locExpansion = expList->GetExp(0);
89 int nDim = locExpansion->GetShapeDimension();
91 ASSERTL0(nDim == 3,
"Preconditioner type only valid in 3D");
119 std::shared_ptr<MultiRegions::ExpList> expList =
120 ((
m_linsys.lock())->GetLocMat()).lock();
125 int nVerts, nEdges, nFaces;
126 int eid, fid, n, cnt, nmodes, nedgemodes, nfacemodes;
130 int vMap1, vMap2, sign1, sign2;
131 int m, v, eMap1, eMap2, fMap1, fMap2;
132 int offset, globalrow, globalcol, nCoeffs;
138 expList->GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
156 int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
157 int nNonDirVerts = asmMap->GetNumNonDirVertexModes();
161 nNonDirVerts, nNonDirVerts, zero, vertstorage);
167 int n_exp = expList->GetNumElmts();
168 int nNonDirEdgeIDs = asmMap->GetNumNonDirEdges();
169 int nNonDirFaceIDs = asmMap->GetNumNonDirFaces();
173 map<int, int> uniqueEdgeMap;
174 map<int, int> uniqueFaceMap;
181 expList->GetBndCondExpansions();
184 &bndConditions = expList->GetBndConditions();
191 for (i = 0; i < extradiredges.size(); ++i)
193 meshEdgeId = extradiredges[i];
194 edgeDirMap.insert(meshEdgeId);
198 for (i = 0; i < bndCondExp.size(); i++)
201 for (j = 0; j < bndCondExp[i]->GetNumElmts(); j++)
204 std::dynamic_pointer_cast<LocalRegions::Expansion2D>(
205 bndCondExp[i]->GetExp(j));
206 if (bndConditions[i]->GetBoundaryConditionType() ==
209 for (k = 0; k < bndCondFaceExp->GetNtraces(); k++)
211 meshEdgeId = bndCondFaceExp->GetGeom()->GetEid(k);
212 if (edgeDirMap.count(meshEdgeId) == 0)
214 edgeDirMap.insert(meshEdgeId);
217 meshFaceId = bndCondFaceExp->GetGeom()->GetGlobalID();
218 faceDirMap.insert(meshFaceId);
226 int nlocalNonDirEdges = 0;
227 int nlocalNonDirFaces = 0;
228 int ntotalentries = 0;
230 map<int, int> EdgeSize;
231 map<int, int> FaceSize;
232 map<int, pair<int, int>> FaceModes;
235 for (n = 0; n < n_exp; ++n)
240 for (j = 0; j < nEdges; ++j)
242 int nEdgeInteriorCoeffs = locExpansion->GetEdgeNcoeffs(j) - 2;
246 if (EdgeSize.count(meshEdgeId) == 0)
248 EdgeSize[meshEdgeId] = nEdgeInteriorCoeffs;
252 EdgeSize[meshEdgeId] =
253 min(EdgeSize[meshEdgeId], nEdgeInteriorCoeffs);
257 nFaces = locExpansion->GetNtraces();
258 for (j = 0; j < nFaces; ++j)
260 int nFaceInteriorCoeffs = locExpansion->GetTraceIntNcoeffs(j);
261 meshFaceId = locExpansion->GetGeom3D()->GetFid(j);
263 if (FaceSize.count(meshFaceId) == 0)
265 FaceSize[meshFaceId] = nFaceInteriorCoeffs;
268 locExpansion->GetTraceNumModes(j, m0, m1,
269 locExpansion->GetTraceOrient(j));
270 FaceModes[meshFaceId] = pair<int, int>(m0, m1);
274 if (nFaceInteriorCoeffs < FaceSize[meshFaceId])
276 FaceSize[meshFaceId] = nFaceInteriorCoeffs;
278 locExpansion->GetTraceNumModes(
279 j, m0, m1, locExpansion->GetTraceOrient(j));
280 FaceModes[meshFaceId] = pair<int, int>(m0, m1);
286 bool verbose = expList->GetSession()->DefinesCmdLineArgument(
"verbose");
290 if (
m_comm->GetSize() > 1)
292 int EdgeSizeLen = EdgeSize.size();
293 int FaceSizeLen = FaceSize.size();
297 map<int, int>::iterator it;
301 for (it = EdgeSize.begin(); it != EdgeSize.end(); ++it, ++cnt)
303 FacetMap[cnt] = it->first;
304 maxid = max(it->first, maxid);
305 FacetLen[cnt] = it->second;
311 for (it = FaceSize.begin(); it != FaceSize.end(); ++it, ++cnt)
313 FacetMap[cnt] = it->first + maxid;
314 FacetLen[cnt] = it->second;
322 for (it = EdgeSize.begin(); it != EdgeSize.end(); ++it, ++cnt)
324 it->second = (int)FacetLen[cnt];
327 for (it = FaceSize.begin(); it != FaceSize.end(); ++it, ++cnt)
329 it->second = (int)FacetLen[cnt];
336 int matrixlocation = 0;
339 for (
auto &pIt : periodicEdges)
341 meshEdgeId = pIt.first;
343 if (edgeDirMap.count(meshEdgeId) == 0)
345 dof = EdgeSize[meshEdgeId];
346 if (uniqueEdgeMap.count(meshEdgeId) == 0 && dof > 0)
348 bool SetUpNewEdge =
true;
350 for (i = 0; i < pIt.second.size(); ++i)
352 if (!pIt.second[i].isLocal)
357 int meshEdgeId2 = pIt.second[i].id;
358 if (edgeDirMap.count(meshEdgeId2) == 0)
360 if (uniqueEdgeMap.count(meshEdgeId2) != 0)
363 uniqueEdgeMap[meshEdgeId] =
364 uniqueEdgeMap[meshEdgeId2];
365 SetUpNewEdge =
false;
370 edgeDirMap.insert(meshEdgeId);
371 SetUpNewEdge =
false;
377 uniqueEdgeMap[meshEdgeId] = matrixlocation;
378 globaloffset[matrixlocation] += ntotalentries;
379 modeoffset[matrixlocation] = dof * dof;
380 ntotalentries += dof * dof;
381 nblks[matrixlocation++] = dof;
387 for (cnt = n = 0; n < n_exp; ++n)
391 for (j = 0; j < locExpansion->GetNedges(); ++j)
393 meshEdgeId = locExpansion->
GetGeom()->GetEid(j);
394 dof = EdgeSize[meshEdgeId];
395 maxEdgeDof = (dof > maxEdgeDof ? dof : maxEdgeDof);
397 if (edgeDirMap.count(meshEdgeId) == 0)
399 if (uniqueEdgeMap.count(meshEdgeId) == 0 && dof > 0)
402 uniqueEdgeMap[meshEdgeId] = matrixlocation;
404 globaloffset[matrixlocation] += ntotalentries;
405 modeoffset[matrixlocation] = dof * dof;
406 ntotalentries += dof * dof;
407 nblks[matrixlocation++] = dof;
409 nlocalNonDirEdges += dof * dof;
417 for (
auto &pIt : periodicFaces)
419 meshFaceId = pIt.first;
421 if (faceDirMap.count(meshFaceId) == 0)
423 dof = FaceSize[meshFaceId];
425 if (uniqueFaceMap.count(meshFaceId) == 0 && dof > 0)
427 bool SetUpNewFace =
true;
429 if (pIt.second[0].isLocal)
431 int meshFaceId2 = pIt.second[0].id;
433 if (faceDirMap.count(meshFaceId2) == 0)
435 if (uniqueFaceMap.count(meshFaceId2) != 0)
438 uniqueFaceMap[meshFaceId] =
439 uniqueFaceMap[meshFaceId2];
440 SetUpNewFace =
false;
445 faceDirMap.insert(meshFaceId);
446 SetUpNewFace =
false;
452 uniqueFaceMap[meshFaceId] = matrixlocation;
454 modeoffset[matrixlocation] = dof * dof;
455 globaloffset[matrixlocation] += ntotalentries;
456 ntotalentries += dof * dof;
457 nblks[matrixlocation++] = dof;
463 for (cnt = n = 0; n < n_exp; ++n)
467 for (j = 0; j < locExpansion->GetNtraces(); ++j)
469 meshFaceId = locExpansion->
GetGeom()->GetFid(j);
471 dof = FaceSize[meshFaceId];
472 maxFaceDof = (dof > maxFaceDof ? dof : maxFaceDof);
474 if (faceDirMap.count(meshFaceId) == 0)
476 if (uniqueFaceMap.count(meshFaceId) == 0 && dof > 0)
478 uniqueFaceMap[meshFaceId] = matrixlocation;
479 modeoffset[matrixlocation] = dof * dof;
480 globaloffset[matrixlocation] += ntotalentries;
481 ntotalentries += dof * dof;
482 nblks[matrixlocation++] = dof;
484 nlocalNonDirFaces += dof * dof;
495 nlocalNonDirEdges + nlocalNonDirFaces, -1);
501 int matrixoffset = 0;
503 int uniEdgeOffset = 0;
507 for (n = 0; n < n_exp; ++n)
509 for (j = 0; j < locExpansion->GetNedges(); ++j)
515 uniEdgeOffset = max(meshEdgeId, uniEdgeOffset);
521 uniEdgeOffset = uniEdgeOffset * maxEdgeDof * maxEdgeDof;
523 for (n = 0; n < n_exp; ++n)
528 for (j = 0; j < locExpansion->GetNedges(); ++j)
531 meshEdgeId = locExpansion->
GetGeom()->GetEid(j);
533 nedgemodes = EdgeSize[meshEdgeId];
535 if (edgeDirMap.count(meshEdgeId) == 0 && nedgemodes > 0)
538 int edgeOffset = globaloffset[uniqueEdgeMap[meshEdgeId]];
541 int uniOffset = meshEdgeId;
542 auto pIt = periodicEdges.find(meshEdgeId);
543 if (pIt != periodicEdges.end())
545 for (
int l = 0; l < pIt->second.size(); ++l)
547 uniOffset = min(uniOffset, pIt->second[l].id);
550 uniOffset = uniOffset * maxEdgeDof * maxEdgeDof;
552 for (k = 0; k < nedgemodes * nedgemodes; ++k)
554 vGlobal = edgeOffset + k;
555 localToGlobalMatrixMap[matrixoffset + k] = vGlobal;
556 BlockToUniversalMap[vGlobal] = uniOffset + k + 1;
558 matrixoffset += nedgemodes * nedgemodes;
565 for (j = 0; j < locExpansion->GetNtraces(); ++j)
568 meshFaceId = locExpansion->GetGeom()->GetFid(j);
570 nfacemodes = FaceSize[meshFaceId];
573 if (faceDirMap.count(meshFaceId) == 0 && nfacemodes > 0)
576 int faceOffset = globaloffset[uniqueFaceMap[meshFaceId]];
578 int uniOffset = meshFaceId;
580 auto pIt = periodicFaces.find(meshFaceId);
581 if (pIt != periodicFaces.end())
583 uniOffset = min(uniOffset, pIt->second[0].id);
585 uniOffset = uniOffset * maxFaceDof * maxFaceDof;
587 for (k = 0; k < nfacemodes * nfacemodes; ++k)
589 vGlobal = faceOffset + k;
591 localToGlobalMatrixMap[matrixoffset + k] = vGlobal;
593 BlockToUniversalMap[vGlobal] =
594 uniOffset + uniEdgeOffset + k + 1;
596 matrixoffset += nfacemodes * nfacemodes;
603 map<int, int>::iterator it;
605 n_blks[0] = nNonDirVerts;
606 for (i = 1, it = nblks.begin(); it != nblks.end(); ++it)
608 n_blks[i++] = it->second;
617 for (cnt = n = 0; n < n_exp; ++n)
629 nVerts = locExpansion->GetGeom()->GetNumVerts();
630 nEdges = locExpansion->GetGeom()->GetNumEdges();
631 nFaces = locExpansion->GetGeom()->GetNumFaces();
634 loc_mat = (
m_linsys.lock())->GetStaticCondBlock(n);
637 bnd_mat = loc_mat->GetBlock(0, 0);
640 offset = bnd_mat->GetRows();
644 DNekMat Sloc(nCoeffs, nCoeffs);
649 for (
int i = 0; i < nCoeffs; ++i)
651 for (
int j = 0; j < nCoeffs; ++j)
655 Sloc.SetValue(i, j, val);
664 for (v = 0; v < nVerts; ++v)
666 vMap1 = locExpansion->GetVertexMap(v);
669 globalrow = asmMap->GetLocalToGlobalBndMap(cnt + vMap1) - nDirBnd;
673 for (m = 0; m < nVerts; ++m)
675 vMap2 = locExpansion->GetVertexMap(m);
680 asmMap->GetLocalToGlobalBndMap(cnt + vMap2) - nDirBnd;
683 if (globalcol == globalrow)
686 sign1 = asmMap->GetLocalToGlobalBndSign(cnt + vMap1);
687 sign2 = asmMap->GetLocalToGlobalBndSign(cnt + vMap2);
689 vertArray[globalrow] +=
690 sign1 * sign2 * RSRT(vMap1, vMap2);
692 meshVertId = locExpansion->GetGeom()->GetVid(v);
694 auto pIt = periodicVerts.find(meshVertId);
695 if (pIt != periodicVerts.end())
697 for (k = 0; k < pIt->second.size(); ++k)
699 meshVertId = min(meshVertId, pIt->second[k].id);
703 VertBlockToUniversalMap[globalrow] = meshVertId + 1;
710 for (eid = 0; eid < nEdges; ++eid)
716 nedgemodes = EdgeSize[meshEdgeId];
722 nedgemodes, nedgemodes, zero, storage);
725 locExpansion->GetEdgeInverseBoundaryMap(eid);
727 if (edgeDirMap.count(meshEdgeId) == 0)
729 for (v = 0; v < nedgemodesloc; ++v)
731 eMap1 = edgemodearray[v];
732 sign1 = asmMap->GetLocalToGlobalBndSign(cnt + eMap1);
739 for (m = 0; m < nedgemodesloc; ++m)
741 eMap2 = edgemodearray[m];
745 asmMap->GetLocalToGlobalBndSign(cnt + eMap2);
748 sign1 * sign2 * RSRT(eMap1, eMap2);
753 BlockArray[matrixoffset + v * nedgemodes + m] =
758 matrixoffset += nedgemodes * nedgemodes;
764 for (fid = 0; fid < nFaces; ++fid)
770 nfacemodes = FaceSize[meshFaceId];
775 nfacemodes, nfacemodes, zero, storage);
777 if (faceDirMap.count(meshFaceId) == 0)
781 locExpansion->GetTraceOrient(fid);
783 auto pIt = periodicFaces.find(meshFaceId);
784 if (pIt != periodicFaces.end())
786 if (meshFaceId == min(meshFaceId, pIt->second[0].id))
789 faceOrient, pIt->second[0].orient);
793 facemodearray = locExpansion->GetTraceInverseBoundaryMap(
794 fid, faceOrient, FaceModes[meshFaceId].first,
795 FaceModes[meshFaceId].second);
797 for (v = 0; v < nfacemodes; ++v)
799 fMap1 = facemodearray[v];
801 sign1 = asmMap->GetLocalToGlobalBndSign(cnt + fMap1);
804 "Something is wrong since we "
805 "shoudl only be extracting modes for "
806 "lowest order expansion");
808 for (m = 0; m < nfacemodes; ++m)
810 fMap2 = facemodearray[m];
814 asmMap->GetLocalToGlobalBndSign(cnt + fMap2);
817 "Something is wrong since "
818 "we shoudl only be extracting modes for "
819 "lowest order expansion");
824 sign1 * sign2 * RSRT(fMap1, fMap2);
828 BlockArray[matrixoffset + v * nfacemodes + m] =
832 matrixoffset += nfacemodes * nfacemodes;
849 Vmath::Assmb(BlockArray.size(), BlockArray, localToGlobalMatrixMap,
858 for (
int i = 0; i < nNonDirVerts; ++i)
860 VertBlk->SetValue(i, i, 1.0 / vertArray[i]);
871 for (
int loc = 0;
loc < n_blks.size() - 1; ++
loc)
873 nmodes = n_blks[1 +
loc];
877 for (v = 0; v < nmodes; ++v)
879 for (m = 0; m < nmodes; ++m)
881 NekDouble Value = GlobalBlock[offset + v * nmodes + m];
882 gmat->SetValue(v, m, Value);
886 offset += modeoffset[
loc];
890 int totblks =
m_BlkMat->GetNumberOfBlockRows();
891 for (i = 1; i < totblks; ++i)
893 unsigned int nmodes =
m_BlkMat->GetNumberOfRowsInBlockRow(i);
917 ASSERTL0(isLocal ==
false,
"PreconditionerLowEnergy is only currently "
918 "set up for Global iteratives sovles");
921 int nNonDir = nGlobal - nDir;
938 std::shared_ptr<MultiRegions::ExpList> expList =
939 ((
m_linsys.lock())->GetLocMat()).lock();
942 map<int, int> EdgeSize;
946 std::map<ShapeType, DNekScalMatSharedPtr> maxRmat;
947 map<ShapeType, LocalRegions::Expansion3DSharedPtr> maxElmt;
948 map<ShapeType, Array<OneD, unsigned int>> vertMapMaxR;
949 map<ShapeType, Array<OneD, Array<OneD, unsigned int>>> edgeMapMaxR;
958 int n_exp = expList->GetNumElmts();
964 nbdry_size, nbdry_size, blkmatStorage);
966 nbdry_size, nbdry_size, blkmatStorage);
976 for (n = 0; n < n_exp; ++n)
982 int nbndcoeffs = locExp->NumBndryCoeffs();
986 rmat =
ExtractLocMat(locExp, maxRmat[eltype], maxElmt[eltype],
987 vertMapMaxR[eltype], edgeMapMaxR[eltype]);
989 m_RBlk->SetBlock(n, n, rmat);
997 m_sameBlock.push_back(pair<int, int>(1, nbndcoeffs));
1006 for (
int i = 0; i < 3; ++i)
1008 if (locExpSav->GetBasis(i) != locExp->GetBasis(i))
1017 m_RBlk->SetBlock(n, n, rmat);
1021 (pair<int, int>(
m_sameBlock[offset].first + 1, nbndcoeffs));
1025 rmat =
ExtractLocMat(locExp, maxRmat[eltype], maxElmt[eltype],
1026 vertMapMaxR[eltype], edgeMapMaxR[eltype]);
1029 m_RBlk->SetBlock(n, n, rmat);
1036 m_sameBlock.push_back(pair<int, int>(1, nbndcoeffs));
1054 int nLocBndDofs =
m_locToGloMap.lock()->GetNumLocalBndCoeffs();
1071 Blas::Dgemm(
'N',
'N', nbndcoeffs, nexp, nbndcoeffs, 1.0,
1072 &(R.GetBlock(cnt1, cnt1)->GetPtr()[0]), nbndcoeffs,
1073 pLocalIn.get() + cnt, nbndcoeffs, 0.0, pInOut.get() + cnt,
1075 cnt += nbndcoeffs * nexp;
1092 int nLocBndDofs =
m_locToGloMap.lock()->GetNumLocalBndCoeffs();
1094 ASSERTL1(pInOut.size() >= nLocBndDofs,
1095 "Output array is not greater than the nLocBndDofs");
1109 Blas::Dgemm(
'T',
'N', nbndcoeffs, nexp, nbndcoeffs, 1.0,
1110 &(R.GetBlock(cnt1, cnt1)->GetPtr()[0]), nbndcoeffs,
1111 pLocalIn.get() + cnt, nbndcoeffs, 0.0, pInOut.get() + cnt,
1113 cnt += nbndcoeffs * nexp;
1130 int nLocBndDofs =
m_locToGloMap.lock()->GetNumLocalBndCoeffs();
1132 ASSERTL1(pInput.size() >= nLocBndDofs,
1133 "Input array is smaller than nLocBndDofs");
1134 ASSERTL1(pOutput.size() >= nLocBndDofs,
1135 "Output array is smaller than nLocBndDofs");
1149 Blas::Dgemm(
'N',
'N', nbndcoeffs, nexp, nbndcoeffs, 1.0,
1150 &(invR.GetBlock(cnt1, cnt1)->GetPtr()[0]), nbndcoeffs,
1151 pLocalIn.get() + cnt, nbndcoeffs, 0.0, pOutput.get() + cnt,
1153 cnt += nbndcoeffs * nexp;
1168 int nLocBndDofs =
m_locToGloMap.lock()->GetNumLocalBndCoeffs();
1170 ASSERTL1(pInput.size() >= nLocBndDofs,
1171 "Input array is less than nLocBndDofs");
1172 ASSERTL1(pOutput.size() >= nLocBndDofs,
1173 "Output array is less than nLocBndDofs");
1187 Blas::Dgemm(
'T',
'N', nbndcoeffs, nexp, nbndcoeffs, 1.0,
1188 &(invR.GetBlock(cnt1, cnt1)->GetPtr()[0]), nbndcoeffs,
1189 pLocalIn.get() + cnt, nbndcoeffs, 0.0, pOutput.get() + cnt,
1191 cnt += nbndcoeffs * nexp;
1204 int n,
int bndoffset,
const std::shared_ptr<DNekScalMat> &loc_mat)
1206 std::shared_ptr<MultiRegions::ExpList> expList =
1207 ((
m_linsys.lock())->GetLocMat()).lock();
1210 locExpansion = expList->GetExp(n);
1211 unsigned int nbnd = locExpansion->NumBndryCoeffs();
1228 for (
int i = 0; i < nbnd; ++i)
1230 for (
int j = 0; j < nbnd; ++j)
1234 Sloc.SetValue(i, j, val);
1255 unsigned int nLocBnd = asmMap->GetNumLocalBndCoeffs();
1266 asmMap->GetLocalToGlobalBndSign();
1268 for (i = 0; i < nLocBnd; ++i)
1285 const int three = 3;
1286 const int nVerts = 6;
1287 const double point[][3] = {
1292 {0, -1,
sqrt(
double(3))},
1293 {0, 1,
sqrt(
double(3))},
1298 for (
int i = 0; i < nVerts; ++i)
1301 three, i, point[i][0], point[i][1], point[i][2]);
1303 const int nEdges = 9;
1304 const int vertexConnectivity[][2] = {{0, 1}, {1, 2}, {3, 2}, {0, 3}, {0, 4},
1305 {1, 4}, {2, 5}, {3, 5}, {4, 5}};
1309 for (
int i = 0; i < nEdges; ++i)
1312 for (
int j = 0; j < 2; ++j)
1314 vertsArray[j] = verts[vertexConnectivity[i][j]];
1317 i, three, vertsArray);
1324 const int nFaces = 5;
1326 const int quadEdgeConnectivity[][4] = {
1327 {0, 1, 2, 3}, {1, 6, 8, 5}, {3, 7, 8, 4}};
1329 const int quadId[] = {0, -1, 1, -1, 2};
1332 const int triEdgeConnectivity[][3] = {{0, 5, 4}, {2, 6, 7}};
1334 const int triId[] = {-1, 0, -1, 1, -1};
1338 for (
int f = 0; f < nFaces; ++f)
1340 if (f == 1 || f == 3)
1344 for (
int j = 0; j < 3; ++j)
1346 edgeArray[j] = edges[triEdgeConnectivity[i][j]];
1356 for (
int j = 0; j < 4; ++j)
1358 edgeArray[j] = edges[quadEdgeConnectivity[i][j]];
1382 const int nVerts = 5;
1383 const double point[][3] = {{-1, -1, 0},
1387 {0, 0,
sqrt(
double(2))}};
1390 const int three = 3;
1392 for (
int i = 0; i < nVerts; ++i)
1395 three, i, point[i][0], point[i][1], point[i][2]);
1397 const int nEdges = 8;
1398 const int vertexConnectivity[][2] = {{0, 1}, {1, 2}, {2, 3}, {3, 0},
1399 {0, 4}, {1, 4}, {2, 4}, {3, 4}};
1403 for (
int i = 0; i < nEdges; ++i)
1406 for (
int j = 0; j < 2; ++j)
1408 vertsArray[j] = verts[vertexConnectivity[i][j]];
1411 i, three, vertsArray);
1418 const int nFaces = 5;
1420 const int quadEdgeConnectivity[][4] = {{0, 1, 2, 3}};
1423 const int triEdgeConnectivity[][3] = {
1424 {0, 5, 4}, {1, 6, 5}, {2, 7, 6}, {3, 4, 7}};
1428 for (
int f = 0; f < nFaces; ++f)
1433 for (
int j = 0; j < 4; ++j)
1435 edgeArray[j] = edges[quadEdgeConnectivity[f][j]];
1446 for (
int j = 0; j < 3; ++j)
1448 edgeArray[j] = edges[triEdgeConnectivity[i][j]];
1473 const int three = 3;
1474 const int nVerts = 4;
1475 const double point[][3] = {{-1, -1 /
sqrt(
double(3)), -1 /
sqrt(
double(6))},
1476 {1, -1 /
sqrt(
double(3)), -1 /
sqrt(
double(6))},
1477 {0, 2 /
sqrt(
double(3)), -1 /
sqrt(
double(6))},
1478 {0, 0, 3 /
sqrt(
double(6))}};
1480 std::shared_ptr<SpatialDomains::PointGeom> verts[4];
1481 for (i = 0; i < nVerts; ++i)
1484 three, i, point[i][0], point[i][1], point[i][2]);
1492 const int nEdges = 6;
1493 const int vertexConnectivity[][2] = {{0, 1}, {1, 2}, {0, 2},
1494 {0, 3}, {1, 3}, {2, 3}};
1498 for (i = 0; i < nEdges; ++i)
1500 std::shared_ptr<SpatialDomains::PointGeom> vertsArray[2];
1501 for (j = 0; j < 2; ++j)
1503 vertsArray[j] = verts[vertexConnectivity[i][j]];
1507 i, three, vertsArray);
1514 const int nFaces = 4;
1515 const int edgeConnectivity[][3] = {
1516 {0, 1, 2}, {0, 4, 3}, {1, 5, 4}, {2, 5, 3}};
1520 for (i = 0; i < nFaces; ++i)
1523 for (j = 0; j < 3; ++j)
1525 edgeArray[j] = edges[edgeConnectivity[i][j]];
1548 const int three = 3;
1550 const int nVerts = 8;
1551 const double point[][3] = {{0, 0, 0}, {1, 0, 0}, {1, 1, 0}, {0, 1, 0},
1552 {0, 0, 1}, {1, 0, 1}, {1, 1, 1}, {0, 1, 1}};
1556 for (
int i = 0; i < nVerts; ++i)
1559 three, i, point[i][0], point[i][1], point[i][2]);
1567 const int nEdges = 12;
1568 const int vertexConnectivity[][2] = {{0, 1}, {1, 2}, {2, 3}, {0, 3},
1569 {0, 4}, {1, 5}, {2, 6}, {3, 7},
1570 {4, 5}, {5, 6}, {6, 7}, {4, 7}};
1574 for (
int i = 0; i < nEdges; ++i)
1577 for (
int j = 0; j < 2; ++j)
1579 vertsArray[j] = verts[vertexConnectivity[i][j]];
1582 i, three, vertsArray);
1589 const int nFaces = 6;
1590 const int edgeConnectivity[][4] = {{0, 1, 2, 3}, {0, 5, 8, 4},
1591 {1, 6, 9, 5}, {2, 7, 10, 6},
1592 {3, 7, 11, 4}, {8, 9, 10, 11}};
1596 for (
int i = 0; i < nFaces; ++i)
1599 for (
int j = 0; j < 4; ++j)
1601 edgeArray[j] = edges[edgeConnectivity[i][j]];
1620 std::map<ShapeType, DNekScalMatSharedPtr> &maxRmat,
1621 map<ShapeType, LocalRegions::Expansion3DSharedPtr> &maxElmt,
1625 std::shared_ptr<MultiRegions::ExpList> expList =
1626 ((
m_linsys.lock())->GetLocMat()).lock();
1632 map<ShapeType, Array<OneD, Array<OneD, unsigned int>>> faceMapMaxR;
1635 int nummodesmax = 0;
1638 for (
int n = 0; n < expList->GetNumElmts(); ++n)
1640 locExp = expList->GetExp(n);
1642 nummodesmax = max(nummodesmax, locExp->GetBasisNumModes(0));
1643 nummodesmax = max(nummodesmax, locExp->GetBasisNumModes(1));
1644 nummodesmax = max(nummodesmax, locExp->GetBasisNumModes(2));
1646 Shapes[locExp->DetShapeType()] = 1;
1693 HexBa, HexBb, HexBc, hexgeom);
1698 HexExp->GetInverseBoundaryMaps(vmap, emap, fmap);
1706 maxRmat[
eHexahedron] = HexExp->GetLocMatrix(HexR);
1716 PointsKey(nummodesmax, eGaussRadauMAlpha1Beta0));
1718 PointsKey(nummodesmax, eGaussRadauMAlpha2Beta0));
1724 TetBa, TetBb, TetBc, tetgeom);
1728 TetExp->GetInverseBoundaryMaps(vmap, emap, fmap);
1741 edgeMapMaxR, faceMapMaxR);
1755 PointsKey(nummodesmax, eGaussRadauMAlpha2Beta0));
1761 PyrBa, PyrBb, PyrBc, pyrgeom);
1766 PyrExp->GetInverseBoundaryMaps(vmap, emap, fmap);
1773 SetUpPyrMaxRMat(nummodesmax, PyrExp, maxRmat, vertMapMaxR, edgeMapMaxR,
1788 PointsKey(nummodesmax, eGaussRadauMAlpha1Beta0));
1794 PrismBa, PrismBb, PrismBc, prismgeom);
1795 maxElmt[
ePrism] = PrismExp;
1798 PrismExp->GetInverseBoundaryMaps(vmap, emap, fmap);
1799 vertMapMaxR[
ePrism] = vmap;
1800 edgeMapMaxR[
ePrism] = emap;
1802 faceMapMaxR[
ePrism] = fmap;
1807 edgeMapMaxR, faceMapMaxR,
false);
1814 maxRmat[
ePrism] = PrismExp->GetLocMatrix(PrismR);
1819 edgeMapMaxR, faceMapMaxR,
true);
1827 std::map<ShapeType, DNekScalMatSharedPtr> &maxRmat,
1832 int nRows = PyrExp->NumBndryCoeffs();
1839 for (
int i = 0; i < nRows; ++i)
1841 newmat->SetValue(i, i, 1.0);
1849 const int nadjedge[] = {3, 3, 3, 3, 4};
1850 const int VEHexVert[][4] = {{0, 0, 0, -1},
1855 const int VEHexEdge[][4] = {{0, 3, 4, -1},
1860 const int VEPyrEdge[][4] = {{0, 3, 4, -1},
1867 for (
int v = 0; v < 5; ++v)
1869 for (
int e = 0; e < nadjedge[v]; ++e)
1871 for (
int i = 0; i < nummodesmax - 2; ++i)
1878 newmat->SetValue(vertMapMaxR[
ePyramid][v],
1879 edgeMapMaxR[
ePyramid][VEPyrEdge[v][e]][i],
1886 nfacemodes = (nummodesmax - 2) * (nummodesmax - 2);
1888 for (
int v = 0; v < 4; ++v)
1890 for (
int i = 0; i < nfacemodes; ++i)
1894 newmat->SetValue(vertMapMaxR[
ePyramid][v],
1899 const int nadjface[] = {2, 2, 2, 2, 4};
1900 const int VFTetVert[][4] = {{0, 0, -1, -1},
1905 const int VFTetFace[][4] = {{1, 3, -1, -1},
1910 const int VFPyrFace[][4] = {{1, 4, -1, -1},
1917 nfacemodes = (nummodesmax - 3) * (nummodesmax - 2) / 2;
1918 for (
int v = 0; v < 5; ++v)
1920 for (
int f = 0; f < nadjface[v]; ++f)
1922 for (
int i = 0; i < nfacemodes; ++i)
1927 newmat->SetValue(vertMapMaxR[
ePyramid][v],
1928 faceMapMaxR[
ePyramid][VFPyrFace[v][f]][i],
1936 nfacemodes = (nummodesmax - 2) * (nummodesmax - 2);
1937 for (
int e = 0; e < 4; ++e)
1939 for (
int i = 0; i < nummodesmax - 2; ++i)
1941 for (
int j = 0; j < nfacemodes; ++j)
1946 val = (*maxRmat[
eHexahedron])(edgemapid, facemapid);
1947 newmat->SetValue(edgeMapMaxR[
ePyramid][e][i],
1953 const int nadjface1[] = {1, 1, 1, 1, 2, 2, 2, 2};
1954 const int EFTetEdge[][2] = {{0, -1}, {1, -1}, {0, -1}, {2, -1},
1955 {3, 3}, {4, 4}, {5, 5}, {3, 5}};
1956 const int EFTetFace[][2] = {{1, -1}, {2, -1}, {1, -1}, {3, -1},
1957 {1, 3}, {1, 2}, {2, 3}, {1, 3}};
1958 const int EFPyrFace[][2] = {{1, -1}, {2, -1}, {3, -1}, {4, -1},
1959 {1, 4}, {1, 2}, {2, 3}, {3, 4}};
1962 nfacemodes = (nummodesmax - 3) * (nummodesmax - 2) / 2;
1963 for (
int e = 0; e < 8; ++e)
1965 for (
int f = 0; f < nadjface1[e]; ++f)
1967 for (
int i = 0; i < nummodesmax - 2; ++i)
1969 for (
int j = 0; j < nfacemodes; ++j)
1977 newmat->SetValue(edgeMapMaxR[
ePyramid][e][i],
1978 faceMapMaxR[
ePyramid][EFPyrFace[e][f]][j],
1992 std::map<ShapeType, DNekScalMatSharedPtr> &maxRmat,
1997 boost::ignore_unused(faceMapMaxR);
1999 int nRows = TetExp->NumBndryCoeffs();
2006 for (
int i = 0; i < nRows; ++i)
2008 for (
int j = 0; j < nRows; ++j)
2011 newmat->SetValue(i, j, val);
2020 const int VEHexVert[][4] = {{0, 0, 0}, {1, 1, 1}, {2, 2, 2}, {4, 5, 6}};
2021 const int VEHexEdge[][4] = {{0, 3, 4}, {0, 1, 5}, {1, 2, 6}, {4, 5, 6}};
2022 const int VETetEdge[][4] = {{0, 2, 3}, {0, 1, 4}, {1, 2, 5}, {3, 4, 5}};
2025 for (
int v = 0; v < 4; ++v)
2027 for (
int e = 0; e < 3; ++e)
2029 for (
int i = 0; i < nummodesmax - 2; ++i)
2051 std::map<ShapeType, DNekScalMatSharedPtr> &maxRmat,
2057 int nRows = PrismExp->NumBndryCoeffs();
2068 for (
int i = 0; i < nRows; ++i)
2070 for (
int j = 0; j < nRows; ++j)
2072 val = (*maxRmat[
ePrism])(i, j);
2073 newmat->SetValue(i, j, val);
2078 const int VETetVert[][2] = {{0, 0}, {1, 1}, {1, 1},
2079 {0, 0}, {3, 3}, {3, 3}};
2080 const int VETetEdge[][2] = {{0, 3}, {0, 4}, {0, 4},
2081 {0, 3}, {3, 4}, {4, 3}};
2082 const int VEPrismEdge[][2] = {{0, 4}, {0, 5}, {2, 6},
2083 {2, 7}, {4, 5}, {6, 7}};
2086 for (
int v = 0; v < 6; ++v)
2088 for (
int e = 0; e < 2; ++e)
2090 for (
int i = 0; i < nummodesmax - 2; ++i)
2097 newmat->SetValue(vertMapMaxR[
ePrism][v],
2098 edgeMapMaxR[
ePrism][VEPrismEdge[v][e]][i],
2107 for (
int i = 0; i < nRows; ++i)
2109 newmat->SetValue(i, i, 1.0);
2119 const int VEHexVert[][3] = {{0, 0, 0}, {1, 1, 1}, {2, 2, 2},
2120 {3, 3, 3}, {4, 5, 5}, {6, 7, 7}};
2121 const int VEHexEdge[][3] = {{0, 3, 4}, {0, 1, 5}, {1, 2, 6},
2122 {2, 3, 7}, {4, 5, 9}, {6, 7, 11}};
2123 const int VEPrismEdge[][3] = {{0, 3, 4}, {0, 1, 5}, {1, 2, 6},
2124 {2, 3, 7}, {4, 5, 8}, {6, 7, 8}};
2127 for (
int v = 0; v < 6; ++v)
2129 for (
int e = 0; e < 3; ++e)
2131 for (
int i = 0; i < nummodesmax - 2; ++i)
2138 newmat->SetValue(vertMapMaxR[
ePrism][v],
2139 edgeMapMaxR[
ePrism][VEPrismEdge[v][e]][i],
2146 const int VFHexVert[][2] = {{0, 0}, {1, 1}, {4, 5},
2147 {2, 2}, {3, 3}, {6, 7}};
2148 const int VFHexFace[][2] = {{0, 4}, {0, 2}, {4, 2},
2149 {0, 2}, {0, 4}, {2, 4}};
2151 const int VQFPrismVert[][2] = {{0, 0}, {1, 1}, {4, 4},
2152 {2, 2}, {3, 3}, {5, 5}};
2153 const int VQFPrismFace[][2] = {{0, 4}, {0, 2}, {4, 2},
2154 {0, 2}, {0, 4}, {2, 4}};
2156 nfacemodes = (nummodesmax - 2) * (nummodesmax - 2);
2158 for (
int v = 0; v < 6; ++v)
2160 for (
int f = 0; f < 2; ++f)
2162 for (
int i = 0; i < nfacemodes; ++i)
2167 newmat->SetValue(vertMapMaxR[
ePrism][VQFPrismVert[v][f]],
2168 faceMapMaxR[
ePrism][VQFPrismFace[v][f]][i],
2175 const int nadjface[] = {1, 2, 1, 2, 1, 1, 1, 1, 2};
2176 const int EFHexEdge[][2] = {{0, -1}, {1, 1}, {2, -1}, {3, 3}, {4, -1},
2177 {5, -1}, {6, -1}, {7, -1}, {9, 11}};
2178 const int EFHexFace[][2] = {{0, -1}, {0, 2}, {0, -1}, {0, 4}, {4, -1},
2179 {2, -1}, {2, -1}, {4, -1}, {2, 4}};
2180 const int EQFPrismEdge[][2] = {{0, -1}, {1, 1}, {2, -1},
2181 {3, 3}, {4, -1}, {5, -1},
2182 {6, -1}, {7, -1}, {8, 8}};
2183 const int EQFPrismFace[][2] = {{0, -1}, {0, 2}, {0, -1},
2184 {0, 4}, {4, -1}, {2, -1},
2185 {2, -1}, {4, -1}, {2, 4}};
2188 nfacemodes = (nummodesmax - 2) * (nummodesmax - 2);
2189 for (
int e = 0; e < 9; ++e)
2191 for (
int f = 0; f < nadjface[e]; ++f)
2193 for (
int i = 0; i < nummodesmax - 2; ++i)
2195 for (
int j = 0; j < nfacemodes; ++j)
2202 val = (*maxRmat[
eHexahedron])(edgemapid, facemapid);
2205 edgeMapMaxR[
ePrism][EQFPrismEdge[e][f]][i];
2207 faceMapMaxR[
ePrism][EQFPrismFace[e][f]][j];
2208 newmat->SetValue(edgemapid1, facemapid1, val);
2215 const int VFTetVert[] = {0, 1, 3, 1, 0, 3};
2216 const int VFTetFace[] = {1, 1, 1, 1, 1, 1};
2217 const int VTFPrismVert[] = {0, 1, 4, 2, 3, 5};
2218 const int VTFPrismFace[] = {1, 1, 1, 3, 3, 3};
2221 nfacemodes = (nummodesmax - 3) * (nummodesmax - 2) / 2;
2222 for (
int v = 0; v < 6; ++v)
2224 for (
int i = 0; i < nfacemodes; ++i)
2230 newmat->SetValue(vertMapMaxR[
ePrism][VTFPrismVert[v]],
2231 faceMapMaxR[
ePrism][VTFPrismFace[v]][i], val);
2236 const int EFTetEdge[] = {0, 3, 4, 0, 4, 3};
2237 const int EFTetFace[] = {1, 1, 1, 1, 1, 1};
2238 const int ETFPrismEdge[] = {0, 4, 5, 2, 6, 7};
2239 const int ETFPrismFace[] = {1, 1, 1, 3, 3, 3};
2243 nfacemodes = (nummodesmax - 3) * (nummodesmax - 2) / 2;
2244 for (
int e = 0; e < 6; ++e)
2246 for (
int i = 0; i < nummodesmax - 2; ++i)
2248 for (
int j = 0; j < nfacemodes; ++j)
2250 int edgemapid = edgeMapMaxR[
eTetrahedron][EFTetEdge[e]][i];
2251 int facemapid = faceMapMaxR[
eTetrahedron][EFTetFace[e]][j];
2254 newmat->SetValue(edgeMapMaxR[
ePrism][ETFPrismEdge[e]][i],
2255 faceMapMaxR[
ePrism][ETFPrismFace[e]][j], val);
2262 maxRmat[
ePrism] = PrismR;
2273 int nRows = locExp->NumBndryCoeffs();
2281 locExp->GetInverseBoundaryMaps(vlocmap, elocmap, flocmap);
2284 for (
int i = 0; i < nRows; ++i)
2287 newmat->SetValue(i, i, val);
2290 int nverts = locExp->GetNverts();
2291 int nedges = locExp->GetNedges();
2292 int nfaces = locExp->GetNtraces();
2295 for (
int e = 0; e < nedges; ++e)
2297 int nEdgeInteriorCoeffs = locExp->GetEdgeNcoeffs(e) - 2;
2299 for (
int v = 0; v < nverts; ++v)
2301 for (
int i = 0; i < nEdgeInteriorCoeffs; ++i)
2303 val = (*maxRmat)(vmap[v], emap[e][i]);
2304 newmat->SetValue(vlocmap[v], elocmap[e][i], val);
2309 for (
int f = 0; f < nfaces; ++f)
2316 int nFaceInteriorCoeffs = locExp->GetTraceIntNcoeffs(f);
2318 locExp->GetTraceNumModes(f, m0, m1, FwdOrient);
2321 maxExp->GetTraceInverseBoundaryMap(f, FwdOrient, m0, m1);
2324 for (
int v = 0; v < nverts; ++v)
2326 for (
int i = 0; i < nFaceInteriorCoeffs; ++i)
2328 val = (*maxRmat)(vmap[v], fmapRmat[i]);
2329 newmat->SetValue(vlocmap[v], flocmap[f][i], val);
2334 for (
int e = 0; e < nedges; ++e)
2336 int nEdgeInteriorCoeffs = locExp->GetEdgeNcoeffs(e) - 2;
2338 for (
int j = 0; j < nEdgeInteriorCoeffs; ++j)
2340 for (
int i = 0; i < nFaceInteriorCoeffs; ++i)
2342 val = (*maxRmat)(emap[e][j], fmapRmat[i]);
2343 newmat->SetValue(elocmap[e][j], flocmap[f][i], val);
#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
Describes the specification for a Basis.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Defines a specification for a set of points.
SpatialDomains::GeometrySharedPtr GetGeom() const
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Describe a linear system.
const StdRegions::ConstFactorMap & GetConstFactors() const
Returns all the constants.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
LibUtilities::CommSharedPtr m_comm
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap
Array< OneD, NekDouble > m_variablePmask
virtual void v_DoTransformCoeffsFromLowEnergy(Array< OneD, NekDouble > &pInOut) override
transform the solution coeffiicents from low energy back to the original coefficient space.
SpatialDomains::HexGeomSharedPtr CreateRefHexGeom(void)
Sets up the reference hexahedral element needed to construct a low energy basis.
void ReSetPrismMaxRMat(int nummodesmax, LocalRegions::PrismExpSharedPtr &PirsmExp, ShapeToDNekMap &maxRmat, ShapeToIntArrayMap &vertMapMaxR, ShapeToIntArrayArrayMap &edgeMapMaxR, ShapeToIntArrayArrayMap &faceMapMaxR, bool UseTetOnly)
std::vector< std::pair< int, int > > m_sameBlock
DNekBlkMatSharedPtr m_InvRBlk
virtual void v_DoTransformCoeffsToLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
Multiply by the block tranposed inverse transformation matrix (R^T)^{-1} which is equivlaent to trans...
virtual DNekScalMatSharedPtr v_TransformedSchurCompl(int n, int offset, const std::shared_ptr< DNekScalMat > &loc_mat) override
Set up the transformed block matrix system.
static std::string className
Name of class.
void SetupBlockTransformationMatrix(void)
DNekBlkMatSharedPtr m_RBlk
DNekMatSharedPtr ExtractLocMat(LocalRegions::Expansion3DSharedPtr &locExp, DNekScalMatSharedPtr &maxRmat, LocalRegions::Expansion3DSharedPtr &expMax, Array< OneD, unsigned int > &vertMapMaxR, Array< OneD, Array< OneD, unsigned int > > &edgeMapMaxR)
virtual void v_BuildPreconditioner() override
Construct the low energy preconditioner from .
void ReSetTetMaxRMat(int nummodesmax, LocalRegions::TetExpSharedPtr &TetExp, ShapeToDNekMap &maxRmat, ShapeToIntArrayMap &vertMapMaxR, ShapeToIntArrayArrayMap &edgeMapMaxR, ShapeToIntArrayArrayMap &faceMapMaxR)
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
PreconditionerLowEnergy(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
SpatialDomains::TetGeomSharedPtr CreateRefTetGeom(void)
Sets up the reference tretrahedral element needed to construct a low energy basis.
virtual void v_InitObject() override
virtual void v_DoTransformBasisFromLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
Multiply by the block inverse transformation matrix This transforms the bassi from Low Energy to orig...
virtual void v_DoTransformBasisToLowEnergy(Array< OneD, NekDouble > &pInOut) override
Transform the basis vector to low energy space.
SpatialDomains::PrismGeomSharedPtr CreateRefPrismGeom(void)
Sets up the reference prismatic element needed to construct a low energy basis.
void SetUpReferenceElements(ShapeToDNekMap &maxRmat, ShapeToExpMap &maxElmt, ShapeToIntArrayMap &vertMapMaxR, ShapeToIntArrayArrayMap &edgeMapMaxR)
Loop expansion and determine different variants of the transformation matrix.
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false) override
void SetUpPyrMaxRMat(int nummodesmax, LocalRegions::PyrExpSharedPtr &PyrExp, ShapeToDNekMap &maxRmat, ShapeToIntArrayMap &vertMapMaxR, ShapeToIntArrayArrayMap &edgeMapMaxR, ShapeToIntArrayArrayMap &faceMapMaxR)
void CreateVariablePMask(void)
DNekBlkMatSharedPtr m_BlkMat
SpatialDomains::PyrGeomSharedPtr CreateRefPyrGeom(void)
Sets up the reference prismatic element needed to construct a low energy basis mapping arrays.
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
int GetNedges() const
return the number of edges in 3D expansion
int NumBndryCoeffs(void) const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
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.
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eModified_B
Principle Modified Functions .
@ eModified_C
Principle Modified Functions .
@ eModifiedPyr_C
Principle Modified Functions.
@ eModified_A
Principle Modified Functions .
std::shared_ptr< PrismExp > PrismExpSharedPtr
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< HexExp > HexExpSharedPtr
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
std::shared_ptr< PyrExp > PyrExpSharedPtr
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
std::shared_ptr< TetExp > TetExpSharedPtr
PreconFactory & GetPreconFactory()
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.
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
std::shared_ptr< PrismGeom > PrismGeomSharedPtr
std::shared_ptr< HexGeom > HexGeomSharedPtr
std::shared_ptr< SegGeom > SegGeomSharedPtr
std::shared_ptr< PyrGeom > PyrGeomSharedPtr
std::shared_ptr< TetGeom > TetGeomSharedPtr
std::shared_ptr< PointGeom > PointGeomSharedPtr
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
std::shared_ptr< TriGeom > TriGeomSharedPtr
@ eDir1FwdDir1_Dir2FwdDir2
std::vector< double > z(NPUPPER)
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
void Assmb(int n, const T *x, const int *y, T *z)
Assemble z[y[i]] += x[i]; z should be zero'd first.
scalarT< T > sqrt(scalarT< T > in)