MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMechanicLargeDisplacement.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Carlos Augusto Teixera Mendes
4 ** All rights reserved.
5 **
6 ** This file is part of the "GeMA" software. It's use should respect
7 ** the terms in the license agreement that can be found together
8 ** with this source code.
9 ** It is provided AS IS, with NO WARRANTY OF ANY KIND,
10 ** INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR
11 ** A PARTICULAR PURPOSE.
12 **
13 ************************************************************************/
14 
24 #ifndef _GEMA_PLUGIN_MECHANICAL_PHYSICS_MECHANIC_LARGEDISPLACEMENT_H_
25 #define _GEMA_PLUGIN_MECHANICAL_PHYSICS_MECHANIC_LARGEDISPLACEMENT_H_
26 
27 
28 #include "gmpMechanicalPhysics.h"
29 #include <gmTrace.h>
30 #include <math.h>
31 #include <gmGaussAccessor.h>
32 #include "gmpMechanicPoint.h"
33 
34 //class GmpMechanicalLargeDisplacementSolid;
36 template <class T> class GmpMechanicLargeDisplacement : public T
37 {
38 public:
39  GmpMechanicLargeDisplacement(const char* pluginType, GmSimulationData* simulation, QString id, QString description,
40  const GmpFemPhysicsCommonMaterialFactory* matFactory, const GmLogCategory& logger)
41  : T(pluginType, simulation, id, description, matFactory, logger)
42  {}
43 
44 protected:
47  virtual double fillBuMatrix(const GmElement* e, const GmShape* shape, const GmVector& coord,
48  const GmMatrix& X, const GmVector& N, const GmMatrix& J, GmMatrix& Bu, const GmMatrix& F)
49  {
50  S_TRACE();
51  assert(e); Q_UNUSED(N); Q_UNUSED(X);
52 
53  // integration point index
54  //int ip = mp->_index;
55  // mesh dimension
56  int d = GmpFemPhysicsCommon::nodeDim();
57  double k, l, m;
58  GmMatrix B;
59 
60  // example to read strain and stress
61  /*
62  // new strain
63  const GmVector e_new(mp->_newStrain->valueAt(e, ip, &coord), 2 * d);
64 
65  // new stresses
66  const GmVector s_new(mp->_newStress->valueAt(e, ip, &coord), 2*d);
67  */
68 
69  double detJ = shape->shapeCartesianPartialsFromJacobian(coord, J, B, true);
70 
71  assert(Bu.n_rows == 2*d && Bu.n_cols == e->numNodes()*d);
72  Bu.zeros();
73  GmMatrix L(d, d);
74  L.zeros();
75 
76  // Checking reference
78 
79  // 2D solid element
80  if (d == 2)
81  {
82  if (modeLD == "totalLagrangian")
83  {
84  for (int i = 0; i < e->numNodes(); i++)
85  {
86  k = 2 * i + 1;
87  l = k - 1;
88 
89  Bu(0, l) = F(0, 0)*B(0, i);
90  Bu(1, l) = F(0, 1)*B(1, i);
91  Bu(3, l) = F(0, 0)*B(1, i) + F(0, 1)*B(0, i);
92 
93  Bu(0, k) = F(1, 0)*B(0, i);
94  Bu(1, k) = F(1, 1)*B(1, i);
95  Bu(3, k) = F(1, 0)*B(1, i) + F(1, 1)*B(0, i);
96  }
97  }
98  else if(modeLD == "updatedLagrangian")
99  {
100  for (int i = 0; i < e->numNodes(); i++)
101  {
102  k = 2 * i + 1;
103  l = 2 * i;
104  Bu(0, l) = B(0, i);
105  Bu(1, k) = B(1, i);
106  Bu(3, l) = B(1, i);
107  Bu(3, k) = B(0, i);
108  }
109  }
110 
111  }
112  // 3D solid element
113  else if (d == 3)
114  {
115  if (modeLD == "totalLagrangian")
116  {
117  for (int i = 0; i < e->numNodes(); i++)
118  {
119  m = 3 * i + 2;
120  l = m - 1;
121  k = m - 2;
122 
123  Bu(0, k) = F(0, 0)*B(0, i);
124  Bu(1, k) = F(0, 1)*B(1, i);
125  Bu(2, k) = F(0, 2)*B(2, i);
126  Bu(3, k) = F(0, 0)*B(1, i) + F(0, 1)*B(0, i);
127  Bu(4, k) = F(0, 2)*B(0, i) + F(0, 0)*B(2, i);
128  Bu(5, k) = F(0, 1)*B(2, i) + F(0, 2)*B(1, i);
129 
130  Bu(0, l) = F(1, 0)*B(0, i);
131  Bu(1, l) = F(1, 1)*B(1, i);
132  Bu(2, l) = F(1, 2)*B(2, i);
133  Bu(3, l) = F(1, 0)*B(1, i) + F(1, 1)*B(0, i);
134  Bu(4, l) = F(1, 2)*B(0, i) + F(1, 0)*B(2, i);
135  Bu(5, l) = F(1, 1)*B(2, i) + F(1, 2)*B(1, i);
136 
137  Bu(0, m) = F(2, 0)*B(0, i);
138  Bu(1, m) = F(2, 1)*B(1, i);
139  Bu(2, m) = F(2, 2)*B(2, i);
140  Bu(3, m) = F(2, 0)*B(1, i) + F(2, 1)*B(0, i);
141  Bu(4, m) = F(2, 2)*B(0, i) + F(2, 0)*B(2, i);
142  Bu(5, m) = F(2, 1)*B(2, i) + F(2, 2)*B(1, i);
143  }
144  }
145  else if (modeLD == "updatedLagrangian")
146  {
147  for (int i = 0; i < e->numNodes(); i++)
148  {
149  m = 3 * i + 2;
150  k = m - 1;
151  l = m - 2;
152  Bu(0, l) = B(0, i);
153  Bu(1, k) = B(1, i);
154  Bu(2, m) = B(2, i);
155  Bu(3, l) = B(1, i);
156  Bu(3, k) = B(0, i);
157  Bu(4, l) = B(2, i);
158  Bu(4, m) = B(0, i);
159  Bu(5, k) = B(2, i);
160  Bu(5, m) = B(1, i);
161  }
162  }
163  }
164  return detJ;
165  }
166 
169  virtual void fillBnlMatrix(const GmElement* e, const GmShape* shape, const GmVector& coord,
170  const GmMatrix& X, const GmVector& N, const GmMatrix& J, GmMatrix& Bnl, const GmMatrix& F)
171  {
172  S_TRACE();
173  assert(e); Q_UNUSED(X); Q_UNUSED(N); Q_UNUSED(F);
174 
175  // integration point index
176  //int ip = mp->_index;
177  // mesh dimension
178  int d = GmpFemPhysicsCommon::nodeDim();
179  double k, l, m;
180  GmMatrix B;
181 
182  // example to read strain and stress
183  /*
184  // new strain
185  const GmVector e_new(mp->_newStrain->valueAt(e, ip, &coord), 2 * d);
186 
187  // new stresses
188  const GmVector s_new(mp->_newStress->valueAt(e, ip, &coord), 2*d);
189  */
190 
191  shape->shapeCartesianPartialsFromJacobian(coord, J, B, true);
192 
193  if (d == 2)
194  {
195  assert(Bnl.n_rows == d * d + 1 && Bnl.n_cols == e->numNodes()*d);
196  Bnl.zeros();
197  for (int i = 0; i < e->numNodes(); i++)
198  {
199  k = 2 * i + 1;
200  l = 2 * i;
201  Bnl(0, l) = B(0, i);
202  Bnl(1, l) = B(1, i);
203  Bnl(2, k) = B(0, i);
204  Bnl(3, k) = B(1, i);
205  }
206  }
207  else if (d == 3)
208  {
209  assert(Bnl.n_rows == d * d && Bnl.n_cols == e->numNodes()*d);
210  Bnl.zeros();
211  for (int i = 0; i < e->numNodes(); i++)
212  {
213  m = 3 * i + 2;
214  k = m - 1;
215  l = m - 2;
216 
217  Bnl(0, l) = B(0, i);
218  Bnl(1, l) = B(1, i);
219  Bnl(2, l) = B(2, i);
220  Bnl(3, k) = B(0, i);
221  Bnl(4, k) = B(1, i);
222  Bnl(5, k) = B(2, i);
223  Bnl(6, m) = B(0, i);
224  Bnl(7, m) = B(1, i);
225  Bnl(8, m) = B(2, i);
226  }
227  }
228  }
229 
232  virtual bool DeformationGradient(const GmElement* e, const GmMatrix& X0, const GmMatrix& Xt, const GmVector& N, const GmMatrix& J0, const GmMatrix& Jt, GmMatrix& F)
233  {
234  S_TRACE();
235  Q_UNUSED(e); Q_UNUSED(X0); Q_UNUSED(Xt); Q_UNUSED(N);
236 
237  // Deformation gradient
238  assert(F.n_rows == 3 && F.n_cols == 3);
239  F.zeros();
240 
241  // mesh dimension
242  int d = GmpFemPhysicsCommon::nodeDim();
243 
244  GmMatrix Ftemp;
245 
246  // Checking the dimension (2D, 3D)
247  switch (d) {
248  case 2:
249  if (e->type() == GM_TRI3 || e->type() == GM_TRI6)
250  {
251  //Obtaining jacobians in relation to a orthogonal set of coordinates
252  /* see Zienkiewicz, Taylor, David (2014) The Finite Element Method for Solid and Structural Mechanics
253  Section A.4*/
254  GmMatrix dXdzeta(2, 2);
255  dXdzeta = J0.submat(0, 1, 1, 2);
256  dXdzeta.each_row() -= J0.submat(2, 1, 2, 2);
257 
258  GmMatrix dxdzeta(2, 2);
259  dxdzeta = Jt.submat(0, 1, 1, 2);
260  dxdzeta.each_row() -= Jt.submat(2, 1, 2, 2);
261 
262  //Deformation gradient
263  Ftemp = arma::solve(dXdzeta, dxdzeta);
264  F(arma::span(0, 1), arma::span(0, 1)) = Ftemp.t();
265  }
266  else
267  {
268  // Deformation gradient
269  Ftemp = arma::solve(J0, Jt);
270  F(arma::span(0, 1), arma::span(0, 1)) = Ftemp.t();
271  }
272 
273  // Plane strain condition
274  F.at(2, 2) = 1;
275 
276  break;
277  case 3:
278  if (e->type() == GM_TET10 || e->type() == GM_TET4)
279  {
280  //Obtaining jacobians in relation to a orthogonal set of coordinates
281  /* see Zienkiewicz, Taylor, David (2014) The Finite Element Method for Solid and Structural Mechanics
282  Section A.5*/
283  GmMatrix dXdzeta(3, 3);
284  dXdzeta = J0.submat(0, 1, 2, 3);
285  dXdzeta.each_row() -= J0.submat(3, 1, 3, 3);
286 
287  GmMatrix dxdzeta(3, 3);
288  dxdzeta = Jt.submat(0, 1, 2, 3);
289  dxdzeta.each_row() -= Jt.submat(3, 1, 3, 3);
290 
291  //Deformation gradient
292  Ftemp = arma::solve(dXdzeta, dxdzeta);
293  F = Ftemp.t();
294  }
295  else
296  {
297  // Deformation gradient
298  Ftemp = arma::solve(J0, Jt);
299  F = Ftemp.t();
300  }
301  break;
302  }
303 
304  // Checking the gradient of deformation validity
305  return arma::det(F) > 0;
306  }
307 };
308 
309 template <class T> class GmpMechanicPlaneStressLargeDisplacement : public T
310 {
311 public:
312  GmpMechanicPlaneStressLargeDisplacement(const char* pluginType, GmSimulationData* simulation, QString id, QString description,
313  const GmpFemPhysicsCommonMaterialFactory* matFactory, const GmLogCategory& logger)
314  : T(pluginType, simulation, id, description, matFactory, logger)
315  {}
316 
317 protected:
318  // Comments on the base class
319  virtual double fillBuMatrix(const GmElement* e, const GmShape* shape, const GmVector& ncoord,
320  const GmMatrix& X, const GmVector& N, const GmMatrix& J, GmMatrix& Bu, const GmMatrix& F)
321  {
322  S_TRACE();
323  assert(e); Q_UNUSED(N); Q_UNUSED(X);
324 
325  double k, l;
326  GmMatrix B;
327  double detJ = shape->shapeCartesianPartialsFromJacobian(ncoord, J, B, true);
328 
329  assert(Bu.n_rows == 2 * nodeDim() && Bu.n_cols == e->numNodes()*nodeDim());
330  Bu.zeros();
331  GmMatrix L(nodeDim(), nodeDim());
332  L.zeros();
333 
334  // 2D solid element
335  assert(nodeDim() == 2);
336 
337  // Checking reference
338  QString modeLD = attribute(GmpMechanicalLargeDisplacementSolid::REFERENCE_LD_ID).toString();
339 
340  if (modeLD == "totalLagrangian")
341  {
342  for (int i = 0; i < e->numNodes(); i++)
343  {
344  k = 2 * i + 1;
345  l = k - 1;
346 
347  Bu(0, l) = F(0, 0)*B(0, i);
348  Bu(1, l) = F(0, 1)*B(1, i);
349  Bu(3, l) = F(0, 0)*B(1, i) + F(0, 1)*B(0, i);
350 
351  Bu(0, k) = F(1, 0)*B(0, i);
352  Bu(1, k) = F(1, 1)*B(1, i);
353  Bu(3, k) = F(1, 0)*B(1, i) + F(1, 1)*B(0, i);
354  }
355  }
356  else if (modeLD == "updatedLagrangian")
357  {
358  for (int i = 0; i < e->numNodes(); i++)
359  {
360  k = 2 * i + 1;
361  l = k - 1;
362  Bu(0, l) = B(0, i);
363  Bu(1, k) = B(1, i);
364  Bu(3, l) = B(1, i);
365  Bu(3, k) = B(0, i);
366  }
367  }
368  return detJ;
369  }
370 
371  virtual void fillBnlMatrix(const GmElement* e, const GmShape* shape, const GmVector& ncoord,
372  const GmMatrix& X, const GmVector& N, const GmMatrix& J, GmMatrix& Bnl, const GmMatrix& F)
373  {
374  S_TRACE();
375  assert(e); Q_UNUSED(X); Q_UNUSED(N); Q_UNUSED(F);
376 
377  double k, l;
378  GmMatrix B;
379  shape->shapeCartesianPartialsFromJacobian(ncoord, J, B, true);
380 
381  // 2D solid element
382  assert(nodeDim() == 2);
383  assert(Bnl.n_rows == nodeDim() * nodeDim() + 1 && Bnl.n_cols == e->numNodes()*nodeDim());
384  Bnl.zeros();
385  for (int i = 0; i < e->numNodes(); i++)
386  {
387  k = 2 * i + 1;
388  l = 2 * i;
389  Bnl(0, l) = B(0, i);
390  Bnl(1, l) = B(1, i);
391  Bnl(2, k) = B(0, i);
392  Bnl(3, k) = B(1, i);
393  }
394  }
395 
396  virtual bool DeformationGradient(const GmElement* e, const GmMatrix& X0, const GmMatrix& Xt, const GmVector& N, const GmMatrix& J0, const GmMatrix& Jt, GmMatrix& F)
397  {
398  S_TRACE();
399  Q_UNUSED(e); Q_UNUSED(X0); Q_UNUSED(Xt); Q_UNUSED(N);
400 
401  // Deformation gradient
402  assert(F.n_rows == 3 && F.n_cols == 3);
403  F.zeros();
404 
405  // Deformation gradient transpose
406  GmMatrix Ftemp;
407 
408  if (e->type() == GM_TRI3 || e->type() == GM_TRI6)
409  {
410  //Obtaining jacobians in relation to a orthogonal set of coordinates
411  /* see Zienkiewicz, Taylor, David (2014) The Finite Element Method for Solid and Structural Mechanics
412  Section A.4*/
413  GmMatrix dXdzeta(2, 2);
414  dXdzeta = J0.submat(0, 1, 1, 2);
415  dXdzeta.each_row() -= J0.submat(2, 1, 2, 2);
416 
417  GmMatrix dxdzeta(2, 2);
418  dxdzeta = Jt.submat(0, 1, 1, 2);
419  dxdzeta.each_row() -= Jt.submat(2, 1, 2, 2);
420 
421  //Deformation gradient
422  Ftemp = arma::solve(dXdzeta, dxdzeta);
423  F(arma::span(0, 1), arma::span(0, 1)) = Ftemp.t();
424  }
425  else
426  {
427  // Deformation gradient
428  Ftemp = arma::solve(J0, Jt);
429  F(arma::span(0, 1), arma::span(0, 1)) = Ftemp.t();
430  }
431 
432  F.at(2, 2) = 1;
433 
434  // Checking the gradient of deformation validity
435  return arma::det(F) > 0;
436  }
437 
438  // Comments on the base class
439  virtual bool isPlaneStress() { return true; }
440 };
441 
443 template <class T> class GmpMechanicAxyLargeDisplacement : public T
444 {
445 public:
446  GmpMechanicAxyLargeDisplacement(const char* pluginType, GmSimulationData* simulation, QString id, QString description,
447  const GmpFemPhysicsCommonMaterialFactory* matFactory, const GmLogCategory& logger)
448  : T(pluginType, simulation, id, description, matFactory, logger)
449  {}
450 
451 protected:
452  // Comments on the base class
453  virtual double fillBuMatrix(const GmElement* e, const GmShape* shape, const GmVector& ncoord,
454  const GmMatrix& X, const GmVector& N, const GmMatrix& J, GmMatrix& Bu, const GmMatrix& F)
455  {
456  S_TRACE();
457  assert(e);
458 
459  double k, l;
460  GmMatrix B;
461  double detJ = shape->shapeCartesianPartialsFromJacobian(ncoord, J, B, true);
462 
463  double rGauss = 0.0;
464  for (int j = 0; j < e->numNodes(); ++j)
465  rGauss += N(j)*X(j, 0);
466 
467  assert(Bu.n_rows == 2 * nodeDim() && Bu.n_cols == e->numNodes()*nodeDim());
468  Bu.zeros();
469  GmMatrix L(nodeDim() + 1, nodeDim() + 1);
470  L.zeros();
471 
472  // 2D solid element
473  assert(nodeDim() == 2);
474 
475  // Checking reference
476  QString modeLD = attribute(GmpMechanicalLargeDisplacementSolid::REFERENCE_LD_ID).toString();
477 
478  if (modeLD == "totalLagrangian")
479  {
480  for (int i = 0; i < e->numNodes(); i++)
481  {
482  k = 2 * i + 1;
483  l = 2 * i;
484 
485  Bu(0, l) = F(0, 0)*B(0, i);
486  Bu(1, l) = F(0, 1)*B(1, i);
487  Bu(2, l) = F(2, 2)*N(i) / rGauss;
488  Bu(3, l) = F(0, 0)*B(1, i) + F(0, 1)*B(0, i);
489 
490  Bu(0, k) = F(1, 0)*B(0, i);
491  Bu(1, k) = F(1, 1)*B(1, i);
492  Bu(3, k) = F(1, 0)*B(1, i) + F(1, 1)*B(0, i);
493 
494  }
495  }
496  else if (modeLD == "updatedLagrangian")
497  {
498  for (int i = 0; i < e->numNodes(); i++)
499  {
500  k = 2 * i + 1;
501  l = 2 * i;
502  Bu(0, l) = B(0, i);
503  Bu(1, k) = B(1, i);
504  Bu(2, l) = N(i) / rGauss;
505  Bu(3, l) = B(1, i);
506  Bu(3, k) = B(0, i);
507  }
508  }
509 
510  return detJ;
511  }
512 
513  virtual void fillBnlMatrix(const GmElement* e, const GmShape* shape, const GmVector& ncoord,
514  const GmMatrix& X, const GmVector& N, const GmMatrix& J, GmMatrix& Bnl, const GmMatrix& F)
515  {
516  S_TRACE();
517  assert(e);
518  Q_UNUSED(F);
519 
520  double k, l;
521  GmMatrix B;
522  shape->shapeCartesianPartialsFromJacobian(ncoord, J, B, true);
523 
524  double rGauss = 0.0;
525  for (int j = 0; j < e->numNodes(); ++j)
526  rGauss += N(j)*X(j, 0);
527 
528  assert(Bnl.n_rows == 2 * nodeDim() && Bnl.n_cols == e->numNodes() * nodeDim());
529  Bnl.zeros();
530  // 2D solid element
531  assert(nodeDim() == 2);
532 
533  for (int i = 0; i < e->numNodes(); i++)
534  {
535  k = 2 * i + 1;
536  l = 2 * i;
537  Bnl(0, l) = B(0, i);
538  Bnl(1, l) = B(1, i);
539  Bnl(2, k) = B(0, i);
540  Bnl(3, k) = B(1, i);
541  Bnl(4, l) = N(i) / rGauss;
542  }
543 
544  //return detJ;
545  }
546 
547  virtual bool DeformationGradient(const GmElement* e, const GmMatrix& X0, const GmMatrix& Xt, const GmVector& N, const GmMatrix& J0, const GmMatrix& Jt, GmMatrix& F)
548  {
549  S_TRACE();
550  assert(e);
551 
552  // Deformation gradient
553  assert(F.n_rows == 3 && F.n_cols == 3);
554  F.zeros();
555 
556  // Radius of gauss point in initial and current configuration
557  double RGauss = 0.0;
558  double rGauss = 0.0;
559  for (int j = 0; j < e->numNodes(); ++j)
560  {
561  RGauss += N(j)*X0(j, 0);
562  rGauss += N(j)*Xt(j, 0);
563  }
564 
565  GmMatrix Ftemp;
566 
567  if (e->type() == GM_TRI3 || e->type() == GM_TRI6)
568  {
569  //Obtaining jacobians in relation to a orthogonal set of coordinates
570  /* see Zienkiewicz, Taylor, David (2014) The Finite Element Method for Solid and Structural Mechanics, Section A.4 */
571  GmMatrix dXdzeta(2, 2);
572  dXdzeta = J0.submat(0, 1, 1, 2);
573  dXdzeta.each_row() -= J0.submat(2, 1, 2, 2);
574 
575  GmMatrix dxdzeta(2, 2);
576  dxdzeta = Jt.submat(0, 1, 1, 2);
577  dxdzeta.each_row() -= Jt.submat(2, 1, 2, 2);
578 
579  // Deformation gradient
580  Ftemp = arma::solve(dXdzeta, dxdzeta);
581  F(arma::span(0, 1), arma::span(0, 1)) = Ftemp.t();
582  }
583  else
584  {
585  // Deformation gradient
586  Ftemp = arma::solve(J0, Jt);
587  F(arma::span(0, 1), arma::span(0, 1)) = Ftemp.t();
588  }
589 
590  // Deformation gradient for axissymetric
591  F.at(2, 2) = rGauss / RGauss;
592 
593  // Checking the gradient of deformation validity
594  return arma::det(F) > 0;
595  }
596 
597  // Comments on the base class
598  virtual double axisSymetricFactor(const GmElement* e, const GmMatrix& X, const GmVector& N)
599  {
600  int n = e->numNodes();
601 
602  double rGauss = 0.0;
603  for (int j = 0; j < n; ++j)
604  rGauss += N(j)*X(j, 0);
605  return 2 * arma::datum::pi *rGauss;
606  }
607 
608 };
609 
610 #endif
Declaration of the GmpMechanicalPhysics class.
#define S_TRACE()
GM_TET4
virtual bool DeformationGradient(const GmElement *e, const GmMatrix &X0, const GmMatrix &Xt, const GmVector &N, const GmMatrix &J0, const GmMatrix &Jt, GmMatrix &F)
Fills the F deformation gradient matrix.
Definition: gmpMechanicLargeDisplacement.h:232
Definition: gmpMechanicLargeDisplacement.h:309
GM_TRI6
GM_TET10
virtual void fillBnlMatrix(const GmElement *e, const GmShape *shape, const GmVector &coord, const GmMatrix &X, const GmVector &N, const GmMatrix &J, GmMatrix &Bnl, const GmMatrix &F)
Fills the Bnl strain- displacement matrix for large displacements.
Definition: gmpMechanicLargeDisplacement.h:169
GM_TRI3
Declaration of the GmpMechanicPoint class.
virtual GmCellType type() const=0
virtual int numNodes() const=0
virtual double shapeCartesianPartialsFromJacobian(const GmVector &ncoord, const GmMatrix &J, GmMatrix &dN, bool transposed=false) const
Id for retrieving the reference large displacement physics attribute.
Definition: gmpMechanicalLargeDisplacementSolid.h:50
Mechanic Axissymetric for Quad and Tri, Elements.
Definition: gmpMechanicLargeDisplacement.h:443
arma::vec GmVector
arma::mat GmMatrix
virtual double fillBuMatrix(const GmElement *e, const GmShape *shape, const GmVector &coord, const GmMatrix &X, const GmVector &N, const GmMatrix &J, GmMatrix &Bu, const GmMatrix &F)
Fills the Bu strain- displacement matrix for large displacements.
Definition: gmpMechanicLargeDisplacement.h:47
Mechanic 2D & 3D for Quad, Tri, Hex and Tet Elements.
Definition: gmpMechanicLargeDisplacement.h:36