GetFEM  5.5
nonlinear_elastostatic.cc
Go to the documentation of this file.
1 /*===========================================================================
2 
3  Copyright (C) 2002-2026 Yves Renard, Julien Pommier.
4 
5  This file is a part of GetFEM
6 
7  GetFEM is free software; you can redistribute it and/or modify it
8  under the terms of the GNU Lesser General Public License as published
9  by the Free Software Foundation; either version 3 of the License, or
10  (at your option) any later version along with the GCC Runtime Library
11  Exception either version 3.1 or (at your option) any later version.
12  This program is distributed in the hope that it will be useful, but
13  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
14  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15  License and GCC Runtime Library Exception for more details.
16  You should have received a copy of the GNU Lesser General Public License
17  along with this program. If not, see https://www.gnu.org/licenses/.
18 
19 ===========================================================================*/
20 
21 /**
22  @file nonlinear_elastostatic.cc
23  @brief Nonlinear Elastostatic problem (large strain).
24 
25  A rubber bar is submitted to a large torsion.
26 
27  This program is used to check that getfem++ is working. This is also
28  a good example of use of GetFEM.
29 */
30 
31 #include "getfem/getfem_export.h"
35 #include "gmm/gmm.h"
36 
37 using std::endl; using std::cout; using std::cerr;
38 using std::ends; using std::cin;
39 
40 /* some GetFEM types that we will be using */
41 using bgeot::base_small_vector; /* special class for small (dim<16) vectors */
42 using bgeot::base_node; /* geometrical nodes(derived from base_small_vector)*/
43 using bgeot::base_vector;
44 using bgeot::scalar_type; /* = double */
45 using bgeot::size_type; /* = unsigned long */
46 using bgeot::dim_type;
47 using bgeot::base_matrix; /* small dense matrix. */
48 
49 /* definition of some matrix/vector types. These ones are built
50  * using the predefined types in Gmm++
51  */
53 typedef getfem::modeling_standard_sparse_matrix sparse_matrix;
54 typedef getfem::modeling_standard_plain_vector plain_vector;
55 
56 /*
57  structure for the elastostatic problem
58 */
59 struct elastostatic_problem {
60 
61  enum { DIRICHLET_BOUNDARY_NUM = 0, NEUMANN_BOUNDARY_NUM = 1};
62  getfem::mesh mesh; /* the mesh */
63  getfem::mesh_im mim; /* the integration methods */
64  getfem::mesh_fem mf_u; /* main mesh_fem, for the elastostatic solution */
65  getfem::mesh_fem mf_p; /* mesh_fem for the pressure. */
66  getfem::mesh_fem mf_rhs; /* mesh_fem for the right hand side (f(x),..) */
67  getfem::mesh_fem mf_coef; /* mesh_fem used to represent pde coefficients */
68  scalar_type p1, p2, p3; /* elastic coefficients. */
69 
70  scalar_type residual; /* max residual for the iterative solvers */
71 
72  std::string datafilename;
73  bgeot::md_param PARAM;
74 
75  bool solve(plain_vector &U);
76  void init(void);
77  elastostatic_problem(void) : mim(mesh), mf_u(mesh), mf_p(mesh), mf_rhs(mesh),
78  mf_coef(mesh) {}
79 };
80 
81 
82 /* Read parameters from the .param file, build the mesh, set finite element
83  and integration methods and selects the boundaries.
84 
85  (this is boilerplate code, not very interesting)
86  */
87 void elastostatic_problem::init(void) {
88  std::string MESH_TYPE = PARAM.string_value("MESH_TYPE","Mesh type ");
89  std::string FEM_TYPE = PARAM.string_value("FEM_TYPE","FEM name");
90  std::string FEM_TYPE_P= PARAM.string_value("FEM_TYPE_P",
91  "FEM name for the pressure");
92  std::string INTEGRATION = PARAM.string_value("INTEGRATION",
93  "Name of integration method");
94  cout << "MESH_TYPE=" << MESH_TYPE << "\n";
95  cout << "FEM_TYPE=" << FEM_TYPE << "\n";
96  cout << "INTEGRATION=" << INTEGRATION << "\n";
97 
98  /* First step : build the mesh */
101  size_type N = pgt->dim();
102  std::vector<size_type> nsubdiv(N);
103  std::fill(nsubdiv.begin(),nsubdiv.end(),
104  PARAM.int_value("NX", "Nomber of space steps "));
105  nsubdiv[1] = PARAM.int_value("NY") ? PARAM.int_value("NY") : nsubdiv[0];
106  if (N>2)
107  nsubdiv[2] = PARAM.int_value("NZ") ? PARAM.int_value("NZ") : nsubdiv[0];
108  getfem::regular_unit_mesh(mesh, nsubdiv, pgt,
109  PARAM.int_value("MESH_NOISED") != 0);
110 
111  bgeot::base_matrix M(N,N);
112  for (size_type i=0; i < N; ++i) {
113  static const char *t[] = {"LX","LY","LZ"};
114  M(i,i) = (i<3) ? PARAM.real_value(t[i],t[i]) : 1.0;
115  }
116  if (N>1) { M(0,1) = PARAM.real_value("INCLINE") * PARAM.real_value("LY"); }
117 
118  /* scale the unit mesh to [LX,LY,..] and incline it */
119  mesh.transformation(M);
120 
121  datafilename = PARAM.string_value("ROOTFILENAME","Base name of data files.");
122  residual = PARAM.real_value("RESIDUAL"); if (residual == 0.) residual = 1e-10;
123 
124  p1 = PARAM.real_value("P1", "First Elastic coefficient");
125  p2 = PARAM.real_value("P2", "Second Elastic coefficient");
126  p3 = PARAM.real_value("P3", "Third Elastic coefficient");
127 
128  mf_u.set_qdim(dim_type(N));
129 
130  /* set the finite element on the mf_u */
131  getfem::pfem pf_u =
132  getfem::fem_descriptor(FEM_TYPE);
133  getfem::pintegration_method ppi =
134  getfem::int_method_descriptor(INTEGRATION);
135 
136  mim.set_integration_method(ppi);
137  mf_u.set_finite_element(pf_u);
138 
139  mf_p.set_finite_element(getfem::fem_descriptor(FEM_TYPE_P));
140 
141  /* set the finite element on mf_rhs (same as mf_u is DATA_FEM_TYPE is
142  not used in the .param file */
143  std::string data_fem_name = PARAM.string_value("DATA_FEM_TYPE");
144  if (data_fem_name.size() == 0) {
145  GMM_ASSERT1(pf_u->is_lagrange(), "You are using a non-lagrange FEM"
146  ". In that case you need to set "
147  << "DATA_FEM_TYPE in the .param file");
148  mf_rhs.set_finite_element(mesh.convex_index(), pf_u);
149  } else {
150  mf_rhs.set_finite_element(mesh.convex_index(),
151  getfem::fem_descriptor(data_fem_name));
152  }
153 
154  /* set the finite element on mf_coef. Here we use a very simple element
155  * since the only function that need to be interpolated on the mesh_fem
156  * is f(x)=1 ... */
157  mf_coef.set_finite_element(mesh.convex_index(),
158  getfem::classical_fem(pgt,0));
159 
160  /* set boundary conditions
161  * (Neuman on the upper face, Dirichlet elsewhere) */
162  cout << "Selecting Neumann and Dirichlet boundaries\n";
163  getfem::mesh_region border_faces;
164  getfem::outer_faces_of_mesh(mesh, border_faces);
165  for (getfem::mr_visitor it(border_faces); !it.finished(); ++it) {
166  assert(it.is_face());
167  base_node un = mesh.normal_of_face_of_convex(it.cv(), it.f());
168  un /= gmm::vect_norm2(un);
169  if (gmm::abs(un[N-1] - 1.0) < 1.0E-7) {
170  mesh.region(DIRICHLET_BOUNDARY_NUM).add(it.cv(),it.f());
171  } else if (gmm::abs(un[N-1] + 1.0) < 1.0E-7) {
172  mesh.region(DIRICHLET_BOUNDARY_NUM).add(it.cv(),it.f());
173  }
174  }
175 }
176 
177 /**************************************************************************/
178 /* Model. */
179 /**************************************************************************/
180 
181 bool elastostatic_problem::solve(plain_vector &U) {
182  size_type nb_dof_rhs = mf_rhs.nb_dof();
183  size_type N = mesh.dim();
184  size_type law_num = PARAM.int_value("LAW");
185  // Linearized elasticity brick.
186  base_vector p(3); p[0] = p1; p[1] = p2; p[2] = p3;
187 
188  /* choose the material law */
189  std::string lawname;
190  switch (law_num) {
191  case 0:
192  case 1: lawname = "Saint_Venant_Kirchhoff"; p.resize(2); break;
193  case 2: lawname = "Ciarlet_Geymonat"; p.resize(3); break;
194  case 3: lawname = "Incompressible_Mooney_Rivlin"; p.resize(2); break;
195  default: GMM_ASSERT1(false, "no such law");
196  }
197 
198  if (N == 2) {
199  cout << "2D plane strain hyper-elasticity\n";
200  lawname = "Plane_Strain_"+lawname;
201  }
202 
203  getfem::model model;
204 
205  // Main unknown of the problem (displacement).
206  model.add_fem_variable("u", mf_u);
207 
208  // Nonlinear elasticity brick
209  model.add_initialized_fixed_size_data("params", p);
210  getfem::add_finite_strain_elasticity_brick(model, mim, lawname, "u","params");
211 
212  // Incompressibility
213  if (law_num == 1 || law_num == 3) {
214  model.add_fem_variable("p", mf_p);
216  }
217 
218  // Defining the volumic source term.
219  base_vector f(N);
220  f[0] = PARAM.real_value("FORCEX","Amplitude of the gravity");
221  f[1] = PARAM.real_value("FORCEY","Amplitude of the gravity");
222  if (N>2)
223  f[2] = PARAM.real_value("FORCEZ","Amplitude of the gravity");
224  plain_vector F(nb_dof_rhs * N);
225  for (size_type i = 0; i < nb_dof_rhs; ++i) {
226  gmm::copy(f, gmm::sub_vector(F, gmm::sub_interval(i*N, N)));
227  }
228  // Volumic source term brick.
229  model.add_initialized_fem_data("VolumicData", mf_rhs, F);
230  getfem::add_source_term_brick(model, mim, "u", "VolumicData");
231 
232  // Dirichlet condition
233  plain_vector F2(nb_dof_rhs * N);
234  model.add_initialized_fem_data("DirichletData", mf_rhs, F2);
235  if (PARAM.int_value("DIRICHLET_VERSION") == 0)
237  (model, mim, "u", mf_u, DIRICHLET_BOUNDARY_NUM, "DirichletData");
238  else
240  (model, mim, "u", 1E15, DIRICHLET_BOUNDARY_NUM, "DirichletData");
241 
242 
243  gmm::iteration iter;
244  gmm::set_traces_level(1); // To avoid some trace messages
245 
246 
247  /* prepare the export routine for OpenDX (all time steps will be exported)
248  (can be viewed with "dx -edit nonlinear_elastostatic.net")
249  */
250  getfem::dx_export exp(datafilename + ".dx",
251  PARAM.int_value("VTK_EXPORT")==1);
253  sl.build(mesh, getfem::slicer_boundary(mesh),8);
254  exp.exporting(sl,true); exp.exporting_mesh_edges();
255  //exp.begin_series("deformationsteps");
256  exp.write_point_data(mf_u, U, "stepinit");
257  exp.serie_add_object("deformationsteps");
258 
259  GMM_ASSERT1(!mf_rhs.is_reduced(), "To be adapted for reduced mesh_fems");
260 
261  int nb_step = int(PARAM.int_value("NBSTEP"));
262  size_type maxit = PARAM.int_value("MAXITER");
263 
264  for (int step = 0; step < nb_step; ++step) {
265  plain_vector DF(F);
266 
267  gmm::copy(gmm::scaled(F, (step+1.)/(scalar_type)nb_step), DF);
268  gmm::copy(DF, model.set_real_variable("VolumicData"));
269 
270  if (N>2) {
271  /* Apply the gradual torsion/extension */
272  scalar_type torsion = PARAM.real_value("TORSION",
273  "Amplitude of the torsion");
274  torsion *= (step+1)/scalar_type(nb_step);
275  scalar_type extension = PARAM.real_value("EXTENSION",
276  "Amplitude of the extension");
277  extension *= (step+1)/scalar_type(nb_step);
278  base_node G(N); G[0] = G[1] = 0.5;
279  for (size_type i = 0; i < nb_dof_rhs; ++i) {
280  const base_node P = mf_rhs.point_of_basic_dof(i) - G;
281  scalar_type r = sqrt(P[0]*P[0]+P[1]*P[1]),
282  theta = atan2(P[1],P[0]);
283  F2[i*N+0] = r*cos(theta + (torsion*P[2])) - P[0];
284  F2[i*N+1] = r*sin(theta + (torsion*P[2])) - P[1];
285  F2[i*N+2] = extension * P[2];
286  }
287  }
288  /* update the imposed displacement */
289  gmm::copy(F2, model.set_real_variable("DirichletData"));
290 
291  cout << "step " << step << ", number of variables : "
292  << model.nb_dof() << endl;
293  iter = gmm::iteration(residual, int(PARAM.int_value("NOISY")),
294  maxit ? maxit : 40000);
295 
296  /* let the non-linear solve (Newton) do its job */
297  getfem::standard_solve(model, iter);
298 
299  gmm::copy(model.real_variable("u"), U);
300  //char s[100]; snprintf(s, 100, "step%d", step+1);
301 
302  /* append the new displacement to the exported opendx file */
303  exp.write_point_data(mf_u, U); //, s);
304  exp.serie_add_object("deformationsteps");
305  }
306 
307  // Solution extraction
308  gmm::copy(model.real_variable("u"), U);
309 
310  return (iter.converged());
311 }
312 
313 /**************************************************************************/
314 /* main program. */
315 /**************************************************************************/
316 
317 int main(int argc, char *argv[]) {
318 
319  GETFEM_MPI_INIT(argc, argv);
320  GMM_SET_EXCEPTION_DEBUG; // Exceptions make a memory fault, to debug.
321  FE_ENABLE_EXCEPT; // Enable floating point exception for Nan.
322 
323  try {
324  elastostatic_problem p;
325  p.PARAM.read_command_line(argc, argv);
326  p.init();
327  p.mesh.write_to_file(p.datafilename + ".mesh");
328  p.mf_u.write_to_file(p.datafilename + ".mf", true);
329  p.mf_rhs.write_to_file(p.datafilename + ".mfd", true);
330  plain_vector U(p.mf_u.nb_dof());
331  GMM_ASSERT1(p.solve(U), "Solve has failed");
332  if (p.PARAM.int_value("VTK_EXPORT")) {
333  gmm::vecsave(p.datafilename + ".U", U);
334  cout << "export to " << p.datafilename + ".vtk" << "..\n";
335  getfem::vtk_export exp(p.datafilename + ".vtk",
336  p.PARAM.int_value("VTK_EXPORT")==1);
337  exp.exporting(p.mf_u);
338  exp.write_point_data(p.mf_u, U, "elastostatic_displacement");
339  cout << "export done, you can view the data file with (for example)\n"
340  "mayavi2 -d " << p.datafilename << ".vtk -f ExtractVectorNorm -f "
341  "WarpVector -m Surface -m Outline\n";
342  }
343  }
344  GMM_STANDARD_CATCH_ERROR;
345 
346  GETFEM_MPI_FINALIZE;
347 
348  return 0;
349 }
A (quite large) class for exportation of data to IBM OpenDX.
Describe a finite element method linked to a mesh.
Describe an integration method linked to a mesh.
"iterator" class for regions.
structure used to hold a set of convexes and/or convex faces.
Describe a mesh (collection of convexes (elements) and points).
Definition: getfem_mesh.h:98
`‘Model’' variables store the variables, the data and the description of a model.
size_type nb_dof(bool with_internal=false) const
Total number of degrees of freedom in the model.
void add_fem_variable(const std::string &name, const mesh_fem &mf, size_type niter=1)
Add a variable being the dofs of a finite element method to the model.
void add_initialized_fixed_size_data(const std::string &name, const VECT &v)
Add a fixed size data (assumed to be a vector) to the model and initialized with v.
void add_initialized_fem_data(const std::string &name, const mesh_fem &mf, const VECT &v)
Add an initialized fixed size data to the model, assumed to be a vector field if the size of the vect...
model_real_plain_vector & set_real_variable(const std::string &name, size_type niter) const
Gives the write access to the vector value of a variable.
const model_real_plain_vector & real_variable(const std::string &name, size_type niter) const
Gives the access to the vector value of a variable.
Extraction of the boundary of a slice.
The output of a getfem::mesh_slicer which has been recorded.
void build(const getfem::mesh &m, const slicer_action &a, size_type nrefine=1)
Build the slice, by applying a slicer_action operation.
VTK/VTU export.
Definition: getfem_export.h:67
The Iteration object calculates whether the solution has reached the desired accuracy,...
Definition: gmm_iter.h:52
sparse vector built upon std::vector.
Definition: gmm_vector.h:995
Export solutions to various formats.
Standard solvers for model bricks.
Non-linear elasticty and incompressibility bricks.
Build regular meshes.
Include common gmm files.
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:976
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
Definition: gmm_blas.h:556
void APIDECL outer_faces_of_mesh(const mesh &m, const dal::bit_vector &cvlst, convex_face_ct &flist)
returns a list of "exterior" faces of a mesh (i.e.
Definition: getfem_mesh.cc:821
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
Definition: getfem_fem.h:243
pfem fem_descriptor(const std::string &name)
get a fem descriptor from its string name.
Definition: getfem_fem.cc:4659
pfem classical_fem(bgeot::pgeometric_trans pgt, short_type k, bool complete=false)
Give a pointer on the structures describing the classical polynomial fem of degree k on a given conve...
Definition: getfem_fem.cc:4566
pgeometric_trans geometric_trans_descriptor(std::string name)
Get the geometric transformation from its string name.
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
size_type APIDECL add_Dirichlet_condition_with_multipliers(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region, const std::string &dataname=std::string())
Add a Dirichlet condition on the variable varname and the mesh region region.
size_type add_finite_strain_elasticity_brick(model &md, const mesh_im &mim, const std::string &lawname, const std::string &varname, const std::string &params, size_type region=size_type(-1))
Add a finite strain elasticity brick to the model with respect to the variable varname (the displacem...
size_type add_finite_strain_incompressibility_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region=size_type(-1))
Add a finite strain incompressibility term (for large strain elasticity) to the model with respect to...
void regular_unit_mesh(mesh &m, std::vector< size_type > nsubdiv, bgeot::pgeometric_trans pgt, bool noised=false)
Build a regular mesh of the unit square/cube/, etc.
pintegration_method int_method_descriptor(std::string name, bool throw_if_not_found=true)
Get an integration method from its name .
size_type APIDECL add_Dirichlet_condition_with_penalization(model &md, const mesh_im &mim, const std::string &varname, scalar_type penalization_coeff, size_type region, const std::string &dataname=std::string(), const mesh_fem *mf_mult=0)
Add a Dirichlet condition on the variable varname and the mesh region region.
size_type APIDECL add_source_term_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region=size_type(-1), const std::string &directdataname=std::string())
Add a source term on the variable varname.
void standard_solve(model &md, gmm::iteration &iter, rmodel_plsolver_type lsolver, abstract_newton_line_search &ls)
A default solver for the model brick system.