13 #include "../general/forall.hpp" 20 mfem_error (
"NonlinearFormIntegrator::GetLocalStateEnergyPA(...)\n" 21 " is not implemented for this class.");
27 mfem_error (
"NonlinearFormIntegrator::AssemblePA(...)\n" 28 " is not implemented for this class.");
34 mfem_error (
"NonlinearFormIntegrator::AssemblePA(...)\n" 35 " is not implemented for this class.");
41 mfem_error (
"NonlinearFormIntegrator::AssembleGradPA(...)\n" 42 " is not implemented for this class.");
47 mfem_error (
"NonlinearFormIntegrator::AddMultPA(...)\n" 48 " is not implemented for this class.");
53 mfem_error (
"NonlinearFormIntegrator::AddMultGradPA(...)\n" 54 " is not implemented for this class.");
59 mfem_error (
"NonlinearFormIntegrator::AssembleGradDiagonalPA(...)\n" 60 " is not implemented for this class.");
65 mfem_error (
"NonlinearFormIntegrator::AssembleMF(...)\n" 66 " is not implemented for this class.");
71 mfem_error (
"NonlinearFormIntegrator::AddMultMF(...)\n" 72 " is not implemented for this class.");
79 mfem_error(
"NonlinearFormIntegrator::AssembleElementVector" 80 " is not overloaded!");
87 mfem_error(
"NonlinearFormIntegrator::AssembleFaceVector" 88 " is not overloaded!");
95 mfem_error(
"NonlinearFormIntegrator::AssembleElementGrad" 96 " is not overloaded!");
104 mfem_error(
"NonlinearFormIntegrator::AssembleFaceGrad" 105 " is not overloaded!");
111 mfem_error(
"NonlinearFormIntegrator::GetElementEnergy" 112 " is not overloaded!");
123 mfem_error(
"BlockNonlinearFormIntegrator::AssembleElementVector" 124 " is not overloaded!");
134 mfem_error(
"BlockNonlinearFormIntegrator::AssembleFaceVector" 135 " is not overloaded!");
144 mfem_error(
"BlockNonlinearFormIntegrator::AssembleElementGrad" 145 " is not overloaded!");
155 mfem_error(
"BlockNonlinearFormIntegrator::AssembleFaceGrad" 156 " is not overloaded!");
164 mfem_error(
"BlockNonlinearFormIntegrator::GetElementEnergy" 165 " is not overloaded!");
174 return 0.5*(
Z*
Z)/J.
Det();
187 for (
int i = 0; i <
dim; i++)
220 for (
int i = 0; i < dof; i++)
221 for (
int j = 0; j <= i; j++)
224 for (
int d = 0; d <
dim; d++)
229 for (
int k = 0; k <
dim; k++)
230 for (
int l = 0; l <= k; l++)
233 A(i+k*dof,j+l*dof) +=
b;
236 A(j+k*dof,i+l*dof) +=
b;
240 A(i+l*dof,j+k*dof) +=
b;
243 A(j+l*dof,i+k*dof) +=
b;
250 for (
int i = 1; i < dof; i++)
251 for (
int j = 0; j < i; j++)
253 for (
int k = 1; k <
dim; k++)
254 for (
int l = 0; l < k; l++)
257 weight*(
C(i,l)*
G(j,k) -
C(i,k)*
G(j,l) +
258 C(j,k)*
G(i,l) -
C(j,l)*
G(i,k) +
259 t*(
G(i,k)*
G(j,l) -
G(i,l)*
G(j,k)));
261 A(i+k*dof,j+l*dof) +=
a;
262 A(j+l*dof,i+k*dof) +=
a;
264 A(i+l*dof,j+k*dof) -=
a;
265 A(j+k*dof,i+l*dof) -=
a;
292 double bI1 = pow(dJ, -2.0/
dim)*(J*J);
294 return 0.5*(
mu*(bI1 -
dim) +
K*(sJ - 1.0)*(sJ - 1.0));
310 double a =
mu*pow(dJ, -2.0/
dim);
311 double b =
K*(dJ/
g - 1.0)/
g -
a*(J*J)/(
dim*dJ);
334 double a =
mu*pow(dJ, -2.0/
dim);
335 double bc =
a*(J*J)/
dim;
336 double b = bc -
K*sJ*(sJ - 1.0);
337 double c = 2.0*bc/
dim +
K*sJ*(2.0*sJ - 1.0);
350 for (
int i = 0; i < dof; i++)
351 for (
int k = 0; k <= i; k++)
354 for (
int d = 0; d <
dim; d++)
356 s += DS(i,d)*DS(k,d);
360 for (
int d = 0; d <
dim; d++)
362 A(i+d*dof,k+d*dof) +=
s;
366 for (
int d = 0; d <
dim; d++)
368 A(k+d*dof,i+d*dof) +=
s;
375 for (
int i = 0; i < dof; i++)
376 for (
int j = 0; j <
dim; j++)
377 for (
int k = 0; k < dof; k++)
378 for (
int l = 0; l <
dim; l++)
380 A(i+j*dof,k+l*dof) +=
381 a*(
C(i,j)*
G(k,l) +
G(i,j)*
C(k,l)) +
382 b*
G(i,l)*
G(k,j) + c*
G(i,j)*
G(k,l);
457 model->
EvalP(Jpt, P);
507 mfem_error(
"IncompressibleNeoHookeanIntegrator::GetElementEnergy" 508 " has incorrect block finite element space size!");
511 int dof_u = el[0]->GetDof();
512 int dim = el[0]->GetDim();
520 int intorder = 2*el[0]->GetOrder() + 3;
532 el[0]->CalcDShape(ip, DSh_u);
552 mfem_error(
"IncompressibleNeoHookeanIntegrator::AssembleElementVector" 553 " has finite element space of incorrect block number");
556 int dof_u = el[0]->GetDof();
557 int dof_p = el[1]->GetDof();
559 int dim = el[0]->GetDim();
564 mfem_error(
"IncompressibleNeoHookeanIntegrator::AssembleElementVector" 565 " is not defined on manifold meshes");
582 int intorder = 2*el[0]->GetOrder() + 3;
594 el[0]->CalcDShape(ip, DSh_u);
595 Mult(DSh_u, J0i, DS_u);
598 el[1]->CalcShape(ip, Sh_p);
600 double pres = Sh_p * *elfun[1];
601 double mu = c_mu->
Eval(Tr, ip);
608 P.
Add(-1.0 * pres * dJ, FinvT);
613 elvec[1]->Add(ip.
weight * Tr.
Weight() * (dJ - 1.0), Sh_p);
624 int dof_u = el[0]->GetDof();
625 int dof_p = el[1]->GetDof();
627 int dim = el[0]->GetDim();
632 elmats(1,1)->
SetSize(dof_p, dof_p);
649 int intorder = 2*el[0]->GetOrder() + 3;
658 el[0]->CalcDShape(ip, DSh_u);
659 Mult(DSh_u, J0i, DS_u);
662 el[1]->CalcShape(ip, Sh_p);
663 double pres = Sh_p * *elfun[1];
664 double mu = c_mu->
Eval(Tr, ip);
671 for (
int i_u = 0; i_u < dof_u; ++i_u)
673 for (
int i_dim = 0; i_dim <
dim; ++i_dim)
675 for (
int j_u = 0; j_u < dof_u; ++j_u)
677 for (
int j_dim = 0; j_dim <
dim; ++j_dim)
683 for (
int n=0; n<
dim; ++n)
685 for (
int l=0; l<
dim; ++l)
687 (*elmats(0,0))(i_u + i_dim*dof_u, j_u + j_dim*dof_u) +=
688 dJ * (
mu * F(i_dim, l) - pres * FinvT(i_dim,l)) *
689 FinvT(j_dim,n) * DS_u(i_u,l) * DS_u(j_u, n) *
692 if (j_dim == i_dim && n==l)
694 (*elmats(0,0))(i_u + i_dim*dof_u, j_u + j_dim*dof_u) +=
695 dJ *
mu * DS_u(i_u, l) * DS_u(j_u,n) *
701 (*elmats(0,0))(i_u + i_dim*dof_u, j_u + j_dim*dof_u) +=
702 dJ * pres * FinvT(i_dim, n) *
703 FinvT(j_dim,l) * DS_u(i_u,l) * DS_u(j_u,n) *
713 for (
int i_p = 0; i_p < dof_p; ++i_p)
715 for (
int j_u = 0; j_u < dof_u; ++j_u)
717 for (
int dim_u = 0; dim_u <
dim; ++dim_u)
719 for (
int l=0; l<
dim; ++l)
721 dJ_FinvT_DS = dJ * FinvT(dim_u,l) * DS_u(j_u, l) * Sh_p(i_p) *
723 (*elmats(1,0))(i_p, j_u + dof_u * dim_u) += dJ_FinvT_DS;
724 (*elmats(0,1))(j_u + dof_u * dim_u, i_p) -= dJ_FinvT_DS;
738 const int order = 2 * fe.
GetOrder() + T.OrderGrad(&fe);
748 const int nd = el.
GetDof();
768 double w = ip.
weight * T.Weight();
769 if (Q) { w *= Q->
Eval(T, ip); }
773 gradEF.
Mult(vec1, vec2);
785 const int nd = el.
GetDof();
806 trans.SetIntPoint(&ip);
811 Mult(dshape,
trans.InverseJacobian(), dshapex);
823 trans.AdjugateJacobian().Mult(vec1, vec2);
826 dshape.
Mult(vec2, vec3);
827 MultVWt(shape, vec3, elmat_comp);
829 for (
int ii = 0; ii <
dim; ii++)
831 elmat.
AddMatrix(elmat_comp, ii * nd, ii * nd);
840 for (
int ii = 0; ii <
dim; ii++)
842 for (
int jj = 0; jj <
dim; jj++)
844 elmat.
AddMatrix(w * gradEF(ii, jj), elmat_comp, ii * nd, jj * nd);
857 const int nd = el.
GetDof();
877 trans.SetIntPoint(&ip);
886 trans.AdjugateJacobian().Mult(vec1, vec2);
889 dshape.
Mult(vec2, vec3);
890 MultVWt(shape, vec3, elmat_comp);
892 for (
int ii = 0; ii <
dim; ii++)
894 elmat.
AddMatrix(elmat_comp, ii * nd, ii * nd);
906 const int nd = el.
GetDof();
929 trans.SetIntPoint(&ip);
934 Mult(dshape,
trans.InverseJacobian(), dshapex);
940 trans.AdjugateJacobian().Mult(vec1, vec2);
943 dshape.
Mult(vec2, vec3);
944 MultVWt(shape, vec3, elmat_comp);
947 for (
int ii = 0; ii <
dim; ii++)
949 elmat.
AddMatrix(.5, elmat_comp, ii * nd, ii * nd);
950 elmat.
AddMatrix(-.5, elmat_comp_T, ii * nd, ii * nd);
Abstract class for all finite elements.
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
void trans(const Vector &u, Vector &x)
int GetNPoints() const
Returns the number of the points in the integration rule.
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
Class for an integration rule - an Array of IntegrationPoint.
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
void SetSize(int s)
Resize the vector to size s.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Data type dense matrix using column-major storage.
virtual double EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt).
virtual double EvalW(const DenseMatrix &J) const
Evaluate the strain energy density function, W = W(Jpt).
virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
virtual double GetElementEnergy(const Array< const FiniteElement *> &el, ElementTransformation &Tr, const Array< const Vector *> &elfun)
Compute the local energy.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const =0
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
double Trace() const
Trace of a square matrix.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
virtual double EvalW(const DenseMatrix &J) const
Evaluate the strain energy density function, W = W(Jpt).
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
void SetSize(int m, int n)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void CalcAdjugateTranspose(const DenseMatrix &a, DenseMatrix &adjat)
Calculate the transposed adjugate of a matrix (for NxN matrices, N=1,2,3)
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
void Mult(const double *x, double *y) const
Matrix vector multiplication.
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun)
Computes the integral of W(Jacobian(Trt)) over a target zone.
Dynamic 2D array using row-major layout.
double * GetData() const
Return a pointer to the beginning of the Vector data.
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
void Transpose()
(*this) = (*this)^t
int GetDim() const
Returns the reference space dimension for the finite element.
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void SetTransformation(ElementTransformation &Ttr_)
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0<=i<A.Height, 0<=j<A.Width.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva)
Calculate the inverse transpose of a matrix (for NxN matrices, N=1,2,3)
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
ElementTransformation * Ttr
virtual void AssembleElementGrad(const Array< const FiniteElement *> &el, ElementTransformation &Tr, const Array< const Vector *> &elfun, const Array2D< DenseMatrix *> &elmats)
Assemble the local gradient matrix.
Class for integration point with weight.
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const =0
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
int Size() const
Return the logical size of the array.
virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
void UseExternalData(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
static const IntegrationRule & GetRule(const FiniteElement &fe, ElementTransformation &T)
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
virtual void AssembleElementVector(const Array< const FiniteElement *> &el, ElementTransformation &Tr, const Array< const Vector *> &elfun, const Array< Vector *> &elvec)
Perform the local action of the NonlinearFormIntegrator.