44 using std::endl;
using std::cout;
using std::cerr;
45 using std::ends;
using std::cin;
51 using bgeot::scalar_type;
53 using bgeot::dim_type;
54 using bgeot::base_matrix;
60 typedef getfem::modeling_standard_sparse_matrix sparse_matrix;
61 typedef getfem::modeling_standard_plain_vector plain_vector;
67 gmm::row_matrix<base_small_vector> sol_K;
68 static scalar_type sol_lambda, sol_mu, alph = 0.3;
71 base_small_vector sol_u(
const base_node &x) {
72 int N = x.size(); base_small_vector res(N);
75 for (
int i = 0; i < N; ++i)
76 res[i] = alph * sin(gmm::vect_sp(sol_K.row(i), x));
80 base_small_vector trans(x.size());
82 base_node y = x - trans;
84 for (
int i = 0; i < N; ++i) res[i] = a;
89 base_small_vector trans(x.size());
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;
100 base_small_vector sol_f(
const base_node &x) {
102 base_small_vector res(N);
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));
115 base_small_vector trans(x.size());
117 base_node y = x - trans;
119 scalar_type r2 = r*r;
120 scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
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;
129 base_small_vector trans(x.size());
131 base_node y = x - trans;
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;
146 base_matrix sol_sigma(
const base_node &x) {
148 base_matrix res(N,N);
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))
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));
165 base_small_vector trans(x.size());
167 base_node y = x - trans;
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;
180 base_small_vector trans(x.size());
182 base_node y = x - trans;
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;
201 struct elastostatic_problem {
203 enum { DIRICHLET_BOUNDARY_NUM = 0, NEUMANN_BOUNDARY_NUM = 1};
210 scalar_type lambda, mu;
212 scalar_type residual;
213 bool mixed_pressure, refine;
216 std::string datafilename;
217 bgeot::md_param PARAM;
219 void select_boundaries(
void);
220 bool solve(plain_vector &U);
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) {}
229 void elastostatic_problem::select_boundaries(
void) {
234 base_node un = mesh.normal_of_face_of_convex(i.cv(), i.f());
236 if (gmm::abs(un[N-1] - 1.0) < 0.5) {
237 mesh.region(NEUMANN_BOUNDARY_NUM).add(i.cv(), i.f());
239 mesh.region(DIRICHLET_BOUNDARY_NUM).add(i.cv(), i.f());
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");
254 cout <<
"MESH_FILE=" << MESH_FILE <<
"\n";
255 cout <<
"FEM_TYPE=" << FEM_TYPE <<
"\n";
256 cout <<
"INTEGRATION=" << INTEGRATION <<
"\n";
258 #if GETFEM_PARA_LEVEL > 1
259 double t_init=MPI_Wtime();
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;
272 GMM_ASSERT1(N == mesh.dim(),
"The mesh has not the right dimension");
274 #if GETFEM_PARA_LEVEL > 1
275 cout<<
"temps creation maillage "<< MPI_Wtime()-t_init<<endl;
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;
288 sol_K(i,j) = (i == j) ? FT : -FT;
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));
301 getfem::pintegration_method ppi =
304 mim.set_integration_method(ppi);
305 mf_u.set_finite_element(pf_u);
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);
311 cout <<
"DIRICHLET_FEM_TYPE=" << dirichlet_fem_name <<
"\n";
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");
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);
336 cout <<
"Selecting Neumann and Dirichlet boundaries\n";
339 #if GETFEM_PARA_LEVEL > 1
342 t_init = MPI_Wtime();
344 mf_u.nb_dof(); mf_rhs.nb_dof(); mf_mult.nb_dof();
346 cout<<
"enumerate dof time "<< MPI_Wtime()-t_init<<endl;
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;
355 void elastostatic_problem::compute_error(plain_vector &U) {
357 std::vector<scalar_type> V(mf_rhs.nb_basic_dof()*N);
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)));
367 mf_rhs.set_qdim(dim_type(N));
371 if (getfem::MPI_IS_MASTER())
372 cout <<
"L2 error = " << l2 << endl
373 <<
"H1 error = " << h1 << endl
389 bool elastostatic_problem::solve(plain_vector &U) {
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;
405 (model, mim,
"u",
"lambda",
"mu");
408 if (mixed_pressure) {
412 (model, mim,
"u",
"p",
size_type(-1),
"incomp_coeff");
416 std::vector<scalar_type> F(mf_rhs.nb_dof()*N);
426 (model, mim,
"u",
"NeumannData", NEUMANN_BOUNDARY_NUM);
433 (model, mim,
"u", mf_u, DIRICHLET_BOUNDARY_NUM,
"DirichletData");
436 #if GETFEM_PARA_LEVEL > 1
437 double t_init=MPI_Wtime();
439 dal::bit_vector cvref;
443 cout <<
"Total number of variables : " << model.
nb_dof() << endl;
467 plain_vector ERR(mesh.convex_index().last_true()+1);
468 getfem::error_estimate(mim, mf_u, U, ERR);
472 scalar_type threshold = 0.0001, min_ = 1e18;
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]);
478 cout <<
"min = " << min_ << endl;
479 cout <<
"Nb elt to be refined : " << cvref.card() << endl;
481 mesh.Bank_refine(cvref);
484 }
while (refine && cvref.card() > 0);
486 #if GETFEM_PARA_LEVEL > 1
487 cout<<
"temps standard solve "<< MPI_Wtime()-t_init<<endl;
491 return (iter.converged());
498 int main(
int argc,
char *argv[]) {
500 GETFEM_MPI_INIT(argc, argv);
502 GMM_SET_EXCEPTION_DEBUG;
509 elastostatic_problem p;
510 p.PARAM.read_command_line(argc, argv);
511 #if GETFEM_PARA_LEVEL > 1
512 double t_ref=MPI_Wtime();
515 #if GETFEM_PARA_LEVEL > 1
516 cout <<
"temps init "<< MPI_Wtime()-t_ref << endl;
518 if (getfem::MPI_IS_MASTER())
519 p.mesh.write_to_file(p.datafilename +
".mesh");
523 #if GETFEM_PARA_LEVEL > 1
525 cout<<
"begining resol"<<endl;
527 if (!p.solve(U)) GMM_ASSERT1(
false,
"Solve has failed");
529 #if GETFEM_PARA_LEVEL > 1
530 cout <<
"temps Resol "<< MPI_Wtime()-t_ref << endl;
534 #if GETFEM_PARA_LEVEL > 1
535 cout <<
"temps error "<< MPI_Wtime()-t_ref << endl;
541 if (p.PARAM.int_value(
"VTK_EXPORT") && getfem::MPI_IS_MASTER()) {
542 cout <<
"export to " << p.datafilename +
".vtk" <<
"..\n";
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";
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).
`‘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.
The Iteration object calculates whether the solution has reached the desired accuracy,...
sparse vector built upon std::vector.
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.
Include common gmm files.
void copy(const L1 &l1, L2 &l2)
*â€/
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
void fill(L &l, typename gmm::linalg_traits< L >::value_type x)
*â€/
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norminf(const V &v)
Infinity norm of a vector.
void resize(V &v, size_type n)
*â€/
strongest_value_type< V1, V2 >::value_type vect_sp(const V1 &v1, const V2 &v2)
*â€/
void add(const L1 &l1, L2 &l2)
*â€/
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.
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
pfem fem_descriptor(const std::string &name)
get a fem descriptor from its string name.
size_t size_type
used as the common size type in the library
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.