73 int nel =
m_f->m_exp[0]->GetExpSize();
77 int nfields =
m_f->m_variables.size();
81 for (
int i = 0; i < nfields + coordim; ++i)
85 vector<Array<OneD, int>> ptsConn;
89 coordim,
m_f->m_variables, pts);
91 m_f->m_fieldPts->SetConnectivity(ptsConn);
95 int coordim =
m_f->m_exp[0]->GetCoordim(0);
96 int shapedim =
m_f->m_exp[0]->GetExp(0)->GetShapeDimension();
97 int npts =
m_f->m_exp[0]->GetTotPoints();
101 bool homogeneous1D =
false;
102 if (
m_f->m_numHomogeneousDir == 1)
106 homogeneous1D =
true;
108 else if (
m_f->m_numHomogeneousDir == 2)
110 ASSERTL0(
false,
"Homegeneous2D case not supported");
115 int newtotpoints = 0;
123 map<int, StdRegions::Orientation> face0orient;
124 set<int> prismorient;
129 vector<Array<OneD, int>> ptsConn;
132 for (
int i = 0; i < nel; ++i)
134 e =
m_f->m_exp[0]->GetExp(i);
138 int fid = e->GetGeom()->GetFid(0);
139 if (face0orient.count(fid))
141 prismorient.insert(i);
161 for (
int i = 0; i < nel; ++i)
163 e =
m_f->m_exp[0]->GetExp(i);
166 int fid = e->GetGeom()->GetFid(2);
168 if (face0orient.count(fid))
183 prismorient.insert(i);
190 prismorient.insert(i);
197 for (
int i = 0; i < nel; ++i)
199 e =
m_f->m_exp[0]->GetExp(i);
202 if (
m_f->m_exp[0]->GetExp(i)->DetShapeType() !=
209 switch (e->DetShapeType())
213 int npoints0 = e->GetBasis(0)->GetNumPoints();
221 int np0 = e->GetBasis(0)->GetNumPoints();
222 int np1 = e->GetBasis(1)->GetNumPoints();
223 int np =
max(np0, np1);
230 int np0 = e->GetBasis(0)->GetNumPoints();
231 int np1 = e->GetBasis(1)->GetNumPoints();
232 int np =
max(np0, np1);
240 int np0 = e->GetBasis(0)->GetNumPoints();
241 int np1 = e->GetBasis(1)->GetNumPoints();
242 int np2 = e->GetBasis(2)->GetNumPoints();
243 int np =
max(np0,
max(np1, np2));
251 int np0 = e->GetBasis(0)->GetNumPoints();
252 int np1 = e->GetBasis(1)->GetNumPoints();
253 int np2 = e->GetBasis(2)->GetNumPoints();
254 int np =
max(np0,
max(np1, np2));
262 int np0 = e->GetBasis(0)->GetNumPoints();
263 int np1 = e->GetBasis(1)->GetNumPoints();
264 int np2 = e->GetBasis(2)->GetNumPoints();
265 int np =
max(np0,
max(np1, np2));
273 int np0 = e->GetBasis(0)->GetNumPoints();
274 int np1 = e->GetBasis(1)->GetNumPoints();
275 int np2 = e->GetBasis(2)->GetNumPoints();
276 int np =
max(np0,
max(np1, np2));
284 ASSERTL0(
false,
"Points not known");
288 ppe.push_back(newpoints);
289 newtotpoints += newpoints;
293 bool standard =
true;
295 if (prismorient.count(i))
300 e->GetSimplexEquiSpacedConnectivity(conn, standard);
305 if ((prevNcoeffs != e->GetNcoeffs()) ||
306 (prevNpoints != e->GetTotPoints()))
308 prevNcoeffs = e->GetNcoeffs();
309 prevNpoints = e->GetTotPoints();
311 e->GetSimplexEquiSpacedConnectivity(conn);
315 for (
int j = 0; j < conn.size(); ++j)
317 newconn[j] = conn[j] + cnt;
320 ptsConn.push_back(newconn);
324 nfields =
m_f->m_variables.size();
328 for (
int i = 0; i < nfields + coordim; ++i)
334 for (
int i = 0; i < coordim; ++i)
339 for (
int i = coordim; i < 3; ++i)
344 m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
348 for (
int n = 0; n < coordim; ++n)
352 for (
int i = 0; i < nel; ++i)
354 m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
355 coords[n] + cnt, tmp = pts[n] + cnt1);
357 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
361 for (
int n = 0; n <
m_f->m_variables.size(); ++n)
366 if (
m_config[
"modalenergy"].m_beenSet &&
370 for (
int i = 0; i < nel; ++i)
374 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
380 for (
int i = 0; i < nel; ++i)
382 m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
383 phys + cnt, tmp = pts[coordim + n] + cnt1);
385 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
391 coordim,
m_f->m_variables, pts);
401 else if (shapedim == 3)
405 m_f->m_fieldPts->SetConnectivity(ptsConn);
412 m_f->m_fieldPts->SetPointsPerElement(ppe);
415 m_f->m_exp = vector<MultiRegions::ExpListSharedPtr>();
421 int nel =
m_f->m_exp[0]->GetPlane(0)->GetExpSize();
422 int nPlanes =
m_f->m_exp[0]->GetZIDs().size();
423 int npts =
m_f->m_fieldPts->GetNpoints();
424 int nptsPerPlane = npts / nPlanes;
427 if (
m_f->m_exp[0]->GetHomogeneousBasis()->GetBasisType() ==
429 m_f->m_exp[0]->GetHomogeneousBasis()->GetPointsType() ==
434 m_f->m_fieldPts->GetPts(pts);
436 for (
int i = 0; i < pts.size(); i++)
443 if (
m_f->m_comm->GetColumnComm()->GetSize() == 1)
445 if (i == (coordim - 1))
449 pts[i], 1, extraPlane, 1);
459 int size =
m_f->m_comm->GetColumnComm()->GetSize();
460 int rank =
m_f->m_comm->GetColumnComm()->GetRank();
461 int fromRank = (rank + 1) % size;
462 int toRank = (rank == 0) ? size - 1 : rank - 1;
465 m_f->m_comm->GetColumnComm()->SendRecv(toRank, send, fromRank,
468 if (i == (coordim - 1) && (rank == size - 1))
471 extraPlane, 1, extraPlane, 1);
475 m_f->m_fieldPts->SetPts(newPts);
477 vector<Array<OneD, int>> oldConn;
479 m_f->m_fieldPts->GetConnectivity(oldConn);
480 vector<Array<OneD, int>> newConn = oldConn;
481 int connPerPlane = oldConn.size() / nPlanes;
482 for (
int i = 0; i < connPerPlane; i++)
485 for (
int j = 0; j < conn.size(); j++)
487 conn[j] = oldConn[i][j] + npts;
489 newConn.push_back(conn);
491 m_f->m_fieldPts->SetConnectivity(newConn);
494 npts += nptsPerPlane;
497 vector<Array<OneD, int>> oldConn;
498 vector<Array<OneD, int>> newConn;
500 m_f->m_fieldPts->GetConnectivity(oldConn);
505 for (
int i = 0; i < nel; ++i)
507 int nLines = oldConn[i].size() / 2;
512 for (
int n = 0; n < nLines; n++)
515 for (
int p = 0; p < nPlanes - 1; p++)
517 conn[cnt++] = oldConn[i + p * nel][2 * n + 0];
518 conn[cnt++] = oldConn[i + p * nel][2 * n + 1];
519 conn[cnt++] = oldConn[i + (p + 1) * nel][2 * n + 0];
521 conn[cnt++] = oldConn[i + p * nel][2 * n + 1];
522 conn[cnt++] = oldConn[i + (p + 1) * nel][2 * n + 0];
523 conn[cnt++] = oldConn[i + (p + 1) * nel][2 * n + 1];
526 newConn.push_back(conn);
533 for (
int i = 0; i < nel; ++i)
535 e =
m_f->m_exp[0]->GetPlane(0)->GetExp(i);
537 int np0 = e->GetBasis(0)->GetNumPoints();
538 int np1 = e->GetBasis(1)->GetNumPoints();
539 int np =
max(np0, np1);
540 maxN =
max(np, maxN);
547 for (
int i = 0; i < nel; ++i)
549 e =
m_f->m_exp[0]->GetPlane(0)->GetExp(i);
550 int np0 = e->GetBasis(0)->GetNumPoints();
551 int np1 = e->GetBasis(1)->GetNumPoints();
552 int np =
max(np0, np1);
553 switch (e->DetShapeType())
562 vId[cnt1] = e->GetGeom()->GetVid(0);
563 vId[cnt1 + np - 1] = e->GetGeom()->GetVid(1);
564 vId[cnt1 + newpoints - 1] = e->GetGeom()->GetVid(2);
570 for (
int n = 1; n < np - 1; n++)
573 edgeOrient = e->GetGeom()->GetEorient(0);
577 4 * nel + maxN * e->GetGeom()->GetEid(0) + n;
581 vId[cnt1 + np - 1 - n] =
582 4 * nel + maxN * e->GetGeom()->GetEid(0) + n;
586 edgeOrient = e->GetGeom()->GetEorient(1);
590 vId[cnt1 + np - 1 + edge1] =
591 4 * nel + maxN * e->GetGeom()->GetEid(1) + n;
596 vId[cnt1 + newpoints - 1 - edge1] =
597 4 * nel + maxN * e->GetGeom()->GetEid(1) + n;
601 edgeOrient = e->GetGeom()->GetEorient(2);
605 vId[cnt1 + newpoints - 1 - edge2] =
606 4 * nel + maxN * e->GetGeom()->GetEid(2) + n;
612 4 * nel + maxN * e->GetGeom()->GetEid(2) + n;
618 for (
int n = 1; n < np - 1; n++)
621 for (
int m = 1; m < np - n - 1; m++)
623 vId[cnt1 + edge2 + m] =
624 4 * nel + maxN * 4 * nel + cnt2;
637 vId[cnt1] = e->GetGeom()->GetVid(0);
638 vId[cnt1 + np - 1] = e->GetGeom()->GetVid(1);
639 vId[cnt1 + np * np - 1] = e->GetGeom()->GetVid(2);
640 vId[cnt1 + np * (np - 1)] = e->GetGeom()->GetVid(3);
644 for (
int n = 1; n < np - 1; n++)
647 edgeOrient = e->GetGeom()->GetEorient(0);
651 4 * nel + maxN * e->GetGeom()->GetEid(0) + n;
655 vId[cnt1 + np - 1 - n] =
656 4 * nel + maxN * e->GetGeom()->GetEid(0) + n;
660 edgeOrient = e->GetGeom()->GetEorient(1);
663 vId[cnt1 + np - 1 + n * np] =
664 4 * nel + maxN * e->GetGeom()->GetEid(1) + n;
668 vId[cnt1 + np * np - 1 - n * np] =
669 4 * nel + maxN * e->GetGeom()->GetEid(1) + n;
673 edgeOrient = e->GetGeom()->GetEorient(2);
676 vId[cnt1 + np * np - 1 - n] =
677 4 * nel + maxN * e->GetGeom()->GetEid(2) + n;
681 vId[cnt1 + np * (np - 1) + n] =
682 4 * nel + maxN * e->GetGeom()->GetEid(2) + n;
686 edgeOrient = e->GetGeom()->GetEorient(3);
689 vId[cnt1 + np * (np - 1) - n * np] =
690 4 * nel + maxN * e->GetGeom()->GetEid(3) + n;
695 4 * nel + maxN * e->GetGeom()->GetEid(3) + n;
700 for (
int n = 1; n < np - 1; n++)
702 for (
int m = 1; m < np - 1; m++)
704 vId[cnt1 + m * np + n] =
705 4 * nel + maxN * 4 * nel + cnt2;
714 ASSERTL0(
false,
"Points not known");
720 for (
int i = 0; i < nel; ++i)
722 int nTris = oldConn[i].size() / 3;
727 for (
int n = 0; n < nTris; n++)
730 int vId0 = vId[oldConn[i][3 * n + 0]];
731 int vId1 = vId[oldConn[i][3 * n + 1]];
732 int vId2 = vId[oldConn[i][3 * n + 2]];
736 if ((vId0 < vId1) && (vId0 < vId2))
750 else if ((vId1 < vId0) && (vId1 < vId2))
780 for (
int p = 0; p < nPlanes - 1; p++)
782 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[0]];
783 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[1]];
784 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[2]];
785 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[2]];
787 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[0]];
788 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[1]];
789 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[2]];
790 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[1]];
792 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[0]];
793 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[1]];
794 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[2]];
795 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[0]];
798 newConn.push_back(conn);
803 ASSERTL0(
false,
"Points type not supported.");
806 m_f->m_fieldPts->SetConnectivity(newConn);