MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
nonlininteg.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.googlecode.com.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "fem.hpp"
13 
14 namespace mfem
15 {
16 
18  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
19  DenseMatrix &elmat)
20 {
21  mfem_error("NonlinearFormIntegrator::AssembleElementGrad"
22  " is not overloaded!");
23 }
24 
26  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
27 {
28  mfem_error("NonlinearFormIntegrator::GetElementEnergy"
29  " is not overloaded!");
30  return 0.0;
31 }
32 
33 
35 {
36  Z.SetSize(J.Width());
38  return 0.5*(Z*Z)/J.Det();
39 }
40 
42 {
43  int dim = J.Width();
44  double t;
45 
46  Z.SetSize(dim);
47  S.SetSize(dim);
49  MultAAt(Z, S);
50  t = 0.5*S.Trace();
51  for (int i = 0; i < dim; i++)
52  S(i,i) -= t;
53  t = J.Det();
54  S *= -1.0/(t*t);
55  Mult(S, Z, P);
56 }
57 
59  const DenseMatrix &J, const DenseMatrix &DS, const double weight,
60  DenseMatrix &A) const
61 {
62  int dof = DS.Height(), dim = DS.Width();
63  double t;
64 
65  Z.SetSize(dim);
66  S.SetSize(dim);
67  G.SetSize(dof, dim);
68  C.SetSize(dof, dim);
69 
71  MultAAt(Z, S);
72 
73  t = 1.0/J.Det();
74  Z *= t; // Z = J^{-t}
75  S *= t; // S = |J| (J.J^t)^{-1}
76  t = 0.5*S.Trace();
77 
78  MultABt(DS, Z, G); // G = DS.J^{-1}
79  Mult(G, S, C);
80 
81  // 1.
82  for (int i = 0; i < dof; i++)
83  for (int j = 0; j <= i; j++)
84  {
85  double a = 0.0;
86  for (int d = 0; d < dim; d++)
87  a += G(i,d)*G(j,d);
88  a *= weight;
89  for (int k = 0; k < dim; k++)
90  for (int l = 0; l <= k; l++)
91  {
92  double b = a*S(k,l);
93  A(i+k*dof,j+l*dof) += b;
94  if (i != j)
95  A(j+k*dof,i+l*dof) += b;
96  if (k != l)
97  {
98  A(i+l*dof,j+k*dof) += b;
99  if (i != j)
100  A(j+l*dof,i+k*dof) += b;
101  }
102  }
103  }
104 
105  // 2.
106  for (int i = 0; i < dof; i++)
107  for (int j = 0; j < i; j++)
108  {
109  for (int k = 0; k < dim; k++)
110  for (int l = 0; l < k; l++)
111  {
112  double a =
113  weight*(C(i,l)*G(j,k) - C(i,k)*G(j,l) +
114  C(j,k)*G(i,l) - C(j,l)*G(i,k) +
115  t*(G(i,k)*G(j,l) - G(i,l)*G(j,k)));
116 
117  A(i+k*dof,j+l*dof) += a;
118  A(j+l*dof,i+k*dof) += a;
119 
120  A(i+l*dof,j+k*dof) -= a;
121  A(j+k*dof,i+l*dof) -= a;
122  }
123  }
124 }
125 
126 
127 inline void NeoHookeanModel::EvalCoeffs() const
128 {
129  mu = c_mu->Eval(*T, T->GetIntPoint());
130  K = c_K->Eval(*T, T->GetIntPoint());
131  if (c_g)
132  g = c_g->Eval(*T, T->GetIntPoint());
133 }
134 
135 double NeoHookeanModel::EvalW(const DenseMatrix &J) const
136 {
137  int dim = J.Width();
138 
139  if (have_coeffs)
140  EvalCoeffs();
141 
142  double dJ = J.Det();
143  double sJ = dJ/g;
144  double bI1 = pow(dJ, -2.0/dim)*(J*J); // \bar{I}_1
145 
146  return 0.5*(mu*(bI1 - dim) + K*(sJ - 1.0)*(sJ - 1.0));
147 }
148 
150 {
151  int dim = J.Width();
152 
153  if (have_coeffs)
154  EvalCoeffs();
155 
156  Z.SetSize(dim);
158 
159  double dJ = J.Det();
160  double a = mu*pow(dJ, -2.0/dim);
161  double b = K*(dJ/g - 1.0)/g - a*(J*J)/(dim*dJ);
162 
163  P = 0.0;
164  P.Add(a, J);
165  P.Add(b, Z);
166 }
167 
169  const double weight, DenseMatrix &A) const
170 {
171  int dof = DS.Height(), dim = DS.Width();
172 
173  if (have_coeffs)
174  EvalCoeffs();
175 
176  Z.SetSize(dim);
177  G.SetSize(dof, dim);
178  C.SetSize(dof, dim);
179 
180  double dJ = J.Det();
181  double sJ = dJ/g;
182  double a = mu*pow(dJ, -2.0/dim);
183  double bc = a*(J*J)/dim;
184  double b = bc - K*sJ*(sJ - 1.0);
185  double c = 2.0*bc/dim + K*sJ*(2.0*sJ - 1.0);
186 
188  Z *= (1.0/dJ); // Z = J^{-t}
189 
190  MultABt(DS, J, C); // C = DS J^t
191  MultABt(DS, Z, G); // G = DS J^{-1}
192 
193  a *= weight;
194  b *= weight;
195  c *= weight;
196 
197  // 1.
198  for (int i = 0; i < dof; i++)
199  for (int k = 0; k <= i; k++)
200  {
201  double s = 0.0;
202  for (int d = 0; d < dim; d++)
203  s += DS(i,d)*DS(k,d);
204  s *= a;
205 
206  for (int d = 0; d < dim; d++)
207  A(i+d*dof,k+d*dof) += s;
208 
209  if (k != i)
210  for (int d = 0; d < dim; d++)
211  A(k+d*dof,i+d*dof) += s;
212  }
213 
214  a *= (-2.0/dim);
215 
216  // 2.
217  for (int i = 0; i < dof; i++)
218  for (int j = 0; j < dim; j++)
219  for (int k = 0; k < dof; k++)
220  for (int l = 0; l < dim; l++)
221  {
222  A(i+j*dof,k+l*dof) +=
223  a*(C(i,j)*G(k,l) + G(i,j)*C(k,l)) +
224  b*G(i,l)*G(k,j) + c*G(i,j)*G(k,l);
225  }
226 }
227 
228 
231  const Vector &elfun)
232 {
233  int dof = el.GetDof(), dim = el.GetDim();
234  double energy;
235 
236  DSh.SetSize(dof, dim);
237  J0i.SetSize(dim);
238  J1.SetSize(dim);
239  J.SetSize(dim);
240  PMatI.UseExternalData(elfun.GetData(), dof, dim);
241 
242  int intorder = 2*el.GetOrder() + 3; // <---
243  const IntegrationRule &ir = IntRules.Get(el.GetGeomType(), intorder);
244 
245  energy = 0.0;
246  model->SetTransformation(Tr);
247  for (int i = 0; i < ir.GetNPoints(); i++)
248  {
249  const IntegrationPoint &ip = ir.IntPoint(i);
250  Tr.SetIntPoint(&ip);
251  CalcInverse(Tr.Jacobian(), J0i);
252 
253  el.CalcDShape(ip, DSh);
254  MultAtB(PMatI, DSh, J1);
255  Mult(J1, J0i, J);
256 
257  energy += ip.weight*Tr.Weight()*model->EvalW(J);
258  }
259 
260  return energy;
261 }
262 
264  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
265  Vector &elvect)
266 {
267  int dof = el.GetDof(), dim = el.GetDim();
268 
269  DSh.SetSize(dof, dim);
270  DS.SetSize(dof, dim);
271  J0i.SetSize(dim);
272  J.SetSize(dim);
273  P.SetSize(dim);
274  PMatI.UseExternalData(elfun.GetData(), dof, dim);
275  elvect.SetSize(dof*dim);
276  PMatO.UseExternalData(elvect.GetData(), dof, dim);
277 
278  int intorder = 2*el.GetOrder() + 3; // <---
279  const IntegrationRule &ir = IntRules.Get(el.GetGeomType(), intorder);
280 
281  elvect = 0.0;
282  model->SetTransformation(Tr);
283  for (int i = 0; i < ir.GetNPoints(); i++)
284  {
285  const IntegrationPoint &ip = ir.IntPoint(i);
286  Tr.SetIntPoint(&ip);
287  CalcInverse(Tr.Jacobian(), J0i);
288 
289  el.CalcDShape(ip, DSh);
290  Mult(DSh, J0i, DS);
291  MultAtB(PMatI, DS, J);
292 
293  model->EvalP(J, P);
294 
295  P *= ip.weight*Tr.Weight();
296  AddMultABt(DS, P, PMatO);
297  }
298 }
299 
301  const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun,
302  DenseMatrix &elmat)
303 {
304  int dof = el.GetDof(), dim = el.GetDim();
305 
306  DSh.SetSize(dof, dim);
307  DS.SetSize(dof, dim);
308  J0i.SetSize(dim);
309  J.SetSize(dim);
310  PMatI.UseExternalData(elfun.GetData(), dof, dim);
311  elmat.SetSize(dof*dim);
312 
313  int intorder = 2*el.GetOrder() + 3; // <---
314  const IntegrationRule &ir = IntRules.Get(el.GetGeomType(), intorder);
315 
316  elmat = 0.0;
317  model->SetTransformation(Tr);
318  for (int i = 0; i < ir.GetNPoints(); i++)
319  {
320  const IntegrationPoint &ip = ir.IntPoint(i);
321  Tr.SetIntPoint(&ip);
322  CalcInverse(Tr.Jacobian(), J0i);
323 
324  el.CalcDShape(ip, DSh);
325  Mult(DSh, J0i, DS);
326  MultAtB(PMatI, DS, J);
327 
328  model->AssembleH(J, DS, ip.weight*Tr.Weight(), elmat);
329  }
330 }
331 
333 {
334  PMatI.ClearExternalData();
335  PMatO.ClearExternalData();
336 }
337 
338 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:202
Abstract class for Finite Elements.
Definition: fe.hpp:42
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
Compute the local energy.
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
Definition: densemat.cpp:2782
virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P=P(J).
virtual double EvalW(const DenseMatrix &J) const
Evaluate the strain energy density function, W=W(J).
int GetDim() const
Returns the space dimension for the finite element.
Definition: fe.hpp:80
Class for integration rule.
Definition: intrules.hpp:63
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:194
void SetSize(int s)
Resizes the vector if the new size is different.
Definition: vector.hpp:248
virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Definition: nonlininteg.cpp:58
double Det() const
Calculates the determinant of the matrix (for 2x2 or 3x3 matrices)
Definition: densemat.cpp:321
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:351
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols.
Definition: operator.hpp:41
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:33
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:89
Data type dense matrix.
Definition: densemat.hpp:22
virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P=P(J).
Definition: nonlininteg.cpp:41
const IntegrationPoint & GetIntPoint()
Definition: eltrans.hpp:35
double * GetData() const
Definition: vector.hpp:80
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:205
void CalcAdjugateTranspose(const DenseMatrix &a, DenseMatrix &adjat)
Calculate the transposed adjugate of a matrix (for NxN matrices, N=1,2,3)
Definition: densemat.cpp:2554
IntegrationRules IntRules(0)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:264
void EvalCoeffs() const
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Definition: densemat.cpp:375
void SetTransformation(ElementTransformation &_T)
An element transformation that can be used to evaluate coefficients.
Definition: nonlininteg.hpp:58
virtual double Weight()=0
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:2588
int GetGeomType() const
Returns the geometry type:
Definition: fe.hpp:83
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: nonlininteg.cpp:17
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
double Trace() const
Trace of a square matrix.
Definition: densemat.cpp:301
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
Definition: densemat.cpp:2859
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:86
void mfem_error(const char *msg)
Definition: error.cpp:23
void ClearExternalData()
Definition: densemat.hpp:57
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
Definition: densemat.cpp:2732
virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS, const double weight, DenseMatrix &A) const =0
virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Class for integration point with weight.
Definition: intrules.hpp:25
virtual double EvalW(const DenseMatrix &J) const
Evaluate the strain energy density function, W=W(J).
Definition: nonlininteg.cpp:34
virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const =0
Evaluate the 1st Piola-Kirchhoff stress tensor, P=P(J).
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
Definition: densemat.cpp:2912
Vector data type.
Definition: vector.hpp:29
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
Compute the local energy.
Definition: nonlininteg.cpp:25
virtual const DenseMatrix & Jacobian()=0
void UseExternalData(double *d, int h, int w)
Definition: densemat.hpp:54
void SetSize(int s)
If the matrix is not a square matrix of size s then recreate it.
Definition: densemat.cpp:73
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
virtual double EvalW(const DenseMatrix &J) const =0
Evaluate the strain energy density function, W=W(J).
ElementTransformation * T
Definition: nonlininteg.hpp:52