38 #include <boost/core/ignore_unused.hpp>
39 #include <boost/math/special_functions/fpclassify.hpp>
53 ModuleKey ProcessEquiSpacedOutput::className =
56 ProcessEquiSpacedOutput::create,
57 "Write data as equi-spaced output using simplices to represent the "
58 "data for connecting points");
64 ConfigOption(
true,
"0",
"Only process tetrahedral elements");
78 int nel =
m_f->m_exp[0]->GetExpSize();
82 int nfields =
m_f->m_variables.size();
86 for (
int i = 0; i < nfields + coordim; ++i)
90 vector<Array<OneD, int> > ptsConn;
94 coordim,
m_f->m_variables, pts);
96 m_f->m_fieldPts->SetConnectivity(ptsConn);
100 int coordim =
m_f->m_exp[0]->GetCoordim(0);
101 int shapedim =
m_f->m_exp[0]->GetExp(0)->GetShapeDimension();
102 int npts =
m_f->m_exp[0]->GetTotPoints();
106 bool homogeneous1D =
false;
107 if (
m_f->m_numHomogeneousDir == 1)
111 homogeneous1D =
true;
113 else if (
m_f->m_numHomogeneousDir == 2)
115 ASSERTL0(
false,
"Homegeneous2D case not supported");
120 int newtotpoints = 0;
128 map<int, StdRegions::Orientation> face0orient;
129 set<int> prismorient;
134 vector<Array<OneD, int> > ptsConn;
137 for (
int i = 0; i < nel; ++i)
139 e =
m_f->m_exp[0]->GetExp(i);
143 int fid = e->GetGeom()->GetFid(0);
144 if (face0orient.count(fid))
146 prismorient.insert(i);
166 for (
int i = 0; i < nel; ++i)
168 e =
m_f->m_exp[0]->GetExp(i);
171 int fid = e->GetGeom()->GetFid(2);
173 if (face0orient.count(fid))
188 prismorient.insert(i);
195 prismorient.insert(i);
202 for (
int i = 0; i < nel; ++i)
204 e =
m_f->m_exp[0]->GetExp(i);
207 if (
m_f->m_exp[0]->GetExp(i)->DetShapeType() !=
214 switch (e->DetShapeType())
218 int npoints0 = e->GetBasis(0)->GetNumPoints();
226 int np0 = e->GetBasis(0)->GetNumPoints();
227 int np1 = e->GetBasis(1)->GetNumPoints();
228 int np = max(np0, np1);
235 int np0 = e->GetBasis(0)->GetNumPoints();
236 int np1 = e->GetBasis(1)->GetNumPoints();
237 int np = max(np0, np1);
245 int np0 = e->GetBasis(0)->GetNumPoints();
246 int np1 = e->GetBasis(1)->GetNumPoints();
247 int np2 = e->GetBasis(2)->GetNumPoints();
248 int np = max(np0, max(np1, np2));
256 int np0 = e->GetBasis(0)->GetNumPoints();
257 int np1 = e->GetBasis(1)->GetNumPoints();
258 int np2 = e->GetBasis(2)->GetNumPoints();
259 int np = max(np0, max(np1, np2));
267 int np0 = e->GetBasis(0)->GetNumPoints();
268 int np1 = e->GetBasis(1)->GetNumPoints();
269 int np2 = e->GetBasis(2)->GetNumPoints();
270 int np = max(np0, max(np1, np2));
278 int np0 = e->GetBasis(0)->GetNumPoints();
279 int np1 = e->GetBasis(1)->GetNumPoints();
280 int np2 = e->GetBasis(2)->GetNumPoints();
281 int np = max(np0, max(np1, np2));
289 ASSERTL0(
false,
"Points not known");
293 ppe.push_back(newpoints);
294 newtotpoints += newpoints;
298 bool standard =
true;
300 if (prismorient.count(i))
305 e->GetSimplexEquiSpacedConnectivity(conn, standard);
310 if ((prevNcoeffs != e->GetNcoeffs()) ||
311 (prevNpoints != e->GetTotPoints()))
313 prevNcoeffs = e->GetNcoeffs();
314 prevNpoints = e->GetTotPoints();
316 e->GetSimplexEquiSpacedConnectivity(conn);
320 for (
int j = 0; j < conn.size(); ++j)
322 newconn[j] = conn[j] + cnt;
325 ptsConn.push_back(newconn);
329 nfields =
m_f->m_variables.size();
333 for (
int i = 0; i < nfields + coordim; ++i)
339 for (
int i = 0; i < coordim; ++i)
344 for (
int i = coordim; i < 3; ++i)
349 m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
353 for (
int n = 0; n < coordim; ++n)
357 for (
int i = 0; i < nel; ++i)
359 m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
360 coords[n] + cnt, tmp = pts[n] + cnt1);
362 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
366 for (
int n = 0; n <
m_f->m_variables.size(); ++n)
371 if (
m_config[
"modalenergy"].as<bool>())
374 for (
int i = 0; i < nel; ++i)
378 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
384 for (
int i = 0; i < nel; ++i)
386 m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
387 phys + cnt, tmp = pts[coordim + n] + cnt1);
389 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
395 coordim,
m_f->m_variables, pts);
405 else if (shapedim == 3)
409 m_f->m_fieldPts->SetConnectivity(ptsConn);
416 m_f->m_exp = vector<MultiRegions::ExpListSharedPtr>();
422 int nel =
m_f->m_exp[0]->GetPlane(0)->GetExpSize();
423 int nPlanes =
m_f->m_exp[0]->GetZIDs().size();
424 int npts =
m_f->m_fieldPts->GetNpoints();
425 int nptsPerPlane = npts / nPlanes;
428 if (
m_f->m_exp[0]->GetHomogeneousBasis()->GetBasisType() ==
430 m_f->m_exp[0]->GetHomogeneousBasis()->GetPointsType() ==
435 m_f->m_fieldPts->GetPts(pts);
437 for (
int i = 0; i < pts.size(); i++)
444 if (
m_f->m_comm->GetColumnComm()->GetSize() == 1)
446 if ( i == (coordim-1))
450 pts[i], 1, extraPlane, 1);
460 int size =
m_f->m_comm->GetColumnComm()->GetSize();
461 int rank =
m_f->m_comm->GetColumnComm()->GetRank();
462 int fromRank = (rank + 1) % size;
463 int toRank = (rank == 0) ? size - 1 : rank - 1;
466 m_f->m_comm->GetColumnComm()->SendRecv(toRank, send, fromRank,
469 if (i == (coordim-1) && (rank == size - 1) )
472 extraPlane, 1, extraPlane, 1);
476 m_f->m_fieldPts->SetPts(newPts);
478 vector<Array<OneD, int> > oldConn;
480 m_f->m_fieldPts->GetConnectivity(oldConn);
481 vector<Array<OneD, int> > newConn = oldConn;
482 int connPerPlane = oldConn.size() / nPlanes;
483 for (
int i = 0; i < connPerPlane; i++)
486 for (
int j = 0; j < conn.size(); j++)
488 conn[j] = oldConn[i][j] + npts;
490 newConn.push_back(conn);
492 m_f->m_fieldPts->SetConnectivity(newConn);
495 npts += nptsPerPlane;
498 vector<Array<OneD, int> > oldConn;
499 vector<Array<OneD, int> > newConn;
501 m_f->m_fieldPts->GetConnectivity(oldConn);
506 for(
int i = 0; i < nel; ++i)
508 int nLines = oldConn[i].size()/2;
513 for(
int n = 0; n<nLines; n++)
516 for (
int p = 0;
p<nPlanes-1;
p++)
518 conn[cnt++] = oldConn[i+
p *nel][2*n+0];
519 conn[cnt++] = oldConn[i+
p *nel][2*n+1];
520 conn[cnt++] = oldConn[i+(
p+1)*nel][2*n+0];
522 conn[cnt++] = oldConn[i+
p *nel][2*n+1];
523 conn[cnt++] = oldConn[i+(
p+1)*nel][2*n+0];
524 conn[cnt++] = oldConn[i+(
p+1)*nel][2*n+1];
527 newConn.push_back(conn);
534 for(
int i = 0; i < nel; ++i)
536 e =
m_f->m_exp[0]->GetPlane(0)->GetExp(i);
538 int np0 = e->GetBasis(0)->GetNumPoints();
539 int np1 = e->GetBasis(1)->GetNumPoints();
540 int np = max(np0,np1);
541 maxN = max(np, maxN);
548 for(
int i = 0; i < nel; ++i)
550 e =
m_f->m_exp[0]->GetPlane(0)->GetExp(i);
551 int np0 = e->GetBasis(0)->GetNumPoints();
552 int np1 = e->GetBasis(1)->GetNumPoints();
553 int np = max(np0,np1);
554 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);
576 vId[cnt1+n] = 4*nel +
577 maxN*e->GetGeom()->GetEid(0) + n;
581 vId[cnt1+np-1-n] = 4*nel +
582 maxN*e->GetGeom()->GetEid(0) + n;
586 edgeOrient = e->GetGeom()->GetEorient(1);
590 vId[cnt1+np-1+edge1] = 4*nel +
591 maxN*e->GetGeom()->GetEid(1) + n;
596 vId[cnt1+newpoints-1-edge1] = 4*nel +
597 maxN*e->GetGeom()->GetEid(1) + n;
601 edgeOrient = e->GetGeom()->GetEorient(2);
605 vId[cnt1+newpoints-1-edge2] = 4*nel +
606 maxN*e->GetGeom()->GetEid(2) + n;
611 vId[cnt1+edge2] = 4*nel +
612 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] = 4*nel + maxN*4*nel + cnt2;
635 vId[cnt1] = e->GetGeom()->GetVid(0);
636 vId[cnt1+np-1] = e->GetGeom()->GetVid(1);
637 vId[cnt1+np*np-1] = e->GetGeom()->GetVid(2);
638 vId[cnt1+np*(np-1)] = e->GetGeom()->GetVid(3);
642 for (
int n = 1; n < np-1; n++)
645 edgeOrient = e->GetGeom()->GetEorient(0);
648 vId[cnt1+n] = 4*nel +
649 maxN*e->GetGeom()->GetEid(0) + n;
653 vId[cnt1+np-1-n] = 4*nel +
654 maxN*e->GetGeom()->GetEid(0) + n;
658 edgeOrient = e->GetGeom()->GetEorient(1);
661 vId[cnt1+np-1+n*np] = 4*nel +
662 maxN*e->GetGeom()->GetEid(1) + n;
666 vId[cnt1+np*np-1-n*np] = 4*nel +
667 maxN*e->GetGeom()->GetEid(1) + n;
671 edgeOrient = e->GetGeom()->GetEorient(2);
674 vId[cnt1+np*np-1-n] = 4*nel +
675 maxN*e->GetGeom()->GetEid(2) + n;
679 vId[cnt1+np*(np-1)+n] = 4*nel +
680 maxN*e->GetGeom()->GetEid(2) + n;
684 edgeOrient = e->GetGeom()->GetEorient(3);
687 vId[cnt1+np*(np-1)-n*np] = 4*nel +
688 maxN*e->GetGeom()->GetEid(3) + n;
692 vId[cnt1+n*np] = 4*nel +
693 maxN*e->GetGeom()->GetEid(3) + n;
698 for (
int n = 1; n < np-1; n++)
700 for (
int m = 1; m < np-1; m++)
702 vId[cnt1+m*np+n] = 4*nel + maxN*4*nel + cnt2;
717 for(
int i = 0; i < nel; ++i)
719 int nTris = oldConn[i].size()/3;
724 for(
int n = 0; n<nTris; n++)
727 int vId0 = vId[oldConn[i][3*n+0]];
728 int vId1 = vId[oldConn[i][3*n+1]];
729 int vId2 = vId[oldConn[i][3*n+2]];
733 if ( (vId0<vId1) && (vId0<vId2))
747 else if ( (vId1<vId0) && (vId1<vId2))
777 for (
int p = 0;
p<nPlanes-1;
p++)
779 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[0]];
780 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[1]];
781 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[2]];
782 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[2]];
784 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[0]];
785 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[1]];
786 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[2]];
787 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[1]];
789 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[0]];
790 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[1]];
791 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[2]];
792 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[0]];
795 newConn.push_back(conn);
800 ASSERTL0(
false,
"Points type not supported.");
803 m_f->m_fieldPts->SetConnectivity(newConn);
812 e =
m_f->m_exp[0]->GetExp(n);
814 switch (e->DetShapeType())
818 int np0 = e->GetBasis(0)->GetNumPoints();
819 int np1 = e->GetBasis(1)->GetNumPoints();
820 int np = max(np0, np1);
834 e->GetBasis(1)->GetBasisKey(), phys, Ba, Bb,
842 int np0 = e->GetBasis(0)->GetNumPoints();
843 int np1 = e->GetBasis(1)->GetNumPoints();
844 int np = max(np0, np1);
856 e->GetBasis(1)->GetBasisKey(), phys, Ba, Bb,
863 ASSERTL0(
false,
"Shape needs setting up");
#define ASSERTL0(condition, msg)
FieldSharedPtr m_f
Field object.
std::map< std::string, ConfigOption > m_config
List of configuration values.
virtual ~ProcessEquiSpacedOutput()
virtual void Process(po::variables_map &vm)
Write mesh to output file.
void SetHomogeneousConnectivity(void)
void GenOrthoModes(int n, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &coeffs)
Abstract base class for processing modules.
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.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
std::shared_ptr< Field > FieldSharedPtr
std::pair< ModuleType, std::string > ModuleKey
ModuleFactory & GetModuleFactory()
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb)
int getNumberOfCoefficients(int Na)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb)
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
@ eOrtho_A
Principle Orthogonal Functions .
@ eOrtho_B
Principle Orthogonal Functions .
@ eFourier
Fourier Expansion .
std::shared_ptr< Expansion > ExpansionSharedPtr
@ eDir1BwdDir2_Dir2BwdDir1
@ eDir1BwdDir1_Dir2BwdDir2
@ eDir1BwdDir2_Dir2FwdDir1
@ eDir1BwdDir1_Dir2FwdDir2
The above copyright notice and this permission notice shall be included.
static Array< OneD, NekDouble > NullNekDouble1DArray
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha - x.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Represents a command-line configuration option.