GetFEM  5.5
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 elastostatic.cc
23  @brief Linear Elastostatic problem. A dummy linear
24  elastotatic problem is solved on a regular mesh, and is compared to
25  the analytical solution.
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  @see laplacian.cc
31  @see nonlinear_elastostatic.cc
32 */
33 
34 #include "getfem/getfem_config.h"
35 #include "getfem/getfem_assembling.h" /* import assembly methods (and norms comp.) */
36 #include "getfem/getfem_export.h" /* export functions (save solution in a file) */
39 #include "gmm/gmm.h"
42 #include "getfem/getfem_import.h"
43 
44 using std::endl; using std::cout; using std::cerr;
45 using std::ends; using std::cin;
46 
47 
48 /* some GetFEM types that we will be using */
49 using bgeot::base_small_vector; /* special class for small (dim<16) vectors */
50 using bgeot::base_node; /* geometrical nodes(derived from base_small_vector)*/
51 using bgeot::scalar_type; /* = double */
52 using bgeot::size_type; /* = unsigned long */
53 using bgeot::dim_type;
54 using bgeot::base_matrix; /* small dense matrix. */
55 
56 /* definition of some matrix/vector types.
57  * default types of getfem_model_solvers.h
58  */
60 typedef getfem::modeling_standard_sparse_matrix sparse_matrix;
61 typedef getfem::modeling_standard_plain_vector plain_vector;
62 
63 /**************************************************************************/
64 /* Exact solution. */
65 /**************************************************************************/
66 
67 gmm::row_matrix<base_small_vector> sol_K;
68 static scalar_type sol_lambda, sol_mu, alph = 0.3;
69 int sol_sing;
70 
71 base_small_vector sol_u(const base_node &x) {
72  int N = x.size(); base_small_vector res(N);
73  switch(sol_sing) {
74  case 0 :
75  for (int i = 0; i < N; ++i)
76  res[i] = alph * sin(gmm::vect_sp(sol_K.row(i), x));
77  break;
78  case 1 :
79  {
80  base_small_vector trans(x.size());
81  gmm::fill(trans, M_PI / 10.0);
82  base_node y = x - trans;
83  scalar_type a = gmm::vect_norm2(y);
84  for (int i = 0; i < N; ++i) res[i] = a;
85  break;
86  }
87  case 2 :
88  {
89  base_small_vector trans(x.size());
90  gmm::fill(trans, M_PI / 10.0);
91  base_node y = x - trans;
92  scalar_type a = gmm::sqrt(gmm::vect_norm2(y));
93  for (int i = 0; i < N; ++i) res[i] = a;
94  break;
95  }
96  }
97  return res;
98 }
99 
100 base_small_vector sol_f(const base_node &x) {
101  int N = x.size();
102  base_small_vector res(N);
103  switch (sol_sing) {
104  case 0 :
105  for (int i = 0; i < N; i++) {
106  res[i] = alph * ( sol_mu * gmm::vect_sp(sol_K.row(i), sol_K.row(i)) )
107  * sin(gmm::vect_sp(sol_K.row(i), x));
108  for (int j = 0; j < N; j++)
109  res[i] += alph * ( (sol_lambda + sol_mu) * sol_K(j,j) * sol_K(j,i))
110  * sin(gmm::vect_sp(sol_K.row(j), x));
111  }
112  break;
113  case 1 :
114  {
115  base_small_vector trans(x.size());
116  gmm::fill(trans, M_PI / 10.0);
117  base_node y = x - trans;
118  scalar_type r = gmm::vect_norm2(y) + 1e-100;
119  scalar_type r2 = r*r;
120  scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
121 
122  for (int i = 0; i < N; i++)
123  res[i] = sol_lambda * (y[i]*tr / r2 - 1.0) / r
124  + sol_mu * (y[i]*tr/r2 - N) / r;
125  }
126  break;
127  case 2 :
128  {
129  base_small_vector trans(x.size());
130  gmm::fill(trans, M_PI / 10.0);
131  base_node y = x - trans;
132 
133  scalar_type r = gmm::vect_norm2(y) + 1e-100;
134  scalar_type a = gmm::sqrt(1.0/r);
135  scalar_type b = a*a*a, c = b*b*a;
136  scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
137  for (int i = 0; i < N; ++i)
138  res[i] -= b * (sol_lambda + sol_mu * (N+1-3.0/2.0)) * 0.5
139  - c * 3.0 * (sol_lambda + sol_mu) * (y[i]*tr) / 4.0;
140  }
141  break;
142  }
143  return res;
144 }
145 
146 base_matrix sol_sigma(const base_node &x) {
147  int N = x.size();
148  base_matrix res(N,N);
149  switch (sol_sing) {
150  case 0 :
151  for (int i = 0; i < N; i++)
152  for (int j = 0; j <= i; j++) {
153  res(j,i) = res(i,j) = alph * sol_mu *
154  ( sol_K(i,j) * cos(gmm::vect_sp(sol_K.row(i), x))
155  + sol_K(j,i) * cos(gmm::vect_sp(sol_K.row(j), x))
156  );
157  if (i == j)
158  for (int k = 0; k < N; k++)
159  res(i,j) += alph * sol_lambda * sol_K(k,k)
160  * cos(gmm::vect_sp(sol_K.row(k), x));
161  }
162  break;
163  case 1 :
164  {
165  base_small_vector trans(x.size());
166  gmm::fill(trans, M_PI / 10.0);
167  base_node y = x - trans;
168  scalar_type r = gmm::vect_norm2(y) + 1e-100;
169  scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
170  for (int i = 0; i < N; i++) {
171  res(i, i) += sol_lambda * tr / r;
172  for (int j = 0; j < N; j++)
173  res(i, j) += sol_mu * (y[i] + y[j]) / r;
174  }
175  }
176 
177  break;
178  case 2 :
179  {
180  base_small_vector trans(x.size());
181  gmm::fill(trans, M_PI / 10.0);
182  base_node y = x - trans;
183  scalar_type r = gmm::vect_norm2(y) + 1e-100;
184  scalar_type a = gmm::sqrt(1.0/r);
185  scalar_type b = a*a*a;
186  scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
187  for (int i = 0; i < N; i++) {
188  res(i, i) += sol_lambda * tr * b * 0.5;
189  for (int j = 0; j < N; j++)
190  res(i, j) += sol_mu * b * (y[i] + y[j]) * 0.5;
191  }
192  }
193  }
194 
195  return res;
196 }
197 
198 /*
199  structure for the elastostatic problem
200 */
201 struct elastostatic_problem {
202 
203  enum { DIRICHLET_BOUNDARY_NUM = 0, NEUMANN_BOUNDARY_NUM = 1};
204  getfem::mesh mesh; /* the mesh */
205  getfem::mesh_im mim; /* the integration methods. */
206  getfem::mesh_fem mf_u; /* main mesh_fem, for the elastostatic solution */
207  getfem::mesh_fem mf_mult; /* mesh_fem for the Dirichlet condition. */
208  getfem::mesh_fem mf_rhs; /* mesh_fem for the right hand side (f(x),..) */
209  getfem::mesh_fem mf_p; /* mesh_fem for the pressure for mixed form */
210  scalar_type lambda, mu; /* Lamé coefficients. */
211 
212  scalar_type residual; /* max residual for iterative solvers */
213  bool mixed_pressure, refine;
214  size_type dirichlet_version;
215 
216  std::string datafilename;
217  bgeot::md_param PARAM;
218 
219  void select_boundaries(void);
220  bool solve(plain_vector &U);
221  void init(void);
222  void compute_error(plain_vector &U);
223  elastostatic_problem(void) : mim(mesh),mf_u(mesh), mf_mult(mesh),
224  mf_rhs(mesh),mf_p(mesh) {}
225 };
226 
227 /* Selects the boundaries */
228 
229 void elastostatic_problem::select_boundaries(void) {
230  size_type N = mesh.dim();
231  getfem::mesh_region border_faces;
232  getfem::outer_faces_of_mesh(mesh, border_faces);
233  for (getfem::mr_visitor i(border_faces); !i.finished(); ++i) {
234  base_node un = mesh.normal_of_face_of_convex(i.cv(), i.f());
235  un /= gmm::vect_norm2(un);
236  if (gmm::abs(un[N-1] - 1.0) < 0.5) { // new Neumann face
237  mesh.region(NEUMANN_BOUNDARY_NUM).add(i.cv(), i.f());
238  } else {
239  mesh.region(DIRICHLET_BOUNDARY_NUM).add(i.cv(), i.f());
240  }
241  }
242 }
243 
244 
245 /* Read parameters from the .param file, build the mesh, set finite element
246  * and integration methods and selects the boundaries.
247  */
248 void elastostatic_problem::init(void) {
249  std::string MESH_FILE = PARAM.string_value("MESH_FILE", "Mesh file");
250  std::string FEM_TYPE = PARAM.string_value("FEM_TYPE","FEM name");
251  std::string INTEGRATION = PARAM.string_value("INTEGRATION",
252  "Name of integration method");
253 
254  cout << "MESH_FILE=" << MESH_FILE << "\n";
255  cout << "FEM_TYPE=" << FEM_TYPE << "\n";
256  cout << "INTEGRATION=" << INTEGRATION << "\n";
257 
258 #if GETFEM_PARA_LEVEL > 1
259  double t_init=MPI_Wtime();
260 #endif
261 
262  size_type NX = PARAM.int_value("NX");
263  size_type N = PARAM.int_value("N");
264  std::stringstream filename; filename << MESH_FILE;
265  if ((MESH_FILE.compare(0,11,"structured:") == 0) && NX > 0) {
266  filename << ";NSUBDIV=[" << NX;
267  for (size_type i = 1; i < N; ++i) filename << "," << NX;
268  filename << "];";
269  }
270  getfem::import_mesh(filename.str(), mesh);
271 
272  GMM_ASSERT1(N == mesh.dim(), "The mesh has not the right dimension");
273 
274 #if GETFEM_PARA_LEVEL > 1
275  cout<<"temps creation maillage "<< MPI_Wtime()-t_init<<endl;
276 #endif
277 
278  dirichlet_version
279  = size_type(PARAM.int_value("DIRICHLET_VERSION",
280  "Dirichlet version"));
281  datafilename = PARAM.string_value("ROOTFILENAME","Base name of data files.");
282  scalar_type FT = PARAM.real_value("FT", "parameter for exact solution");
283  residual = PARAM.real_value("RESIDUAL");
284  if (residual == 0.) residual = 1e-10;
285  gmm::resize(sol_K, N, N);
286  for (size_type i = 0; i < N; i++)
287  for (size_type j = 0; j < N; j++)
288  sol_K(i,j) = (i == j) ? FT : -FT;
289 
290  mu = PARAM.real_value("MU", "Lamé coefficient mu");
291  lambda = PARAM.real_value("LAMBDA", "Lamé coefficient lambda");
292  sol_sing = int(PARAM.int_value("SOL_SING", "Optional singular solution"));
293  refine = (PARAM.int_value("REFINE", "Optional refinement") != 0);
294  sol_lambda = lambda; sol_mu = mu;
295  mf_u.set_qdim(dim_type(N));
296  mf_mult.set_qdim(dim_type(N));
297 
298  /* set the finite element on the mf_u */
299  getfem::pfem pf_u =
300  getfem::fem_descriptor(FEM_TYPE);
301  getfem::pintegration_method ppi =
302  getfem::int_method_descriptor(INTEGRATION);
303 
304  mim.set_integration_method(ppi);
305  mf_u.set_finite_element(pf_u);
306 
307  std::string dirichlet_fem_name = PARAM.string_value("DIRICHLET_FEM_TYPE");
308  if (dirichlet_fem_name.size() == 0)
309  mf_mult.set_finite_element(pf_u);
310  else {
311  cout << "DIRICHLET_FEM_TYPE=" << dirichlet_fem_name << "\n";
312  mf_mult.set_finite_element(getfem::fem_descriptor(dirichlet_fem_name));
313  }
314 
315  mixed_pressure =
316  (PARAM.int_value("MIXED_PRESSURE","Mixed version or not.") != 0);
317  if (mixed_pressure) {
318  std::string FEM_TYPE_P = PARAM.string_value("FEM_TYPE_P","FEM name P");
319  mf_p.set_finite_element(getfem::fem_descriptor(FEM_TYPE_P));
320  }
321 
322  /* set the finite element on mf_rhs (same as mf_u is DATA_FEM_TYPE is
323  not used in the .param file */
324  std::string data_fem_name = PARAM.string_value("DATA_FEM_TYPE");
325  if (data_fem_name.size() == 0) {
326  GMM_ASSERT1(pf_u->is_lagrange(), "You are using a non-lagrange FEM. "
327  << "In that case you need to set "
328  << "DATA_FEM_TYPE in the .param file");
329  mf_rhs.set_finite_element(pf_u);
330  } else {
331  mf_rhs.set_finite_element(getfem::fem_descriptor(data_fem_name));
332  }
333 
334  /* set boundary conditions
335  * (Neuman on the upper face, Dirichlet elsewhere) */
336  cout << "Selecting Neumann and Dirichlet boundaries\n";
337  select_boundaries();
338 
339 #if GETFEM_PARA_LEVEL > 1
340 
341 
342  t_init = MPI_Wtime();
343 
344  mf_u.nb_dof(); mf_rhs.nb_dof(); mf_mult.nb_dof();
345 
346  cout<<"enumerate dof time "<< MPI_Wtime()-t_init<<endl;
347 #else
348  double t_init = gmm::uclock_sec();
349  mf_u.nb_dof(); mf_rhs.nb_dof(); mf_mult.nb_dof();
350  cout << "enumerate dof time " << gmm::uclock_sec() - t_init << endl;
351 #endif
352 }
353 
354 /* compute the error with respect to the exact solution */
355 void elastostatic_problem::compute_error(plain_vector &U) {
356  size_type N = mesh.dim();
357  std::vector<scalar_type> V(mf_rhs.nb_basic_dof()*N);
358  getfem::interpolation(mf_u, mf_rhs, U, V);
359  for (size_type i = 0; i < mf_rhs.nb_basic_dof(); ++i) {
360  gmm::add(gmm::scaled(sol_u(mf_rhs.point_of_basic_dof(i)), -1.0),
361  gmm::sub_vector(V, gmm::sub_interval(i*N, N)));
362  }
363 
364 
365 
366  cout.precision(16);
367  mf_rhs.set_qdim(dim_type(N));
368  scalar_type l2 = getfem::asm_L2_norm(mim, mf_rhs, V);
369  scalar_type h1 = getfem::asm_H1_norm(mim, mf_rhs, V);
370 
371  if (getfem::MPI_IS_MASTER())
372  cout << "L2 error = " << l2 << endl
373  << "H1 error = " << h1 << endl
374  << "Linfty error = " << gmm::vect_norminf(V) << endl;
375 
376 /* getfem::vtk_export exp(datafilename + "_err.vtk", */
377 /* PARAM.int_value("VTK_EXPORT")==1); */
378 /* exp.exporting(mf_rhs); */
379 /* exp.write_point_data(mf_rhs, V, "elastostatic_displacement"); */
380 
381  mf_rhs.set_qdim(1);
382 }
383 
384 /**************************************************************************/
385 /* Model. */
386 /**************************************************************************/
387 
388 
389 bool elastostatic_problem::solve(plain_vector &U) {
390 
391  size_type N = mesh.dim();
392 
393  if (mixed_pressure) cout << "Number of dof for P: " << mf_p.nb_dof() << endl;
394  cout << "Number of dof for u: " << mf_u.nb_dof() << endl;
395 
396  getfem::model model;
397 
398  // Main unknown of the problem.
399  model.add_fem_variable("u", mf_u);
400 
401  // Linearized elasticity brick.
402  model.add_initialized_scalar_data("lambda", mixed_pressure ? 0.0 : lambda);
403  model.add_initialized_scalar_data("mu", mu);
405  (model, mim, "u", "lambda", "mu");
406 
407  // Linearized incompressibility condition brick.
408  if (mixed_pressure) {
409  model.add_initialized_scalar_data("incomp_coeff", 1.0/lambda);
410  model.add_fem_variable("p", mf_p); // Adding the pressure as a variable
412  (model, mim, "u", "p", size_type(-1), "incomp_coeff");
413  }
414 
415  // Volumic source term.
416  std::vector<scalar_type> F(mf_rhs.nb_dof()*N);
417  getfem::interpolation_function(mf_rhs, F, sol_f);
418  model.add_initialized_fem_data("VolumicData", mf_rhs, F);
419  getfem::add_source_term_brick(model, mim, "u", "VolumicData");
420 
421  // Neumann condition.
422  gmm::resize(F, mf_rhs.nb_dof()*N*N);
423  getfem::interpolation_function(mf_rhs, F, sol_sigma, NEUMANN_BOUNDARY_NUM);
424  model.add_initialized_fem_data("NeumannData", mf_rhs, F);
426  (model, mim, "u", "NeumannData", NEUMANN_BOUNDARY_NUM);
427 
428  // Dirichlet condition.
429  gmm::resize(F, mf_rhs.nb_dof()*N);
430  getfem::interpolation_function(mf_rhs, F, sol_u);
431  model.add_initialized_fem_data("DirichletData", mf_rhs, F);
433  (model, mim, "u", mf_u, DIRICHLET_BOUNDARY_NUM, "DirichletData");
434 
435  gmm::iteration iter(residual, 1, 40000);
436 #if GETFEM_PARA_LEVEL > 1
437  double t_init=MPI_Wtime();
438 #endif
439  dal::bit_vector cvref;
440 
441  do { // solve with optional refinement
442 
443  cout << "Total number of variables : " << model.nb_dof() << endl;
444 
445  // Defining the volumic source term.
446  size_type nb_dof_rhs = mf_rhs.nb_dof();
447  gmm::resize(F, nb_dof_rhs * N);
448  getfem::interpolation_function(mf_rhs, F, sol_f);
449  gmm::copy(F, model.set_real_variable("VolumicData"));
450 
451  // Defining the Neumann source term.
452  gmm::resize(F, nb_dof_rhs * N * N);
453  getfem::interpolation_function(mf_rhs, F, sol_sigma, NEUMANN_BOUNDARY_NUM);
454  gmm::copy(F, model.set_real_variable("NeumannData"));
455 
456  // Defining the Dirichlet condition value.
457  gmm::resize(F, nb_dof_rhs * N);
458  getfem::interpolation_function(mf_rhs, F, sol_u, DIRICHLET_BOUNDARY_NUM);
459  gmm::copy(F, model.set_real_variable("DirichletData"));
460 
461  iter.init();
462  getfem::standard_solve(model, iter);
463  gmm::resize(U, mf_u.nb_dof());
464  gmm::copy(model.real_variable("u"), U);
465 
466  if (refine) {
467  plain_vector ERR(mesh.convex_index().last_true()+1);
468  getfem::error_estimate(mim, mf_u, U, ERR);
469 
470  cout << "max = " << gmm::vect_norminf(ERR) << endl;
471  // scalar_type threshold = gmm::vect_norminf(ERR) * 0.7;
472  scalar_type threshold = 0.0001, min_ = 1e18;
473  cvref.clear();
474  for (dal::bv_visitor i(mesh.convex_index()); !i.finished(); ++i) {
475  if (ERR[i] > threshold) cvref.add(i);
476  min_ = std::min(min_, ERR[i]);
477  }
478  cout << "min = " << min_ << endl;
479  cout << "Nb elt to be refined : " << cvref.card() << endl;
480 
481  mesh.Bank_refine(cvref);
482  }
483 
484  } while (refine && cvref.card() > 0);
485 
486 #if GETFEM_PARA_LEVEL > 1
487  cout<<"temps standard solve "<< MPI_Wtime()-t_init<<endl;
488 #endif
489 
490 
491  return (iter.converged());
492 }
493 
494 /**************************************************************************/
495 /* main program. */
496 /**************************************************************************/
497 
498 int main(int argc, char *argv[]) {
499 
500  GETFEM_MPI_INIT(argc, argv); // For parallelized version
501 
502  GMM_SET_EXCEPTION_DEBUG; // Exceptions make a memory fault, to debug.
503  FE_ENABLE_EXCEPT; // Enable floating point exception for Nan.
504 
505  // try {
506 
507 
508 
509  elastostatic_problem p;
510  p.PARAM.read_command_line(argc, argv);
511 #if GETFEM_PARA_LEVEL > 1
512  double t_ref=MPI_Wtime();
513 #endif
514  p.init();
515 #if GETFEM_PARA_LEVEL > 1
516  cout << "temps init "<< MPI_Wtime()-t_ref << endl;
517 #endif
518  if (getfem::MPI_IS_MASTER())
519  p.mesh.write_to_file(p.datafilename + ".mesh");
520 
521  plain_vector U;
522 
523 #if GETFEM_PARA_LEVEL > 1
524  t_ref=MPI_Wtime();
525  cout<<"begining resol"<<endl;
526 #endif
527  if (!p.solve(U)) GMM_ASSERT1(false, "Solve has failed");
528 
529 #if GETFEM_PARA_LEVEL > 1
530  cout << "temps Resol "<< MPI_Wtime()-t_ref << endl;
531  t_ref = MPI_Wtime();
532 #endif
533  p.compute_error(U);
534 #if GETFEM_PARA_LEVEL > 1
535  cout << "temps error "<< MPI_Wtime()-t_ref << endl;
536  t_ref = MPI_Wtime();
537 #endif
538 
539  // if (getfem::MPI_IS_MASTER()) { p.mesh.write_to_file("toto.mesh"); }
540 
541  if (p.PARAM.int_value("VTK_EXPORT") && getfem::MPI_IS_MASTER()) {
542  cout << "export to " << p.datafilename + ".vtk" << "..\n";
543  getfem::vtk_export exp(p.datafilename + ".vtk",
544  p.PARAM.int_value("VTK_EXPORT")==1);
545  exp.exporting(p.mf_u);
546  exp.write_point_data(p.mf_u, U, "elastostatic_displacement");
547  cout << "export done, you can view the data file with (for example)\n"
548  "mayavi2 -d " << p.datafilename << ".vtk -f ExtractVectorNorm -f "
549  "WarpVector -m Surface -m Outline\n";
550  }
551 
552  // } GMM_STANDARD_CATCH_ERROR;
553 
554  GETFEM_MPI_FINALIZE;
555 
556  return 0;
557 }
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_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...
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.
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....
defines and typedefs for namespace getfem
Definition of a posteriori error estimates.
Export solutions to various formats.
Import mesh files from various formats.
Interpolation of fields from a mesh_fem onto another.
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
void fill(L &l, typename gmm::linalg_traits< L >::value_type x)
*â€/
Definition: gmm_blas.h:103
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 resize(V &v, size_type n)
*â€/
Definition: gmm_blas.h:209
strongest_value_type< V1, V2 >::value_type vect_sp(const V1 &v1, const V2 &v2)
*â€/
Definition: gmm_blas.h:263
void add(const L1 &l1, L2 &l2)
*â€/
Definition: gmm_blas.h:1275
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.
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
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
size_type APIDECL add_normal_source_term_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region)
Add a source term on the variable varname on a boundary region.
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.
void interpolation_function(mesh_fem &mf_target, const VECT &VV, F &f, mesh_region rg=mesh_region::all_convexes())
interpolation of a function f on mf_target.
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 import_mesh(const std::string &filename, const std::string &format, mesh &m)
imports a mesh file.
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 APIDECL add_linear_incompressibility(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname_pressure, size_type region=size_type(-1), const std::string &dataexpr_penal_coeff=std::string())
Mixed linear incompressibility condition brick.
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.