GetFEM  5.5
plate.cc
Go to the documentation of this file.
1 /*===========================================================================
2 
3  Copyright (C) 2002-2026 Yves Renard, Michel Salaün.
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 plate.cc
23  @brief Linear Elastostatic plate problem.
24 
25  This program is used to check that getfem++ is working. This is
26  also a good example of use of GetFEM.
27 */
28 
29 #include "getfem/getfem_assembling.h" /* import assembly methods (and norms comp.) */
31 #include "getfem/getfem_export.h" /* export functions (save solution in a file) */
34 #include "gmm/gmm.h"
35 
36 using std::endl; using std::cout; using std::cerr;
37 using std::ends; using std::cin;
38 
39 /* some GetFEM types that we will be using */
40 using bgeot::base_small_vector; /* special class for small (dim<16) vectors */
41 using bgeot::base_node; /* geometrical nodes(derived from base_small_vector)*/
42 using bgeot::scalar_type; /* = double */
43 using bgeot::size_type; /* = unsigned long */
44 using bgeot::dim_type;
45 using bgeot::base_matrix; /* small dense matrix. */
46 
47 /* definition of some matrix/vector types. These ones are built
48  * using the predefined types in Gmm++
49  */
51 typedef getfem::modeling_standard_sparse_matrix sparse_matrix;
52 typedef getfem::modeling_standard_plain_vector plain_vector;
53 
54 /*
55  structure for the elastostatic problem
56 */
57 struct plate_problem {
58 
59  enum { SIMPLY_FIXED_BOUNDARY_NUM = 0 };
60  getfem::mesh mesh; /* the mesh */
61  getfem::mesh_im mim, mim_subint;
62  getfem::mesh_fem mf_ut;
63  getfem::mesh_fem mf_u3;
64  getfem::mesh_fem mf_theta;
65  getfem::mesh_fem mf_rhs; /* mesh_fem for the right hand side (f(x),..) */
66  getfem::mesh_fem mf_coef; /* mesh_fem used to represent pde coefficients */
67  scalar_type lambda, mu; /* Lamé coefficients. */
68  scalar_type E, nu; /* Lamé coefficients. */
69  scalar_type epsilon; /* thickness of the plate. */
70  scalar_type pressure;
71  scalar_type residual; /* max residual for the iterative solvers */
72  scalar_type LX , LY ; // default : LX = LY = 1
73  bool mitc;
74  int sol_ref; // sol_ref = 0 : simple support on the vertical edges
75  // sol_ref = 1 : homogeneous on the vertical edges
76  // sol_ref = 2 : homogeneous on the 4 vertical
77  // edges with solution u3 = sin²(x)*sin²(y)
78  scalar_type eta; // useful only if sol_ref == 2 :
79  // eta = 0 => Kirchhoff-Love
80  // eta = small => Mindlin
81  size_type N_Four ;
82  base_matrix theta1_Four, theta2_Four, u3_Four ;
83 
84  int study_flag; // if studyflag = 1, then the loadings applied are chosen
85  // in order to have a maximal vertical displacement equal to one.
86  // Nothing is done if study_flag has another value.
87 
88  std::string datafilename;
89  bgeot::md_param PARAM;
90 
91  base_small_vector theta_exact(base_node P);
92  scalar_type u3_exact(base_node P);
93 
94  bool solve(plain_vector &Ut, plain_vector &U3, plain_vector &THETA);
95  void init(void);
96  void compute_error(plain_vector &Ut, plain_vector &U3, plain_vector &THETA);
97  plate_problem(void) : mim(mesh), mim_subint(mesh), mf_ut(mesh), mf_u3(mesh),
98  mf_theta(mesh), mf_rhs(mesh), mf_coef(mesh) {}
99 };
100 
101 /* Read parameters from the .param file, build the mesh, set finite element
102  * and integration methods and selects the boundaries.
103  */
104 void plate_problem::init(void) {
105  std::string MESH_FILE = PARAM.string_value("MESH_FILE");
106  std::string MESH_TYPE = PARAM.string_value("MESH_TYPE","Mesh type ");
107  std::string FEM_TYPE_UT = PARAM.string_value("FEM_TYPE_UT","FEM name");
108  std::string FEM_TYPE_U3 = PARAM.string_value("FEM_TYPE_U3","FEM name");
109  std::string FEM_TYPE_THETA = PARAM.string_value("FEM_TYPE_THETA","FEM name");
110  std::string INTEGRATION = PARAM.string_value("INTEGRATION",
111  "Name of integration method");
112  std::string INTEGRATION_CT = PARAM.string_value("INTEGRATION_CT",
113  "Name of integration method");
114  cout << "MESH_TYPE=" << MESH_TYPE << "\n";
115  cout << "FEM_TYPE_UT=" << FEM_TYPE_UT << "\n";
116  cout << "INTEGRATION=" << INTEGRATION << "\n";
117  cout << "INTEGRATION_CT=" << INTEGRATION_CT << "\n";
118 
119  /* First step : build the mesh */
120  size_type N;
123  if (!MESH_FILE.empty()) {
124  cout << "MESH_FILE=" << MESH_FILE << "\n";
125  mesh.read_from_file(MESH_FILE);
127  (mesh.trans_of_convex(mesh.convex_index().first_true()));
128  cout << "MESH_TYPE=" << MESH_TYPE << "\n";
129  N = mesh.dim();
130  } else {
131  N = pgt->dim();
132  GMM_ASSERT1(N == 2, "For a plate problem, N should be 2");
133  std::vector<size_type> nsubdiv(N);
134  std::fill(nsubdiv.begin(),nsubdiv.end(),
135  PARAM.int_value("NX", "Number of space steps "));
136  getfem::regular_unit_mesh(mesh, nsubdiv, pgt,
137  PARAM.int_value("MESH_NOISED") != 0);
138  }
139 
140  datafilename = PARAM.string_value("ROOTFILENAME","Base name of data files.");
141  residual = PARAM.real_value("RESIDUAL"); if (residual == 0.) residual = 1e-10;
142  mitc = (PARAM.int_value("MITC", "Mitc version ?") != 0);
143 
144 
145  cout << "MITC = " ;
146  if (mitc) cout << "true \n" ; else cout << "false \n" ;
147  sol_ref = int(PARAM.int_value("SOL_REF") ) ;
148  study_flag = int(PARAM.int_value("STUDY_FLAG") ) ;
149  eta = (PARAM.real_value("ETA") );
150  N_Four = (PARAM.int_value("N_Four") ) ;
151 
152 
153  LX = PARAM.real_value("LX");
154  LY = PARAM.real_value("LY");
155  mu = PARAM.real_value("MU", "Lamé coefficient mu");
156  lambda = PARAM.real_value("LAMBDA", "Lamé coefficient lambda");
157  epsilon = PARAM.real_value("EPSILON", "thickness of the plate");
158  pressure = PARAM.real_value("PRESSURE",
159  "pressure on the top surface of the plate.");
160 
161 
162  cout << "SOL_REF = " ;
163  if (sol_ref==0) cout << "appui simple aux 2 bords verticaux\n" ;
164  if (sol_ref==1) cout << "encastrement aux 2 bords verticaux\n" ;
165  if (sol_ref==2) {
166  cout << "encastrement aux 4 bords verticaux, solution en sin(x)^2*sin(y)^2\n" ;
167  cout << "eta = " << eta <<"\n";
168  }
169  if (sol_ref==4) {
170  cout << "bord en appuis simple\n" ;
171  cout << "nombre de terme pour calcul sol exacte : " << N_Four << " \n" ;
172  // Calcul des coeeficients de Fourier de la solution exacte :
173  // Cas où le chargement est seulement vertical (pas de moment appliqué)
174  gmm::resize( theta1_Four, N_Four, N_Four) ;
175  gmm::resize( theta2_Four, N_Four, N_Four) ;
176  gmm::resize( u3_Four, N_Four, N_Four) ;
177  base_matrix Jmn(3, 3) ;
178  base_small_vector Bmn(3), Xmn(3) ;
179  scalar_type /*det_Jmn, */ A, B, e2, Pmn ;
180  E = 4.*mu*(mu+lambda) / (2. * mu + lambda);
181  nu = lambda / (2. * mu + lambda);
182  e2 = epsilon * epsilon / 4. ;
183  for(size_type i = 0 ; i < N_Four ; i++) {
184  for(size_type j = 0 ; j < N_Four ; j++) {
185  A = scalar_type(j + 1) * M_PI / LX ;
186  B = scalar_type(i + 1) * M_PI / LY ;
187  Jmn(0, 0) = 2. * A * A / (1. - nu) + B * B + 3. / e2 ;
188  Jmn(0, 1) = A * B * (1. +nu) / (1. - nu) ;
189  Jmn(0, 2) = A * 3. / e2 ;
190  Jmn(1, 0) = A * B * (1. +nu) / (1. - nu) ;
191  Jmn(1, 1) = 2. * B * B / (1. - nu) + A * A + 3. / e2 ;
192  Jmn(1, 2) = B * 3. / e2 ;
193  Jmn(2, 0) = - A ;
194  Jmn(2, 1) = - B ;
195  Jmn(2, 2) = A * A + B * B ;
196  gmm::scale(Jmn, - E*(epsilon/2.) / (1. + nu) ) ;
197 
198  // calcul du développement de Fourrier du chargement :
199  if ( ( (i + 1) % 2 == 1 ) && ( (j + 1) % 2 == 1) ) {
200  Pmn = 16. * pressure / ( scalar_type(i + 1) * scalar_type(j + 1) * M_PI * M_PI) ; }
201  else {
202  Pmn = 0. ; }
203  Bmn[0] = 0. ;
204  Bmn[1] = 0. ;
205  Bmn[2] = Pmn ;
206  gmm::lu_solve(Jmn, Xmn, Bmn) ;
207  theta1_Four(i, j) = Xmn[0] ;
208  theta1_Four(i, j) = Xmn[1] ;
209  u3_Four(i, j) = Xmn[2] ;
210  }
211  }
212  }
213 
214  mf_ut.set_qdim(dim_type(N));
215  mf_theta.set_qdim(dim_type(N));
216 
217  /* set the finite element on the mf_u */
218  getfem::pfem pf_ut = getfem::fem_descriptor(FEM_TYPE_UT);
219  getfem::pfem pf_u3 = getfem::fem_descriptor(FEM_TYPE_U3);
220  getfem::pfem pf_theta = getfem::fem_descriptor(FEM_TYPE_THETA);
221 
222  getfem::pintegration_method ppi =
223  getfem::int_method_descriptor(INTEGRATION);
224  getfem::pintegration_method ppi_ct =
225  getfem::int_method_descriptor(INTEGRATION_CT);
226 
227  mim.set_integration_method(mesh.convex_index(), ppi);
228  mim_subint.set_integration_method(mesh.convex_index(), ppi_ct);
229  mf_ut.set_finite_element(mesh.convex_index(), pf_ut);
230  mf_u3.set_finite_element(mesh.convex_index(), pf_u3);
231  mf_theta.set_finite_element(mesh.convex_index(), pf_theta);
232 
233  /* set the finite element on mf_rhs (same as mf_u is DATA_FEM_TYPE is
234  not used in the .param file */
235  std::string data_fem_name = PARAM.string_value("DATA_FEM_TYPE");
236  if (data_fem_name.size() == 0) {
237  GMM_ASSERT1(pf_ut->is_lagrange(), "You are using a non-lagrange FEM. "
238  << "In that case you need to set "
239  << "DATA_FEM_TYPE in the .param file");
240  mf_rhs.set_finite_element(mesh.convex_index(), pf_ut);
241  } else {
242  mf_rhs.set_finite_element(mesh.convex_index(),
243  getfem::fem_descriptor(data_fem_name));
244  }
245 
246  /* set the finite element on mf_coef. Here we use a very simple element
247  * since the only function that need to be interpolated on the mesh_fem
248  * is f(x)=1 ... */
249  mf_coef.set_finite_element(mesh.convex_index(),
250  getfem::classical_fem(pgt,0));
251 
252  /* set boundary conditions
253  * (Neuman on the upper face, Dirichlet elsewhere) */
254  cout << "Selecting Neumann and Dirichlet boundaries\n";
255  getfem::mesh_region border_faces;
256  getfem::outer_faces_of_mesh(mesh, border_faces);
257  for (getfem::mr_visitor i(border_faces); !i.finished(); ++i) {
258  assert(i.is_face());
259  base_node un = mesh.normal_of_face_of_convex(i.cv(), i.f());
260  un /= gmm::vect_norm2(un);
261  switch(sol_ref){
262  case 0 :
263  if (gmm::abs(un[1]) <= 1.0E-7) // new Neumann face
264  mesh.region(SIMPLY_FIXED_BOUNDARY_NUM).add(i.cv(), i.f());
265  break ;
266  case 1 :
267  if (gmm::abs(un[1]) <= 1.0E-7) // new Neumann face
268  mesh.region(SIMPLY_FIXED_BOUNDARY_NUM).add(i.cv(), i.f());
269  break ;
270  case 2 :
271  if ( (gmm::abs(un[0]) <= 1.0E-7) || (gmm::abs(un[1]) <= 1.0E-7) )
272  mesh.region(SIMPLY_FIXED_BOUNDARY_NUM).add(i.cv(), i.f());
273  break ;
274  case 3 :
275  if (un[0] <= (- 1. + 1.0E-7)) // new Neumann face
276  mesh.region(SIMPLY_FIXED_BOUNDARY_NUM).add(i.cv(), i.f());
277  break ;
278  case 4 :
279  if ( (gmm::abs(un[0]) <= 1.0E-7) || (gmm::abs(un[1]) <= 1.0E-7) )
280  mesh.region(SIMPLY_FIXED_BOUNDARY_NUM).add(i.cv(), i.f());
281  break ;
282  default :
283  GMM_ASSERT1(false, "SOL_REF parameter is undefined");
284  break ;
285  }
286  }
287 }
288 
289 base_small_vector plate_problem::theta_exact(base_node P) {
290  base_small_vector theta(2);
291  if (sol_ref == 0) { // appui simple aux 2 bords
292  theta[0] = - (-pressure / (32. * mu * epsilon * epsilon * epsilon / 8.))
293  * (4. * pow(P[0] - .5, 3.) - 3 * (P[0] - .5));
294  theta[1] = 0.;
295  }
296  if (sol_ref == 1) { // encastrement aux 2 bords
297  theta[0] = - (-pressure / (16. * mu * epsilon * epsilon * epsilon / 8.))
298  * P[0] * ( 2.*P[0]*P[0] - 3.* P[0] + 1. ) ;
299  theta[1] = 0.;
300  }
301  if (sol_ref == 2) { // encastrement aux 4 bords et sols vraiment 2D, non polynomiale
302  theta[0] = (1. + eta) * M_PI * sin(2.*M_PI*P[0]) * sin(M_PI*P[1])* sin(M_PI*P[1]) ;
303  theta[1] = (1. + eta) * M_PI * sin(2.*M_PI*P[1]) * sin(M_PI*P[0])* sin(M_PI*P[0]) ;
304  }
305  if (sol_ref == 3) { // plaque cantilever
306  theta[0] = - (- 3. * pressure / (8. * mu * epsilon * epsilon * epsilon / 8.))
307  * P[0] * ( 0.25 * P[0] * P[0] - P[0] + 1. ) ;
308  theta[1] = 0.;
309  }
310  if (sol_ref == 4) { // bord entier en appui simple
311  theta[0] = 0. ;
312  theta[1] = 0. ;
313  for(size_type i = 0 ; i < N_Four ; i ++) {
314  for(size_type j = 0 ; j < N_Four ; j ++) {
315  theta[0] -= theta1_Four(i, j) * cos( scalar_type(j + 1) * M_PI * P[0] / LX ) * sin( scalar_type(i + 1) * M_PI * P[1] / LY ) ;
316  theta[0] -= theta2_Four(i, j) * sin( scalar_type(j + 1) * M_PI * P[0] / LX ) * cos( scalar_type(i + 1) * M_PI * P[1] / LY ) ;
317  }
318  }
319  }
320  return theta;
321 }
322 
323 scalar_type plate_problem::u3_exact(base_node P) {
324  switch(sol_ref) {
325  case 0 : return (pressure / (32. * mu * epsilon * epsilon * epsilon / 8.))
326  * P[0] * (P[0] - 1.)
327  * (gmm::sqr(P[0] - .5) -1.25-(8.* epsilon*epsilon / 4.));
328  break ;
329  case 1 : return (pressure /(32.* mu * epsilon * epsilon * epsilon / 8.))
330  * P[0] * (P[0] - 1.)
331  * ( P[0] * P[0] - P[0] - 8. * epsilon *epsilon / 4.) ;
332  break ;
333  case 2 : return gmm::sqr(sin(M_PI*P[0])) * gmm::sqr(sin(M_PI*P[1]));
334  break ;
335  case 3 : return (3. * pressure / (4. * mu * epsilon * epsilon * epsilon / 8. ))
336  * P[0] * ( P[0] * P[0] * P[0] / 24. - P[0] * P[0] / 6. + P[0] / 4.
337  - (epsilon * epsilon / 4.) * P[0] / 3.
338  + 2. * (epsilon * epsilon / 4.) / 3.) ;
339  break ;
340  case 4 :
341  scalar_type u3_local ;
342  u3_local = 0. ;
343  for(size_type i = 0 ; i < N_Four ; i ++) {
344  for(size_type j = 0 ; j < N_Four ; j ++)
345  u3_local += u3_Four(i, j) * sin( scalar_type(j + 1) * M_PI * P[0] / LX ) * sin( scalar_type(i + 1) * M_PI * P[1] / LY ) ;
346  }
347  return (u3_local) ;
348  break ;
349  default : GMM_ASSERT1(false, "indice de solution de référence incorrect");
350  }
351 }
352 
353 
354 /* compute the error with respect to the exact solution */
355 void plate_problem::compute_error(plain_vector &Ut, plain_vector &U3, plain_vector &THETA) {
356  cout.precision(16);
357  if (PARAM.int_value("SOL_EXACTE") == 1) {
358  gmm::clear(Ut); gmm::clear(U3); gmm::clear(THETA);
359  }
360 
361  std::vector<scalar_type> V(mf_rhs.nb_dof()*2);
362 
363  getfem::interpolation(mf_ut, mf_rhs, Ut, V);
364  mf_rhs.set_qdim(2);
365  scalar_type l2 = gmm::sqr(getfem::asm_L2_norm(mim, mf_rhs, V));
366  scalar_type h1 = gmm::sqr(getfem::asm_H1_norm(mim, mf_rhs, V));
367  scalar_type linf = gmm::vect_norminf(V);
368  mf_rhs.set_qdim(1);
369  cout << "L2 error = " << sqrt(l2) << endl
370  << "H1 error = " << sqrt(h1) << endl
371  << "Linfty error = " << linf << endl;
372 
373  getfem::interpolation(mf_theta, mf_rhs, THETA, V);
374  GMM_ASSERT1(!mf_rhs.is_reduced(),
375  "To be adapted, use interpolation_function");
376  for (size_type i = 0; i < mf_rhs.nb_dof(); ++i) {
377  gmm::add(gmm::scaled(theta_exact(mf_rhs.point_of_basic_dof(i)), -1.0),
378  gmm::sub_vector(V, gmm::sub_interval(i*2, 2)));
379  }
380  mf_rhs.set_qdim(2);
381  l2 += gmm::sqr(getfem::asm_L2_norm(mim, mf_rhs, V));
382  h1 += gmm::sqr(getfem::asm_H1_semi_norm(mim, mf_rhs, V));
383  linf = std::max(linf, gmm::vect_norminf(V));
384  mf_rhs.set_qdim(1);
385  cout << "L2 error theta:" << sqrt(l2) << endl
386  << "H1 error theta:" << sqrt(h1) << endl
387  << "Linfty error = " << linf << endl;
388 
389  gmm::resize(V, mf_rhs.nb_dof());
390  getfem::interpolation(mf_u3, mf_rhs, U3, V);
391 
392  for (size_type i = 0; i < mf_rhs.nb_dof(); ++i)
393  V[i] -= u3_exact(mf_rhs.point_of_basic_dof(i));
394 
395  l2 = gmm::sqr(getfem::asm_L2_norm(mim, mf_rhs, V));
396  h1 = gmm::sqr(getfem::asm_H1_semi_norm(mim, mf_rhs, V));
397  linf = std::max(linf, gmm::vect_norminf(V));
398 
399  cout.precision(16);
400  cout << "L2 error u3:" << sqrt(l2) << endl
401  << "H1 error u3:" << sqrt(h1) << endl
402  << "Linfty error = " << linf << endl;
403 
404  // stockage de l'erreur H1 :
405  if (PARAM.int_value("SAUV")){
406  std::ofstream f_out("errH1.don");
407  if (!f_out) throw std :: runtime_error("Impossible to open file") ;
408  f_out << sqrt(h1) << "\n" ;
409  }
410 
411 
412 }
413 
414 /**************************************************************************/
415 /* Model. */
416 /**************************************************************************/
417 
418 bool plate_problem::solve(plain_vector &Ut, plain_vector &U3, plain_vector &THETA) {
419  size_type nb_dof_rhs = mf_rhs.nb_dof();
420 
421  cout << "Number of dof for ut: " << mf_ut.nb_dof() << endl;
422  cout << "Number of dof for u3: " << mf_u3.nb_dof() << endl;
423  cout << "Number of dof for theta: " << mf_theta.nb_dof() << endl;
424 
425  E = 4.*mu*(mu+lambda) / (2. * mu + lambda);
426  nu = lambda / (2. * mu + lambda);
427  scalar_type kappa = 5./6.;
428 
429  getfem::model md;
430  md.add_fem_variable("ut", mf_ut);
431  md.add_fem_variable("u3", mf_u3);
432  md.add_fem_variable("theta", mf_theta);
433 
434  // Linearized plate brick.
435  md.add_initialized_scalar_data("E", E);
436  md.add_initialized_scalar_data("nu", nu);
437  md.add_initialized_scalar_data("lambda", lambda);
438  md.add_initialized_scalar_data("mu", mu);
439  md.add_initialized_scalar_data("epsilon", epsilon);
440  md.add_initialized_scalar_data("kappa", kappa);
441  getfem::add_Mindlin_Reissner_plate_brick(md, mim, mim_subint, "u3", "theta",
442  "E", "nu", "epsilon", "kappa",
443  (mitc) ? 2 : 1);
444  getfem::add_isotropic_linearized_elasticity_brick(md, mim, "ut", "lambda", "mu");
445 
446  // Defining the surface source term.
447  if (study_flag == 1 ){
448  cout << "Attention : l'intensité de la pression verticale " ;
449  cout << "a été choisie pour que le déplacement maximal soit unitaire." ;
450  cout << "Pour annuler cette option, faire STUDY_FLAG = 0\n" ;
451  switch(sol_ref) {
452  case 0 :
453  pressure = 128. * mu * epsilon * epsilon * epsilon / 8. ;
454  pressure /= 1.25 + 8. * epsilon * epsilon / 4. ;
455  break ;
456  case 1 :
457  pressure = 128. * mu * epsilon * epsilon * epsilon / 8. ;
458  pressure /= 0.25 + 8. * epsilon * epsilon / 4. ;
459  break ;
460  case 3 :
461  pressure = 32. * mu * epsilon * epsilon * epsilon / 8.;
462  pressure /= 3. + 8. * epsilon * epsilon / 4.;
463  default :
464  break ;
465  }
466  }
467  plain_vector F(nb_dof_rhs);
468  plain_vector M(nb_dof_rhs * 2);
469  if (sol_ref == 2) {
470  base_small_vector P(2) ;
471  scalar_type sx, sy, cx, cy, s2x, s2y, c2x, c2y ;
472  E = 4.*mu*(mu+lambda) / (2. * mu + lambda);
473  nu = lambda / (2. * mu + lambda);
474  for (size_type i = 0; i < nb_dof_rhs; ++i) {
475  P = mf_rhs.point_of_basic_dof(i);
476  sx = sin(M_PI*P[0]) ;
477  cx = cos(M_PI*P[0]) ;
478  sy = sin(M_PI*P[1]) ;
479  cy = cos(M_PI*P[1]) ;
480  c2x = cos(2.*M_PI*P[0]) ;
481  c2y = cos(2.*M_PI*P[1]) ;
482  s2x = sin(2.*M_PI*P[0]) ;
483  s2y = sin(2.*M_PI*P[1]) ;
484  F[i] = 2. * (epsilon / 2.) * E * M_PI * M_PI * eta *
485  ( sy * sy * c2x + sx * sx * c2y ) / ( 1. + nu ) ;
486  M[2*i] = -((epsilon * epsilon * epsilon / 8.) * E * M_PI * s2x / 3. / (1. + nu))
487  * ( (4. * M_PI * M_PI * (1. + eta) * (2. * c2y - 1.) / (1.- nu))
488  - 3. * eta * sy * sy / (epsilon/2.) / (epsilon/2.) ) ;
489  M[2*i+1] = -((epsilon * epsilon * epsilon/8.) * E * M_PI * s2y / 3. / (1. + nu))
490  * ( (4. * M_PI * M_PI * (1. + eta) * (2. * c2x - 1.) / (1.- nu))
491  - 3. * eta * sx * sx / (epsilon/2.) / (epsilon/2.) ) ;
492  }
493  }
494  else // sol_ref = 0 ou 1 ou 3 ou 4: pression verticale uniforme
495  for (size_type i = 0; i < nb_dof_rhs; ++i) {
496  F[i] = pressure;
497  }
498 
499  md.add_initialized_fem_data("VF", mf_rhs, F);
500  getfem::add_source_term_brick(md, mim, "u3", "VF");
501  md.add_initialized_fem_data("VM", mf_rhs, M);
502  getfem::add_source_term_brick(md, mim, "theta", "VM");
503 
505  (md, mim, "u3", mf_u3, SIMPLY_FIXED_BOUNDARY_NUM);
507  (md, mim, "ut", mf_ut, SIMPLY_FIXED_BOUNDARY_NUM);
508 
509  if (sol_ref == 1 || sol_ref == 2 || sol_ref == 3)
511  (md, mim, "theta", mf_ut, SIMPLY_FIXED_BOUNDARY_NUM);
512 
513 
514  // Generic solve.
515  gmm::iteration iter(residual, 1, 40000);
516  getfem::standard_solve(md, iter);
517 
518 
519  gmm::resize(Ut, mf_ut.nb_dof());
520  gmm::copy(md.real_variable("ut"), Ut);
521  gmm::resize(U3, mf_u3.nb_dof());
522  gmm::copy(md.real_variable("u3"), U3);
523  gmm::resize(THETA, mf_theta.nb_dof());
524  gmm::copy(md.real_variable("theta"), THETA);
525 
526  if (PARAM.int_value("VTK_EXPORT")) {
527  cout << "export to " << datafilename + ".vtk" << "..\n";
528  getfem::vtk_export exp(datafilename + ".vtk",
529  PARAM.int_value("VTK_EXPORT")==1);
530  exp.exporting(mf_u3);
531  exp.write_point_data(mf_u3, U3, "plate_normal_displacement");
532  cout << "export done, you can view the data file with (for example)\n"
533  "mayavi2 -d " << datafilename << ".vtk -f "
534  "WarpScalar -m Surface -m Outline\n";
535 // cout << "export done, you can view the data file with (for example)\n"
536 // "mayavi -d " << datafilename << ".vtk -f ExtractVectorNorm -f "
537 // "WarpVector -m BandedSurfaceMap -m Outline\n";
538  }
539  if (PARAM.int_value("DX_EXPORT")) {
540  cout << "export to " << datafilename + ".dx" << ".\n";
541  getfem::dx_export exp(datafilename + ".dx",
542  PARAM.int_value("DX_EXPORT")==1);
543  exp.exporting(mf_u3);
544  exp.write_point_data(mf_u3, U3, "plate_normal_displacement");
545  }
546 
547 
548  return (iter.converged());
549 }
550 
551 /**************************************************************************/
552 /* main program. */
553 /**************************************************************************/
554 
555 int main(int argc, char *argv[]) {
556 
557  GETFEM_MPI_INIT(argc, argv);
558  GMM_SET_EXCEPTION_DEBUG; // Exceptions make a memory fault, to debug.
559  FE_ENABLE_EXCEPT; // Enable floating point exception for Nan.
560 
561  try {
562  plate_problem p;
563  p.PARAM.read_command_line(argc, argv);
564  p.init();
565  if ((p.study_flag != 1)&&((p.sol_ref == 0) || (p.sol_ref ==1)))
566  p.pressure *= p.epsilon * p.epsilon * p.epsilon / 8.;
567  p.mesh.write_to_file(p.datafilename + ".mesh");
568  plain_vector Ut, U3, THETA;
569  bool ok = p.solve(Ut, U3, THETA);
570  p.compute_error(Ut, U3, THETA);
571  GMM_ASSERT1(ok, "Solve has failed");
572 
573  }
574  GMM_STANDARD_CATCH_ERROR;
575 
576  GETFEM_MPI_FINALIZE;
577 
578  return 0;
579 }
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.
void add_initialized_scalar_data(const std::string &name, T e)
Add a scalar data (i.e.
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_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...
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.
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
Miscelleanous assembly routines for common terms. Use the low-level generic assembly....
Export solutions to various formats.
Reissner-Mindlin plate model brick.
Standard solvers for model 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
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norminf(const V &v)
Infinity norm of a vector.
Definition: gmm_blas.h:692
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:58
void resize(V &v, size_type n)
*‍/
Definition: gmm_blas.h:209
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1275
void lu_solve(const DenseMatrix &LU, const Pvector &pvector, VectorX &x, const VectorB &b)
LU Solve : Solve equation Ax=b, given an LU factored matrix.
Definition: gmm_dense_lu.h:129
scalar_type asm_L2_norm(const mesh_im &mim, const mesh_fem &mf, const VEC &U, const mesh_region &rg=mesh_region::all_convexes())
compute , U might be real or complex
scalar_type asm_H1_norm(const mesh_im &mim, const mesh_fem &mf, const VEC &U, const mesh_region &rg=mesh_region::all_convexes())
compute the H1 norm of U.
scalar_type asm_H1_semi_norm(const mesh_im &mim, const mesh_fem &mf, const VEC &U, const mesh_region &rg=mesh_region::all_convexes())
compute , U might be real or complex
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
std::string name_of_geometric_trans(pgeometric_trans p)
Get the string name of a geometric transformation.
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 APIDECL add_isotropic_linearized_elasticity_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataname_lambda, const std::string &dataname_mu, size_type region=size_type(-1), const std::string &dataname_preconstraint=std::string())
Linear elasticity brick ( ).
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.
void interpolation(const mesh_fem &mf_source, const mesh_fem &mf_target, const VECTU &U, VECTV &V, int extrapolation=0, double EPS=1E-10, mesh_region rg_source=mesh_region::all_convexes(), mesh_region rg_target=mesh_region::all_convexes())
interpolation/extrapolation of (mf_source, U) on mf_target.
size_type add_Mindlin_Reissner_plate_brick(model &md, const mesh_im &mim, const mesh_im &mim_reduced, const std::string &u3, const std::string &Theta, const std::string &param_E, const std::string &param_nu, const std::string &param_epsilon, const std::string &param_kappa, size_type variant=size_type(2), size_type region=size_type(-1))
Add a term corresponding to the classical Reissner-Mindlin plate model for which u3 is the transverse...
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_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.