Modified the ASM-Integrand interface by introducing the class FiniteElement (with subclass MxFiniteElement for mixed problems) incapsulating the various finite element quantities evaluated at the integration point. The number of evalInt interface are thereby also reduced, making the code (hopefully) more readable, etc.

git-svn-id: http://svn.sintef.no/trondheim/IFEM/trunk@880 e10b68d5-8a6e-419e-a041-bce267b0401d
This commit is contained in:
kmo
2011-03-29 17:27:23 +00:00
committed by Knut Morten Okstad
parent 087c8bce35
commit 1b4e69e5db
31 changed files with 755 additions and 789 deletions

View File

@@ -17,6 +17,7 @@
#include "ASMs3D.h"
#include "TimeDomain.h"
#include "FiniteElement.h"
#include "GlobalIntegral.h"
#include "IntegrandBase.h"
#include "CoordinateMapping.h"
@@ -1132,10 +1133,9 @@ bool ASMs3D::integrate (Integrand& integrand,
const int nel1 = n1 - p1 + 1;
const int nel2 = n2 - p2 + 1;
Vector N(p1*p2*p3), Navg;
Matrix dNdu, dNdX, Xnod, Jac;
Matrix3D d2Ndu2, d2NdX2, Hess;
double h = 0.0;
FiniteElement fe(p1*p2*p3);
Matrix dNdu, Xnod, Jac;
Matrix3D d2Ndu2, Hess;
Vec4 X;
@@ -1157,13 +1157,13 @@ bool ASMs3D::integrate (Integrand& integrand,
// Compute characteristic element length, if needed
if (integrand.getIntegrandType() == 2)
h = getElmSize(p1,p2,p3,Xnod);
fe.h = getElmSize(p1,p2,p3,Xnod);
else if (integrand.getIntegrandType() == 3)
{
// --- Compute average value of basis functions over the element -----
Navg.resize(p1*p2*p3,true);
fe.Navg.resize(p1*p2*p3,true);
double vol = 0.0;
int ip = (((i3-p3)*nGauss*nel2 + i2-p2)*nGauss*nel1 + i1-p1)*nGauss;
for (int k = 0; k < nGauss; k++, ip += nGauss*(nel2-1)*nGauss*nel1)
@@ -1171,20 +1171,20 @@ bool ASMs3D::integrate (Integrand& integrand,
for (int i = 0; i < nGauss; i++, ip++)
{
// Fetch basis function derivatives at current integration point
extractBasis(spline[ip],N,dNdu);
extractBasis(spline[ip],fe.N,dNdu);
// Compute Jacobian determinant of coordinate mapping
// and multiply by weight of current integration point
double detJac = utl::Jacobian(Jac,fe.dNdX,Xnod,dNdu,false);
double weight = 0.125*dV*wg[i]*wg[j]*wg[k];
double detJxW = utl::Jacobian(Jac,dNdX,Xnod,dNdu,false)*weight;
// Numerical quadrature
Navg += N*detJxW;
vol += detJxW;
fe.Navg.add(fe.N,detJac*weight);
vol += detJac*weight;
}
// Divide by element volume
Navg /= vol;
fe.Navg /= vol;
}
else if (integrand.getIntegrandType() == 4)
@@ -1218,50 +1218,39 @@ bool ASMs3D::integrate (Integrand& integrand,
for (int j = 0; j < nGauss; j++, ip += nGauss*(nel1-1))
for (int i = 0; i < nGauss; i++, ip++)
{
// Weight of current integration point
double weight = 0.125*dV*wg[i]*wg[j]*wg[k];
// Parameter values of current integration point
fe.u = gpar[0](i+1,i1-p1+1);
fe.v = gpar[1](j+1,i2-p2+1);
fe.w = gpar[2](k+1,i3-p3+1);
// Fetch basis function derivatives at current integration point
if (integrand.getIntegrandType() == 2)
extractBasis(spline2[ip],N,dNdu,d2Ndu2);
extractBasis(spline2[ip],fe.N,dNdu,d2Ndu2);
else
extractBasis(spline[ip],N,dNdu);
extractBasis(spline[ip],fe.N,dNdu);
// Compute Jacobian inverse of coordinate mapping and derivatives
double detJ = utl::Jacobian(Jac,dNdX,Xnod,dNdu);
if (detJ == 0.0) continue; // skip singular points
fe.detJxW = utl::Jacobian(Jac,fe.dNdX,Xnod,dNdu);
if (fe.detJxW == 0.0) continue; // skip singular points
// Compute Hessian of coordinate mapping and 2nd order derivatives
if (integrand.getIntegrandType() == 2)
if (!utl::Hessian(Hess,d2NdX2,Jac,Xnod,d2Ndu2,dNdu))
if (!utl::Hessian(Hess,fe.d2NdX2,Jac,Xnod,d2Ndu2,dNdu))
return false;
#if SP_DEBUG > 4
std::cout <<"\niel, ip = "<< iel <<" "<< ip
<<"\nN ="<< N <<"dNdX ="<< dNdX << std::endl;
<<"\nN ="<< fe.N <<"dNdX ="<< fe.dNdX << std::endl;
#endif
// Cartesian coordinates of current integration point
X = Xnod * N;
X = Xnod * fe.N;
X.t = time.t;
// Evaluate the integrand and accumulate element contributions
bool ok = false;
switch (integrand.getIntegrandType())
{
case 2:
ok = integrand.evalInt(elmInt,time,detJ*weight,
N,dNdX,d2NdX2,X,h);
break;
case 3:
ok = integrand.evalInt(elmInt,time,detJ*weight,
N,dNdX,Navg,X);
break;
default:
ok = integrand.evalInt(elmInt,time,detJ*weight,
N,dNdX,X);
}
if (!ok) return false;
fe.detJxW *= 0.125*dV*wg[i]*wg[j]*wg[k];
if (!integrand.evalInt(elmInt,fe,time,X))
return false;
}
// Finalize the element quantities
@@ -1344,8 +1333,8 @@ bool ASMs3D::integrate (Integrand& integrand, int lIndex,
const int nel1 = n1 - p1 + 1;
const int nel2 = n2 - p2 + 1;
Vector N(p1*p2*p3);
Matrix dNdu, dNdX, Xnod, Jac;
FiniteElement fe(p1*p2*p3);
Matrix dNdu, Xnod, Jac;
Vec4 X;
Vec3 normal;
@@ -1405,24 +1394,22 @@ bool ASMs3D::integrate (Integrand& integrand, int lIndex,
for (int j = 0; j < nGauss; j++, ip += nGauss*(nf1-1))
for (int i = 0; i < nGauss; i++, ip++)
{
// Weight of current integration point
double weight = 0.25*dA*wg[i]*wg[j];
// Fetch basis function derivatives at current integration point
extractBasis(spline[ip],N,dNdu);
extractBasis(spline[ip],fe.N,dNdu);
// Compute basis function derivatives and the face normal
double dS = utl::Jacobian(Jac,normal,dNdX,Xnod,dNdu,t1,t2);
if (dS == 0.0) continue; // skip singular points
fe.detJxW = utl::Jacobian(Jac,normal,fe.dNdX,Xnod,dNdu,t1,t2);
if (fe.detJxW == 0.0) continue; // skip singular points
if (faceDir < 0) normal *= -1.0;
// Cartesian coordinates of current integration point
X = Xnod * N;
X = Xnod * fe.N;
X.t = time.t;
// Evaluate the integrand and accumulate element contributions
if (!integrand.evalBou(elmInt,time,dS*weight,N,dNdX,X,normal))
fe.detJxW *= 0.25*dA*wg[i]*wg[j];
if (!integrand.evalBou(elmInt,fe,time,X,normal))
return false;
}
@@ -1496,10 +1483,10 @@ bool ASMs3D::integrateEdge (Integrand& integrand, int lEdge,
const int p2 = svol->order(1);
const int p3 = svol->order(2);
Vector N(p1*p2*p3);
Matrix dNdu, dNdX, Xnod, Jac;
FiniteElement fe(p1*p2*p3);
Matrix dNdu, Xnod, Jac;
Vec4 X;
Vec3 tangent;
Vec3 tang;
// === Assembly loop over all elements on the patch edge =====================
@@ -1564,22 +1551,20 @@ bool ASMs3D::integrateEdge (Integrand& integrand, int lEdge,
LocalIntegral* elmInt = 0;
for (int i = 0; i < nGauss; i++, ip++)
{
// Weight of current integration point
double weight = 0.5*dS*wg[i];
// Fetch basis function derivatives at current integration point
extractBasis(spline[ip],N,dNdu);
extractBasis(spline[ip],fe.N,dNdu);
// Compute basis function derivatives and the edge tangent
double dT = utl::Jacobian(Jac,tangent,dNdX,Xnod,dNdu,1+(lEdge-1)/4);
if (dT == 0.0) continue; // skip singular points
// Compute basis function derivatives and the edge tang
fe.detJxW = utl::Jacobian(Jac,tang,fe.dNdX,Xnod,dNdu,1+(lEdge-1)/4);
if (fe.detJxW == 0.0) continue; // skip singular points
// Cartesian coordinates of current integration point
X = Xnod * N;
X = Xnod * fe.N;
X.t = time.t;
// Evaluate the integrand and accumulate element contributions
if (!integrand.evalBou(elmInt,time,dT*weight,N,dNdX,X,tangent))
fe.detJxW *= 0.5*dS*wg[i];
if (!integrand.evalBou(elmInt,fe,time,X,tang))
return false;
}
@@ -1604,12 +1589,14 @@ bool ASMs3D::getGridParameters (RealArray& prm, int dir, int nSegPerSpan) const
}
DoubleIter uit = svol->basis(dir).begin();
double ucurr, uprev = *(uit++);
double ucurr = 0.0, uprev = *(uit++);
while (uit != svol->basis(dir).end())
{
ucurr = *(uit++);
if (ucurr > uprev)
for (int i = 0; i < nSegPerSpan; i++)
if (nSegPerSpan == 1)
prm.push_back(uprev);
else for (int i = 0; i < nSegPerSpan; i++)
{
double xg = (double)(2*i-nSegPerSpan)/(double)nSegPerSpan;
prm.push_back(0.5*(ucurr*(1.0+xg) + uprev*(1.0-xg)));
@@ -1617,7 +1604,8 @@ bool ASMs3D::getGridParameters (RealArray& prm, int dir, int nSegPerSpan) const
uprev = ucurr;
}
prm.push_back(svol->basis(dir).endparam());
if (ucurr > prm.back())
prm.push_back(ucurr);
return true;
}
@@ -1745,7 +1733,7 @@ bool ASMs3D::evalSolution (Matrix& sField, const Integrand& integrand,
if (npe)
{
// Compute parameter values of the result sampling points
DoubleVec gpar[3];
RealArray gpar[3];
for (int dir = 0; dir < 3 && retVal; dir++)
retVal = this->getGridParameters(gpar[dir],dir,npe[dir]-1);