1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241
|
/*===========================================================================
Copyright (C) 2002-2020 Yves Renard, Julien Pommier.
This file is a part of GetFEM
GetFEM is free software; you can redistribute it and/or modify it
under the terms of the GNU Lesser General Public License as published
by the Free Software Foundation; either version 3 of the License, or
(at your option) any later version along with the GCC Runtime Library
Exception either version 3.1 or (at your option) any later version.
This program is distributed in the hope that it will be useful, but
WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
License and GCC Runtime Library Exception for more details.
You should have received a copy of the GNU Lesser General Public License
along with this program; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
===========================================================================*/
/**
@file stokes.cc
@brief Solve the stokes problem (incompressible viscuous fluid).
This program is used to check that getfem++ is working. This is also
a good example of use of GetFEM.
This file is almost identical to @link elastostatic.cc
tests/elastostatic.cc@endlink, except than a
linear incompressibility brick is inserted.
*/
#include "getfem/getfem_assembling.h" /* import assembly methods (and norms comp.) */
#include "getfem/getfem_export.h" /* export functions (save solution in a file) */
#include "getfem/getfem_regular_meshes.h"
#include "getfem/getfem_model_solvers.h"
#include "gmm/gmm.h"
using std::endl; using std::cout; using std::cerr;
using std::ends; using std::cin;
/* some GetFEM types that we will be using */
using bgeot::base_small_vector; /* special class for small (dim<16) vectors */
using bgeot::base_node; /* geometrical nodes(derived from base_small_vector)*/
using bgeot::scalar_type; /* = double */
using bgeot::size_type; /* = unsigned long */
using bgeot::base_matrix; /* small dense matrix. */
/* definition of some matrix/vector types. These ones are built
* using the predefined types in Gmm++
*/
typedef getfem::modeling_standard_sparse_vector sparse_vector;
typedef getfem::modeling_standard_sparse_matrix sparse_matrix;
typedef getfem::modeling_standard_plain_vector plain_vector;
/*
* structure for the stokes problem
*/
struct stokes_problem {
enum { DIRICHLET_BOUNDARY_NUM = 0, NEUMANN_BOUNDARY_NUM = 1};
getfem::mesh mesh; /* the mesh */
getfem::mesh_im mim; /* integration methods. */
getfem::mesh_fem mf_u; /* main mesh_fem, for the velocity */
getfem::mesh_fem mf_p; /* mesh_fem for the pressure */
getfem::mesh_fem mf_rhs; /* mesh_fem for the right hand side (f(x),..) */
scalar_type nu; /* Lam coefficients. */
scalar_type residual; /* max residual for the iterative solvers */
std::string datafilename;
bgeot::md_param PARAM;
bool solve(plain_vector &U);
void init(void);
stokes_problem(void) : mim(mesh), mf_u(mesh), mf_p(mesh),
mf_rhs(mesh) {}
};
/* Read parameters from the .param file, build the mesh, set finite element
* and integration methods and selects the boundaries.
*/
void stokes_problem::init(void) {
std::string MESH_TYPE = PARAM.string_value("MESH_TYPE","Mesh type ");
std::string FEM_TYPE = PARAM.string_value("FEM_TYPE","FEM name");
std::string FEM_TYPE_P = PARAM.string_value("FEM_TYPE_P","FEM name P");
std::string INTEGRATION = PARAM.string_value("INTEGRATION",
"Name of integration method");
cout << "MESH_TYPE=" << MESH_TYPE << "\n";
cout << "FEM_TYPE=" << FEM_TYPE << "\n";
cout << "INTEGRATION=" << INTEGRATION << "\n";
/* First step : build the mesh */
bgeot::pgeometric_trans pgt =
bgeot::geometric_trans_descriptor(MESH_TYPE);
size_type N = pgt->dim();
std::vector<size_type> nsubdiv(N);
std::fill(nsubdiv.begin(),nsubdiv.end(),
PARAM.int_value("NX", "Nomber of space steps "));
getfem::regular_unit_mesh(mesh, nsubdiv, pgt,
PARAM.int_value("MESH_NOISED") != 0);
/* scale the unit mesh to [LX,LY,..] and incline it */
bgeot::base_matrix M(N,N);
for (size_type i=0; i < N; ++i) {
static const char *t[] = {"LX","LY","LZ"};
M(i,i) = (i<3) ? PARAM.real_value(t[i],t[i]) : 1.0;
}
mesh.transformation(M);
datafilename = PARAM.string_value("ROOTFILENAME","Base name of data files.");
residual = PARAM.real_value("RESIDUAL"); if (residual == 0.) residual = 1e-10;
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;
}
|