35#include <boost/core/ignore_unused.hpp>
49 std::bind(&
Expansion::CreateIndexMap, this, std::placeholders::_1),
50 std::string(
"ExpansionIndexMap")),
51 m_geom(pGeom), m_metricinfo(m_geom->GetGeomFactors()),
52 m_elementTraceLeft(-1), m_elementTraceRight(-1)
62 string type =
"regular";
69 err << nDim <<
"D " << type <<
" Jacobian not positive "
70 <<
"(element ID = " <<
m_geom->GetGlobalID() <<
") "
71 <<
"(first vertex ID = " <<
m_geom->GetVid(0) <<
")";
79 : StdExpansion(pSrc), m_indexMapManager(pSrc.m_indexMapManager),
80 m_geom(pSrc.m_geom), m_metricinfo(pSrc.m_metricinfo)
111 const NekDouble *data,
const std::vector<unsigned int> &nummodes,
112 const int nmodes_offset,
NekDouble *coeffs,
113 std::vector<LibUtilities::BasisType> &fromType)
119 const int edge,
const std::shared_ptr<Expansion> &EdgeExp,
127 const int edge,
const std::shared_ptr<Expansion> &EdgeExp,
134 const int face,
const std::shared_ptr<Expansion> &FaceExp,
146 v_DGDeriv(dir, inarray, EdgeExp, coeffs, outarray);
212 ASSERTL0(
false,
"Boundary Index Map not implemented yet.");
218 ASSERTL0(
false,
"Boundary Index Map not implemented yet.");
224 ASSERTL0(
false,
"Boundary Index Map not implemented yet.");
229 ASSERTL0(
false,
"Vertex Index Map not implemented yet.");
234 ASSERTL0(
false,
"The Index Map you are requiring "
235 "is not between the possible options.");
241 for (
int i = 0; i < map.size(); i++)
243 (*returnval)[i].index = map[i];
244 (*returnval)[i].sign =
sign[i];
257 std::map<int, NormalVector>::const_iterator x;
272 boost::ignore_unused(mkey);
282 "Geometric information is not set up");
286 unsigned int nint = (
unsigned int)(
m_ncoeffs - nbdry);
287 unsigned int exp_size[] = {nbdry, nint};
288 unsigned int nblks = 2;
290 nblks, nblks, exp_size, exp_size);
303 goto UseLocRegionsMatrix;
308 goto UseStdRegionsMatrix;
314 goto UseLocRegionsMatrix;
327 factor, Asubmat = mat->GetBlock(0, 0)));
331 one, Asubmat = mat->GetBlock(0, 1)));
335 factor, Asubmat = mat->GetBlock(1, 0)));
339 invfactor, Asubmat = mat->GetBlock(1, 1)));
362 for (i = 0; i < nbdry; ++i)
364 for (j = 0; j < nbdry; ++j)
366 (*A)(i, j) = mat(bmap[i], bmap[j]);
369 for (j = 0; j < nint; ++j)
371 (*B)(i, j) = mat(bmap[i], imap[j]);
375 for (i = 0; i < nint; ++i)
377 for (j = 0; j < nbdry; ++j)
379 (*C)(i, j) = mat(imap[i], bmap[j]);
382 for (j = 0; j < nint; ++j)
384 (*D)(i, j) = mat(imap[i], imap[j]);
393 (*A) = (*A) - (*B) * (*C);
420 boost::ignore_unused(mkey);
491 for (
int i = 0; i < ntraces; ++i)
503 for (
int i = 0; i < ndir; ++i)
520 for (
int j = 0; j < ntraces; ++j)
523 for (
int k = 0; k < nTracePts; ++k)
525 for (
int d = 0;
d < ndir; ++
d)
527 (*DerivMat[
d])(i, cnt + k) = Deriv[
d][traceids[j][k]];
544 const int expDim =
m_base.size();
551 for (
int i = 0; i < expDim; ++i)
553 CBasis[i] =
m_geom->GetXmap()->GetBasis(i);
554 nqGeom *= CBasis[i]->GetNumPoints();
556 m_base[i]->GetBasisKey().SamePoints(CBasis[i]->GetBasisKey());
565 for (
int i = 0; i <
m_geom->GetCoordim(); ++i)
567 m_geom->GetXmap()->BwdTrans(
m_geom->GetCoeffs(i), tmp[i]);
572 for (
int i = 0; i <
m_geom->GetCoordim(); ++i)
575 m_geom->GetXmap()->BwdTrans(
m_geom->GetCoeffs(i), tmpGeom);
582 CBasis[0]->GetPointsKey(), &tmpGeom[0],
583 m_base[0]->GetPointsKey(), &tmp[i][0]);
589 CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
590 &tmpGeom[0],
m_base[0]->GetPointsKey(),
591 m_base[1]->GetPointsKey(), &tmp[i][0]);
597 CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
598 CBasis[2]->GetPointsKey(), &tmpGeom[0],
600 m_base[2]->GetPointsKey(), &tmp[i][0]);
612 int shapedim = dfdir.size();
614 int nqtot = direction.size() / coordim;
616 for (
int j = 0; j < shapedim; j++)
619 for (
int k = 0; k < coordim; k++)
624 &direction[k * nqtot], 1, &dfdir[j][0], 1,
630 &direction[k * nqtot], 1, &dfdir[j][0], 1,
643 int nquad0, nquad1, nquad2;
649 nquad0 =
m_base[0]->GetNumPoints();
650 nquad1 =
m_base[1]->GetNumPoints();
651 nqtot = nquad0 * nquad1;
656 nquad0 =
m_base[0]->GetNumPoints();
657 nquad1 =
m_base[1]->GetNumPoints();
658 nquad2 =
m_base[2]->GetNumPoints();
659 nqtot = nquad0 * nquad1 * nquad2;
678 for (
int k = 0; k < coordim; k++)
680 StdRegions::VarCoeffMap::const_iterator MFdir =
681 varcoeffs.find(MMFCoeffs[5 * dir + k]);
682 tmp = MFdir->second.GetValue();
706 StdRegions::VarCoeffMap::const_iterator MFdir =
707 varcoeffs.find(MMFCoeffs[5 * dir + indxDiv]);
729 StdRegions::VarCoeffMap::const_iterator MFdir =
730 varcoeffs.find(MMFCoeffs[5 * dir + indxMag]);
739 boost::ignore_unused(r_bnd, matrixType);
747 boost::ignore_unused(r_bnd);
753 const NekDouble *data,
const std::vector<unsigned int> &nummodes,
754 const int nmodes_offset,
NekDouble *coeffs,
755 std::vector<LibUtilities::BasisType> &fromType)
757 boost::ignore_unused(data, nummodes, nmodes_offset, coeffs, fromType);
762 const int edge,
const std::shared_ptr<Expansion> &EdgeExp,
766 boost::ignore_unused(edge, EdgeExp, Fx, Fy, outarray);
771 const int edge,
const std::shared_ptr<Expansion> &EdgeExp,
774 boost::ignore_unused(edge, EdgeExp, Fn, outarray);
779 const int face,
const std::shared_ptr<Expansion> &FaceExp,
782 boost::ignore_unused(face, FaceExp, Fn, outarray);
792 boost::ignore_unused(dir, inarray, EdgeExp, coeffs, outarray);
799 boost::ignore_unused(vec);
801 "This function is only valid for "
802 "shape expansion in LocalRegions, not parant class");
811 boost::ignore_unused(
factors, d0factors, d1factors);
813 "This function is only valid for "
814 "shape expansion in LocalRegions, not parant class");
819 boost::ignore_unused(trace);
827 boost::ignore_unused(dir, inarray, outarray);
834 boost::ignore_unused(trace, outarray);
836 "Method does not exist for this shape or library");
844 boost::ignore_unused(trace, TraceExp, inarray, outarray, orient);
846 "Method does not exist for this shape or library");
851 boost::ignore_unused(edge, outarray);
853 "Method does not exist for this shape or library");
860 boost::ignore_unused(orient, idmap, nq0, nq1);
862 "Method does not exist for this shape or library");
867 boost::ignore_unused(traceid, exp);
869 "Method does not exist for this shape or library");
874 boost::ignore_unused(
id);
875 ASSERTL0(
false,
"Cannot compute trace normal for this expansion.");
886 boost::ignore_unused(normal);
892 boost::ignore_unused(edge);
900 boost::ignore_unused(trace, primCoeffs, inoutmat);
908 boost::ignore_unused(traceid, primCoeffs, incoeffs, coeffs);
914 boost::ignore_unused(traceid, h,
p);
919 const int nbnd)
const
923 "m_elmtBndNormDirElmtLen normal not computed.");
931 boost::ignore_unused(dir, inarray, outarray);
#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
std::map< int, NormalVector > m_traceNormals
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
virtual void v_SetPhysNormals(Array< OneD, const NekDouble > &normal)
SpatialDomains::GeometrySharedPtr GetGeom() const
virtual void v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
SpatialDomains::GeometrySharedPtr m_geom
void GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
void DropLocMatrix(const LocalRegions::MatrixKey &mkey)
const SpatialDomains::GeomFactorsSharedPtr & GetMetricInfo() const
DNekMatSharedPtr BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
virtual void v_DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp)
void ComputeLaplacianMetric()
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_AddRobinTraceContribution(const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
virtual void v_GetTraceQFactors(const int trace, Array< OneD, NekDouble > &outarray)
virtual void v_ComputeTraceNormal(const int id)
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Array< OneD, NekDouble > GetMFMag(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
virtual DNekMatSharedPtr v_BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
virtual DNekScalMatSharedPtr v_GetLocMatrix(const LocalRegions::MatrixKey &mkey)
virtual void v_AddFaceNormBoundaryInt(const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_DGDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
virtual DNekMatSharedPtr v_BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen(const int nbnd) const
virtual void v_AddRobinMassMatrix(const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
std::map< int, ExpansionWeakPtr > m_traceExp
virtual void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
void ComputeQuadratureMetric()
void AddEdgeNormBoundaryInt(const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddFaceNormBoundaryInt(const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_DropLocMatrix(const LocalRegions::MatrixKey &mkey)
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec)
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
virtual StdRegions::Orientation v_GetTraceOrient(int trace)
void DGDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
void StdDerivBaseOnTraceMat(Array< OneD, DNekMatSharedPtr > &DerivMat)
void ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
virtual void v_MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1=-1)
void NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
virtual void v_ComputeLaplacianMetric()
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals()
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
virtual void v_AddEdgeNormBoundaryInt(const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
virtual void v_SetUpPhysNormals(const int id)
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Array< OneD, NekDouble > GetMFDiv(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
const NormalVector & GetTraceNormal(const int id)
Array< OneD, NekDouble > GetMF(const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
virtual void v_GetTracePhysVals(const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
DNekMatSharedPtr BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec)
StdRegions::Orientation GetIndexOrientation() const
IndexMapType GetIndexMapType() const
int GetIndexEntity() const
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
int GetTraceNumPoints(const int i) const
This function returns the number of quadrature points belonging to the i-th trace.
int NumBndryCoeffs(void) const
const LibUtilities::PointsKeyVector GetPointsKeys() const
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
int GetNtraces() const
Returns the number of trace elements connected to this element.
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
void StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis,...
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,...
std::vector< PointsKey > PointsKeyVector
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eNoGeomType
No type defined.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< Geometry > GeometrySharedPtr
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
std::vector< double > d(NPUPPER *NPUPPER)
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
static DNekMatSharedPtr NullDNekMatSharedPtr
static Array< OneD, NekDouble > NullNekDouble1DArray
static DNekScalMatSharedPtr NullDNekScalMatSharedPtr
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 Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
void Vdiv(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 Zero(int n, T *x, const int incx)
Zero vector.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)