21 mfem_error(
"NonlinearFormIntegrator::AssembleElementVector"
22 " is not overloaded!");
29 mfem_error(
"NonlinearFormIntegrator::AssembleFaceVector"
30 " is not overloaded!");
37 mfem_error(
"NonlinearFormIntegrator::AssembleElementGrad"
38 " is not overloaded!");
46 mfem_error(
"NonlinearFormIntegrator::AssembleFaceGrad"
47 " is not overloaded!");
53 mfem_error(
"NonlinearFormIntegrator::GetElementEnergy"
54 " is not overloaded!");
65 mfem_error(
"BlockNonlinearFormIntegrator::AssembleElementVector"
66 " is not overloaded!");
76 mfem_error(
"BlockNonlinearFormIntegrator::AssembleFaceVector"
77 " is not overloaded!");
86 mfem_error(
"BlockNonlinearFormIntegrator::AssembleElementGrad"
87 " is not overloaded!");
97 mfem_error(
"BlockNonlinearFormIntegrator::AssembleFaceGrad"
98 " is not overloaded!");
106 mfem_error(
"BlockNonlinearFormIntegrator::GetElementEnergy"
107 " is not overloaded!");
116 return 0.5*(
Z*
Z)/J.
Det();
129 for (
int i = 0; i <
dim; i++)
162 for (
int i = 0; i < dof; i++)
163 for (
int j = 0; j <= i; j++)
166 for (
int d = 0; d <
dim; d++)
171 for (
int k = 0; k <
dim; k++)
172 for (
int l = 0; l <= k; l++)
175 A(i+k*dof,j+l*dof) += b;
178 A(j+k*dof,i+l*dof) += b;
182 A(i+l*dof,j+k*dof) += b;
185 A(j+l*dof,i+k*dof) += b;
192 for (
int i = 1; i < dof; i++)
193 for (
int j = 0; j < i; j++)
195 for (
int k = 1; k <
dim; k++)
196 for (
int l = 0; l < k; l++)
199 weight*(
C(i,l)*
G(j,k) -
C(i,k)*
G(j,l) +
200 C(j,k)*
G(i,l) -
C(j,l)*
G(i,k) +
201 t*(
G(i,k)*
G(j,l) -
G(i,l)*
G(j,k)));
203 A(i+k*dof,j+l*dof) += a;
204 A(j+l*dof,i+k*dof) += a;
206 A(i+l*dof,j+k*dof) -= a;
207 A(j+k*dof,i+l*dof) -= a;
234 double bI1 = pow(dJ, -2.0/dim)*(J*J);
236 return 0.5*(
mu*(bI1 -
dim) +
K*(sJ - 1.0)*(sJ - 1.0));
252 double a =
mu*pow(dJ, -2.0/dim);
253 double b =
K*(dJ/
g - 1.0)/
g - a*(J*J)/(dim*dJ);
276 double a =
mu*pow(dJ, -2.0/
dim);
277 double bc = a*(J*J)/
dim;
278 double b = bc -
K*sJ*(sJ - 1.0);
279 double c = 2.0*bc/
dim +
K*sJ*(2.0*sJ - 1.0);
292 for (
int i = 0; i < dof; i++)
293 for (
int k = 0; k <= i; k++)
296 for (
int d = 0; d <
dim; d++)
298 s += DS(i,d)*DS(k,d);
302 for (
int d = 0; d <
dim; d++)
304 A(i+d*dof,k+d*dof) += s;
308 for (
int d = 0; d <
dim; d++)
310 A(k+d*dof,i+d*dof) += s;
317 for (
int i = 0; i < dof; i++)
318 for (
int j = 0; j <
dim; j++)
319 for (
int k = 0; k < dof; k++)
320 for (
int l = 0; l <
dim; l++)
322 A(i+j*dof,k+l*dof) +=
323 a*(
C(i,j)*
G(k,l) +
G(i,j)*
C(k,l)) +
324 b*
G(i,l)*
G(k,j) + c*
G(i,j)*
G(k,l);
399 model->
EvalP(Jpt, P);
449 mfem_error(
"IncompressibleNeoHookeanIntegrator::GetElementEnergy"
450 " has incorrect block finite element space size!");
453 int dof_u = el[0]->GetDof();
454 int dim = el[0]->GetDim();
462 int intorder = 2*el[0]->GetOrder() + 3;
474 el[0]->CalcDShape(ip, DSh_u);
478 mu = c_mu->
Eval(Tr, ip);
494 mfem_error(
"IncompressibleNeoHookeanIntegrator::AssembleElementVector"
495 " has finite element space of incorrect block number");
498 int dof_u = el[0]->GetDof();
499 int dof_p = el[1]->GetDof();
501 int dim = el[0]->GetDim();
506 mfem_error(
"IncompressibleNeoHookeanIntegrator::AssembleElementVector"
507 " is not defined on manifold meshes");
524 int intorder = 2*el[0]->GetOrder() + 3;
536 el[0]->CalcDShape(ip, DSh_u);
537 Mult(DSh_u, J0i, DS_u);
540 el[1]->CalcShape(ip, Sh_p);
542 double pres = Sh_p * *elfun[1];
543 double mu = c_mu->
Eval(Tr, ip);
550 P.
Add(-1.0 * pres * dJ, FinvT);
555 elvec[1]->Add(ip.
weight * Tr.
Weight() * (dJ - 1.0), Sh_p);
566 int dof_u = el[0]->GetDof();
567 int dof_p = el[1]->GetDof();
569 int dim = el[0]->GetDim();
571 elmats(0,0)->
SetSize(dof_u*dim, dof_u*dim);
572 elmats(0,1)->
SetSize(dof_u*dim, dof_p);
573 elmats(1,0)->
SetSize(dof_p, dof_u*dim);
574 elmats(1,1)->
SetSize(dof_p, dof_p);
591 int intorder = 2*el[0]->GetOrder() + 3;
600 el[0]->CalcDShape(ip, DSh_u);
601 Mult(DSh_u, J0i, DS_u);
604 el[1]->CalcShape(ip, Sh_p);
605 double pres = Sh_p * *elfun[1];
606 double mu = c_mu->
Eval(Tr, ip);
613 for (
int i_u = 0; i_u < dof_u; ++i_u)
615 for (
int i_dim = 0; i_dim <
dim; ++i_dim)
617 for (
int j_u = 0; j_u < dof_u; ++j_u)
619 for (
int j_dim = 0; j_dim <
dim; ++j_dim)
625 for (
int n=0; n<
dim; ++n)
627 for (
int l=0; l<
dim; ++l)
629 (*elmats(0,0))(i_u + i_dim*dof_u, j_u + j_dim*dof_u) +=
630 dJ * (mu * F(i_dim, l) - pres * FinvT(i_dim,l)) *
631 FinvT(j_dim,n) * DS_u(i_u,l) * DS_u(j_u, n) *
634 if (j_dim == i_dim && n==l)
636 (*elmats(0,0))(i_u + i_dim*dof_u, j_u + j_dim*dof_u) +=
637 dJ * mu * DS_u(i_u, l) * DS_u(j_u,n) *
643 (*elmats(0,0))(i_u + i_dim*dof_u, j_u + j_dim*dof_u) +=
644 dJ * pres * FinvT(i_dim, n) *
645 FinvT(j_dim,l) * DS_u(i_u,l) * DS_u(j_u,n) *
655 for (
int i_p = 0; i_p < dof_p; ++i_p)
657 for (
int j_u = 0; j_u < dof_u; ++j_u)
659 for (
int dim_u = 0; dim_u <
dim; ++dim_u)
661 for (
int l=0; l<
dim; ++l)
663 dJ_FinvT_DS = dJ * FinvT(dim_u,l) * DS_u(j_u, l) * Sh_p(i_p) *
665 (*elmats(1,0))(i_p, j_u + dof_u * dim_u) += dJ_FinvT_DS;
666 (*elmats(0,1))(j_u + dof_u * dim_u, i_p) -= dJ_FinvT_DS;
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for 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.
int Size() const
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 EvalW(const DenseMatrix &J) const
Evaluate the strain energy density function, W = W(Jpt).
int GetDim() const
Returns the reference space dimension for the finite element.
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 an integration rule - an Array of IntegrationPoint.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void SetSize(int s)
Resize the vector to size 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...
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().
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
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 void EvalP(const DenseMatrix &J, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
double * GetData() const
Return a pointer to the beginning of the Vector data.
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const =0
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(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)
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the 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 double GetElementEnergy(const Array< const FiniteElement * > &el, ElementTransformation &Tr, const Array< const Vector * > &elfun)
Compute the local energy.
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.
void SetTransformation(ElementTransformation &_Ttr)
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
double Trace() const
Trace of a square matrix.
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
void CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva)
Calculate the inverse transpose of a matrix (for NxN matrices, N=1,2,3)
ElementTransformation * Ttr
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...
Class for integration point with weight.
virtual double EvalW(const DenseMatrix &J) const
Evaluate the strain energy density function, W = W(Jpt).
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 AssembleElementVector(const Array< const FiniteElement * > &el, ElementTransformation &Tr, const Array< const Vector * > &elfun, const Array< Vector * > &elvec)
Perform the local action of the NonlinearFormIntegrator.
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 MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
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 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...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)