52 namespace MultiRegions
64 DisContField::DisContField()
65 :
ExpList(), m_bndConditions(), m_bndCondExpansions(),
84 const std::string &variable,
const bool SetUpJustDG,
85 const bool DeclareCoeffPhysArrays,
87 const std::string bcvariable)
88 :
ExpList(pSession, graph, DeclareCoeffPhysArrays, variable, ImpType),
92 if (bcvariable ==
"NotSet")
101 if (bcvar.compare(
"DefaultVar") != 0)
106 DeclareCoeffPhysArrays);
107 if (DeclareCoeffPhysArrays)
124 for (
int e = 0; e < locExpList->GetExpSize(); ++e)
126 (*m_exp)[ElmtID[cnt + e]]->SetTraceExp(TraceID[cnt + e],
127 locExpList->GetExp(e));
128 locExpList->GetExp(e)->SetAdjacentElementExp(
129 TraceID[cnt + e], (*
m_exp)[ElmtID[cnt + e]]);
140 if (
m_session->DefinesSolverInfo(
"PROJECTION"))
142 std::string ProjectStr =
m_session->GetSolverInfo(
"PROJECTION");
143 if ((ProjectStr ==
"MixedCGDG") ||
144 (ProjectStr ==
"Mixed_CG_Discontinuous"))
178 m_comm,
true,
"DefaultVar", ImpType);
188 if (
m_session->DefinesCmdLineArgument(
"verbose"))
201 for (
int i = 0; i <
m_exp->size(); ++i)
203 for (
int j = 0; j < (*m_exp)[i]->GetNtraces(); ++j)
208 exp->SetAdjacentElementExp(j, (*
m_exp)[i]);
209 (*m_exp)[i]->SetTraceExp(j, exp);
231 m_traceMap->GetBndCondIDToGlobalTraceID(cnt + v));
238 std::unordered_map<int, pair<int, int>> perTraceToExpMap;
239 for (cnt = n = 0; n <
m_exp->size(); ++n)
241 for (
int v = 0; v < (*m_exp)[n]->GetNtraces(); ++v, ++cnt)
243 auto it = periodicTraces.find((*
m_exp)[n]->GetGeom()->GetTid(v));
245 if (it != periodicTraces.end())
247 perTraceToExpMap[it->first] = make_pair(n, v);
256 for (cnt = n = 0; n <
m_exp->size(); ++n)
258 for (
int v = 0; v < (*m_exp)[n]->GetNtraces(); ++v, ++cnt)
265 for (cnt = n = 0; n <
m_exp->size(); ++n)
267 for (
int v = 0; v < (*m_exp)[n]->GetNtraces(); ++v, ++cnt)
269 int GeomId = (*m_exp)[n]->GetGeom()->GetTid(v);
272 auto it = periodicTraces.find(GeomId);
274 if (it != periodicTraces.end())
277 auto it2 = perTraceToExpMap.find(ent.
id);
279 if (it2 == perTraceToExpMap.end())
281 if (
m_session->GetComm()->GetRowComm()->GetSize() > 1 &&
288 ASSERTL1(
false,
"Periodic trace not found!");
293 "Periodic trace in non-forward space?");
295 int offset =
m_trace->GetPhys_Offset(
296 (
m_traceMap->GetElmtToTrace())[n][v]->GetElmtId());
298 int offset2 =
m_trace->GetPhys_Offset(
299 (
m_traceMap->GetElmtToTrace())[it2->second.first]
315 int nquad = elmtToTrace[n][v]->GetNumPoints(0);
317 vector<int> tmpBwd(nquad);
318 vector<int> tmpFwd(nquad);
322 for (
int i = 0; i < nquad; ++i)
324 tmpBwd[i] = offset2 + i;
325 tmpFwd[i] = offset + i;
330 for (
int i = 0; i < nquad; ++i)
332 tmpBwd[i] = offset2 + i;
333 tmpFwd[i] = offset + nquad - i - 1;
337 for (
int i = 0; i < nquad; ++i)
348 int nquad1 = elmtToTrace[n][v]->GetNumPoints(0);
349 int nquad2 = elmtToTrace[n][v]->GetNumPoints(1);
351 vector<int> tmpBwd(nquad1 * nquad2);
352 vector<int> tmpFwd(nquad1 * nquad2);
362 for (
int i = 0; i < nquad2; ++i)
364 for (
int j = 0; j < nquad1; ++j)
366 tmpBwd[i * nquad1 + j] =
367 offset2 + i * nquad1 + j;
368 tmpFwd[i * nquad1 + j] =
369 offset + j * nquad2 + i;
375 for (
int i = 0; i < nquad2; ++i)
377 for (
int j = 0; j < nquad1; ++j)
379 tmpBwd[i * nquad1 + j] =
380 offset2 + i * nquad1 + j;
381 tmpFwd[i * nquad1 + j] =
382 offset + i * nquad1 + j;
396 for (
int i = 0; i < nquad2; ++i)
398 for (
int j = 0; j < nquad1 / 2; ++j)
400 swap(tmpFwd[i * nquad1 + j],
401 tmpFwd[i * nquad1 + nquad1 - j - 1]);
415 for (
int j = 0; j < nquad1; ++j)
417 for (
int i = 0; i < nquad2 / 2; ++i)
419 swap(tmpFwd[i * nquad1 + j],
420 tmpFwd[(nquad2 - i - 1) * nquad1 + j]);
425 for (
int i = 0; i < nquad1 * nquad2; ++i)
488 if (traceEl->GetLeftAdjacentElementTrace() == -1 ||
489 traceEl->GetRightAdjacentElementTrace() == -1)
502 else if (traceEl->GetLeftAdjacentElementTrace() != -1 &&
503 traceEl->GetRightAdjacentElementTrace() != -1)
506 fwd = (traceEl->GetLeftAdjacentElementExp().get()) == (*
m_exp)[n].get();
510 ASSERTL2(
false,
"Unconnected trace element!");
522 const std::string &variable)
524 boost::ignore_unused(variable);
530 map<int, int> GeometryToRegionsMap;
538 for (
auto &it : bregions)
540 for (
auto &bregionIt : *it.second)
544 int id = bregionIt.second->m_geomVec[0]->GetGlobalID();
545 GeometryToRegionsMap[id] = it.first;
549 map<int, SpatialDomains::GeometrySharedPtr> EndOfDomain;
552 for (
auto &domIt : domain)
555 for (
int i = 0; i < geomvector->m_geomVec.size(); ++i)
557 for (
int j = 0; j < 2; ++j)
559 int vid = geomvector->m_geomVec[i]->GetVid(j);
560 if (EndOfDomain.count(vid) == 0)
562 EndOfDomain[vid] = geomvector->m_geomVec[i]->GetVertex(j);
566 EndOfDomain.erase(vid);
571 ASSERTL1(EndOfDomain.size() == 2,
"Did not find two ends of domain");
573 for (
auto ®It : EndOfDomain)
575 if (GeometryToRegionsMap.count(regIt.first) != 0)
578 auto iter = GeometryToRegionsMap.find(regIt.first);
579 ASSERTL1(iter != GeometryToRegionsMap.end(),
580 "Failied to find GeometryToRegionMap");
582 int regionId = iter->second;
583 auto bregionsIter = bregions.find(regionId);
584 ASSERTL1(bregionsIter != bregions.end(),
585 "Failed to find boundary region");
588 returnval->AddBoundaryRegions(regionId, breg);
590 auto bconditionsIter = bconditions.find(regionId);
591 ASSERTL1(bconditionsIter != bconditions.end(),
592 "Failed to find boundary collection");
595 bconditionsIter->second;
596 returnval->AddBoundaryConditions(regionId, bcond);
614 const std::string &variable,
616 bool SetToOneSpaceDimension,
618 :
ExpList(pSession, domain, graph1D, true, variable, SetToOneSpaceDimension,
621 if (variable.compare(
"DefaultVar") != 0)
640 const bool DeclareCoeffPhysArrays)
641 :
ExpList(In, DeclareCoeffPhysArrays), m_bndConditions(In.m_bndConditions),
642 m_bndCondExpansions(In.m_bndCondExpansions),
643 m_globalBndMat(In.m_globalBndMat), m_traceMap(In.m_traceMap),
644 m_boundaryTraces(In.m_boundaryTraces),
645 m_periodicVerts(In.m_periodicVerts),
646 m_periodicFwdCopy(In.m_periodicFwdCopy),
647 m_periodicBwdCopy(In.m_periodicBwdCopy),
648 m_leftAdjacentTraces(In.m_leftAdjacentTraces),
649 m_locTraceToTraceMap(In.m_locTraceToTraceMap)
654 *In.
m_trace, DeclareCoeffPhysArrays);
664 const std::string &variable,
const bool SetUpJustDG,
665 const bool DeclareCoeffPhysArrays)
666 :
ExpList(In, DeclareCoeffPhysArrays)
673 if (variable.compare(
"DefaultVar") != 0)
678 if (DeclareCoeffPhysArrays)
706 for (e = 0; e < locExpList->GetExpSize(); ++e)
708 (*m_exp)[ElmtID[cnt + e]]->SetTraceExp(
709 TraceID[cnt + e], locExpList->GetExp(e));
710 locExpList->GetExp(e)->SetAdjacentElementExp(
711 TraceID[cnt + e], (*
m_exp)[ElmtID[cnt + e]]);
717 if (
m_session->DefinesSolverInfo(
"PROJECTION"))
719 std::string ProjectStr =
722 if ((ProjectStr ==
"MixedCGDG") ||
723 (ProjectStr ==
"Mixed_CG_Discontinuous"))
767 for (e = 0; e < locExpList->GetExpSize(); ++e)
769 (*m_exp)[ElmtID[cnt + e]]->SetTraceExp(
770 TraceID[cnt + e], locExpList->GetExp(e));
771 locExpList->GetExp(e)->SetAdjacentElementExp(
772 TraceID[cnt + e], (*
m_exp)[ElmtID[cnt + e]]);
818 const bool DeclareCoeffPhysArrays)
836 for (
auto &it : bregions)
841 m_session, *(it.second), graph, DeclareCoeffPhysArrays, variable,
842 false, bc->GetComm());
853 boost::iequals(type,
"I") || boost::iequals(type,
"CalcBC"))
882 int i, region1ID, region2ID;
885 map<int, int> BregionToVertMap;
889 for (
auto &it : bregions)
894 if (locBCond->GetBoundaryConditionType() !=
900 it.second->begin()->second->m_geomVec[0]->GetGlobalID();
902 BregionToVertMap[it.first] = id;
907 int n = vComm->GetSize();
908 int p = vComm->GetRank();
911 nregions[
p] = BregionToVertMap.size();
918 for (i = 1; i < n; ++i)
920 regOffset[i] = regOffset[i - 1] + nregions[i - 1];
927 for (
auto &iit : BregionToVertMap)
929 bregid[i] = iit.first;
930 bregmap[i++] = iit.second;
931 islocal.insert(iit.first);
937 for (
int i = 0; i < totRegions; ++i)
939 BregionToVertMap[bregid[i]] = bregmap[i];
943 for (
auto &it : bregions)
948 if (locBCond->GetBoundaryConditionType() !=
955 region1ID = it.first;
957 std::static_pointer_cast<
959 ->m_connectedBoundaryRegion;
961 ASSERTL0(BregionToVertMap.count(region1ID) != 0,
962 "Cannot determine vertex of region1ID");
964 ASSERTL0(BregionToVertMap.count(region2ID) != 0,
965 "Cannot determine vertex of region2ID");
969 islocal.count(region2ID) != 0);
977 int region1ID, region2ID, i, j, k, cnt;
981 m_graph->GetCompositeOrdering();
983 m_graph->GetBndRegionOrdering();
990 map<int, int> perComps;
991 map<int, vector<int>> allVerts;
993 map<int, pair<int, StdRegions::Orientation>> allEdges;
996 for (i = 0; i < (*m_exp).size(); ++i)
998 for (j = 0; j < (*m_exp)[i]->GetNverts(); ++j)
1000 int id = (*m_exp)[i]->GetGeom()->GetVid(j);
1001 locVerts.insert(
id);
1006 for (
auto &it : bregions)
1011 if (locBCond->GetBoundaryConditionType() !=
1018 region1ID = it.first;
1020 std::static_pointer_cast<
1022 ->m_connectedBoundaryRegion;
1027 if (vComm->GetSize() == 1)
1029 cId1 = it.second->begin()->first;
1030 cId2 = bregions.find(region2ID)->second->begin()->first;
1034 cId1 = bndRegOrder.find(region1ID)->second[0];
1035 cId2 = bndRegOrder.find(region2ID)->second[0];
1039 "Boundary region " +
1040 boost::lexical_cast<string>(region1ID) +
1041 " should only contain 1 composite.");
1043 vector<unsigned int> tmpOrder;
1048 it.second->begin()->second;
1050 for (i = 0; i < c->m_geomVec.size(); ++i)
1053 std::dynamic_pointer_cast<SpatialDomains::SegGeom>(
1055 ASSERTL0(segGeom,
"Unable to cast to shared ptr");
1058 m_graph->GetElementsFromEdge(segGeom);
1060 "The periodic boundaries belong to "
1061 "more than one element of the mesh");
1064 std::dynamic_pointer_cast<SpatialDomains::Geometry2D>(
1067 allEdges[c->m_geomVec[i]->GetGlobalID()] =
1068 make_pair(elmt->at(0).second,
1069 geom->GetEorient(elmt->at(0).second));
1073 if (vComm->GetSize() == 1)
1075 tmpOrder.push_back(c->m_geomVec[i]->GetGlobalID());
1078 vector<int> vertList(2);
1079 vertList[0] = segGeom->GetVid(0);
1080 vertList[1] = segGeom->GetVid(1);
1081 allVerts[c->m_geomVec[i]->GetGlobalID()] = vertList;
1084 if (vComm->GetSize() == 1)
1086 compOrder[it.second->begin()->first] = tmpOrder;
1091 if (perComps.count(cId1) == 0)
1093 if (perComps.count(cId2) == 0)
1095 perComps[cId1] = cId2;
1099 std::stringstream ss;
1100 ss <<
"Boundary region " << cId2 <<
" should be "
1101 <<
"periodic with " << perComps[cId2] <<
" but "
1102 <<
"found " << cId1 <<
" instead!";
1103 ASSERTL0(perComps[cId2] == cId1, ss.str());
1108 std::stringstream ss;
1109 ss <<
"Boundary region " << cId1 <<
" should be "
1110 <<
"periodic with " << perComps[cId1] <<
" but "
1111 <<
"found " << cId2 <<
" instead!";
1112 ASSERTL0(perComps[cId1] == cId1, ss.str());
1119 for (
auto &cIt : perComps)
1121 idmax = max(idmax, cIt.first);
1122 idmax = max(idmax, cIt.second);
1127 for (
auto &cIt : perComps)
1129 perid[cIt.first] = cIt.second;
1133 for (
int i = 0; i < idmax; ++i)
1139 if (perComps.count(perid[i]))
1143 perComps[i] = perid[i];
1149 int n = vComm->GetSize();
1150 int p = vComm->GetRank();
1156 edgecounts[
p] = allEdges.size();
1160 for (i = 1; i < n; ++i)
1162 edgeoffset[i] = edgeoffset[i - 1] + edgecounts[i - 1];
1171 auto sIt = allEdges.begin();
1173 for (i = 0; sIt != allEdges.end(); ++sIt)
1175 edgeIds[edgeoffset[
p] + i] = sIt->first;
1176 edgeIdx[edgeoffset[
p] + i] = sIt->second.first;
1177 edgeOrient[edgeoffset[
p] + i] = sIt->second.second;
1178 edgeVerts[edgeoffset[
p] + i++] = allVerts[sIt->first].size();
1201 for (i = 0; i < n; ++i)
1203 if (edgecounts[i] > 0)
1206 edgeVerts + edgeoffset[i], 1);
1215 for (i = 1; i < n; ++i)
1217 vertoffset[i] = vertoffset[i - 1] + procVerts[i - 1];
1221 for (i = 0, sIt = allEdges.begin(); sIt != allEdges.end(); ++sIt)
1223 for (j = 0; j < allVerts[sIt->first].size(); ++j)
1225 traceIds[vertoffset[
p] + i++] = allVerts[sIt->first][j];
1232 map<int, StdRegions::Orientation> orientMap;
1233 map<int, vector<int>> vertMap;
1235 for (cnt = i = 0; i < totEdges; ++i)
1237 ASSERTL0(orientMap.count(edgeIds[i]) == 0,
1238 "Already found edge in orientation map!");
1258 orientMap[edgeIds[i]] = o;
1260 vector<int> verts(edgeVerts[i]);
1262 for (j = 0; j < edgeVerts[i]; ++j)
1264 verts[j] = traceIds[cnt++];
1266 vertMap[edgeIds[i]] = verts;
1273 map<int, int> allCompPairs;
1282 for (
auto &cIt : perComps)
1285 const int id1 = cIt.first;
1286 const int id2 = cIt.second;
1287 std::string id1s = boost::lexical_cast<string>(id1);
1288 std::string id2s = boost::lexical_cast<string>(id2);
1290 if (compMap.count(id1) > 0)
1292 c[0] = compMap[id1];
1295 if (compMap.count(id2) > 0)
1297 c[1] = compMap[id2];
1303 map<int, int> compPairs;
1306 "Unable to find composite " + id1s +
" in order map.");
1308 "Unable to find composite " + id2s +
" in order map.");
1309 ASSERTL0(compOrder[id1].size() == compOrder[id2].size(),
1310 "Periodic composites " + id1s +
" and " + id2s +
1311 " should have the same number of elements.");
1312 ASSERTL0(compOrder[id1].size() > 0,
"Periodic composites " +
1313 id1s +
" and " + id2s +
1317 for (i = 0; i < compOrder[id1].size(); ++i)
1319 int eId1 = compOrder[id1][i];
1320 int eId2 = compOrder[id2][i];
1322 ASSERTL0(compPairs.count(eId1) == 0,
"Already paired.");
1324 if (compPairs.count(eId2) != 0)
1326 ASSERTL0(compPairs[eId2] == eId1,
"Pairing incorrect");
1328 compPairs[eId1] = eId2;
1334 for (i = 0; i < 2; ++i)
1341 if (c[i]->m_geomVec.size() > 0)
1343 for (j = 0; j < c[i]->m_geomVec.size(); ++j)
1345 locEdges.insert(c[i]->m_geomVec[j]->GetGlobalID());
1351 for (
auto &pIt : compPairs)
1353 int ids[2] = {pIt.first, pIt.second};
1354 bool local[2] = {locEdges.count(pIt.first) > 0,
1355 locEdges.count(pIt.second) > 0};
1357 ASSERTL0(orientMap.count(ids[0]) > 0 &&
1358 orientMap.count(ids[1]) > 0,
1359 "Unable to find edge in orientation map.");
1361 allCompPairs[pIt.first] = pIt.second;
1362 allCompPairs[pIt.second] = pIt.first;
1364 for (i = 0; i < 2; ++i)
1371 int other = (i + 1) % 2;
1374 orientMap[ids[i]] == orientMap[ids[other]]
1382 for (i = 0; i < 2; ++i)
1384 int other = (i + 1) % 2;
1387 orientMap[ids[i]] == orientMap[ids[other]]
1392 vector<int> perVerts1 = vertMap[ids[i]];
1393 vector<int> perVerts2 = vertMap[ids[other]];
1395 map<int, pair<int, bool>> tmpMap;
1398 tmpMap[perVerts1[0]] = make_pair(
1399 perVerts2[0], locVerts.count(perVerts2[0]) > 0);
1400 tmpMap[perVerts1[1]] = make_pair(
1401 perVerts2[1], locVerts.count(perVerts2[1]) > 0);
1405 tmpMap[perVerts1[0]] = make_pair(
1406 perVerts2[1], locVerts.count(perVerts2[1]) > 0);
1407 tmpMap[perVerts1[1]] = make_pair(
1408 perVerts2[0], locVerts.count(perVerts2[0]) > 0);
1411 for (
auto &mIt : tmpMap)
1417 auto perIt = periodicVerts.find(mIt.first);
1419 if (perIt == periodicVerts.end())
1421 periodicVerts[mIt.first].push_back(ent2);
1422 perIt = periodicVerts.find(mIt.first);
1427 for (j = 0; j < perIt->second.size(); ++j)
1429 if (perIt->second[j].id == mIt.second.first)
1438 perIt->second.push_back(ent2);
1449 for (cnt = i = 0; i < totEdges; ++i)
1451 int edgeId = edgeIds[i];
1453 ASSERTL0(allCompPairs.count(edgeId) > 0,
1454 "Unable to find matching periodic edge.");
1456 int perEdgeId = allCompPairs[edgeId];
1458 for (j = 0; j < edgeVerts[i]; ++j, ++cnt)
1460 int vId = traceIds[cnt];
1462 auto perId = periodicVerts.find(vId);
1464 if (perId == periodicVerts.end())
1473 orientMap[edgeId] == orientMap[perEdgeId]
1474 ? vertMap[perEdgeId][(j + 1) % 2]
1475 : vertMap[perEdgeId][j];
1479 locVerts.count(perVertexId) > 0);
1481 periodicVerts[vId].push_back(ent);
1488 for (
auto &perIt : periodicVerts)
1491 for (i = 0; i < perIt.second.size(); ++i)
1493 auto perIt2 = periodicVerts.find(perIt.second[i].id);
1494 ASSERTL0(perIt2 != periodicVerts.end(),
1495 "Couldn't find periodic vertex.");
1497 for (j = 0; j < perIt2->second.size(); ++j)
1499 if (perIt2->second[j].id == perIt.first)
1505 for (k = 0; k < perIt.second.size(); ++k)
1507 if (perIt2->second[j].id == perIt.second[k].id)
1516 perIt.second.push_back(perIt2->second[j]);
1524 for (
auto &perIt : periodicVerts)
1526 if (locVerts.count(perIt.first) > 0)
1536 m_graph->GetCompositeOrdering();
1538 m_graph->GetBndRegionOrdering();
1553 map<int, RotPeriodicInfo> rotComp;
1554 map<int, int> perComps;
1555 map<int, vector<int>> allVerts;
1556 map<int, SpatialDomains::PointGeomVector> allCoord;
1557 map<int, vector<int>> allEdges;
1558 map<int, vector<StdRegions::Orientation>> allOrient;
1563 int region1ID, region2ID, i, j, k, cnt;
1567 for (i = 0; i < (*m_exp).size(); ++i)
1569 for (j = 0; j < (*m_exp)[i]->GetNverts(); ++j)
1571 int id = (*m_exp)[i]->GetGeom()->GetVid(j);
1572 locVerts.insert(
id);
1575 for (j = 0; j < (*m_exp)[i]->GetGeom()->GetNumEdges(); ++j)
1577 int id = (*m_exp)[i]->GetGeom()->GetEid(j);
1578 locEdges.insert(
id);
1585 for (
auto &it : bregions)
1591 if (locBCond->GetBoundaryConditionType() !=
1598 region1ID = it.first;
1600 std::static_pointer_cast<
1602 ->m_connectedBoundaryRegion;
1606 "Boundary region " +
1607 boost::lexical_cast<string>(region1ID) +
1608 " should only contain 1 composite.");
1616 if (vComm->GetSize() == 1)
1618 cId1 = it.second->begin()->first;
1619 cId2 = bregions.find(region2ID)->second->begin()->first;
1623 cId1 = bndRegOrder.find(region1ID)->second[0];
1624 cId2 = bndRegOrder.find(region2ID)->second[0];
1628 if (boost::icontains(locBCond->GetUserDefined(),
"Rotated"))
1630 vector<string> tmpstr;
1632 boost::split(tmpstr, locBCond->GetUserDefined(),
1633 boost::is_any_of(
":"));
1635 if (boost::iequals(tmpstr[0],
"Rotated"))
1638 "Expected Rotated user defined string to "
1639 "contain direction and rotation angle "
1640 "and optionally a tolerance, "
1641 "i.e. Rotated:dir:PI/2:1e-6");
1643 ASSERTL1((tmpstr[1] ==
"x") || (tmpstr[1] ==
"y") ||
1646 "not specified as x,y or z");
1649 RotInfo.
m_dir = (tmpstr[1] ==
"x") ? 0
1650 : (tmpstr[1] ==
"y") ? 1
1657 if (tmpstr.size() == 4)
1662 boost::lexical_cast<NekDouble>(tmpstr[3]);
1667 "failed to cast tolerance input "
1668 "to a double value in Rotated"
1669 "boundary information");
1674 RotInfo.
m_tol = 1e-8;
1676 rotComp[cId1] = RotInfo;
1681 it.second->begin()->second;
1683 vector<unsigned int> tmpOrder;
1691 for (i = 0; i < c->m_geomVec.size(); ++i)
1694 std::dynamic_pointer_cast<SpatialDomains::Geometry2D>(
1696 ASSERTL1(faceGeom,
"Unable to cast to shared ptr");
1699 int faceId = c->m_geomVec[i]->GetGlobalID();
1700 locFaces.insert(faceId);
1704 if (vComm->GetSize() == 1)
1706 tmpOrder.push_back(c->m_geomVec[i]->GetGlobalID());
1711 vector<int> vertList, edgeList;
1713 vector<StdRegions::Orientation> orientVec;
1714 for (j = 0; j < faceGeom->GetNumVerts(); ++j)
1716 vertList.push_back(faceGeom->GetVid(j));
1717 edgeList.push_back(faceGeom->GetEid(j));
1718 coordVec.push_back(faceGeom->GetVertex(j));
1719 orientVec.push_back(faceGeom->GetEorient(j));
1722 allVerts[faceId] = vertList;
1723 allEdges[faceId] = edgeList;
1724 allCoord[faceId] = coordVec;
1725 allOrient[faceId] = orientVec;
1730 if (vComm->GetSize() == 1)
1732 compOrder[it.second->begin()->first] = tmpOrder;
1738 if (perComps.count(cId1) == 0)
1740 if (perComps.count(cId2) == 0)
1742 perComps[cId1] = cId2;
1746 std::stringstream ss;
1747 ss <<
"Boundary region " << cId2 <<
" should be "
1748 <<
"periodic with " << perComps[cId2] <<
" but "
1749 <<
"found " << cId1 <<
" instead!";
1750 ASSERTL0(perComps[cId2] == cId1, ss.str());
1755 std::stringstream ss;
1756 ss <<
"Boundary region " << cId1 <<
" should be "
1757 <<
"periodic with " << perComps[cId1] <<
" but "
1758 <<
"found " << cId2 <<
" instead!";
1759 ASSERTL0(perComps[cId1] == cId1, ss.str());
1766 for (
auto &cIt : perComps)
1768 idmax = max(idmax, cIt.first);
1769 idmax = max(idmax, cIt.second);
1774 for (
auto &cIt : perComps)
1776 perid[cIt.first] = cIt.second;
1780 for (
int i = 0; i < idmax; ++i)
1786 if (perComps.count(perid[i]))
1790 perComps[i] = perid[i];
1797 int n = vComm->GetSize();
1798 int p = vComm->GetRank();
1808 rotcounts[
p] = rotComp.size();
1814 for (i = 1; i < n; ++i)
1816 rotoffset[i] = rotoffset[i - 1] + rotcounts[i - 1];
1825 auto rIt = rotComp.begin();
1827 for (i = 0; rIt != rotComp.end(); ++rIt)
1829 compid[rotoffset[
p] + i] = rIt->first;
1830 rotdir[rotoffset[
p] + i] = rIt->second.m_dir;
1831 rotangle[rotoffset[
p] + i] = rIt->second.m_angle;
1832 rottol[rotoffset[
p] + i++] = rIt->second.m_tol;
1841 for (i = 0; i < totrot; ++i)
1845 rotComp[compid[i]] = rinfo;
1850 facecounts[
p] = locFaces.size();
1856 for (i = 1; i < n; ++i)
1858 faceoffset[i] = faceoffset[i - 1] + facecounts[i - 1];
1872 auto sIt = locFaces.begin();
1873 for (i = 0; sIt != locFaces.end(); ++sIt)
1875 faceIds[faceoffset[
p] + i] = *sIt;
1876 faceVerts[faceoffset[
p] + i++] = allVerts[*sIt].size();
1899 for (i = 0; i < n; ++i)
1901 if (facecounts[i] > 0)
1904 faceVerts + faceoffset[i], 1);
1915 for (i = 1; i < n; ++i)
1917 vertoffset[i] = vertoffset[i - 1] + procVerts[i - 1];
1931 for (cnt = 0, sIt = locFaces.begin(); sIt != locFaces.end(); ++sIt)
1933 for (j = 0; j < allVerts[*sIt].size(); ++j)
1935 int vertId = allVerts[*sIt][j];
1936 vertIds[vertoffset[
p] + cnt] = vertId;
1937 vertX[vertoffset[
p] + cnt] = (*allCoord[*sIt][j])(0);
1938 vertY[vertoffset[
p] + cnt] = (*allCoord[*sIt][j])(1);
1939 vertZ[vertoffset[
p] + cnt] = (*allCoord[*sIt][j])(2);
1940 edgeIds[vertoffset[
p] + cnt] = allEdges[*sIt][j];
1941 edgeOrt[vertoffset[
p] + cnt++] = allOrient[*sIt][j];
1956 map<int, vector<int>> vertMap;
1957 map<int, vector<int>> edgeMap;
1958 map<int, SpatialDomains::PointGeomVector> coordMap;
1964 map<int, SpatialDomains::PointGeomSharedPtr> vCoMap;
1965 map<int, pair<int, int>> eIdMap;
1967 for (cnt = i = 0; i < totFaces; ++i)
1969 vector<int> edges(faceVerts[i]);
1970 vector<int> verts(faceVerts[i]);
1976 for (j = 0; j < faceVerts[i]; ++j, ++cnt)
1978 edges[j] = edgeIds[cnt];
1979 verts[j] = vertIds[cnt];
1983 vCoMap[vertIds[cnt]] = coord[j];
1986 auto testIns = eIdMap.insert(make_pair(
1988 make_pair(vertIds[tmp + j],
1989 vertIds[tmp + ((j + 1) % faceVerts[i])])));
1991 if (testIns.second ==
false)
2014 swap(testIns.first->second.first,
2015 testIns.first->second.second);
2019 vertMap[faceIds[i]] = verts;
2020 edgeMap[faceIds[i]] = edges;
2021 coordMap[faceIds[i]] = coord;
2039 map<int, map<StdRegions::Orientation, vector<int>>> vmap;
2040 map<int, map<StdRegions::Orientation, vector<int>>> emap;
2042 map<StdRegions::Orientation, vector<int>> quadVertMap;
2052 map<StdRegions::Orientation, vector<int>> quadEdgeMap;
2062 map<StdRegions::Orientation, vector<int>> triVertMap;
2066 map<StdRegions::Orientation, vector<int>> triEdgeMap;
2070 vmap[3] = triVertMap;
2071 vmap[4] = quadVertMap;
2072 emap[3] = triEdgeMap;
2073 emap[4] = quadEdgeMap;
2075 map<int, int> allCompPairs;
2079 map<int, int> fIdToCompId;
2084 for (
auto &cIt : perComps)
2087 const int id1 = cIt.first;
2088 const int id2 = cIt.second;
2089 std::string id1s = boost::lexical_cast<string>(id1);
2090 std::string id2s = boost::lexical_cast<string>(id2);
2092 if (compMap.count(id1) > 0)
2094 c[0] = compMap[id1];
2097 if (compMap.count(id2) > 0)
2099 c[1] = compMap[id2];
2105 map<int, int> compPairs;
2108 "Unable to find composite " + id1s +
" in order map.");
2110 "Unable to find composite " + id2s +
" in order map.");
2111 ASSERTL0(compOrder[id1].size() == compOrder[id2].size(),
2112 "Periodic composites " + id1s +
" and " + id2s +
2113 " should have the same number of elements.");
2114 ASSERTL0(compOrder[id1].size() > 0,
"Periodic composites " +
2115 id1s +
" and " + id2s +
2119 for (i = 0; i < compOrder[id1].size(); ++i)
2121 int eId1 = compOrder[id1][i];
2122 int eId2 = compOrder[id2][i];
2124 ASSERTL0(compPairs.count(eId1) == 0,
"Already paired.");
2127 if (compPairs.count(eId2) != 0)
2129 ASSERTL0(compPairs[eId2] == eId1,
"Pairing incorrect");
2131 compPairs[eId1] = eId2;
2134 fIdToCompId[eId1] = id1;
2135 fIdToCompId[eId2] = id2;
2141 for (
auto &pIt : compPairs)
2143 int ids[2] = {pIt.first, pIt.second};
2144 bool local[2] = {locFaces.count(pIt.first) > 0,
2145 locFaces.count(pIt.second) > 0};
2147 ASSERTL0(coordMap.count(ids[0]) > 0 &&
2148 coordMap.count(ids[1]) > 0,
2149 "Unable to find face in coordinate map");
2151 allCompPairs[pIt.first] = pIt.second;
2152 allCompPairs[pIt.second] = pIt.first;
2157 coordMap[ids[0]], coordMap[ids[1]]};
2159 ASSERTL0(tmpVec[0].size() == tmpVec[1].size(),
2160 "Two periodic faces have different number "
2169 bool rotbnd =
false;
2176 if (rotComp.count(fIdToCompId[pIt.first]))
2179 dir = rotComp[fIdToCompId[pIt.first]].m_dir;
2180 angle = rotComp[fIdToCompId[pIt.first]].m_angle;
2181 tol = rotComp[fIdToCompId[pIt.first]].m_tol;
2185 for (i = 0; i < 2; ++i)
2193 int other = (i + 1) % 2;
2196 sign = (i == 0) ? 1.0 : -1.0;
2199 if (tmpVec[0].size() == 3)
2202 tmpVec[i], tmpVec[other], rotbnd, dir,
2208 tmpVec[i], tmpVec[other], rotbnd, dir,
2218 int nFaceVerts = vertMap[ids[0]].size();
2221 for (i = 0; i < 2; ++i)
2223 int other = (i + 1) % 2;
2226 sign = (i == 0) ? 1.0 : -1.0;
2229 if (tmpVec[0].size() == 3)
2232 tmpVec[i], tmpVec[other], rotbnd, dir,
2238 tmpVec[i], tmpVec[other], rotbnd, dir,
2242 if (nFaceVerts == 3)
2247 "Unsupported face orientation for face " +
2248 boost::lexical_cast<string>(ids[i]));
2252 vector<int> per1 = vertMap[ids[i]];
2253 vector<int> per2 = vertMap[ids[other]];
2257 map<int, pair<int, bool>> tmpMap;
2261 for (j = 0; j < nFaceVerts; ++j)
2263 int v = vmap[nFaceVerts][o][j];
2265 make_pair(per2[v], locVerts.count(per2[v]) > 0);
2269 for (
auto &mIt : tmpMap)
2276 auto perIt = periodicVerts.find(mIt.first);
2278 if (perIt == periodicVerts.end())
2282 periodicVerts[mIt.first].push_back(ent2);
2283 perIt = periodicVerts.find(mIt.first);
2290 for (k = 0; k < perIt->second.size(); ++k)
2292 if (perIt->second[k].id == mIt.second.first)
2298 if (k == perIt->second.size())
2300 perIt->second.push_back(ent2);
2308 for (i = 0; i < 2; ++i)
2310 int other = (i + 1) % 2;
2313 sign = (i == 0) ? 1.0 : -1.0;
2315 if (tmpVec[0].size() == 3)
2318 tmpVec[i], tmpVec[other], rotbnd, dir,
2324 tmpVec[i], tmpVec[other], rotbnd, dir,
2328 vector<int> per1 = edgeMap[ids[i]];
2329 vector<int> per2 = edgeMap[ids[other]];
2331 map<int, pair<int, bool>> tmpMap;
2333 for (j = 0; j < nFaceVerts; ++j)
2335 int e = emap[nFaceVerts][o][j];
2337 make_pair(per2[e], locEdges.count(per2[e]) > 0);
2340 for (
auto &mIt : tmpMap)
2347 auto perIt = periodicEdges.find(mIt.first);
2349 if (perIt == periodicEdges.end())
2351 periodicEdges[mIt.first].push_back(ent2);
2352 perIt = periodicEdges.find(mIt.first);
2356 for (k = 0; k < perIt->second.size(); ++k)
2358 if (perIt->second[k].id == mIt.second.first)
2364 if (k == perIt->second.size())
2366 perIt->second.push_back(ent2);
2376 ASSERTL1(allCompPairs.size() == fIdToCompId.size(),
2377 "At this point the size of allCompPairs "
2378 "should have been the same as fIdToCompId");
2382 map<int, int> eIdToCompId;
2388 for (cnt = i = 0; i < totFaces; ++i)
2390 bool rotbnd =
false;
2395 int faceId = faceIds[i];
2397 ASSERTL0(allCompPairs.count(faceId) > 0,
2398 "Unable to find matching periodic face. faceId = " +
2399 boost::lexical_cast<string>(faceId));
2401 int perFaceId = allCompPairs[faceId];
2404 ASSERTL1((rotComp.size() == 0) || fIdToCompId.count(faceId) > 0,
2405 "Face " + boost::lexical_cast<string>(faceId) +
2406 " not found in fIdtoCompId map");
2407 if (rotComp.count(fIdToCompId[faceId]))
2410 dir = rotComp[fIdToCompId[faceId]].m_dir;
2411 angle = rotComp[fIdToCompId[faceId]].m_angle;
2412 tol = rotComp[fIdToCompId[faceId]].m_tol;
2415 for (j = 0; j < faceVerts[i]; ++j, ++cnt)
2417 int vId = vertIds[cnt];
2419 auto perId = periodicVerts.find(vId);
2421 if (perId == periodicVerts.end())
2430 coordMap[faceId], coordMap[perFaceId]};
2432 int nFaceVerts = tmpVec[0].size();
2436 tmpVec[0], tmpVec[1], rotbnd, dir, angle,
2439 tmpVec[0], tmpVec[1], rotbnd, dir, angle,
2445 vertMap[perFaceId][vmap[nFaceVerts][o][j]];
2449 locVerts.count(perVertexId) > 0);
2451 periodicVerts[vId].push_back(ent);
2454 int eId = edgeIds[cnt];
2456 perId = periodicEdges.find(eId);
2462 eIdToCompId[eId] = fIdToCompId[faceId];
2465 if (perId == periodicEdges.end())
2473 coordMap[faceId], coordMap[perFaceId]};
2475 int nFaceEdges = tmpVec[0].size();
2479 tmpVec[0], tmpVec[1], rotbnd, dir, angle,
2482 tmpVec[0], tmpVec[1], rotbnd, dir, angle,
2488 edgeMap[perFaceId][emap[nFaceEdges][o][j]];
2491 locEdges.count(perEdgeId) > 0);
2493 periodicEdges[eId].push_back(ent);
2499 eIdToCompId[perEdgeId] = fIdToCompId[perFaceId];
2507 for (
auto &perIt : periodicVerts)
2510 for (i = 0; i < perIt.second.size(); ++i)
2513 auto perIt2 = periodicVerts.find(perIt.second[i].id);
2514 ASSERTL0(perIt2 != periodicVerts.end(),
2515 "Couldn't find periodic vertex.");
2520 for (j = 0; j < perIt2->second.size(); ++j)
2522 if (perIt2->second[j].id == perIt.first)
2527 for (k = 0; k < perIt.second.size(); ++k)
2529 if (perIt2->second[j].id == perIt.second[k].id)
2535 if (k == perIt.second.size())
2537 perIt.second.push_back(perIt2->second[j]);
2543 for (
auto &perIt : periodicEdges)
2545 for (i = 0; i < perIt.second.size(); ++i)
2547 auto perIt2 = periodicEdges.find(perIt.second[i].id);
2548 ASSERTL0(perIt2 != periodicEdges.end(),
2549 "Couldn't find periodic edge.");
2551 for (j = 0; j < perIt2->second.size(); ++j)
2553 if (perIt2->second[j].id == perIt.first)
2558 for (k = 0; k < perIt.second.size(); ++k)
2560 if (perIt2->second[j].id == perIt.second[k].id)
2566 if (k == perIt.second.size())
2568 perIt.second.push_back(perIt2->second[j]);
2575 for (
auto &perIt : periodicEdges)
2577 bool rotbnd =
false;
2583 auto eIt = eIdMap.find(perIt.first);
2585 *vCoMap[eIt->second.second]};
2588 if (rotComp.count(eIdToCompId[perIt.first]))
2591 dir = rotComp[eIdToCompId[perIt.first]].m_dir;
2592 angle = rotComp[eIdToCompId[perIt.first]].m_angle;
2593 tol = rotComp[eIdToCompId[perIt.first]].m_tol;
2599 for (i = 0; i < perIt.second.size(); ++i)
2601 eIt = eIdMap.find(perIt.second[i].id);
2604 *vCoMap[eIt->second.first],
2605 *vCoMap[eIt->second.second]};
2607 int vMap[2] = {-1, -1};
2613 r.
Rotate(v[0], dir, angle);
2615 if (r.
dist(w[0]) < tol)
2621 r.
Rotate(v[1], dir, angle);
2622 if (r.
dist(w[0]) < tol)
2629 "Unable to align rotationally "
2630 "periodic edge vertex");
2637 0.5 * (w[0](0) - v[0](0) + w[1](0) - v[1](0));
2639 0.5 * (w[0](1) - v[0](1) + w[1](1) - v[1](1));
2641 0.5 * (w[0](2) - v[0](2) + w[1](2) - v[1](2));
2643 for (j = 0; j < 2; ++j)
2648 for (k = 0; k < 2; ++k)
2654 if (
sqrt((x1 - x) * (x1 - x) +
2655 (y1 - y) * (y1 - y) +
2656 (z1 - z) * (z1 - z)) < 1e-8)
2665 ASSERTL0(vMap[0] >= 0 && vMap[1] >= 0,
2666 "Unable to align periodic edge vertex.");
2667 ASSERTL0((vMap[0] == 0 || vMap[0] == 1) &&
2668 (vMap[1] == 0 || vMap[1] == 1) &&
2669 (vMap[0] != vMap[1]),
2670 "Unable to align periodic edge vertex.");
2684 for (
auto &perIt : periodicVerts)
2686 if (locVerts.count(perIt.first) > 0)
2692 for (
auto &perIt : periodicEdges)
2694 if (locEdges.count(perIt.first) > 0)
2702 ASSERTL1(
false,
"Not setup for this expansion");
2714 "Routine currently only tested for HybridDGHelmholtz");
2717 "Full matrix global systems are not supported for HDG "
2721 "The local to global map is not set up for the requested "
2730 (*m_globalBndMat)[mkey] = glo_matrix;
2734 glo_matrix = matrixIter->second;
2749 bool returnval =
true;
2768 int vSame = (returnval ? 1 : 0);
2771 return (vSame == 1);
2786 for (
int v = 0; v < 2; ++v)
2791 if (vertExp->GetLeftAdjacentElementExp()
2794 (*m_exp)[i]->GetGeom()->GetGlobalID())
2918 bool PutFwdInBwdOnBCs,
bool DoExchange)
2950 for (
int n = 0, cnt = 0; n < nexp; ++n)
2955 for (
int e = 0; e < exp->GetNtraces(); ++e, ++cnt)
2958 m_trace->GetPhys_Offset(elmtToTrace[n][e]->GetElmtId());
2963 exp->GetTracePhysVals(e, elmtToTrace[n][e], field + phys_offset,
2979 m_traceMap->GetAssemblyCommDG()->PerformExchange(Fwd, Bwd);
2990 bool PutFwdInBwdOnBCs)
2993 if (PutFwdInBwdOnBCs)
3002 for (
int e = 0; e < ne; ++e)
3006 auto id2 =
m_trace->GetPhys_Offset(
3007 m_traceMap->GetBndCondIDToGlobalTraceID(cnt + e));
3019 for (
int e = 0; e < ne; ++e)
3023 auto id2 =
m_trace->GetPhys_Offset(
3024 m_traceMap->GetBndCondIDToGlobalTraceID(cnt + e));
3032 ASSERTL0(
false,
"Method not set up for this "
3033 "boundary condition.");
3045 for (
int e = 0; e < ne; ++e)
3050 auto id2 =
m_trace->GetPhys_Offset(
3051 m_traceMap->GetBndCondIDToGlobalTraceID(cnt + e));
3065 for (
int e = 0; e < ne; ++e)
3072 "Method not set up for non-zero "
3073 "Neumann boundary condition");
3074 auto id2 =
m_trace->GetPhys_Offset(
3075 m_traceMap->GetBndCondIDToGlobalTraceID(cnt + e));
3091 "Method not set up for this boundary "
3131 "v_AddTraceQuadPhysToField not coded for eGauss_Lagrange");
3166 m_traceMap->GetAssemblyCommDG()->PerformExchange(outarray, outarray);
3172 int n,
p, offset, phys_offset;
3177 "input array is of insufficient length");
3178 for (n = 0; n < nexp; ++n)
3181 for (
p = 0;
p < (*m_exp)[n]->GetNtraces(); ++
p)
3184 m_trace->GetPhys_Offset(elmtToTrace[n][
p]->GetElmtId());
3185 (*m_exp)[n]->GetTracePhysVals(
p, elmtToTrace[n][
p],
3186 inarray + phys_offset,
3187 t_tmp = outarray + offset);
3215 int n, offset, t_offset;
3222 int e_ncoeffs = (*m_exp)[n]->GetNcoeffs();
3225 if ((*
m_exp)[n]->GetBasis(0)->GetBasisType() !=
3229 GetTrace()->GetCoeff_Offset(elmtToTrace[n][0]->GetElmtId());
3230 if (negatedFluxNormal[2 * n])
3232 outarray[offset] -= Fn[t_offset];
3236 outarray[offset] += Fn[t_offset];
3239 GetTrace()->GetCoeff_Offset(elmtToTrace[n][1]->GetElmtId());
3240 if (negatedFluxNormal[2 * n + 1])
3242 outarray[offset + (*m_exp)[n]->GetVertexMap(1)] -=
3247 outarray[offset + (*m_exp)[n]->GetVertexMap(1)] +=
3255 static int sav_ncoeffs = 0;
3256 if (!m_Ixm.get() || e_ncoeffs != sav_ncoeffs)
3266 m_Ixm = BASE->GetI(coords);
3268 m_Ixp = BASE->GetI(coords);
3269 sav_ncoeffs = e_ncoeffs;
3272 GetTrace()->GetCoeff_Offset(elmtToTrace[n][0]->GetElmtId());
3273 if (negatedFluxNormal[2 * n])
3275 for (j = 0; j < e_ncoeffs; j++)
3277 outarray[offset + j] -=
3278 (m_Ixm->GetPtr())[j] * Fn[t_offset];
3283 for (j = 0; j < e_ncoeffs; j++)
3285 outarray[offset + j] +=
3286 (m_Ixm->GetPtr())[j] * Fn[t_offset];
3290 GetTrace()->GetCoeff_Offset(elmtToTrace[n][1]->GetElmtId());
3291 if (negatedFluxNormal[2 * n + 1])
3293 for (j = 0; j < e_ncoeffs; j++)
3295 outarray[offset + j] -=
3296 (m_Ixp->GetPtr())[j] * Fn[t_offset];
3301 for (j = 0; j < e_ncoeffs; j++)
3303 outarray[offset + j] +=
3304 (m_Ixp->GetPtr())[j] * Fn[t_offset];
3318 m_trace->IProductWRTBase(Fn, Fcoeffs);
3324 int e, n, offset, t_offset;
3333 for (e = 0; e < (*m_exp)[n]->GetNtraces(); ++e)
3335 t_offset =
GetTrace()->GetPhys_Offset(
3336 elmtToTrace[n][e]->GetElmtId());
3337 (*m_exp)[n]->AddEdgeNormBoundaryInt(
3338 e, elmtToTrace[n][e], Fn + t_offset,
3339 e_outarray = outarray + offset);
3348 for (e = 0; e < (*m_exp)[n]->GetNtraces(); ++e)
3350 t_offset =
GetTrace()->GetPhys_Offset(
3351 elmtToTrace[n][e]->GetElmtId());
3352 (*m_exp)[n]->AddFaceNormBoundaryInt(
3353 e, elmtToTrace[n][e], Fn + t_offset,
3354 e_outarray = outarray + offset);
3386 "tested for 1D expansion");
3388 m_trace->IProductWRTBase(Fwd, Coeffs);
3390 m_trace->IProductWRTBase(Bwd, Coeffs);
3401 boost::ignore_unused(varfactors, dirForcing);
3402 int i, n, cnt, nbndry;
3410 if (PhysSpaceForcing)
3422 int GloBndDofs =
m_traceMap->GetNumGlobalBndCoeffs();
3423 int NumDirBCs =
m_traceMap->GetNumLocalDirBndCoeffs();
3430 int LocBndCoeffs =
m_traceMap->GetNumLocalBndCoeffs();
3439 for (cnt = n = 0; n < nexp; ++n)
3441 nbndry = (*m_exp)[n]->NumDGBndryCoeffs();
3442 e_ncoeffs = (*m_exp)[n]->GetNcoeffs();
3448 Floc =
Transpose(*(HDGLamToU->GetBlock(n, n))) * ElmtFce;
3452 m_traceMap->GetBndCondCoeffsToLocalTraceMap();
3468 locid = bndCondMap[cnt + j];
3469 loclambda[locid] = Sign[locid] * bndcoeffs[j];
3480 locid = bndCondMap[cnt + j];
3481 bndrhs[locid] += Sign[locid] * bndcoeffs[j];
3498 if (GloBndDofs - NumDirBCs > 0)
3503 LinSys->Solve(bndrhs, loclambda,
m_traceMap);
3521 out = (*InvHDGHelm) * F + (*HDGLamToU) * LocLambda;
3563 const std::string varName,
3579 npoints = locExpList->GetNpoints();
3584 locExpList->GetCoords(x0, x1, x2);
3608 0, (std::static_pointer_cast<
3611 ->m_dirichletCondition)
3612 .Evaluate(x0[0], x1[0], x2[0], time));
3620 0, (std::static_pointer_cast<
3623 ->m_neumannCondition)
3624 .Evaluate(x0[0], x1[0], x2[0], time));
3630 0, (std::static_pointer_cast<
3634 .Evaluate(x0[0], x1[0], x2[0], time));
3648 "This type of BC not implemented yet");
3657 std::static_pointer_cast<
3662 valuesExp(npoints, 1.0);
3664 string bcfilename = bcPtr->m_filename;
3665 string exprbcs = bcPtr->m_expr;
3667 if (bcfilename !=
"")
3669 locExpList->ExtractCoeffsFromFile(
3670 bcfilename, bcPtr->GetComm(), varName,
3671 locExpList->UpdateCoeffs());
3672 locExpList->BwdTrans(locExpList->GetCoeffs(),
3673 locExpList->UpdatePhys());
3674 valuesFile = locExpList->GetPhys();
3680 std::static_pointer_cast<
3683 ->m_dirichletCondition;
3685 condition.
Evaluate(x0, x1, x2, time, valuesExp);
3689 locExpList->UpdatePhys(), 1);
3690 locExpList->FwdTransBndConstrained(
3691 locExpList->GetPhys(), locExpList->UpdateCoeffs());
3697 std::static_pointer_cast<
3702 if (bcfilename !=
"")
3704 locExpList->ExtractCoeffsFromFile(
3705 bcfilename, bcPtr->GetComm(), varName,
3706 locExpList->UpdateCoeffs());
3707 locExpList->BwdTrans(locExpList->GetCoeffs(),
3708 locExpList->UpdatePhys());
3713 std::static_pointer_cast<
3716 ->m_neumannCondition;
3717 condition.
Evaluate(x0, x1, x2, time,
3718 locExpList->UpdatePhys());
3721 locExpList->IProductWRTBase(locExpList->GetPhys(),
3722 locExpList->UpdateCoeffs());
3728 std::static_pointer_cast<
3734 if (bcfilename !=
"")
3736 locExpList->ExtractCoeffsFromFile(
3737 bcfilename, bcPtr->GetComm(), varName,
3738 locExpList->UpdateCoeffs());
3739 locExpList->BwdTrans(locExpList->GetCoeffs(),
3740 locExpList->UpdatePhys());
3745 std::static_pointer_cast<
3749 condition.
Evaluate(x0, x1, x2, time,
3750 locExpList->UpdatePhys());
3753 locExpList->IProductWRTBase(locExpList->GetPhys(),
3754 locExpList->UpdateCoeffs());
3764 "This type of BC not implemented yet");
3793 id2 =
m_trace->GetPhys_Offset(
3794 m_traceMap->GetBndCondIDToGlobalTraceID(cnt + e));
3809 id2 =
m_trace->GetPhys_Offset(
3810 m_traceMap->GetBndCondIDToGlobalTraceID(cnt + e));
3821 "Method not set up for this boundary condition.");
3838 map<int, int> VertGID;
3863 VertGID[Vid] = cnt++;
3872 for (i = 0; i < exp->GetNverts(); ++i)
3874 id = exp->GetGeom()->GetVid(i);
3876 if (VertGID.count(
id) > 0)
3878 bid = VertGID.find(
id)->second;
3880 "Edge already set");
3888 "Failed to visit all boundary condtiions");
3893 map<int, int> globalIdMap;
3902 globalIdMap[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
3928 m_graph->GetElementsFromEdge(exp1d->GetGeom1D());
3931 globalIdMap[(*tmp)[0].first->GetGlobalID()];
3939 map<int, int> globalIdMap;
3949 globalIdMap[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
3973 m_graph->GetElementsFromFace(exp2d->GetGeom2D());
3975 globalIdMap[tmp->at(0).first->GetGlobalID()];
3982 ASSERTL1(
false,
"Not setup for this expansion");
3992 std::shared_ptr<ExpList> &result,
3993 const bool DeclareCoeffPhysArrays)
3996 std::vector<unsigned int> eIDs;
4002 for (cnt = n = 0; n < i; ++n)
4010 eIDs.push_back(ElmtID[cnt + n]);
4017 if (DeclareCoeffPhysArrays)
4020 int offsetOld, offsetNew;
4022 for (n = 0; n < result->GetExpSize(); ++n)
4024 nq =
GetExp(ElmtID[cnt + n])->GetTotPoints();
4026 offsetNew = result->GetPhys_Offset(n);
4028 tmp2 = result->UpdatePhys() + offsetNew, 1);
4029 nq =
GetExp(ElmtID[cnt + n])->GetNcoeffs();
4031 offsetNew = result->GetCoeff_Offset(n);
4033 tmp2 = result->UpdateCoeffs() + offsetNew, 1);
4069 map<int, RobinBCInfoSharedPtr> returnval;
4085 int npoints = locExpList->GetNpoints();
4091 locExpList->GetCoords(x0, x1, x2);
4094 std::static_pointer_cast<
4096 ->m_robinPrimitiveCoeff;
4099 coeffeqn.
Evaluate(x0, x1, x2, 0.0, coeffphys);
4101 for (e = 0; e < locExpList->GetExpSize(); ++e)
4106 Array_tmp = coeffphys + locExpList->GetPhys_Offset(e));
4108 elmtid = ElmtID[cnt + e];
4110 if (returnval.count(elmtid) != 0)
4112 rInfo->next = returnval.find(elmtid)->second;
4114 returnval[elmtid] = rInfo;
4137 int i, e, ncoeff_edge;
4147 int LocBndCoeffs =
m_traceMap->GetNumLocalBndCoeffs();
4151 edge_lambda = loc_lambda;
4159 int nEdges = (*m_exp)[i]->GetGeom()->GetNumEdges();
4162 for (e = 0; e < nEdges; ++e)
4164 edgedir = (*m_exp)[i]->GetTraceOrient(e);
4165 ncoeff_edge = elmtToTrace[i][e]->GetNcoeffs();
4167 Vmath::Vcopy(ncoeff_edge, edge_lambda, 1, edgeCoeffs[e], 1);
4168 elmtToTrace[i][e]->SetCoeffsToOrientation(edgedir, edgeCoeffs[e],
4170 edge_lambda = edge_lambda + ncoeff_edge;
4173 (*m_exp)[i]->DGDeriv(dir, tmp_coeffs = coeffs +
m_coeff_offset[i],
4174 elmtToTrace[i], edgeCoeffs, out_tmp = out_d + cnt);
4175 cnt += (*m_exp)[i]->GetNcoeffs();
4208 int i, cnt, e, ncoeff_trace;
4215 int nq_elmt, nm_elmt;
4216 int LocBndCoeffs =
m_traceMap->GetNumLocalBndCoeffs();
4221 trace_lambda = loc_lambda;
4225 int num_points[] = {0, 0, 0};
4226 int num_modes[] = {0, 0, 0};
4231 nq_elmt = (*m_exp)[i]->GetTotPoints();
4232 nm_elmt = (*m_exp)[i]->GetNcoeffs();
4236 out_tmp = force + nm_elmt;
4239 for (
int j = 0; j < dim; ++j)
4241 num_points[j] = (*m_exp)[i]->GetBasis(j)->GetNumPoints();
4242 num_modes[j] = (*m_exp)[i]->GetBasis(j)->GetNumModes();
4248 int nTraces = (*m_exp)[i]->GetNtraces();
4251 for (e = 0; e < (*m_exp)[i]->GetNtraces(); ++e)
4253 edgedir = (*m_exp)[i]->GetTraceOrient(e);
4254 ncoeff_trace = elmtToTrace[i][e]->GetNcoeffs();
4256 Vmath::Vcopy(ncoeff_trace, trace_lambda, 1, traceCoeffs[e], 1);
4259 elmtToTrace[i][e]->SetCoeffsToOrientation(
4260 edgedir, traceCoeffs[e], traceCoeffs[e]);
4266 ->SetFaceToGeomOrientation(e, traceCoeffs[e]);
4268 trace_lambda = trace_lambda + ncoeff_trace;
4283 num_modes[0], PkeyQ1);
4285 num_modes[1], PkeyQ2);
4287 std::dynamic_pointer_cast<SpatialDomains::QuadGeom>(
4288 (*
m_exp)[i]->GetGeom());
4290 BkeyQ1, BkeyQ2, qGeom);
4298 num_points[1], LibUtilities::eGaussRadauMAlpha1Beta0);
4300 num_modes[0], PkeyT1);
4302 num_modes[1], PkeyT2);
4304 std::dynamic_pointer_cast<SpatialDomains::TriGeom>(
4305 (*
m_exp)[i]->GetGeom());
4307 BkeyT1, BkeyT2, tGeom);
4319 num_modes[0], PkeyH1);
4321 num_modes[1], PkeyH2);
4323 num_modes[2], PkeyH3);
4325 std::dynamic_pointer_cast<SpatialDomains::HexGeom>(
4326 (*
m_exp)[i]->GetGeom());
4328 BkeyH1, BkeyH2, BkeyH3, hGeom);
4336 num_points[1], LibUtilities::eGaussRadauMAlpha1Beta0);
4338 num_points[2], LibUtilities::eGaussRadauMAlpha2Beta0);
4340 num_modes[0], PkeyT1);
4342 num_modes[1], PkeyT2);
4344 num_modes[2], PkeyT3);
4346 std::dynamic_pointer_cast<SpatialDomains::TetGeom>(
4347 (*
m_exp)[i]->GetGeom());
4349 BkeyT1, BkeyT2, BkeyT3, tGeom);
4354 "Wrong shape type, HDG postprocessing is not "
4361 elmtToTrace[i], traceCoeffs, out_tmp);
4362 (*m_exp)[i]->BwdTrans(out_tmp, qrhs);
4363 ppExp->IProductWRTDerivBase(0, qrhs, force);
4367 elmtToTrace[i], traceCoeffs, out_tmp);
4369 (*m_exp)[i]->BwdTrans(out_tmp, qrhs);
4370 ppExp->IProductWRTDerivBase(1, qrhs, out_tmp);
4372 Vmath::Vadd(nm_elmt, force, 1, out_tmp, 1, force, 1);
4375 (*m_exp)[i]->BwdTrans(tmp_coeffs = coeffs +
m_coeff_offset[i], qrhs);
4376 force[0] = (*m_exp)[i]->Integral(qrhs);
4381 ppExp->DetShapeType(), *ppExp);
4387 out = (*lapsys) * in;
4391 ppExp->BwdTrans(out.
GetPtr(), work);
4392 (*m_exp)[i]->FwdTrans(work, tmp_coeffs = outarray +
m_coeff_offset[i]);
4421 m_trace->IProductWRTBase(FwdFlux, FCoeffs);
4423 m_trace->IProductWRTBase(BwdFlux, FCoeffs);
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
#define sign(a, b)
return the sign(b)*a
Describes the specification for a Basis.
NekDouble Evaluate() const
Interpreter class for the evaluation of mathematical expressions.
int DefineFunction(const std::string &vlist, const std::string &expr)
Defines a function for the purposes of evaluation.
NekDouble Evaluate(const int id)
Evaluate a function which depends only on constants and/or parameters.
Defines a specification for a set of points.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
This class is the abstractio n of a global discontinuous two- dimensional spectral/hp element expansi...
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions() override
std::vector< int > m_periodicBwdCopy
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray) override
Add trace contributions into elemental coefficient spaces.
PeriodicMap m_periodicEdges
A map which identifies pairs of periodic edges.
void SetUpDG(const std::string="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
Set up all DG member variables and maps.
PeriodicMap m_periodicFaces
A map which identifies pairs of periodic faces.
std::set< int > m_boundaryTraces
A set storing the global IDs of any boundary Verts.
LocTraceToTraceMapSharedPtr m_locTraceToTraceMap
virtual ~DisContField()
Destructor.
std::vector< bool > m_negatedFluxNormal
SpatialDomains::BoundaryConditionsSharedPtr GetDomainBCs(const SpatialDomains::CompositeMap &domain, const SpatialDomains::BoundaryConditions &Allbcs, const std::string &variable)
AssemblyMapDGSharedPtr m_traceMap
Local to global DG mapping for trace space.
Array< OneD, int > m_BCtoTraceMap
void EvaluateHDGPostProcessing(const Array< OneD, const NekDouble > &coeffs, Array< OneD, NekDouble > &outarray)
Evaluate HDG post-processing to increase polynomial order of solution.
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces) override
Obtain a copy of the periodic edges and vertices for this field.
GlobalLinSysSharedPtr GetGlobalBndLinSys(const GlobalLinSysKey &mkey)
For a given key, returns the associated global linear system.
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions() override
virtual MultiRegions::ExpListSharedPtr & v_UpdateBndCondExpansion(int i) override
virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value) override
virtual void v_PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd) override
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &TraceID) override
virtual void v_ExtractTracePhys(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
This method extracts the trace (verts in 1D) from the field inarray and puts the values in outarray.
InterfaceMapDGSharedPtr m_interfaceMap
Interfaces mapping for trace space.
virtual void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays) override
virtual void v_FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp) override
Fill the weight with m_bndCondBndWeight.
virtual const Array< OneD, const NekDouble > & v_GetBndCondBwdWeight() override
virtual const Array< OneD, const MultiRegions::ExpListSharedPtr > & v_GetBndCondExpansions() override
void FindPeriodicTraces(const SpatialDomains::BoundaryConditions &bcs, const std::string variable)
Generate a associative map of periodic vertices in a mesh.
Array< OneD, NekDouble > m_bndCondBndWeight
bool SameTypeOfBoundaryConditions(const DisContField &In)
Check to see if expansion has the same BCs as In.
virtual ExpListSharedPtr & v_GetTrace() override
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition structure definition on the diff...
virtual void v_EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble) override
Evaluate all boundary conditions at a given time..
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray) override
Add trace contributions into elemental coefficient spaces.
virtual const LocTraceToTraceMapSharedPtr & v_GetLocTraceToTraceMap(void) const override
virtual void v_AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray) override
virtual void v_FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs) override
GlobalLinSysMapShPtr m_globalBndMat
Global boundary matrix.
PeriodicMap m_periodicVerts
A map which identifies groups of periodic vertices.
ExpListSharedPtr m_trace
Trace space storage for points between elements.
virtual void v_AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field) override
virtual void v_Reset() override
Reset this field, so that geometry information can be updated.
std::vector< bool > & GetNegatedFluxNormal(void)
DisContField()
Default constructor.
virtual std::vector< bool > & v_GetLeftAdjacentTraces(void) override
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
virtual GlobalLinSysKey v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing) override
Solve the Helmholtz equation.
std::vector< int > m_periodicFwdCopy
A vector indicating degress of freedom which need to be copied from forwards to backwards space in ca...
bool IsLeftAdjacentTrace(const int n, const int e)
This routine determines if an element is to the "left" of the adjacent trace, which arises from the i...
std::vector< bool > m_leftAdjacentTraces
NekDouble L2_DGDeriv(const int dir, const Array< OneD, const NekDouble > &coeffs, const Array< OneD, const NekDouble > &soln)
Calculate the error of the derivative using the consistent DG evaluation of .
virtual void v_GetLocTraceFromTracePts(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd) override
Array< OneD, int > m_BCtoElmMap
virtual AssemblyMapDGSharedPtr & v_GetTraceMap(void) override
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo() override
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd) override
void GenerateBoundaryConditionExpansion(const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::BoundaryConditions &bcs, const std::string variable, const bool DeclareCoeffPhysArrays=true)
Discretises the boundary conditions.
Base class for all multi-elemental spectral/hp expansions.
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to all local expansion modes .
std::shared_ptr< ExpList > & GetTrace()
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
const Array< OneD, const NekDouble > & GetPhys() const
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
int GetCoeff_Offset(int n) const
Get the start offset position for a local contiguous list of coeffs correspoinding to element n.
int GetExpSize(void)
This function returns the number of elements in the expansion.
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
std::shared_ptr< GlobalLinSys > GenGlobalBndLinSys(const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocTo...
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
bool m_physState
The state of the array m_phys.
LibUtilities::CommSharedPtr m_comm
Communicator.
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
NekDouble L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error of the global This function calculates the error with respect to...
void ApplyGeomInfo()
Apply geometry information to each expansion.
LibUtilities::SessionReaderSharedPtr m_session
Session.
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the backward transformation of the global spectral/hp element exp...
int GetPhys_Offset(int n) const
Get the start offset position for a local contiguous list of quadrature points in a full array corres...
ExpansionType m_expType
Exapnsion type.
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Describe a linear system.
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
Describes a matrix with ordering defined by a local to global map.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
Array< OneD, DataType > & GetPtr()
const BoundaryConditionCollection & GetBoundaryConditions(void) const
const BoundaryRegionCollection & GetBoundaryRegions(void) const
void Rotate(PointGeom &a, int dir, NekDouble angle)
_this = rotation of a by angle 'angle' around axis dir
NekDouble dist(PointGeom &a)
return distance between this and input a
static StdRegions::Orientation GetFaceOrientation(const QuadGeom &face1, const QuadGeom &face2, bool doRot=false, int dir=0, NekDouble angle=0.0, NekDouble tol=1e-8)
Get the orientation of face1.
static StdRegions::Orientation GetFaceOrientation(const TriGeom &face1, const TriGeom &face2, bool doRot, int dir, NekDouble angle, NekDouble tol)
BasisManagerT & BasisManager(void)
std::shared_ptr< Basis > BasisSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
@ eOrtho_A
Principle Orthogonal Functions .
@ eOrtho_C
Principle Orthogonal Functions .
@ eOrtho_B
Principle Orthogonal Functions .
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
std::shared_ptr< Expansion0D > Expansion0DSharedPtr
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
static AssemblyMapSharedPtr NullAssemblyMapSharedPtr
static ExpListSharedPtr NullExpListSharedPtr
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
std::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static GlobalLinSysKey NullGlobalLinSysKey(StdRegions::eNoMatrixType)
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
static const NekDouble kNekUnsetDouble
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
std::map< int, std::vector< unsigned int > > CompositeOrdering
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
std::shared_ptr< BoundaryRegion > BoundaryRegionShPtr
std::shared_ptr< BoundaryConditions > BoundaryConditionsSharedPtr
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
std::shared_ptr< std::vector< std::pair< GeometrySharedPtr, int > > > GeometryLinkSharedPtr
std::map< int, std::vector< unsigned int > > BndRegionOrdering
std::shared_ptr< DirichletBoundaryCondition > DirichletBCShPtr
std::shared_ptr< Composite > CompositeSharedPtr
std::vector< PointGeomSharedPtr > PointGeomVector
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
std::shared_ptr< HexGeom > HexGeomSharedPtr
std::shared_ptr< SegGeom > SegGeomSharedPtr
std::shared_ptr< TetGeom > TetGeomSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::shared_ptr< NeumannBoundaryCondition > NeumannBCShPtr
std::shared_ptr< BoundaryConditionMap > BoundaryConditionMapShPtr
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
std::shared_ptr< TriGeom > TriGeomSharedPtr
std::shared_ptr< RobinBoundaryCondition > RobinBCShPtr
std::map< int, CompositeSharedPtr > CompositeMap
@ eInvLaplacianWithUnityMean
std::map< ConstFactorType, NekDouble > ConstFactorMap
@ eDir1BwdDir2_Dir2BwdDir1
@ eDir1FwdDir1_Dir2FwdDir2
@ eDir1BwdDir1_Dir2BwdDir2
@ eDir1BwdDir2_Dir2FwdDir1
@ eDir1FwdDir1_Dir2BwdDir2
@ eDir1BwdDir1_Dir2FwdDir2
@ eDir1FwdDir2_Dir2FwdDir1
@ eDir1FwdDir2_Dir2BwdDir1
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
static Array< OneD, NekDouble > NullNekDouble1DArray
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 Neg(int n, T *x, const int incx)
Negate x = -x.
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
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.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
scalarT< T > sqrt(scalarT< T > in)
int id
Geometry ID of entity.
StdRegions::Orientation orient
Orientation of entity within higher dimensional entity.
bool isLocal
Flag specifying if this entity is local to this partition.
NekDouble m_tol
Tolerance to rotation is considered identical.
int m_dir
Axis of rotation. 0 = 'x', 1 = 'y', 2 = 'z'.
NekDouble m_angle
Angle of rotation in radians.