45 using std::endl;
using std::cout;
using std::cerr;
46 using std::ends;
using std::cin;
52 using bgeot::scalar_type;
54 using bgeot::dim_type;
55 using bgeot::base_matrix;
61 typedef getfem::modeling_standard_sparse_matrix sparse_matrix;
62 typedef getfem::modeling_standard_plain_vector plain_vector;
68 gmm::row_matrix<base_small_vector> sol_K;
69 static scalar_type sol_lambda, sol_mu, alph = 0.3;
72 base_small_vector sol_u(
const base_node &x) {
73 int N = x.size(); base_small_vector res(N);
76 for (
int i = 0; i < N; ++i)
77 res[i] = alph * sin(gmm::vect_sp(sol_K.row(i), x));
81 base_small_vector trans(x.size());
82 gmm::fill(trans, M_PI / 10.0);
83 base_node y = x - trans;
85 for (
int i = 0; i < N; ++i) res[i] = a;
90 base_small_vector trans(x.size());
91 gmm::fill(trans, M_PI / 10.0);
92 base_node y = x - trans;
93 scalar_type a = gmm::sqrt(gmm::vect_norm2(y));
94 for (
int i = 0; i < N; ++i) res[i] = a;
101 base_small_vector sol_f(
const base_node &x) {
103 base_small_vector res(N);
106 for (
int i = 0; i < N; i++) {
107 res[i] = alph * ( sol_mu * gmm::vect_sp(sol_K.row(i), sol_K.row(i)) )
108 * sin(gmm::vect_sp(sol_K.row(i), x));
109 for (
int j = 0; j < N; j++)
110 res[i] += alph * ( (sol_lambda + sol_mu) * sol_K(j,j) * sol_K(j,i))
111 * sin(gmm::vect_sp(sol_K.row(j), x));
116 base_small_vector trans(x.size());
117 gmm::fill(trans, M_PI / 10.0);
118 base_node y = x - trans;
120 scalar_type r2 = r*r;
121 scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
123 for (
int i = 0; i < N; i++)
124 res[i] = sol_lambda * (y[i]*tr / r2 - 1.0) / r
125 + sol_mu * (y[i]*tr/r2 - N) / r;
130 base_small_vector trans(x.size());
131 gmm::fill(trans, M_PI / 10.0);
132 base_node y = x - trans;
135 scalar_type a = gmm::sqrt(1.0/r);
136 scalar_type b = a*a*a, c = b*b*a;
137 scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
138 for (
int i = 0; i < N; ++i)
139 res[i] -= b * (sol_lambda + sol_mu * (N+1-3.0/2.0)) * 0.5
140 - c * 3.0 * (sol_lambda + sol_mu) * (y[i]*tr) / 4.0;
147 base_matrix sol_sigma(
const base_node &x) {
149 base_matrix res(N,N);
152 for (
int i = 0; i < N; i++)
153 for (
int j = 0; j <= i; j++) {
154 res(j,i) = res(i,j) = alph * sol_mu *
155 ( sol_K(i,j) * cos(gmm::vect_sp(sol_K.row(i), x))
156 + sol_K(j,i) * cos(gmm::vect_sp(sol_K.row(j), x))
159 for (
int k = 0; k < N; k++)
160 res(i,j) += alph * sol_lambda * sol_K(k,k)
161 * cos(gmm::vect_sp(sol_K.row(k), x));
166 base_small_vector trans(x.size());
167 gmm::fill(trans, M_PI / 10.0);
168 base_node y = x - trans;
170 scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
171 for (
int i = 0; i < N; i++) {
172 res(i, i) += sol_lambda * tr / r;
173 for (
int j = 0; j < N; j++)
174 res(i, j) += sol_mu * (y[i] + y[j]) / r;
181 base_small_vector trans(x.size());
182 gmm::fill(trans, M_PI / 10.0);
183 base_node y = x - trans;
185 scalar_type a = gmm::sqrt(1.0/r);
186 scalar_type b = a*a*a;
187 scalar_type tr(0); tr = std::accumulate(y.begin(), y.end(), tr);
188 for (
int i = 0; i < N; i++) {
189 res(i, i) += sol_lambda * tr * b * 0.5;
190 for (
int j = 0; j < N; j++)
191 res(i, j) += sol_mu * b * (y[i] + y[j]) * 0.5;
202 struct elastostatic_problem {
204 enum { DIRICHLET_BOUNDARY_NUM = 0, NEUMANN_BOUNDARY_NUM = 1};
211 scalar_type lambda, mu;
213 scalar_type residual;
214 bool mixed_pressure, refine;
217 std::string datafilename;
218 bgeot::md_param PARAM;
220 void select_boundaries(
void);
221 bool solve(plain_vector &U);
223 void compute_error(plain_vector &U);
224 elastostatic_problem(
void) : mim(mesh),mf_u(mesh), mf_mult(mesh),
225 mf_rhs(mesh),mf_p(mesh) {}
230 void elastostatic_problem::select_boundaries(
void) {
235 base_node un = mesh.normal_of_face_of_convex(i.cv(), i.f());
237 if (gmm::abs(un[N-1] - 1.0) < 0.5) {
238 mesh.region(NEUMANN_BOUNDARY_NUM).add(i.cv(), i.f());
240 mesh.region(DIRICHLET_BOUNDARY_NUM).add(i.cv(), i.f());
249 void elastostatic_problem::init(
void) {
250 std::string MESH_FILE = PARAM.string_value(
"MESH_FILE",
"Mesh file");
251 std::string FEM_TYPE = PARAM.string_value(
"FEM_TYPE",
"FEM name");
252 std::string INTEGRATION = PARAM.string_value(
"INTEGRATION",
253 "Name of integration method");
255 cout <<
"MESH_FILE=" << MESH_FILE <<
"\n";
256 cout <<
"FEM_TYPE=" << FEM_TYPE <<
"\n";
257 cout <<
"INTEGRATION=" << INTEGRATION <<
"\n";
259 #if GETFEM_PARA_LEVEL > 1
260 double t_init=MPI_Wtime();
265 std::stringstream filename; filename << MESH_FILE;
266 if ((MESH_FILE.compare(0,11,
"structured:") == 0) && NX > 0) {
267 filename <<
";NSUBDIV=[" << NX;
268 for (
size_type i = 1; i < N; ++i) filename <<
"," << NX;
273 GMM_ASSERT1(N == mesh.dim(),
"The mesh has not the right dimension");
275 #if GETFEM_PARA_LEVEL > 1
276 cout<<
"temps creation maillage "<< MPI_Wtime()-t_init<<endl;
280 =
size_type(PARAM.int_value(
"DIRICHLET_VERSION",
281 "Dirichlet version"));
282 datafilename = PARAM.string_value(
"ROOTFILENAME",
"Base name of data files.");
283 scalar_type FT = PARAM.real_value(
"FT",
"parameter for exact solution");
284 residual = PARAM.real_value(
"RESIDUAL");
285 if (residual == 0.) residual = 1e-10;
289 sol_K(i,j) = (i == j) ? FT : -FT;
291 mu = PARAM.real_value(
"MU",
"Lamé coefficient mu");
292 lambda = PARAM.real_value(
"LAMBDA",
"Lamé coefficient lambda");
293 sol_sing = int(PARAM.int_value(
"SOL_SING",
"Optional singular solution"));
294 refine = (PARAM.int_value(
"REFINE",
"Optional refinement") != 0);
295 sol_lambda = lambda; sol_mu = mu;
296 mf_u.set_qdim(dim_type(N));
297 mf_mult.set_qdim(dim_type(N));
302 getfem::pintegration_method ppi =
305 mim.set_integration_method(ppi);
306 mf_u.set_finite_element(pf_u);
308 std::string dirichlet_fem_name = PARAM.string_value(
"DIRICHLET_FEM_TYPE");
309 if (dirichlet_fem_name.size() == 0)
310 mf_mult.set_finite_element(pf_u);
312 cout <<
"DIRICHLET_FEM_TYPE=" << dirichlet_fem_name <<
"\n";
317 (PARAM.int_value(
"MIXED_PRESSURE",
"Mixed version or not.") != 0);
318 if (mixed_pressure) {
319 std::string FEM_TYPE_P = PARAM.string_value(
"FEM_TYPE_P",
"FEM name P");
325 std::string data_fem_name = PARAM.string_value(
"DATA_FEM_TYPE");
326 if (data_fem_name.size() == 0) {
327 GMM_ASSERT1(pf_u->is_lagrange(),
"You are using a non-lagrange FEM. "
328 <<
"In that case you need to set "
329 <<
"DATA_FEM_TYPE in the .param file");
330 mf_rhs.set_finite_element(pf_u);
337 cout <<
"Selecting Neumann and Dirichlet boundaries\n";
340 #if GETFEM_PARA_LEVEL > 1
343 t_init = MPI_Wtime();
345 mf_u.nb_dof(); mf_rhs.nb_dof(); mf_mult.nb_dof();
347 cout<<
"enumerate dof time "<< MPI_Wtime()-t_init<<endl;
349 double t_init = gmm::uclock_sec();
350 mf_u.nb_dof(); mf_rhs.nb_dof(); mf_mult.nb_dof();
351 cout <<
"enumerate dof time " << gmm::uclock_sec() - t_init << endl;
356 void elastostatic_problem::compute_error(plain_vector &U) {
358 std::vector<scalar_type> V(mf_rhs.nb_basic_dof()*N);
360 for (
size_type i = 0; i < mf_rhs.nb_basic_dof(); ++i) {
361 gmm::add(gmm::scaled(sol_u(mf_rhs.point_of_basic_dof(i)), -1.0),
362 gmm::sub_vector(V, gmm::sub_interval(i*N, N)));
368 mf_rhs.set_qdim(dim_type(N));
372 if (getfem::MPI_IS_MASTER())
373 cout <<
"L2 error = " << l2 << endl
374 <<
"H1 error = " << h1 << endl
390 bool elastostatic_problem::solve(plain_vector &U) {
394 if (mixed_pressure) cout <<
"Number of dof for P: " << mf_p.nb_dof() << endl;
395 cout <<
"Number of dof for u: " << mf_u.nb_dof() << endl;
406 (model, mim,
"u",
"lambda",
"mu");
409 if (mixed_pressure) {
413 (model, mim,
"u",
"p",
size_type(-1),
"incomp_coeff");
417 std::vector<scalar_type> F(mf_rhs.nb_dof()*N);
427 (model, mim,
"u",
"NeumannData", NEUMANN_BOUNDARY_NUM);
434 (model, mim,
"u", mf_u, DIRICHLET_BOUNDARY_NUM,
"DirichletData");
437 #if GETFEM_PARA_LEVEL > 1
438 double t_init=MPI_Wtime();
440 dal::bit_vector cvref;
444 cout <<
"Total number of variables : " << model.
nb_dof() << endl;
468 plain_vector ERR(mesh.convex_index().last_true()+1);
469 getfem::error_estimate(mim, mf_u, U, ERR);
473 scalar_type threshold = 0.0001, min_ = 1e18;
475 for (dal::bv_visitor i(mesh.convex_index()); !i.finished(); ++i) {
476 if (ERR[i] > threshold) cvref.add(i);
477 min_ = std::min(min_, ERR[i]);
479 cout <<
"min = " << min_ << endl;
480 cout <<
"Nb elt to be refined : " << cvref.card() << endl;
482 mesh.Bank_refine(cvref);
485 }
while (refine && cvref.card() > 0);
487 #if GETFEM_PARA_LEVEL > 1
488 cout<<
"temps standard solve "<< MPI_Wtime()-t_init<<endl;
492 return (iter.converged());
499 int main(
int argc,
char *argv[]) {
501 GETFEM_MPI_INIT(argc, argv);
503 GMM_SET_EXCEPTION_DEBUG;
510 elastostatic_problem p;
511 p.PARAM.read_command_line(argc, argv);
512 #if GETFEM_PARA_LEVEL > 1
513 double t_ref=MPI_Wtime();
516 #if GETFEM_PARA_LEVEL > 1
517 cout <<
"temps init "<< MPI_Wtime()-t_ref << endl;
519 if (getfem::MPI_IS_MASTER())
520 p.mesh.write_to_file(p.datafilename +
".mesh");
524 #if GETFEM_PARA_LEVEL > 1
526 cout<<
"begining resol"<<endl;
528 if (!p.solve(U)) GMM_ASSERT1(
false,
"Solve has failed");
530 #if GETFEM_PARA_LEVEL > 1
531 cout <<
"temps Resol "<< MPI_Wtime()-t_ref << endl;
535 #if GETFEM_PARA_LEVEL > 1
536 cout <<
"temps error "<< MPI_Wtime()-t_ref << endl;
542 if (p.PARAM.int_value(
"VTK_EXPORT") && getfem::MPI_IS_MASTER()) {
543 cout <<
"export to " << p.datafilename +
".vtk" <<
"..\n";
545 p.PARAM.int_value(
"VTK_EXPORT")==1);
546 exp.exporting(p.mf_u);
547 exp.write_point_data(p.mf_u, U,
"elastostatic_displacement");
548 cout <<
"export done, you can view the data file with (for example)\n"
549 "mayavi2 -d " << p.datafilename <<
".vtk -f ExtractVectorNorm -f "
550 "WarpVector -m Surface -m Outline\n";