multi point element transformation in nginterface_v2

This commit is contained in:
Joachim Schoeberl 2009-07-06 08:16:02 +00:00
parent 5764ae448c
commit 654914c3e1
6 changed files with 532 additions and 292 deletions

View File

@ -24,7 +24,20 @@ public:
template <int DIM> DLL_HEADER int Ng_GetNElements (); template <int DIM>
DLL_HEADER int Ng_GetNElements ();
template <int DIM> DLL_HEADER Ng_Element Ng_GetElement (int nr); template <int DIM>
DLL_HEADER Ng_Element Ng_GetElement (int nr);
/// Curved Elements:
/// xi..... DIM_EL local coordinates
/// sxi ... step xi
/// x ..... DIM_SPACE global coordinates
/// dxdxi...DIM_SPACE x DIM_EL Jacobian matrix (row major storage)
template <int DIM_EL, int DIM_SPACE>
DLL_HEADER void Ng_MultiElementTransformation (int elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi);

View File

@ -1269,7 +1269,7 @@ inline void ScaledJacobiPolynomial (int n, S x, St t, double alpha, double beta,
Vector shapes; Vector shapes;
DenseMatrix dshapes; MatrixFixWidth<2> dshapes;
Array<Vec<3> > coefs; Array<Vec<3> > coefs;
const Element2d & el = mesh[elnr]; const Element2d & el = mesh[elnr];
@ -1359,7 +1359,7 @@ inline void ScaledJacobiPolynomial (int n, S x, St t, double alpha, double beta,
const Element2d & el = mesh[info.elnr]; const Element2d & el = mesh[info.elnr];
shapes.SetSize(info.ndof); shapes.SetSize(info.ndof);
shapes = 0; // shapes = 0;
if (rational && info.order >= 2) if (rational && info.order >= 2)
{ {
@ -1494,15 +1494,15 @@ inline void ScaledJacobiPolynomial (int n, S x, St t, double alpha, double beta,
void CurvedElements :: void CurvedElements ::
CalcElementDShapes (SurfaceElementInfo & info, const Point<2> & xi, DenseMatrix & dshapes) const CalcElementDShapes (SurfaceElementInfo & info, const Point<2> & xi, MatrixFixWidth<2> & dshapes) const
{ {
const Element2d & el = mesh[info.elnr]; const Element2d & el = mesh[info.elnr];
ELEMENT_TYPE type = el.GetType(); ELEMENT_TYPE type = el.GetType();
double lami[4]; double lami[4];
dshapes.SetSize(info.ndof,2); dshapes.SetSize(info.ndof);
dshapes = 0; // dshapes = 0;
// *testout << "calcelementdshapes, info.ndof = " << info.ndof << endl; // *testout << "calcelementdshapes, info.ndof = " << info.ndof << endl;
@ -1555,6 +1555,8 @@ inline void ScaledJacobiPolynomial (int n, S x, St t, double alpha, double beta,
case TRIG: case TRIG:
{ {
dshapes(0,0) = 1; dshapes(0,0) = 1;
dshapes(0,1) = 0.0;
dshapes(1,0) = 0.0;
dshapes(1,1) = 1; dshapes(1,1) = 1;
dshapes(2,0) = -1; dshapes(2,0) = -1;
dshapes(2,1) = -1; dshapes(2,1) = -1;
@ -1756,16 +1758,20 @@ inline void ScaledJacobiPolynomial (int n, S x, St t, double alpha, double beta,
} }
template <int DIM_SPACE>
void CurvedElements :: void CurvedElements ::
GetCoefficients (SurfaceElementInfo & info, Array<Vec<3> > & coefs) const GetCoefficients (SurfaceElementInfo & info, Array<Vec<DIM_SPACE> > & coefs) const
{ {
const Element2d & el = mesh[info.elnr]; const Element2d & el = mesh[info.elnr];
coefs.SetSize (info.ndof); coefs.SetSize (info.ndof);
// coefs = Vec<3> (0,0,0); // coefs = Vec<3> (0,0,0);
for (int i = 0; i < info.nv; i++) for (int i = 0; i < info.nv; i++)
coefs[i] = Vec<3> (mesh[el[i]]); {
Vec<3> hv(mesh[el[i]]);
for (int j = 0; j < DIM_SPACE; j++)
coefs[i](j) = hv(j);
}
if (info.order == 1) return; if (info.order == 1) return;
@ -1776,16 +1782,23 @@ inline void ScaledJacobiPolynomial (int n, S x, St t, double alpha, double beta,
int first = edgecoeffsindex[info.edgenrs[i]]; int first = edgecoeffsindex[info.edgenrs[i]];
int next = edgecoeffsindex[info.edgenrs[i]+1]; int next = edgecoeffsindex[info.edgenrs[i]+1];
for (int j = first; j < next; j++, ii++) for (int j = first; j < next; j++, ii++)
coefs[ii] = edgecoeffs[j]; for (int k = 0; k < DIM_SPACE; k++)
coefs[ii](k) = edgecoeffs[j](k);
} }
int first = facecoeffsindex[info.facenr]; int first = facecoeffsindex[info.facenr];
int next = facecoeffsindex[info.facenr+1]; int next = facecoeffsindex[info.facenr+1];
for (int j = first; j < next; j++, ii++) for (int j = first; j < next; j++, ii++)
coefs[ii] = facecoeffs[j]; for (int k = 0; k < DIM_SPACE; k++)
coefs[ii](k) = facecoeffs[j](k);
} }
template void CurvedElements ::
GetCoefficients<2> (SurfaceElementInfo & info, Array<Vec<2> > & coefs) const;
template void CurvedElements ::
GetCoefficients<3> (SurfaceElementInfo & info, Array<Vec<3> > & coefs) const;
@ -2799,7 +2812,7 @@ inline void ScaledJacobiPolynomial (int n, S x, St t, double alpha, double beta,
Vector shapes; Vector shapes;
DenseMatrix dshapes; MatrixFixWidth<2> dshapes;
Array<Vec<3> > coefs; Array<Vec<3> > coefs;
@ -2871,6 +2884,181 @@ inline void ScaledJacobiPolynomial (int n, S x, St t, double alpha, double beta,
} }
} }
template <int DIM_SPACE>
void CurvedElements ::
CalcMultiPointSurfaceTransformation (SurfaceElementIndex elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi)
{
if (mesh.coarsemesh)
{
const HPRefElement & hpref_el =
(*mesh.hpelements) [mesh[elnr].hp_elnr];
// xi umrechnen
double lami[4];
FlatVector vlami(4, lami);
ArrayMem<Point<2>, 50> coarse_xi (npts);
for (int pi = 0; pi < npts; pi++)
{
vlami = 0;
Point<2> hxi(xi[pi*sxi], xi[pi*sxi+1]);
mesh[elnr].GetShapeNew ( hxi, vlami);
Point<2> cxi(0,0);
for (int i = 0; i < hpref_el.np; i++)
for (int j = 0; j < 2; j++)
cxi(j) += hpref_el.param[i][j] * lami[i];
coarse_xi[pi] = cxi;
}
mesh.coarsemesh->GetCurvedElements().
CalcMultiPointSurfaceTransformation<DIM_SPACE> (hpref_el.coarse_elnr, npts,
&coarse_xi[0](0), &coarse_xi[1](0)-&coarse_xi[0](0),
x, sx, dxdxi, sdxdxi);
Mat<2,2> trans;
Mat<3,2> dxdxic;
if (dxdxi)
{
MatrixFixWidth<2> dlami(4);
dlami = 0;
for (int pi = 0; pi < npts; pi++)
{
Point<2> hxi(xi[pi*sxi], xi[pi*sxi+1]);
mesh[elnr].GetDShapeNew ( hxi, dlami);
trans = 0;
for (int k = 0; k < 2; k++)
for (int l = 0; l < 2; l++)
for (int i = 0; i < hpref_el.np; i++)
trans(l,k) += hpref_el.param[i][l] * dlami(i, k);
Mat<DIM_SPACE,2> hdxdxic, hdxdxi;
for (int k = 0; k < 2*DIM_SPACE; k++)
hdxdxic(k) = dxdxi[pi*sdxdxi+k];
hdxdxi = hdxdxic * trans;
for (int k = 0; k < 2*DIM_SPACE; k++)
dxdxi[pi*sdxdxi+k] = hdxdxi(k);
// dxdxic = (*dxdxi)[pi];
// (*dxdxi)[pi] = dxdxic * trans;
}
}
return;
}
Vector shapes;
MatrixFixWidth<2> dshapes;
Array<Vec<DIM_SPACE> > coefs;
const Element2d & el = mesh[elnr];
ELEMENT_TYPE type = el.GetType();
SurfaceElementInfo info;
info.elnr = elnr;
info.order = order;
switch (type)
{
case TRIG : info.nv = 3; break;
case QUAD : info.nv = 4; break;
case TRIG6: info.nv = 6; break;
default:
cerr << "undef element in CalcMultPointSurfaceTrao" << endl;
}
info.ndof = info.nv;
if (info.order > 1)
{
const MeshTopology & top = mesh.GetTopology();
top.GetSurfaceElementEdges (elnr+1, info.edgenrs);
for (int i = 0; i < info.edgenrs.Size(); i++)
info.edgenrs[i]--;
info.facenr = top.GetSurfaceElementFace (elnr+1)-1;
for (int i = 0; i < info.edgenrs.Size(); i++)
info.ndof += edgecoeffsindex[info.edgenrs[i]+1] - edgecoeffsindex[info.edgenrs[i]];
info.ndof += facecoeffsindex[info.facenr+1] - facecoeffsindex[info.facenr];
}
GetCoefficients (info, coefs);
if (x)
{
for (int j = 0; j < npts; j++)
{
Point<2> vxi(xi[j*sxi], xi[j*sxi+1]);
CalcElementShapes (info, vxi, shapes);
Point<DIM_SPACE> val = 0.0;
for (int i = 0; i < coefs.Size(); i++)
val += shapes(i) * coefs[i];
for (int k = 0; k < DIM_SPACE; k++)
x[j*sx+k] = val(k);
}
}
if (dxdxi)
{
for (int j = 0; j < npts; j++)
{
Point<2> vxi(xi[j*sxi], xi[j*sxi+1]);
CalcElementDShapes (info, vxi, dshapes);
Mat<DIM_SPACE,2> ds;
ds = 0.0;
for (int i = 0; i < coefs.Size(); i++)
for (int j = 0; j < DIM_SPACE; j++)
for (int k = 0; k < 2; k++)
ds(j,k) += dshapes(i,k) * coefs[i](j);
// (*dxdxi)[ip] = ds;
for (int k = 0; k < 2*DIM_SPACE; k++)
dxdxi[j*sdxdxi+k] = ds(k);
}
}
}
template void CurvedElements ::
CalcMultiPointSurfaceTransformation<2> (SurfaceElementIndex elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi);
template void CurvedElements ::
CalcMultiPointSurfaceTransformation<3> (SurfaceElementIndex elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi);
void CurvedElements :: void CurvedElements ::
CalcMultiPointElementTransformation (Array< Point<3> > * xi, ElementIndex elnr, CalcMultiPointElementTransformation (Array< Point<3> > * xi, ElementIndex elnr,
Array< Point<3> > * x, Array< Point<3> > * x,

View File

@ -123,6 +123,13 @@ public:
Array< Point<3> > * x, Array< Point<3> > * x,
Array< Mat<3,2> > * dxdxi); Array< Mat<3,2> > * dxdxi);
template <int DIM_SPACE>
void CalcMultiPointSurfaceTransformation (SurfaceElementIndex elnr, int n,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi);
void CalcMultiPointElementTransformation (Array< Point<3> > * xi, ElementIndex elnr, void CalcMultiPointElementTransformation (Array< Point<3> > * xi, ElementIndex elnr,
Array< Point<3> > * x, Array< Point<3> > * x,
Array< Mat<3,3> > * dxdxi); Array< Mat<3,3> > * dxdxi);
@ -200,8 +207,9 @@ private:
}; };
void CalcElementShapes (SurfaceElementInfo & elinfo, const Point<2> & xi, Vector & shapes) const; void CalcElementShapes (SurfaceElementInfo & elinfo, const Point<2> & xi, Vector & shapes) const;
void GetCoefficients (SurfaceElementInfo & elinfo, Array<Vec<3> > & coefs) const; template <int DIM_SPACE>
void CalcElementDShapes (SurfaceElementInfo & elinfo, const Point<2> & xi, DenseMatrix & dshapes) const; void GetCoefficients (SurfaceElementInfo & elinfo, Array<Vec<DIM_SPACE> > & coefs) const;
void CalcElementDShapes (SurfaceElementInfo & elinfo, const Point<2> & xi, MatrixFixWidth<2> & dshapes) const;
}; };

View File

@ -10,7 +10,6 @@
namespace netgen namespace netgen
{ {

View File

@ -929,38 +929,10 @@ void Ng_GetMultiElementTransformation (int ei, int n,
double * dxdxi, int sdxdxi) double * dxdxi, int sdxdxi)
{ {
if (mesh->GetDimension() == 2) if (mesh->GetDimension() == 2)
{ mesh->GetCurvedElements().CalcMultiPointSurfaceTransformation<2> (ei-1, n, xi, sxi, x, sx, dxdxi, sdxdxi);
for (int i = 0; i < n; i++)
{
Point<2> xl(xi[i*sxi], xi[i*sxi+1]);
Point<3> xg;
Mat<3,2> dx;
mesh->GetCurvedElements().CalcSurfaceTransformation (xl, ei-1, xg, dx);
if (x)
{
x[i*sx ] = xg(0);
x[i*sx+1] = xg(1);
}
if (dxdxi)
{
dxdxi[i*sdxdxi ] = dx(0,0);
dxdxi[i*sdxdxi+1] = dx(0,1);
dxdxi[i*sdxdxi+2] = dx(1,0);
dxdxi[i*sdxdxi+3] = dx(1,1);
}
}
}
else else
{
mesh->GetCurvedElements().CalcMultiPointElementTransformation (ei-1, n, xi, sxi, x, sx, dxdxi, sdxdxi); mesh->GetCurvedElements().CalcMultiPointElementTransformation (ei-1, n, xi, sxi, x, sx, dxdxi, sdxdxi);
} }
}

View File

@ -124,3 +124,63 @@ template <> DLL_HEADER Ng_Element Ng_GetElement<3> (int nr)
} }
template <>
DLL_HEADER void Ng_MultiElementTransformation<3,3> (int elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi)
{
mesh->GetCurvedElements().CalcMultiPointElementTransformation (elnr, npts, xi, sxi, x, sx, dxdxi, sdxdxi);
}
template <>
DLL_HEADER void Ng_MultiElementTransformation<2,2> (int elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi)
{
mesh->GetCurvedElements().CalcMultiPointSurfaceTransformation<2> (elnr, npts, xi, sxi, x, sx, dxdxi, sdxdxi);
}
template <>
DLL_HEADER void Ng_MultiElementTransformation<2,3> (int elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi)
{
mesh->GetCurvedElements().CalcMultiPointSurfaceTransformation<3> (elnr, npts, xi, sxi, x, sx, dxdxi, sdxdxi);
}
template <>
DLL_HEADER void Ng_MultiElementTransformation<1,2> (int elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi)
{
for (int ip = 0; ip < npts; ip++)
{
Point<3> xg;
Vec<3> dx;
mesh->GetCurvedElements().CalcSegmentTransformation (xi[ip*sxi], elnr, xg, dx);
if (x)
for (int i = 0; i < 2; i++)
x[ip*sx+i] = xg(i);
if (dxdxi)
for (int i=0; i<2; i++)
dxdxi[ip*sdxdxi+i] = dx(i);
}
}
template <>
DLL_HEADER void Ng_MultiElementTransformation<1,1> (int elnr, int npts,
const double * xi, int sxi,
double * x, int sx,
double * dxdxi, int sdxdxi)
{
cout << "1D not supported" << endl;
}