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:
@@ -13,6 +13,7 @@
|
||||
|
||||
#include "NonlinearElasticity.h"
|
||||
#include "MaterialBase.h"
|
||||
#include "FiniteElement.h"
|
||||
#include "Utilities.h"
|
||||
#include "Profiler.h"
|
||||
|
||||
@@ -75,15 +76,15 @@ void NonlinearElasticity::setMode (SIM::SolutionMode mode)
|
||||
}
|
||||
|
||||
|
||||
bool NonlinearElasticity::evalInt (LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
bool NonlinearElasticity::evalInt (LocalIntegral*& elmInt,
|
||||
const FiniteElement& fe,
|
||||
const Vec3& X) const
|
||||
{
|
||||
PROFILE3("NonlinearEl::evalInt");
|
||||
|
||||
// Evaluate the kinematic quantities, F and E, at this point
|
||||
Tensor F(nsd);
|
||||
if (!this->kinematics(dNdX,F,E))
|
||||
if (!this->kinematics(fe.dNdX,F,E))
|
||||
return false;
|
||||
|
||||
// Evaluate current tangent at this point, that is
|
||||
@@ -106,30 +107,30 @@ bool NonlinearElasticity::evalInt (LocalIntegral*& elmInt, double detJW,
|
||||
PROFILE4("stiff_TL_");
|
||||
if (nsd == 3)
|
||||
if (fullCmat)
|
||||
stiff_tl3d_(N.size(),detJW,dNdX.ptr(),F.ptr(),
|
||||
stiff_tl3d_(fe.N.size(),fe.detJxW,fe.dNdX.ptr(),F.ptr(),
|
||||
Cmat.ptr(),eKm->ptr());
|
||||
else
|
||||
stiff_tl3d_isoel_(N.size(),detJW,dNdX.ptr(),F.ptr(),
|
||||
stiff_tl3d_isoel_(fe.N.size(),fe.detJxW,fe.dNdX.ptr(),F.ptr(),
|
||||
Cmat(1,1),Cmat(1,2),Cmat(4,4),eKm->ptr());
|
||||
else if (nsd == 2)
|
||||
if (fullCmat)
|
||||
stiff_tl2d_(N.size(),detJW,dNdX.ptr(),F.ptr(),
|
||||
stiff_tl2d_(fe.N.size(),fe.detJxW,fe.dNdX.ptr(),F.ptr(),
|
||||
Cmat.ptr(),eKm->ptr());
|
||||
else
|
||||
stiff_tl2d_isoel_(N.size(),detJW,dNdX.ptr(),F.ptr(),
|
||||
stiff_tl2d_isoel_(fe.N.size(),fe.detJxW,fe.dNdX.ptr(),F.ptr(),
|
||||
Cmat(1,1),Cmat(1,2),Cmat(3,3),eKm->ptr());
|
||||
else if (nsd == 1)
|
||||
for (a = 1; a <= N.size(); a++)
|
||||
for (b = 1; b <= N.size(); b++)
|
||||
(*eKm)(a,b) += dNdX(a,1)*F(1,1)*Cmat(1,1)*F(1,1)*dNdX(b,1);
|
||||
for (a = 1; a <= fe.N.size(); a++)
|
||||
for (b = 1; b <= fe.N.size(); b++)
|
||||
(*eKm)(a,b) += fe.dNdX(a,1)*F(1,1)*Cmat(1,1)*F(1,1)*fe.dNdX(b,1);
|
||||
#else
|
||||
// This is too costly, but is basically what is done in the fortran routines
|
||||
PROFILE4("dNdX^t*F^t*C*F*dNdX");
|
||||
unsigned short int l, m, n;
|
||||
SymmTensor4 D(Cmat,nsd); // fourth-order material tensor
|
||||
Matrix& EM = *eKm;
|
||||
for (a = 1; a <= N.size(); a++)
|
||||
for (b = 1; b <= N.size(); b++)
|
||||
for (a = 1; a <= fe.N.size(); a++)
|
||||
for (b = 1; b <= fe.N.size(); b++)
|
||||
for (m = 1; m <= nsd; m++)
|
||||
for (n = 1; n <= nsd; n++)
|
||||
{
|
||||
@@ -138,39 +139,39 @@ bool NonlinearElasticity::evalInt (LocalIntegral*& elmInt, double detJW,
|
||||
for (j = 1; j <= nsd; j++)
|
||||
for (k = 1; k <= nsd; k++)
|
||||
for (l = 1; l <= nsd; l++)
|
||||
km += dNdX(a,i)*F(m,j)*D(i,j,k,l)*F(n,k)*dNdX(b,l);
|
||||
km += fe.dNdX(a,i)*F(m,j)*D(i,j,k,l)*F(n,k)*fe.dNdX(b,l);
|
||||
|
||||
EM(nsd*(a-1)+m,nsd*(b-1)+n) += km*detJW;
|
||||
EM(nsd*(a-1)+m,nsd*(b-1)+n) += km*fe.detJxW;
|
||||
}
|
||||
#endif
|
||||
}
|
||||
|
||||
if (eKg && haveStrains)
|
||||
// Integrate the geometric stiffness matrix
|
||||
this->formKG(*eKg,dNdX,S,detJW);
|
||||
this->formKG(*eKg,fe.dNdX,S,fe.detJxW);
|
||||
|
||||
if (eM)
|
||||
// Integrate the mass matrix
|
||||
this->formMassMatrix(*eM,N,X,detJW);
|
||||
this->formMassMatrix(*eM,fe.N,X,fe.detJxW);
|
||||
|
||||
if (iS && haveStrains)
|
||||
{
|
||||
// Integrate the internal forces
|
||||
Vector& ES = *iS;
|
||||
for (a = 1; a <= N.size(); a++)
|
||||
for (a = 1; a <= fe.N.size(); a++)
|
||||
for (k = 1; k <= nsd; k++)
|
||||
{
|
||||
double f = 0.0;
|
||||
for (i = 1; i <= nsd; i++)
|
||||
for (j = 1; j <= nsd; j++)
|
||||
f -= dNdX(a,i)*F(k,j)*S(i,j);
|
||||
ES(nsd*(a-1)+k) += f*detJW;
|
||||
f -= fe.dNdX(a,i)*F(k,j)*S(i,j);
|
||||
ES(nsd*(a-1)+k) += f*fe.detJxW;
|
||||
}
|
||||
}
|
||||
|
||||
if (eS)
|
||||
// Integrate the load vector due to gravitation and other body forces
|
||||
this->formBodyForce(*eS,N,X,detJW);
|
||||
this->formBodyForce(*eS,fe.N,X,fe.detJxW);
|
||||
|
||||
return this->getIntegralResult(elmInt);
|
||||
}
|
||||
|
||||
@@ -49,12 +49,9 @@ public:
|
||||
|
||||
//! \brief Evaluates the integrand at an interior point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N Basis function values
|
||||
//! \param[in] dNdX Basis function gradients
|
||||
//! \param[in] fe Finite element data of current integration point
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
virtual bool evalInt(LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
virtual bool evalInt(LocalIntegral*& elmInt, const FiniteElement& fe,
|
||||
const Vec3& X) const;
|
||||
|
||||
//! \brief Evaluates the secondary solution at a result point.
|
||||
|
||||
@@ -13,6 +13,7 @@
|
||||
|
||||
#include "NonlinearElasticity.h"
|
||||
#include "MaterialBase.h"
|
||||
#include "FiniteElement.h"
|
||||
#include "ElmMats.h"
|
||||
#include "Tensor.h"
|
||||
#include "Vec3Oper.h"
|
||||
@@ -82,15 +83,15 @@ void NonlinearElasticityTL::setMode (SIM::SolutionMode mode)
|
||||
}
|
||||
|
||||
|
||||
bool NonlinearElasticityTL::evalInt (LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
bool NonlinearElasticityTL::evalInt (LocalIntegral*& elmInt,
|
||||
const FiniteElement& fe,
|
||||
const Vec3& X) const
|
||||
{
|
||||
// Evaluate the deformation gradient, F, and the Green-Lagrange strains, E,
|
||||
// and compute the nonlinear strain-displacement matrix, B, from dNdX and F
|
||||
Tensor F(nsd);
|
||||
SymmTensor E(nsd), S(nsd);
|
||||
if (!this->kinematics(dNdX,F,E))
|
||||
if (!this->kinematics(fe.dNdX,F,E))
|
||||
return false;
|
||||
|
||||
// Evaluate the constitutive relation
|
||||
@@ -105,36 +106,36 @@ bool NonlinearElasticityTL::evalInt (LocalIntegral*& elmInt, double detJW,
|
||||
if (eKm)
|
||||
{
|
||||
// Integrate the material stiffness matrix
|
||||
CB.multiply(Cmat,Bmat).multiply(detJW); // CB = C*B*|J|*w
|
||||
eKm->multiply(Bmat,CB,true,false,true); // EK += B^T * CB
|
||||
CB.multiply(Cmat,Bmat).multiply(fe.detJxW); // CB = C*B*|J|*w
|
||||
eKm->multiply(Bmat,CB,true,false,true); // EK += B^T * CB
|
||||
}
|
||||
|
||||
if (eKg && lHaveStrains)
|
||||
// Integrate the geometric stiffness matrix
|
||||
this->formKG(*eKg,dNdX,S,detJW);
|
||||
this->formKG(*eKg,fe.dNdX,S,fe.detJxW);
|
||||
|
||||
if (eM)
|
||||
// Integrate the mass matrix
|
||||
this->formMassMatrix(*eM,N,X,detJW);
|
||||
this->formMassMatrix(*eM,fe.N,X,fe.detJxW);
|
||||
|
||||
if (iS && lHaveStrains)
|
||||
{
|
||||
// Integrate the internal forces
|
||||
S *= -detJW;
|
||||
S *= -fe.detJxW;
|
||||
if (!Bmat.multiply(S,*iS,true,true)) // ES -= B^T*S
|
||||
return false;
|
||||
}
|
||||
|
||||
if (eS)
|
||||
// Integrate the load vector due to gravitation and other body forces
|
||||
this->formBodyForce(*eS,N,X,detJW);
|
||||
this->formBodyForce(*eS,fe.N,X,fe.detJxW);
|
||||
|
||||
return this->getIntegralResult(elmInt);
|
||||
}
|
||||
|
||||
|
||||
bool NonlinearElasticityTL::evalBou (LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
bool NonlinearElasticityTL::evalBou (LocalIntegral*& elmInt,
|
||||
const FiniteElement& fe,
|
||||
const Vec3& X, const Vec3& normal) const
|
||||
{
|
||||
if (!tracFld)
|
||||
@@ -163,7 +164,7 @@ bool NonlinearElasticityTL::evalBou (LocalIntegral*& elmInt, double detJW,
|
||||
// Compute the deformation gradient, F
|
||||
Tensor F(nsd);
|
||||
SymmTensor dummy(0);
|
||||
if (!this->kinematics(dNdX,F,dummy)) return false;
|
||||
if (!this->kinematics(fe.dNdX,F,dummy)) return false;
|
||||
|
||||
// Compute its inverse and determinant, J
|
||||
double J = F.inverse();
|
||||
@@ -180,9 +181,9 @@ bool NonlinearElasticityTL::evalBou (LocalIntegral*& elmInt, double detJW,
|
||||
|
||||
// Integrate the force vector
|
||||
Vector& ES = *eS;
|
||||
for (size_t a = 1; a <= N.size(); a++)
|
||||
for (size_t a = 1; a <= fe.N.size(); a++)
|
||||
for (i = 1; i <= nsd; i++)
|
||||
ES(nsd*(a-1)+i) += T[i-1]*N(a)*detJW;
|
||||
ES(nsd*(a-1)+i) += T[i-1]*fe.N(a)*fe.detJxW;
|
||||
|
||||
return this->getIntegralResult(elmInt);
|
||||
}
|
||||
|
||||
@@ -44,26 +44,20 @@ public:
|
||||
|
||||
//! \brief Evaluates the integrand at an interior point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N Basis function values
|
||||
//! \param[in] dNdX Basis function gradients
|
||||
//! \param[in] fe Finite element data of current integration point
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
virtual bool evalInt(LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
virtual bool evalInt(LocalIntegral*& elmInt, const FiniteElement& fe,
|
||||
const Vec3& X) const;
|
||||
|
||||
//! \brief Evaluates the integrand at a boundary point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N Basis function values
|
||||
//! \param[in] dNdX Basis function gradients
|
||||
//! \param[in] fe Finite element data of current integration point
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
//! \param[in] normal Boundary normal vector at current integration point
|
||||
//!
|
||||
//! \details This method is reimplemented in this class to account for
|
||||
//! possibly with-rotated traction fields in the Total-Lagrangian setting.
|
||||
virtual bool evalBou(LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
virtual bool evalBou(LocalIntegral*& elmInt, const FiniteElement& fe,
|
||||
const Vec3& X, const Vec3& normal) const;
|
||||
|
||||
protected:
|
||||
|
||||
@@ -13,6 +13,7 @@
|
||||
|
||||
#include "NonlinearElasticityUL.h"
|
||||
#include "MaterialBase.h"
|
||||
#include "FiniteElement.h"
|
||||
#include "ElmMats.h"
|
||||
#include "ElmNorm.h"
|
||||
#include "TimeDomain.h"
|
||||
@@ -101,14 +102,14 @@ void NonlinearElasticityUL::initResultPoints ()
|
||||
|
||||
|
||||
bool NonlinearElasticityUL::evalInt (LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
const FiniteElement& fe,
|
||||
const TimeDomain& prm,
|
||||
const Vec3& X) const
|
||||
{
|
||||
// Evaluate the deformation gradient, F, and the Green-Lagrange strains, E
|
||||
Tensor F(nsd);
|
||||
SymmTensor E(nsd);
|
||||
if (!this->kinematics(dNdX,F,E))
|
||||
if (!this->kinematics(fe.dNdX,F,E))
|
||||
return false;
|
||||
|
||||
bool lHaveStrains = !E.isZero();
|
||||
@@ -121,12 +122,12 @@ bool NonlinearElasticityUL::evalInt (LocalIntegral*& elmInt,
|
||||
if (J == 0.0) return false;
|
||||
|
||||
// Scale with J=|F| since we are integrating over current configuration
|
||||
detJW *= J;
|
||||
const_cast<double&>(fe.detJxW) *= J;
|
||||
|
||||
if (eKm || iS)
|
||||
{
|
||||
// Push-forward the basis function gradients to current configuration
|
||||
dNdx.multiply(dNdX,Fi); // dNdx = dNdX * F^-1
|
||||
dNdx.multiply(fe.dNdX,Fi); // dNdx = dNdX * F^-1
|
||||
// Compute the small-deformation strain-displacement matrix B from dNdx
|
||||
if (!this->formBmatrix(dNdx)) return false;
|
||||
#ifdef INT_DEBUG
|
||||
@@ -137,7 +138,7 @@ bool NonlinearElasticityUL::evalInt (LocalIntegral*& elmInt,
|
||||
}
|
||||
else if (eKm || iS)
|
||||
// Initial state, no deformation yet
|
||||
if (!this->formBmatrix(dNdX)) return false;
|
||||
if (!this->formBmatrix(fe.dNdX)) return false;
|
||||
|
||||
// Evaluate the constitutive relation
|
||||
SymmTensor sigma(nsd);
|
||||
@@ -151,36 +152,36 @@ bool NonlinearElasticityUL::evalInt (LocalIntegral*& elmInt,
|
||||
if (eKm)
|
||||
{
|
||||
// Integrate the material stiffness matrix
|
||||
CB.multiply(Cmat,Bmat).multiply(detJW); // CB = C*B*|J|*w
|
||||
eKm->multiply(Bmat,CB,true,false,true); // EK += B^T * CB
|
||||
CB.multiply(Cmat,Bmat).multiply(fe.detJxW); // CB = C*B*|J|*w
|
||||
eKm->multiply(Bmat,CB,true,false,true); // EK += B^T * CB
|
||||
}
|
||||
|
||||
if (eKg && lHaveStrains)
|
||||
// Integrate the geometric stiffness matrix
|
||||
this->formKG(*eKg,dNdx,sigma,detJW);
|
||||
this->formKG(*eKg,dNdx,sigma,fe.detJxW);
|
||||
|
||||
if (eM)
|
||||
// Integrate the mass matrix
|
||||
this->formMassMatrix(*eM,N,X,detJW);
|
||||
this->formMassMatrix(*eM,fe.N,X,fe.detJxW);
|
||||
|
||||
if (iS && lHaveStrains)
|
||||
{
|
||||
// Integrate the internal forces
|
||||
sigma *= -detJW;
|
||||
sigma *= -fe.detJxW;
|
||||
if (!Bmat.multiply(sigma,*iS,true,true)) // ES -= B^T*sigma
|
||||
return false;
|
||||
}
|
||||
|
||||
if (eS)
|
||||
// Integrate the load vector due to gravitation and other body forces
|
||||
this->formBodyForce(*eS,N,X,detJW);
|
||||
this->formBodyForce(*eS,fe.N,X,fe.detJxW);
|
||||
|
||||
return this->getIntegralResult(elmInt);
|
||||
}
|
||||
|
||||
|
||||
bool NonlinearElasticityUL::evalBou (LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
bool NonlinearElasticityUL::evalBou (LocalIntegral*& elmInt,
|
||||
const FiniteElement& fe,
|
||||
const Vec3& X, const Vec3& normal) const
|
||||
{
|
||||
if (!tracFld)
|
||||
@@ -207,7 +208,7 @@ bool NonlinearElasticityUL::evalBou (LocalIntegral*& elmInt, double detJW,
|
||||
{
|
||||
// Compute the deformation gradient, F
|
||||
Tensor F(nsd);
|
||||
if (!this->formDefGradient(dNdX,F)) return false;
|
||||
if (!this->formDefGradient(fe.dNdX,F)) return false;
|
||||
|
||||
// Check for with-rotated pressure load
|
||||
if (tracFld->isNormalPressure())
|
||||
@@ -226,14 +227,14 @@ bool NonlinearElasticityUL::evalBou (LocalIntegral*& elmInt, double detJW,
|
||||
}
|
||||
else
|
||||
// Scale with J=|F| since we are integrating over current configuration
|
||||
detJW *= F.det();
|
||||
const_cast<double&>(fe.detJxW) *= F.det();
|
||||
}
|
||||
|
||||
// Integrate the force vector
|
||||
Vector& ES = *eS;
|
||||
for (size_t a = 1; a <= N.size(); a++)
|
||||
for (size_t a = 1; a <= fe.N.size(); a++)
|
||||
for (i = 1; i <= nsd; i++)
|
||||
ES(nsd*(a-1)+i) += T[i-1]*N(a)*detJW;
|
||||
ES(nsd*(a-1)+i) += T[i-1]*fe.N(a)*fe.detJxW;
|
||||
|
||||
return this->getIntegralResult(elmInt);
|
||||
}
|
||||
@@ -323,8 +324,8 @@ void ElasticityNormUL::initIntegration (const TimeDomain& prm)
|
||||
|
||||
|
||||
bool ElasticityNormUL::evalInt (LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector&, const Matrix& dNdX,
|
||||
const FiniteElement& fe,
|
||||
const TimeDomain& prm,
|
||||
const Vec3& X) const
|
||||
{
|
||||
ElmNorm& pnorm = ElasticityNorm::getElmNormBuffer(elmInt);
|
||||
@@ -332,9 +333,9 @@ bool ElasticityNormUL::evalInt (LocalIntegral*& elmInt,
|
||||
NonlinearElasticityUL* ulp = static_cast<NonlinearElasticityUL*>(&problem);
|
||||
|
||||
// Evaluate the deformation gradient, F, and the Green-Lagrange strains, E
|
||||
Tensor F(dNdX.cols());
|
||||
SymmTensor E(dNdX.cols());
|
||||
if (!ulp->kinematics(dNdX,F,E))
|
||||
Tensor F(fe.dNdX.cols());
|
||||
SymmTensor E(fe.dNdX.cols());
|
||||
if (!ulp->kinematics(fe.dNdX,F,E))
|
||||
return false;
|
||||
|
||||
// Evaluate the 2nd Piola Kirchhoff stresses, S
|
||||
@@ -346,7 +347,24 @@ bool ElasticityNormUL::evalInt (LocalIntegral*& elmInt,
|
||||
|
||||
// Integrate the energy norm a(u^h,u^h) = Int_Omega0 (S:E) dV0
|
||||
if (U == 0.0) U = S.innerProd(E);
|
||||
pnorm[0] += U*detJW;
|
||||
pnorm[0] += U*fe.detJxW;
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
|
||||
bool ElasticityNormUL::evalIntMx (LocalIntegral*& elmInt,
|
||||
const MxFiniteElement& fe,
|
||||
const TimeDomain& prm,
|
||||
const Vec3& X) const
|
||||
{
|
||||
return this->evalInt(elmInt,fe,prm,X);
|
||||
}
|
||||
|
||||
|
||||
bool ElasticityNormUL::evalBouMx (LocalIntegral*& elmInt,
|
||||
const MxFiniteElement& fe,
|
||||
const Vec3& X, const Vec3& normal) const
|
||||
{
|
||||
return this->evalBou(elmInt,fe,X,normal);
|
||||
}
|
||||
|
||||
@@ -51,21 +51,15 @@ public:
|
||||
|
||||
//! \brief Evaluates the integrand at an interior point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] fe Finite element data of current integration point
|
||||
//! \param[in] prm Nonlinear solution algorithm parameters
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N Basis function values
|
||||
//! \param[in] dNdX Basis function gradients
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
virtual bool evalInt(LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
const Vec3& X) const;
|
||||
virtual bool evalInt(LocalIntegral*& elmInt, const FiniteElement& fe,
|
||||
const TimeDomain& prm, const Vec3& X) const;
|
||||
|
||||
//! \brief Evaluates the integrand at a boundary point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N Basis function values
|
||||
//! \param[in] dNdX Basis function gradients
|
||||
//! \param[in] fe Finite element data of current integration point
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
//! \param[in] normal Boundary normal vector at current integration point
|
||||
//!
|
||||
@@ -73,8 +67,7 @@ public:
|
||||
//! possibly with-rotated traction fields (non-conservative loads).
|
||||
//! For uni-directional (conservative) loads, it is similar to the
|
||||
//! \a LinearElasticity::evalBou method.
|
||||
virtual bool evalBou(LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
virtual bool evalBou(LocalIntegral*& elmInt, const FiniteElement& fe,
|
||||
const Vec3& X, const Vec3& normal) const;
|
||||
|
||||
//! \brief Returns a pointer to an Integrand for solution norm evaluation.
|
||||
@@ -144,45 +137,26 @@ public:
|
||||
|
||||
//! \brief Evaluates the integrand at an interior point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] fe Finite element data of current integration point
|
||||
//! \param[in] prm Nonlinear solution algorithm parameters
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N Basis function values
|
||||
//! \param[in] dNdX Basis function gradients
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
virtual bool evalInt(LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
const Vec3& X) const;
|
||||
virtual bool evalInt(LocalIntegral*& elmInt, const FiniteElement& fe,
|
||||
const TimeDomain& prm, const Vec3& X) const;
|
||||
//! \brief Evaluates the integrand at an interior point (mixed).
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] fe Mixed finite element data of current integration point
|
||||
//! \param[in] prm Nonlinear solution algorithm parameters
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N1 Basis function values for the displacement field
|
||||
//! \param[in] dN1dX Basis function gradients for the displacement field
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
virtual bool evalInt(LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N1, const Vector&,
|
||||
const Matrix& dN1dX, const Matrix&,
|
||||
const Vec3& X) const
|
||||
{
|
||||
return this->evalInt(elmInt,prm,detJW,N1,dN1dX,X);
|
||||
}
|
||||
virtual bool evalIntMx(LocalIntegral*& elmInt, const MxFiniteElement& fe,
|
||||
const TimeDomain& prm, const Vec3& X) const;
|
||||
|
||||
//! \brief Evaluates the integrand at a boundary point (mixed).
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N1 Basis function values for the displacement field
|
||||
//! \param[in] dN1dX Basis function gradients for the displacement field
|
||||
//! \param[in] fe Mixed finite element data of current integration point
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
//! \param[in] normal Boundary normal vector at current integration point
|
||||
virtual bool evalBou(LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N1, const Vector&,
|
||||
const Matrix& dN1dX, const Matrix&,
|
||||
const Vec3& X, const Vec3& normal) const
|
||||
{
|
||||
return this->ElasticityNorm::evalBou(elmInt,detJW,N1,dN1dX,X,normal);
|
||||
}
|
||||
virtual bool evalBouMx(LocalIntegral*& elmInt, const MxFiniteElement& fe,
|
||||
const Vec3& X, const Vec3& normal) const;
|
||||
};
|
||||
|
||||
#endif
|
||||
|
||||
@@ -13,6 +13,7 @@
|
||||
|
||||
#include "NonlinearElasticityULMX.h"
|
||||
#include "MaterialBase.h"
|
||||
#include "FiniteElement.h"
|
||||
#include "ElmMats.h"
|
||||
#include "ElmNorm.h"
|
||||
#include "TimeDomain.h"
|
||||
@@ -140,26 +141,25 @@ bool NonlinearElasticityULMX::initElement (const std::vector<int>& MNPC,
|
||||
|
||||
|
||||
bool NonlinearElasticityULMX::evalInt (LocalIntegral*& elmInt,
|
||||
const TimeDomain&, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
const Vec3& X) const
|
||||
const FiniteElement& fe,
|
||||
const TimeDomain&, const Vec3& X) const
|
||||
{
|
||||
if (myMats->b.size() < 3) return false;
|
||||
|
||||
ItgPtData& ptData = myData[iP++];
|
||||
|
||||
#if INT_DEBUG > 0
|
||||
std::cout <<"NonlinearElasticityUL::dNdX ="<< dNdX;
|
||||
std::cout <<"NonlinearElasticityUL::dNdX ="<< fe.dNdX;
|
||||
#endif
|
||||
|
||||
// Evaluate the deformation gradient, Fp, at previous configuration
|
||||
const_cast<NonlinearElasticityULMX*>(this)->eV = &myMats->b[2];
|
||||
if (!this->formDefGradient(dNdX,ptData.Fp))
|
||||
if (!this->formDefGradient(fe.dNdX,ptData.Fp))
|
||||
return false;
|
||||
|
||||
// Evaluate the deformation gradient, F, at current configuration
|
||||
const_cast<NonlinearElasticityULMX*>(this)->eV = &myMats->b[1];
|
||||
if (!this->formDefGradient(dNdX,ptData.F))
|
||||
if (!this->formDefGradient(fe.dNdX,ptData.F))
|
||||
return false;
|
||||
|
||||
if (nsd == 2) // In 2D we always assume plane strain so set F(3,3)=1
|
||||
@@ -174,10 +174,10 @@ bool NonlinearElasticityULMX::evalInt (LocalIntegral*& elmInt,
|
||||
if (J == 0.0) return false;
|
||||
|
||||
// Push-forward the basis function gradients to current configuration
|
||||
ptData.dNdx.multiply(dNdX,Fi); // dNdx = dNdX * F^-1
|
||||
ptData.dNdx.multiply(fe.dNdX,Fi); // dNdx = dNdX * F^-1
|
||||
|
||||
ptData.X.assign(X);
|
||||
ptData.detJW = detJW;
|
||||
ptData.detJW = fe.detJxW;
|
||||
|
||||
// Evaluate the pressure modes (generalized coordinates)
|
||||
Vec3 Xg = X-X0;
|
||||
@@ -190,15 +190,15 @@ bool NonlinearElasticityULMX::evalInt (LocalIntegral*& elmInt,
|
||||
|
||||
if (Hh)
|
||||
// Integrate the Hh-matrix
|
||||
Hh->outer_product(ptData.Phi,ptData.Phi*detJW,true);
|
||||
Hh->outer_product(ptData.Phi,ptData.Phi*fe.detJxW,true);
|
||||
|
||||
if (eM)
|
||||
// Integrate the mass matrix
|
||||
this->formMassMatrix(*eM,N,X,J*detJW);
|
||||
this->formMassMatrix(*eM,fe.N,X,J*fe.detJxW);
|
||||
|
||||
if (eS)
|
||||
// Integrate the load vector due to gravitation and other body forces
|
||||
this->formBodyForce(*eS,N,X,J*detJW);
|
||||
this->formBodyForce(*eS,fe.N,X,J*fe.detJxW);
|
||||
|
||||
return this->getIntegralResult(elmInt);
|
||||
}
|
||||
@@ -470,12 +470,10 @@ bool ElasticityNormULMX::initElement (const std::vector<int>& MNPC,
|
||||
|
||||
|
||||
bool ElasticityNormULMX::evalInt (LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
const Vec3& X) const
|
||||
const FiniteElement& fe,
|
||||
const TimeDomain& td, const Vec3& X) const
|
||||
{
|
||||
NonlinearElasticityULMX& elp = static_cast<NonlinearElasticityULMX&>(problem);
|
||||
return elp.evalInt(elmInt,prm,detJW,N,dNdX,X);
|
||||
return static_cast<NonlinearElasticityULMX&>(problem).evalInt(elmInt,fe,td,X);
|
||||
}
|
||||
|
||||
|
||||
|
||||
@@ -63,20 +63,16 @@ public:
|
||||
|
||||
//! \brief Evaluates the integrand at an interior point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] fe Finite element data of current integration point
|
||||
//! \param[in] prm Nonlinear solution algorithm parameters
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N Basis function values
|
||||
//! \param[in] dNdX Basis function gradients
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
//!
|
||||
//! \details This method mainly stores the integration point quantities
|
||||
//! depending on the basis function values in the internal member \a myData.
|
||||
//! The actual numerical integration of the tangent stiffness and internal
|
||||
//! forces is then performed by the \a finalizeElement method.
|
||||
virtual bool evalInt(LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
const Vec3& X) const;
|
||||
virtual bool evalInt(LocalIntegral*& elmInt, const FiniteElement& fe,
|
||||
const TimeDomain& prm, const Vec3& X) const;
|
||||
|
||||
//! \brief Finalizes the element matrices after the numerical integration.
|
||||
//! \details This method is used to implement the multiple integration point
|
||||
@@ -131,20 +127,16 @@ public:
|
||||
|
||||
//! \brief Evaluates the integrand at an interior point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] fe Finite element data of current integration point
|
||||
//! \param[in] prm Nonlinear solution algorithm parameters
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N Basis function values
|
||||
//! \param[in] dNdX Basis function gradients
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
//!
|
||||
//! \details This method just forwards to the \a evalInt method of the
|
||||
//! corresponding NonlinearElasticityULMX object, for which this class is
|
||||
//! declared as friend such that it can access the data members. The actual
|
||||
//! norm integration us then performed by the \a finalizeElement method.
|
||||
virtual bool evalInt(LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N, const Matrix& dNdX,
|
||||
const Vec3& X) const;
|
||||
virtual bool evalInt(LocalIntegral*& elmInt, const FiniteElement& fe,
|
||||
const TimeDomain& prm, const Vec3& X) const;
|
||||
|
||||
//! \brief Finalizes the element norms after the numerical integration.
|
||||
//! \details This method is used to implement the multiple integration point
|
||||
|
||||
@@ -13,6 +13,7 @@
|
||||
|
||||
#include "NonlinearElasticityULMixed.h"
|
||||
#include "MaterialBase.h"
|
||||
#include "FiniteElement.h"
|
||||
#include "ElmMats.h"
|
||||
#include "Utilities.h"
|
||||
|
||||
@@ -290,25 +291,22 @@ bool NonlinearElasticityULMixed::initElementBou (const std::vector<int>& MNPC1,
|
||||
}
|
||||
|
||||
|
||||
bool NonlinearElasticityULMixed::evalInt (LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N1,
|
||||
const Vector& N2,
|
||||
const Matrix& dNdX1,
|
||||
const Matrix& dNdX2,
|
||||
const Vec3& X) const
|
||||
bool NonlinearElasticityULMixed::evalIntMx (LocalIntegral*& elmInt,
|
||||
const MxFiniteElement& fe,
|
||||
const TimeDomain& prm,
|
||||
const Vec3& X) const
|
||||
{
|
||||
// Evaluate the deformation gradient, F, and the Green-Lagrange strains, E
|
||||
Tensor F(nsd);
|
||||
SymmTensor E(nsd);
|
||||
if (!this->kinematics(dNdX1,F,E))
|
||||
if (!this->kinematics(fe.dN1dX,F,E))
|
||||
return false;
|
||||
|
||||
bool lHaveStrains = !E.isZero();
|
||||
|
||||
// Evaluate the volumetric change and pressure fields
|
||||
double Theta = N2.dot(myMats->b[T]) + 1.0;
|
||||
double Press = N2.dot(myMats->b[P]);
|
||||
double Theta = fe.N2.dot(myMats->b[T]) + 1.0;
|
||||
double Press = fe.N2.dot(myMats->b[P]);
|
||||
#if INT_DEBUG > 0
|
||||
std::cout <<"NonlinearElasticityULMixed::b_theta ="<< myMats->b[T];
|
||||
std::cout <<"NonlinearElasticityULMixed::b_p ="<< myMats->b[P];
|
||||
@@ -316,7 +314,7 @@ bool NonlinearElasticityULMixed::evalInt (LocalIntegral*& elmInt,
|
||||
#endif
|
||||
|
||||
// Compute the mixed integration point volume
|
||||
double dVol = Theta*detJW;
|
||||
double dVol = Theta*fe.detJxW;
|
||||
|
||||
// Compute the mixed model deformation gradient, F_bar
|
||||
Fbar = F; // notice that F_bar always has dimension 3
|
||||
@@ -335,20 +333,20 @@ bool NonlinearElasticityULMixed::evalInt (LocalIntegral*& elmInt,
|
||||
if (J == 0.0) return false;
|
||||
|
||||
// Push-forward the basis function gradients to current configuration
|
||||
dNdx.multiply(dNdX1,Fi); // dNdx = dNdX * F^-1
|
||||
dNdx.multiply(fe.dN1dX,Fi); // dNdx = dNdX * F^-1
|
||||
#if INT_DEBUG > 0
|
||||
std::cout <<"NonlinearElasticityULMixed::dNdX ="<< dNdX1;
|
||||
std::cout <<"NonlinearElasticityULMixed::dNdX ="<< fe.dN1dX;
|
||||
std::cout <<"NonlinearElasticityULMixed::dNdx ="<< dNdx;
|
||||
std::cout <<"NonlinearElasticityULMixed::Fbar ="<< Fbar;
|
||||
#endif
|
||||
|
||||
if (eM)
|
||||
// Integrate the mass matrix
|
||||
this->formMassMatrix(*eM,N1,X,J*detJW);
|
||||
this->formMassMatrix(*eM,fe.N1,X,J*fe.detJxW);
|
||||
|
||||
if (eS)
|
||||
// Integrate the load vector due to gravitation and other body forces
|
||||
this->formBodyForce(*eS,N1,X,J*detJW);
|
||||
this->formBodyForce(*eS,fe.N1,X,J*fe.detJxW);
|
||||
|
||||
// Evaluate the constitutive relation
|
||||
SymmTensor Eps(3), Sig(3), Sigma(nsd);
|
||||
@@ -384,34 +382,34 @@ bool NonlinearElasticityULMixed::evalInt (LocalIntegral*& elmInt,
|
||||
// Integrate the material stiffness matrix
|
||||
Dmat *= dVol;
|
||||
if (nsd == 2)
|
||||
acckm2d_(N1.size(),dNdx.ptr(),Dmat.ptr(),eKm->ptr());
|
||||
acckm2d_(fe.N1.size(),dNdx.ptr(),Dmat.ptr(),eKm->ptr());
|
||||
else
|
||||
acckm3d_(N1.size(),dNdx.ptr(),Dmat.ptr(),eKm->ptr());
|
||||
acckm3d_(fe.N1.size(),dNdx.ptr(),Dmat.ptr(),eKm->ptr());
|
||||
#endif
|
||||
|
||||
// Integrate the volumetric change and pressure tangent terms
|
||||
size_t a, b;
|
||||
unsigned short int i, j, k;
|
||||
for (a = 1; a <= N1.size(); a++)
|
||||
for (b = 1; b <= N2.size(); b++)
|
||||
for (a = 1; a <= fe.N1.size(); a++)
|
||||
for (b = 1; b <= fe.N2.size(); b++)
|
||||
for (i = 1; i <= nsd; i++)
|
||||
{
|
||||
myMats->A[Kut](nsd*(a-1)+i,b) += dNdx(a,i)*N2(b) * Dmat(i,7);
|
||||
myMats->A[Kut](nsd*(a-1)+i,b) += dNdx(a,i)*fe.N2(b) * Dmat(i,7);
|
||||
if (nsd == 2)
|
||||
myMats->A[Kut](nsd*(a-1)+i,b) += dNdx(a,3-i)*N2(b) * Dmat(4,7);
|
||||
myMats->A[Kut](nsd*(a-1)+i,b) += dNdx(a,3-i)*fe.N2(b) * Dmat(4,7);
|
||||
else if (i < 3)
|
||||
{
|
||||
j = i + 1;
|
||||
k = j%3 + 1;
|
||||
myMats->A[Kut](nsd*(a-1)+i,b) += dNdx(a,3-i)*N2(b) * Dmat(4,7);
|
||||
myMats->A[Kut](nsd*(a-1)+j,b) += dNdx(a,5-j)*N2(b) * Dmat(5,7);
|
||||
myMats->A[Kut](nsd*(a-1)+k,b) += dNdx(a,4-k)*N2(b) * Dmat(6,7);
|
||||
myMats->A[Kut](nsd*(a-1)+i,b) += dNdx(a,3-i)*fe.N2(b) * Dmat(4,7);
|
||||
myMats->A[Kut](nsd*(a-1)+j,b) += dNdx(a,5-j)*fe.N2(b) * Dmat(5,7);
|
||||
myMats->A[Kut](nsd*(a-1)+k,b) += dNdx(a,4-k)*fe.N2(b) * Dmat(6,7);
|
||||
}
|
||||
myMats->A[Kup](nsd*(a-1)+i,b) += dNdx(a,i)*N2(b) * dVol;
|
||||
myMats->A[Kup](nsd*(a-1)+i,b) += dNdx(a,i)*fe.N2(b) * dVol;
|
||||
}
|
||||
|
||||
myMats->A[Ktt].outer_product(N2,N2*Dmat(7,7),true); // Ktt += N2*N2^T * D77
|
||||
myMats->A[Ktp].outer_product(N2,N2*(-detJW),true); // Ktp -= N2*N2^T * detJW
|
||||
myMats->A[Ktt].outer_product(fe.N2,fe.N2*Dmat(7,7),true); // += N2*N2^T*D77
|
||||
myMats->A[Ktp].outer_product(fe.N2,fe.N2*(-fe.detJxW),true); // -= N2*N2^T*|J|
|
||||
|
||||
if (lHaveStrains)
|
||||
{
|
||||
@@ -420,13 +418,28 @@ bool NonlinearElasticityULMixed::evalInt (LocalIntegral*& elmInt,
|
||||
|
||||
// Integrate the internal forces
|
||||
Sigma *= -dVol;
|
||||
if (!Bmat.multiply(Sigma,myMats->b[Ru],true,true)) // Ru -= B^T * sigma
|
||||
if (!Bmat.multiply(Sigma,myMats->b[Ru],true,true)) // -= B^T*sigma
|
||||
return false;
|
||||
|
||||
// Integrate the volumetric change and pressure forces
|
||||
myMats->b[Rt].add(N2,(Press - Bpres)*detJW); // Rt += N2 * (p-pBar)
|
||||
myMats->b[Rp].add(N2,(Theta - J)*detJW); // Rp += N2 * (Theta-J)
|
||||
myMats->b[Rt].add(fe.N2,(Press - Bpres)*fe.detJxW); // += N2*(p-pBar)*|J|
|
||||
myMats->b[Rp].add(fe.N2,(Theta - J)*fe.detJxW); // += N2*(Theta-J)*|J|
|
||||
}
|
||||
|
||||
return this->getIntegralResult(elmInt);
|
||||
}
|
||||
|
||||
|
||||
/*!
|
||||
The boundary integral is the same as that of the parent class. It does not
|
||||
depend on the pressure and volumetric change fields. Thus, this call is
|
||||
forwarded to the corresponding single-field method of the parent class.
|
||||
*/
|
||||
|
||||
bool NonlinearElasticityULMixed::evalBouMx (LocalIntegral*& elmInt,
|
||||
const MxFiniteElement& fe,
|
||||
const Vec3& X,
|
||||
const Vec3& normal) const
|
||||
{
|
||||
return this->evalBou(elmInt,fe,X,normal);
|
||||
}
|
||||
|
||||
@@ -70,37 +70,19 @@ public:
|
||||
|
||||
//! \brief Evaluates the mixed field problem integrand at an interior point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] fe Mixed finite element data of current integration point
|
||||
//! \param[in] prm Nonlinear solution algorithm parameters
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N1 Basis function values, field 1 (displacements)
|
||||
//! \param[in] N2 Basis function values, field 2 (pressure and volume change)
|
||||
//! \param[in] dN1dX Basis function gradients, field 1
|
||||
//! \param[in] dN2dX Basis function gradients, field 2
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
virtual bool evalInt(LocalIntegral*& elmInt,
|
||||
const TimeDomain& prm, double detJW,
|
||||
const Vector& N1, const Vector& N2,
|
||||
const Matrix& dN1dX, const Matrix& dN2dX,
|
||||
const Vec3& X) const;
|
||||
virtual bool evalIntMx(LocalIntegral*& elmInt, const MxFiniteElement& fe,
|
||||
const TimeDomain& prm, const Vec3& X) const;
|
||||
|
||||
//! \brief Evaluates the integrand at a boundary point.
|
||||
//! \param elmInt The local integral object to receive the contributions
|
||||
//! \param[in] detJW Jacobian determinant times integration point weight
|
||||
//! \param[in] N1 Basis function values, field 1
|
||||
//! \param[in] dN1dX Basis function gradients, field 1
|
||||
//! \param[in] fe Mixed finite element data of current integration point
|
||||
//! \param[in] X Cartesian coordinates of current integration point
|
||||
//! \param[in] normal Boundary normal vector at current integration point
|
||||
//!
|
||||
//! \details The boundary integral is the same as that of the parent class.
|
||||
//! It does not depend on the pressure and volumetric change fields.
|
||||
//! Thus, this call is forwarded to the single-field parent class method.
|
||||
virtual bool evalBou(LocalIntegral*& elmInt, double detJW,
|
||||
const Vector& N1, const Vector&,
|
||||
const Matrix& dN1dX, const Matrix&,
|
||||
const Vec3& X, const Vec3& normal) const
|
||||
{
|
||||
return this->NonlinearElasticityUL::evalBou(elmInt,detJW,N1,dN1dX,X,normal);
|
||||
}
|
||||
virtual bool evalBouMx(LocalIntegral*& elmInt, const MxFiniteElement& fe,
|
||||
const Vec3& X, const Vec3& normal) const;
|
||||
|
||||
protected:
|
||||
mutable Tensor Fbar; //!< Mixed model deformation gradient
|
||||
|
||||
Reference in New Issue
Block a user