GetFEM  5.5
stokes.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 stokes.cc
23  @brief Solve the stokes problem (incompressible viscuous fluid).
24 
25  This program is used to check that getfem++ is working. This is also
26  a good example of use of GetFEM.
27 
28  This file is almost identical to @link elastostatic.cc
29  tests/elastostatic.cc@endlink, except than a
30  linear incompressibility brick is inserted.
31 */
32 
33 #include "getfem/getfem_assembling.h" /* import assembly methods (and norms comp.) */
34 #include "getfem/getfem_export.h" /* export functions (save solution in a file) */
37 #include "gmm/gmm.h"
38 using std::endl; using std::cout; using std::cerr;
39 using std::ends; using std::cin;
40 
41 /* some GetFEM types that we will be using */
42 using bgeot::base_small_vector; /* special class for small (dim<16) vectors */
43 using bgeot::base_node; /* geometrical nodes(derived from base_small_vector)*/
44 using bgeot::scalar_type; /* = double */
45 using bgeot::size_type; /* = unsigned long */
46 using bgeot::base_matrix; /* small dense matrix. */
47 
48 /* definition of some matrix/vector types. These ones are built
49  * using the predefined types in Gmm++
50  */
52 typedef getfem::modeling_standard_sparse_matrix sparse_matrix;
53 typedef getfem::modeling_standard_plain_vector plain_vector;
54 
55 /*
56  * structure for the stokes problem
57  */
58 struct stokes_problem {
59 
60  enum { DIRICHLET_BOUNDARY_NUM = 0, NEUMANN_BOUNDARY_NUM = 1};
61  getfem::mesh mesh; /* the mesh */
62  getfem::mesh_im mim; /* integration methods. */
63  getfem::mesh_fem mf_u; /* main mesh_fem, for the velocity */
64  getfem::mesh_fem mf_p; /* mesh_fem for the pressure */
65  getfem::mesh_fem mf_rhs; /* mesh_fem for the right hand side (f(x),..) */
66  scalar_type nu; /* Lam� coefficients. */
67 
68  scalar_type residual; /* max residual for the iterative solvers */
69 
70  std::string datafilename;
71  bgeot::md_param PARAM;
72 
73  bool solve(plain_vector &U);
74  void init(void);
75  stokes_problem(void) : mim(mesh), mf_u(mesh), mf_p(mesh),
76  mf_rhs(mesh) {}
77 };
78 
79 /* Read parameters from the .param file, build the mesh, set finite element
80  * and integration methods and selects the boundaries.
81  */
82 void stokes_problem::init(void) {
83  std::string MESH_TYPE = PARAM.string_value("MESH_TYPE","Mesh type ");
84  std::string FEM_TYPE = PARAM.string_value("FEM_TYPE","FEM name");
85  std::string FEM_TYPE_P = PARAM.string_value("FEM_TYPE_P","FEM name P");
86  std::string INTEGRATION = PARAM.string_value("INTEGRATION",
87  "Name of integration method");
88  cout << "MESH_TYPE=" << MESH_TYPE << "\n";
89  cout << "FEM_TYPE=" << FEM_TYPE << "\n";
90  cout << "INTEGRATION=" << INTEGRATION << "\n";
91 
92  /* First step : build the mesh */
95  size_type N = pgt->dim();
96  std::vector<size_type> nsubdiv(N);
97  std::fill(nsubdiv.begin(),nsubdiv.end(),
98  PARAM.int_value("NX", "Nomber of space steps "));
99  getfem::regular_unit_mesh(mesh, nsubdiv, pgt,
100  PARAM.int_value("MESH_NOISED") != 0);
101 
102  /* scale the unit mesh to [LX,LY,..] and incline it */
103  bgeot::base_matrix M(N,N);
104  for (size_type i=0; i < N; ++i) {
105  static const char *t[] = {"LX","LY","LZ"};
106  M(i,i) = (i<3) ? PARAM.real_value(t[i],t[i]) : 1.0;
107  }
108  mesh.transformation(M);
109 
110  datafilename = PARAM.string_value("ROOTFILENAME","Base name of data files.");
111  residual = PARAM.real_value("RESIDUAL"); if (residual == 0.) residual = 1e-10;
112 
113  nu = PARAM.real_value("NU", "Viscosit�); mf_u.set_qdim(bgeot::dim_type(N)); /* set the finite element on the mf_u */ getfem::pfem pf_u = getfem::fem_descriptor(FEM_TYPE); getfem::pintegration_method ppi = getfem::int_method_descriptor(INTEGRATION); mim.set_integration_method(mesh.convex_index(), ppi); mf_u.set_finite_element(mesh.convex_index(), pf_u); mf_p.set_finite_element(mesh.convex_index(), getfem::fem_descriptor(FEM_TYPE_P)); /* set the finite element on mf_rhs (same as mf_u is DATA_FEM_TYPE is not used in the .param file */ std::string data_fem_name = PARAM.string_value("DATA_FEM_TYPE"); if (data_fem_name.size() == 0) { GMM_ASSERT1(pf_u->is_lagrange(), "You are using a non-lagrange FEM " << data_fem_name << ". In that case you need to set " << "DATA_FEM_TYPE in the .param file"); mf_rhs.set_finite_element(mesh.convex_index(), pf_u); } else { mf_rhs.set_finite_element(mesh.convex_index(), getfem::fem_descriptor(data_fem_name)); } /* set boundary conditions * (Neuman on the upper face, Dirichlet elsewhere) */ cout << "Selecting Neumann and Dirichlet boundaries\n"; getfem::mesh_region border_faces; getfem::outer_faces_of_mesh(mesh, border_faces); for (getfem::mr_visitor it(border_faces); !it.finished(); ++it) { assert(it.is_face()); base_node un = mesh.normal_of_face_of_convex(it.cv(), it.f()); un /= gmm::vect_norm2(un); if (gmm::abs(un[N-1] - 1.0) < 1.0E-7) { // new Neumann face mesh.region(NEUMANN_BOUNDARY_NUM).add(it.cv(), it.f()); } else { mesh.region(DIRICHLET_BOUNDARY_NUM).add(it.cv(), it.f()); } } } /**************************************************************************/ /* Model. */ /**************************************************************************/ base_small_vector sol_f(const base_small_vector &P) { base_small_vector res(P.size()); res[P.size()-1] = -1.0; return res; } bool stokes_problem::solve(plain_vector &U) { size_type N = mesh.dim(); getfem::model model; // Main unknown of the problem. model.add_fem_variable("u", mf_u); // Linearized elasticity brick. model.add_initialized_fixed_size_data ("lambda", plain_vector(1, 0.0)); model.add_initialized_fixed_size_data("nu", plain_vector(1, nu)); getfem::add_isotropic_linearized_elasticity_brick (model, mim, "u", "lambda", "nu"); // Linearized incompressibility condition brick. model.add_fem_variable("p", mf_p); // Adding the pressure as a variable add_linear_incompressibility(model, mim, "u", "p"); // Volumic source term. std::vector<scalar_type> F(mf_rhs.nb_dof()*N); getfem::interpolation_function(mf_rhs, F, sol_f); model.add_initialized_fem_data("VolumicData", mf_rhs, F); getfem::add_source_term_brick(model, mim, "u", "VolumicData"); // Dirichlet condition. gmm::clear(F); model.add_initialized_fem_data("DirichletData", mf_rhs, F); getfem::add_Dirichlet_condition_with_multipliers (model, mim, "u", mf_u, DIRICHLET_BOUNDARY_NUM, "DirichletData"); gmm::iteration iter(residual, 1, 40000); getfem::standard_solve(model, iter); // Solution extraction gmm::copy(model.real_variable("u"), U); return (iter.converged()); } /**************************************************************************/ /* main program. */ /**************************************************************************/ int main(int argc, char *argv[]) { GETFEM_MPI_INIT(argc, argv); GMM_SET_EXCEPTION_DEBUG; // Exceptions make a memory fault, to debug. FE_ENABLE_EXCEPT; // Enable floating point exception for Nan. try { stokes_problem p; p.PARAM.read_command_line(argc, argv); p.init(); // p.mesh.write_to_file(p.datafilename + ".mesh"); plain_vector U(p.mf_u.nb_dof()); if (!p.solve(U)) GMM_ASSERT1(false, "Solve has failed"); if (p.PARAM.int_value("VTK_EXPORT")) { cout << "export to " << p.datafilename + ".vtk" << "..\n"; getfem::vtk_export exp(p.datafilename + ".vtk", p.PARAM.int_value("VTK_EXPORT")==1); exp.exporting(p.mf_u); exp.write_point_data(p.mf_u, U, "stokes_displacement"); cout << "export done, you can view the data file with (for example)\n" "mayavi2 -d " << p.datafilename << ".vtk -f ExtractVectorNorm -f " "WarpVector -m Surface -m Outline\n"; } } GMM_STANDARD_CATCH_ERROR; GETFEM_MPI_FINALIZE; return 0; } ");
114  mf_u.set_qdim(bgeot::dim_type(N));
115 
116  /* set the finite element on the mf_u */
117  getfem::pfem pf_u =
118  getfem::fem_descriptor(FEM_TYPE);
119  getfem::pintegration_method ppi =
120  getfem::int_method_descriptor(INTEGRATION);
121 
122  mim.set_integration_method(mesh.convex_index(), ppi);
123  mf_u.set_finite_element(mesh.convex_index(), pf_u);
124  mf_p.set_finite_element(mesh.convex_index(),
125  getfem::fem_descriptor(FEM_TYPE_P));
126 
127  /* set the finite element on mf_rhs (same as mf_u is DATA_FEM_TYPE is
128  not used in the .param file */
129  std::string data_fem_name = PARAM.string_value("DATA_FEM_TYPE");
130  if (data_fem_name.size() == 0) {
131  GMM_ASSERT1(pf_u->is_lagrange(), "You are using a non-lagrange FEM "
132  << data_fem_name << ". In that case you need to set "
133  << "DATA_FEM_TYPE in the .param file");
134  mf_rhs.set_finite_element(mesh.convex_index(), pf_u);
135  } else {
136  mf_rhs.set_finite_element(mesh.convex_index(),
137  getfem::fem_descriptor(data_fem_name));
138  }
139 
140  /* set boundary conditions
141  * (Neuman on the upper face, Dirichlet elsewhere) */
142  cout << "Selecting Neumann and Dirichlet boundaries\n";
143  getfem::mesh_region border_faces;
144  getfem::outer_faces_of_mesh(mesh, border_faces);
145  for (getfem::mr_visitor it(border_faces); !it.finished(); ++it) {
146  assert(it.is_face());
147  base_node un = mesh.normal_of_face_of_convex(it.cv(), it.f());
148  un /= gmm::vect_norm2(un);
149  if (gmm::abs(un[N-1] - 1.0) < 1.0E-7) { // new Neumann face
150  mesh.region(NEUMANN_BOUNDARY_NUM).add(it.cv(), it.f());
151  } else {
152  mesh.region(DIRICHLET_BOUNDARY_NUM).add(it.cv(), it.f());
153  }
154  }
155 }
156 
157 /**************************************************************************/
158 /* Model. */
159 /**************************************************************************/
160 
161 base_small_vector sol_f(const base_small_vector &P) {
162  base_small_vector res(P.size());
163  res[P.size()-1] = -1.0;
164  return res;
165 }
166 
167 
168 bool stokes_problem::solve(plain_vector &U) {
169  size_type N = mesh.dim();
170 
171  getfem::model model;
172 
173  // Main unknown of the problem.
174  model.add_fem_variable("u", mf_u);
175 
176  // Linearized elasticity brick.
178  ("lambda", plain_vector(1, 0.0));
179  model.add_initialized_fixed_size_data("nu", plain_vector(1, nu));
181  (model, mim, "u", "lambda", "nu");
182 
183  // Linearized incompressibility condition brick.
184  model.add_fem_variable("p", mf_p); // Adding the pressure as a variable
185  add_linear_incompressibility(model, mim, "u", "p");
186 
187  // Volumic source term.
188  std::vector<scalar_type> F(mf_rhs.nb_dof()*N);
189  getfem::interpolation_function(mf_rhs, F, sol_f);
190  model.add_initialized_fem_data("VolumicData", mf_rhs, F);
191  getfem::add_source_term_brick(model, mim, "u", "VolumicData");
192 
193  // Dirichlet condition.
194  gmm::clear(F);
195  model.add_initialized_fem_data("DirichletData", mf_rhs, F);
197  (model, mim, "u", mf_u, DIRICHLET_BOUNDARY_NUM, "DirichletData");
198 
199  gmm::iteration iter(residual, 1, 40000);
200  getfem::standard_solve(model, iter);
201 
202  // Solution extraction
203  gmm::copy(model.real_variable("u"), U);
204 
205  return (iter.converged());
206 }
207 
208 /**************************************************************************/
209 /* main program. */
210 /**************************************************************************/
211 
212 int main(int argc, char *argv[]) {
213 
214  GETFEM_MPI_INIT(argc, argv);
215  GMM_SET_EXCEPTION_DEBUG; // Exceptions make a memory fault, to debug.
216  FE_ENABLE_EXCEPT; // Enable floating point exception for Nan.
217 
218  try {
219  stokes_problem p;
220  p.PARAM.read_command_line(argc, argv);
221  p.init();
222  // p.mesh.write_to_file(p.datafilename + ".mesh");
223  plain_vector U(p.mf_u.nb_dof());
224  if (!p.solve(U)) GMM_ASSERT1(false, "Solve has failed");
225  if (p.PARAM.int_value("VTK_EXPORT")) {
226  cout << "export to " << p.datafilename + ".vtk" << "..\n";
227  getfem::vtk_export exp(p.datafilename + ".vtk",
228  p.PARAM.int_value("VTK_EXPORT")==1);
229  exp.exporting(p.mf_u);
230  exp.write_point_data(p.mf_u, U, "stokes_displacement");
231  cout << "export done, you can view the data file with (for example)\n"
232  "mayavi2 -d " << p.datafilename << ".vtk -f ExtractVectorNorm -f "
233  "WarpVector -m Surface -m Outline\n";
234  }
235  }
236  GMM_STANDARD_CATCH_ERROR;
237  GETFEM_MPI_FINALIZE;
238 
239  return 0;
240 }
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_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...
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.
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 clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:58
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
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.
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 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.
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.