46 #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");
63 f->m_setUpEquiSpacedFields =
true;
66 ConfigOption(
true,
"NotSet",
"Only process tetrahedral elements");
69 ConfigOption(
true,
"NotSet",
"Write output as modal energy");
85 if (
m_f->m_comm->TreatAsRankZero())
87 cout <<
"Interpolating fields to equispaced..." << endl;
90 int nel =
m_f->m_exp[0]->GetExpSize();
97 int coordim =
m_f->m_exp[0]->GetCoordim(0);
98 int shapedim =
m_f->m_exp[0]->GetExp(0)->GetShapeDimension();
99 int npts =
m_f->m_exp[0]->GetTotPoints();
103 bool homogeneous1D =
false;
104 if (
m_f->m_fielddef.size())
106 if (
m_f->m_fielddef[0]->m_numHomogeneousDir == 1)
110 homogeneous1D =
true;
112 else if (
m_f->m_fielddef[0]->m_numHomogeneousDir == 2)
114 ASSERTL0(
false,
"Homegeneous2D case not supported");
119 if (
m_f->m_session->DefinesSolverInfo(
"HOMOGENEOUS"))
121 std::string HomoStr =
m_f->m_session->GetSolverInfo(
"HOMOGENEOUS");
123 if ((HomoStr ==
"HOMOGENEOUS1D") || (HomoStr ==
"Homogeneous1D") ||
124 (HomoStr ==
"1D") || (HomoStr ==
"Homo1D"))
128 homogeneous1D =
true;
130 if ((HomoStr ==
"HOMOGENEOUS2D") || (HomoStr ==
"Homogeneous2D") ||
131 (HomoStr ==
"2D") || (HomoStr ==
"Homo2D"))
133 ASSERTL0(
false,
"Homegeneous2D case not supported");
140 int newtotpoints = 0;
148 map<int, StdRegions::Orientation> face0orient;
149 set<int> prismorient;
153 vector<std::string> fieldNames;
155 vector<Array<OneD, int> > ptsConn;
158 for (
int i = 0; i < nel; ++i)
160 e =
m_f->m_exp[0]->GetExp(i);
164 int fid = e->GetGeom()->GetFid(0);
165 if (face0orient.count(fid))
167 prismorient.insert(i);
187 for (
int i = 0; i < nel; ++i)
189 e =
m_f->m_exp[0]->GetExp(i);
192 int fid = e->GetGeom()->GetFid(2);
194 if (face0orient.count(fid))
209 prismorient.insert(i);
216 prismorient.insert(i);
223 for (
int i = 0; i < nel; ++i)
225 e =
m_f->m_exp[0]->GetExp(i);
228 if (
m_f->m_exp[0]->GetExp(i)->DetShapeType() !=
235 switch (e->DetShapeType())
239 int npoints0 = e->GetBasis(0)->GetNumPoints();
247 int np0 = e->GetBasis(0)->GetNumPoints();
248 int np1 = e->GetBasis(1)->GetNumPoints();
249 int np = max(np0, np1);
256 int np0 = e->GetBasis(0)->GetNumPoints();
257 int np1 = e->GetBasis(1)->GetNumPoints();
258 int np = max(np0, np1);
266 int np0 = e->GetBasis(0)->GetNumPoints();
267 int np1 = e->GetBasis(1)->GetNumPoints();
268 int np2 = e->GetBasis(2)->GetNumPoints();
269 int np = max(np0, max(np1, np2));
277 int np0 = e->GetBasis(0)->GetNumPoints();
278 int np1 = e->GetBasis(1)->GetNumPoints();
279 int np2 = e->GetBasis(2)->GetNumPoints();
280 int np = max(np0, max(np1, np2));
288 int np0 = e->GetBasis(0)->GetNumPoints();
289 int np1 = e->GetBasis(1)->GetNumPoints();
290 int np2 = e->GetBasis(2)->GetNumPoints();
291 int np = max(np0, max(np1, np2));
299 int np0 = e->GetBasis(0)->GetNumPoints();
300 int np1 = e->GetBasis(1)->GetNumPoints();
301 int np2 = e->GetBasis(2)->GetNumPoints();
302 int np = max(np0, max(np1, np2));
310 ASSERTL0(
false,
"Points not known");
314 ppe.push_back(newpoints);
315 newtotpoints += newpoints;
319 bool standard =
true;
321 if (prismorient.count(i))
326 e->GetSimplexEquiSpacedConnectivity(conn, standard);
331 if ((prevNcoeffs != e->GetNcoeffs()) ||
332 (prevNpoints != e->GetTotPoints()))
334 prevNcoeffs = e->GetNcoeffs();
335 prevNpoints = e->GetTotPoints();
337 e->GetSimplexEquiSpacedConnectivity(conn);
341 for (
int j = 0; j < conn.num_elements(); ++j)
343 newconn[j] = conn[j] + cnt;
346 ptsConn.push_back(newconn);
350 if (
m_f->m_fielddef.size())
352 nfields =
m_f->m_exp.size();
361 for (
int i = 0; i < nfields + coordim; ++i)
367 for (
int i = 0; i < coordim; ++i)
372 for (
int i = coordim; i < 3; ++i)
377 m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
379 int nq1 =
m_f->m_exp[0]->GetTotPoints();
385 m_f->m_exp[0]->GetCoords(x1, y1, z1);
389 for (
int n = 0; n < coordim; ++n)
393 for (
int i = 0; i < nel; ++i)
395 m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
396 coords[n] + cnt, tmp = pts[n] + cnt1);
398 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
402 if (
m_f->m_fielddef.size())
404 ASSERTL0(
m_f->m_fielddef[0]->m_fields.size() ==
m_f->m_exp.size(),
405 "More expansion defined than fields");
407 for (
int n = 0; n <
m_f->m_exp.size(); ++n)
412 if (
m_config[
"modalenergy"].m_beenSet)
415 for (
int i = 0; i < nel; ++i)
419 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
425 for (
int i = 0; i < nel; ++i)
427 m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
428 phys + cnt, tmp = pts[coordim + n] + cnt1);
430 cnt +=
m_f->m_exp[0]->GetExp(i)->GetTotPoints();
435 fieldNames.push_back(
m_f->m_fielddef[0]->m_fields[n]);
440 coordim, fieldNames, pts);
450 else if (shapedim == 3)
454 m_f->m_fieldPts->SetConnectivity(ptsConn);
464 int nel =
m_f->m_exp[0]->GetPlane(0)->GetExpSize();
465 int nPlanes =
m_f->m_exp[0]->GetZIDs().num_elements();
466 int npts =
m_f->m_fieldPts->GetNpoints();
467 int nptsPerPlane = npts / nPlanes;
470 if (
m_f->m_exp[0]->GetHomogeneousBasis()->GetBasisType() ==
472 m_f->m_exp[0]->GetHomogeneousBasis()->GetPointsType() ==
477 m_f->m_fieldPts->GetPts(pts);
479 for (
int i = 0; i < pts.num_elements(); i++)
486 if (
m_f->m_comm->GetColumnComm()->GetSize() == 1)
488 if ( i == (coordim-1))
492 pts[i], 1, extraPlane, 1);
502 int size =
m_f->m_comm->GetColumnComm()->GetSize();
503 int rank =
m_f->m_comm->GetColumnComm()->GetRank();
504 int fromRank = (rank + 1) % size;
505 int toRank = (rank == 0) ? size - 1 : rank - 1;
508 m_f->m_comm->GetColumnComm()->SendRecv(toRank, send, fromRank,
511 if (i == (coordim-1) && (rank == size - 1) )
514 extraPlane, 1, extraPlane, 1);
518 m_f->m_fieldPts->SetPts(newPts);
520 vector<Array<OneD, int> > oldConn;
522 m_f->m_fieldPts->GetConnectivity(oldConn);
523 vector<Array<OneD, int> > newConn = oldConn;
524 int connPerPlane = oldConn.size() / nPlanes;
525 for (
int i = 0; i < connPerPlane; i++)
528 for (
int j = 0; j < conn.num_elements(); j++)
530 conn[j] = oldConn[i][j] +
npts;
532 newConn.push_back(conn);
534 m_f->m_fieldPts->SetConnectivity(newConn);
537 npts += nptsPerPlane;
540 vector<Array<OneD, int> > oldConn;
541 vector<Array<OneD, int> > newConn;
543 m_f->m_fieldPts->GetConnectivity(oldConn);
548 for(
int i = 0; i < nel; ++i)
550 int nLines = oldConn[i].num_elements()/2;
555 for(
int n = 0; n<nLines; n++)
558 for (
int p = 0;
p<nPlanes-1;
p++)
560 conn[cnt++] = oldConn[i+
p *nel][2*n+0];
561 conn[cnt++] = oldConn[i+
p *nel][2*n+1];
562 conn[cnt++] = oldConn[i+(
p+1)*nel][2*n+0];
564 conn[cnt++] = oldConn[i+
p *nel][2*n+1];
565 conn[cnt++] = oldConn[i+(
p+1)*nel][2*n+0];
566 conn[cnt++] = oldConn[i+(
p+1)*nel][2*n+1];
569 newConn.push_back(conn);
576 for(
int i = 0; i < nel; ++i)
578 e =
m_f->m_exp[0]->GetPlane(0)->GetExp(i);
580 int np0 = e->GetBasis(0)->GetNumPoints();
581 int np1 = e->GetBasis(1)->GetNumPoints();
582 int np = max(np0,np1);
583 maxN = max(np, maxN);
590 for(
int i = 0; i < nel; ++i)
592 e =
m_f->m_exp[0]->GetPlane(0)->GetExp(i);
593 int np0 = e->GetBasis(0)->GetNumPoints();
594 int np1 = e->GetBasis(1)->GetNumPoints();
595 int np = max(np0,np1);
596 switch(e->DetShapeType())
604 vId[cnt1] = e->GetGeom()->GetVid(0);
605 vId[cnt1+np-1] = e->GetGeom()->GetVid(1);
606 vId[cnt1+newpoints-1] = e->GetGeom()->GetVid(2);
612 for (
int n = 1; n < np-1; n++)
615 edgeOrient = e->GetGeom()->GetEorient(0);
618 vId[cnt1+n] = 4*nel +
619 maxN*e->GetGeom()->GetEid(0) + n;
623 vId[cnt1+np-1-n] = 4*nel +
624 maxN*e->GetGeom()->GetEid(0) + n;
628 edgeOrient = e->GetGeom()->GetEorient(1);
632 vId[cnt1+np-1+edge1] = 4*nel +
633 maxN*e->GetGeom()->GetEid(1) + n;
638 vId[cnt1+newpoints-1-edge1] = 4*nel +
639 maxN*e->GetGeom()->GetEid(1) + n;
643 edgeOrient = e->GetGeom()->GetEorient(2);
647 vId[cnt1+newpoints-1-edge2] = 4*nel +
648 maxN*e->GetGeom()->GetEid(2) + n;
653 vId[cnt1+edge2] = 4*nel +
654 maxN*e->GetGeom()->GetEid(2) + n;
660 for (
int n = 1; n < np-1; n++)
662 for (
int m = 1; m < np-n-1; m++)
664 vId[cnt1+mode] = 4*nel + maxN*4*nel + cnt2;
677 vId[cnt1] = e->GetGeom()->GetVid(0);
678 vId[cnt1+np-1] = e->GetGeom()->GetVid(1);
679 vId[cnt1+np*np-1] = e->GetGeom()->GetVid(2);
680 vId[cnt1+np*(np-1)] = e->GetGeom()->GetVid(3);
684 for (
int n = 1; n < np-1; n++)
687 edgeOrient = e->GetGeom()->GetEorient(0);
690 vId[cnt1+n] = 4*nel +
691 maxN*e->GetGeom()->GetEid(0) + n;
695 vId[cnt1+np-1-n] = 4*nel +
696 maxN*e->GetGeom()->GetEid(0) + n;
700 edgeOrient = e->GetGeom()->GetEorient(1);
703 vId[cnt1+np-1+n*np] = 4*nel +
704 maxN*e->GetGeom()->GetEid(1) + n;
708 vId[cnt1+np*np-1-n*np] = 4*nel +
709 maxN*e->GetGeom()->GetEid(1) + n;
713 edgeOrient = e->GetGeom()->GetEorient(2);
716 vId[cnt1+np*np-1-n] = 4*nel +
717 maxN*e->GetGeom()->GetEid(2) + n;
721 vId[cnt1+np*(np-1)+n] = 4*nel +
722 maxN*e->GetGeom()->GetEid(2) + n;
726 edgeOrient = e->GetGeom()->GetEorient(3);
729 vId[cnt1+np*(np-1)-n*np] = 4*nel +
730 maxN*e->GetGeom()->GetEid(3) + n;
734 vId[cnt1+n*np] = 4*nel +
735 maxN*e->GetGeom()->GetEid(3) + n;
740 for (
int n = 1; n < np-1; n++)
742 for (
int m = 1; m < np-1; m++)
744 vId[cnt1+m*np+n] = 4*nel + maxN*4*nel + cnt2;
759 for(
int i = 0; i < nel; ++i)
761 int nTris = oldConn[i].num_elements()/3;
766 for(
int n = 0; n<nTris; n++)
769 int vId0 = vId[oldConn[i][3*n+0]];
770 int vId1 = vId[oldConn[i][3*n+1]];
771 int vId2 = vId[oldConn[i][3*n+2]];
775 if ( (vId0<vId1) && (vId0<vId2))
789 else if ( (vId1<vId0) && (vId1<vId2))
819 for (
int p = 0;
p<nPlanes-1;
p++)
821 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[0]];
822 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[1]];
823 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[2]];
824 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[2]];
826 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[0]];
827 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[1]];
828 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[2]];
829 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[1]];
831 conn[cnt2++] = oldConn[i+
p *nel][3*n+rot[0]];
832 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[1]];
833 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[2]];
834 conn[cnt2++] = oldConn[i+(
p+1)*nel][3*n+rot[0]];
837 newConn.push_back(conn);
842 ASSERTL0(
false,
"Points type not supported.");
845 m_f->m_fieldPts->SetConnectivity(newConn);
854 e =
m_f->m_exp[0]->GetExp(n);
856 switch (e->DetShapeType())
860 int np0 = e->GetBasis(0)->GetNumPoints();
861 int np1 = e->GetBasis(1)->GetNumPoints();
862 int np = max(np0, np1);
876 e->GetBasis(1)->GetBasisKey(), phys, Ba, Bb,
884 int np0 = e->GetBasis(0)->GetNumPoints();
885 int np1 = e->GetBasis(1)->GetNumPoints();
886 int np = max(np0, np1);
898 e->GetBasis(1)->GetBasisKey(), phys, Ba, Bb,
905 ASSERTL0(
false,
"Shape needs setting up");
map< string, ConfigOption > m_config
List of configuration values.
#define ASSERTL0(condition, msg)
void SetHomogeneousConnectivity(void)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
static Array< OneD, NekDouble > NullNekDouble1DArray
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Represents a command-line configuration option.
virtual void Process(po::variables_map &vm)
Write mesh to output file.
pair< ModuleType, string > ModuleKey
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Principle Orthogonal Functions .
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...
1D Evenly-spaced points using Fourier Fit
int getNumberOfCoefficients(int Na)
boost::shared_ptr< Field > FieldSharedPtr
int getNumberOfCoefficients(int Na, int Nb)
Principle Orthogonal Functions .
Defines a specification for a set of points.
void SetupEquiSpacedField(void)
boost::shared_ptr< Expansion > ExpansionSharedPtr
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
void GenOrthoModes(int n, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &coeffs)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
static PtsFieldSharedPtr NullPtsField
int getNumberOfCoefficients(int Na, int Nb)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Abstract base class for processing modules.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
virtual ~ProcessEquiSpacedOutput()
Describes the specification for a Basis.
ModuleFactory & GetModuleFactory()
FieldSharedPtr m_f
Field object.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space...