42 #include "libmesh/libmesh_config.h"
43 #include "libmesh/libmesh.h"
44 #include "libmesh/mesh.h"
45 #include "libmesh/mesh_generation.h"
46 #include "libmesh/exodusII_io.h"
47 #include "libmesh/gnuplot_io.h"
48 #include "libmesh/linear_implicit_system.h"
49 #include "libmesh/equation_systems.h"
50 #include "libmesh/fe.h"
51 #include "libmesh/quadrature_gauss.h"
52 #include "libmesh/dof_map.h"
53 #include "libmesh/sparse_matrix.h"
54 #include "libmesh/numeric_vector.h"
55 #include "libmesh/dense_matrix.h"
56 #include "libmesh/dense_submatrix.h"
57 #include "libmesh/dense_vector.h"
58 #include "libmesh/dense_subvector.h"
59 #include "libmesh/perf_log.h"
60 #include "libmesh/elem.h"
61 #include "libmesh/boundary_info.h"
62 #include "libmesh/zero_function.h"
63 #include "libmesh/dirichlet_boundaries.h"
64 #include "libmesh/string_to_enum.h"
65 #include "libmesh/getpot.h"
66 #include "libmesh/solver_configuration.h"
67 #include "libmesh/petsc_linear_solver.h"
68 #include "libmesh/petsc_macro.h"
69 #include "libmesh/enum_solver_package.h"
74 #define BOUNDARY_ID_MIN_Z 0
75 #define BOUNDARY_ID_MIN_Y 1
76 #define BOUNDARY_ID_MAX_X 2
77 #define BOUNDARY_ID_MAX_Y 3
78 #define BOUNDARY_ID_MIN_X 4
79 #define BOUNDARY_ID_MAX_Z 5
80 #define NODE_BOUNDARY_ID 10
81 #define EDGE_BOUNDARY_ID 20
86 #ifdef LIBMESH_HAVE_PETSC
94 _petsc_linear_solver(petsc_linear_solver)
100 PetscErrorCode
ierr = 0;
101 ierr = KSPSetType (_petsc_linear_solver.ksp(), const_cast<KSPType>(KSPCG));
102 CHKERRABORT(_petsc_linear_solver.comm().get(),
ierr);
104 ierr = PCSetType (_petsc_linear_solver.pc(), const_cast<PCType>(PCBJACOBI));
105 CHKERRABORT(_petsc_linear_solver.comm().get(),
ierr);
131 return i == j ? 1. : 0.;
143 const Real poisson_ratio = 0.3;
144 const Real young_modulus = 1.;
147 const Real lambda_1 = (young_modulus*poisson_ratio)/((1.+poisson_ratio)*(1.-2.*poisson_ratio));
148 const Real lambda_2 = young_modulus/(2.*(1.+poisson_ratio));
171 fe->attach_quadrature_rule (&qrule);
175 fe_face->attach_quadrature_rule (&qface);
177 const std::vector<Real> & JxW = fe->get_JxW();
178 const std::vector<std::vector<Real>> & phi = fe->get_phi();
179 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
196 std::vector<dof_id_type> dof_indices;
197 std::vector<std::vector<dof_id_type>> dof_indices_var(3);
202 for (
unsigned int var=0; var<3; var++)
203 dof_map.
dof_indices (elem, dof_indices_var[var], var);
205 const unsigned int n_dofs = dof_indices.size();
206 const unsigned int n_var_dofs = dof_indices_var[0].size();
210 Ke.
resize (n_dofs, n_dofs);
211 for (
unsigned int var_i=0; var_i<3; var_i++)
212 for (
unsigned int var_j=0; var_j<3; var_j++)
213 Ke_var[var_i][var_j].reposition (var_i*n_var_dofs, var_j*n_var_dofs, n_var_dofs, n_var_dofs);
216 for (
unsigned int var=0; var<3; var++)
217 Fe_var[var].reposition (var*n_var_dofs, n_var_dofs);
219 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
222 for (
unsigned int dof_i=0; dof_i<n_var_dofs; dof_i++)
223 for (
unsigned int dof_j=0; dof_j<n_var_dofs; dof_j++)
224 for (
unsigned int i=0; i<3; i++)
225 for (
unsigned int j=0; j<3; j++)
226 for (
unsigned int k=0; k<3; k++)
227 for (
unsigned int l=0; l<3; l++)
228 Ke_var[i][k](dof_i,dof_j) +=
236 for (
unsigned int dof_i=0; dof_i<n_var_dofs; dof_i++)
237 for (
unsigned int i=0; i<3; i++)
238 Fe_var[i](dof_i) += JxW[qp] * (f_vec(i) * phi[dof_i][qp]);
247 for (
auto side : elem->side_index_range())
248 if (elem->neighbor_ptr(side) ==
nullptr)
250 const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
251 const std::vector<Real> & JxW_face = fe_face->get_JxW();
253 fe_face->reinit(elem, side);
256 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
258 for (
unsigned int dof_i=0; dof_i<n_var_dofs; dof_i++)
259 for (
unsigned int i=0; i<3; i++)
260 Fe_var[i](dof_i) += JxW_face[qp] * (g_vec(i) * phi_face[dof_i][qp]);
279 unsigned int displacement_vars[3];
289 fe->attach_quadrature_rule (&qrule);
291 const std::vector<Real> & JxW = fe->get_JxW();
292 const std::vector<std::vector<Real>> & phi = fe->get_phi();
293 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
298 unsigned int sigma_vars[6];
305 unsigned int vonMises_var = stress_system.
variable_number (
"vonMises");
308 std::vector<std::vector<dof_id_type>> dof_indices_var(system.
n_vars());
309 std::vector<dof_id_type> stress_dof_indices_var;
310 std::vector<dof_id_type> vonmises_dof_indices_var;
314 for (
unsigned int var=0; var<3; var++)
315 dof_map.
dof_indices (elem, dof_indices_var[var], displacement_vars[var]);
317 const unsigned int n_var_dofs = dof_indices_var[0].size();
321 std::vector<DenseMatrix<Number>> stress_tensor_qp(qrule.n_points());
322 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
324 stress_tensor_qp[qp].resize(3,3);
328 for (
unsigned int var_i=0; var_i<3; var_i++)
329 for (
unsigned int var_j=0; var_j<3; var_j++)
330 for (
unsigned int j=0; j<n_var_dofs; j++)
331 grad_u(var_i,var_j) += dphi[j][qp](var_j) * system.
current_solution(dof_indices_var[var_i][j]);
333 for (
unsigned int var_i=0; var_i<3; var_i++)
334 for (
unsigned int var_j=0; var_j<3; var_j++)
335 for (
unsigned int k=0; k<3; k++)
336 for (
unsigned int l=0; l<3; l++)
337 stress_tensor_qp[qp](var_i,var_j) +=
elasticity_tensor(var_i,var_j,k,l) * grad_u(k,l);
340 stress_dof_map.dof_indices (elem, vonmises_dof_indices_var, vonMises_var);
341 std::vector<DenseMatrix<Number>> elem_sigma_vec(vonmises_dof_indices_var.size());
342 for (std::size_t index=0; index<elem_sigma_vec.size(); index++)
343 elem_sigma_vec[index].resize(3,3);
348 unsigned int stress_var_index = 0;
349 for (
unsigned int var_i=0; var_i<3; var_i++)
350 for (
unsigned int var_j=var_i; var_j<3; var_j++)
352 stress_dof_map.dof_indices (elem, stress_dof_indices_var, sigma_vars[stress_var_index]);
354 const unsigned int n_proj_dofs = stress_dof_indices_var.size();
357 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
359 for(
unsigned int i=0; i<n_proj_dofs; i++)
360 for(
unsigned int j=0; j<n_proj_dofs; j++)
362 Me(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp]);
367 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
368 for(
unsigned int i=0; i<n_proj_dofs; i++)
370 Fe(i) += JxW[qp] * stress_tensor_qp[qp](var_i,var_j) * phi[i][qp];
376 for(
unsigned int index=0; index<n_proj_dofs; index++)
378 dof_id_type dof_index = stress_dof_indices_var[index];
379 if ((stress_system.
solution->first_local_index() <= dof_index) &&
380 (dof_index < stress_system.solution->last_local_index()))
381 stress_system.
solution->set(dof_index, projected_data(index));
383 elem_sigma_vec[index](var_i,var_j) = projected_data(index);
389 for (std::size_t index=0; index<elem_sigma_vec.size(); index++)
391 elem_sigma_vec[index](1,0) = elem_sigma_vec[index](0,1);
392 elem_sigma_vec[index](2,0) = elem_sigma_vec[index](0,2);
393 elem_sigma_vec[index](2,1) = elem_sigma_vec[index](1,2);
396 Number vonMises_value =
std::sqrt(0.5*(
pow(elem_sigma_vec[index](0,0) - elem_sigma_vec[index](1,1), 2.) +
397 pow(elem_sigma_vec[index](1,1) - elem_sigma_vec[index](2,2), 2.) +
398 pow(elem_sigma_vec[index](2,2) - elem_sigma_vec[index](0,0), 2.) +
399 6.*(
pow(elem_sigma_vec[index](0,1), 2.) +
400 pow(elem_sigma_vec[index](1,2), 2.) +
401 pow(elem_sigma_vec[index](2,0), 2.))));
403 dof_id_type dof_index = vonmises_dof_indices_var[index];
405 if ((stress_system.
solution->first_local_index() <= dof_index) &&
406 (dof_index < stress_system.solution->last_local_index()))
407 stress_system.
solution->set(dof_index, vonMises_value);
419 int main (
int argc,
char ** argv)
426 "--enable-petsc, --enable-trilinos, or --enable-eigen");
429 const unsigned int dim = 3;
432 libmesh_example_requires(
dim == LIBMESH_DIM,
"3D support");
435 #ifndef LIBMESH_ENABLE_DIRICHLET
436 libmesh_example_requires(
false,
"--enable-dirichlet");
459 side_max_x = 0, side_min_y = 0,
460 side_max_y = 0, side_max_z = 0;
463 found_side_max_x =
false, found_side_max_y =
false,
464 found_side_min_y =
false, found_side_max_z =
false;
466 for (
auto side : elem->side_index_range())
471 found_side_max_x =
true;
477 found_side_min_y =
true;
483 found_side_max_y =
true;
489 found_side_max_z =
true;
496 if (found_side_max_x && found_side_max_y && found_side_max_z)
497 for (
auto n : elem->node_index_range())
498 if (elem->is_node_on_side(n, side_max_x) &&
499 elem->is_node_on_side(n, side_max_y) &&
500 elem->is_node_on_side(n, side_max_z))
507 if (found_side_max_x && found_side_min_y)
508 for (
auto e : elem->edge_index_range())
509 if (elem->is_edge_on_side(e, side_max_x) &&
510 elem->is_edge_on_side(e, side_min_y))
522 #ifdef LIBMESH_HAVE_PETSC
528 petsc_linear_solver->set_solver_configuration(petsc_solver_config);
534 #ifdef LIBMESH_ENABLE_DIRICHLET
540 std::set<boundary_id_type> boundary_ids;
541 boundary_ids.insert(BOUNDARY_ID_MIN_X);
542 boundary_ids.insert(NODE_BOUNDARY_ID);
543 boundary_ids.insert(EDGE_BOUNDARY_ID);
546 std::vector<unsigned int> variables;
547 variables.push_back(u_var);
548 variables.push_back(v_var);
549 variables.push_back(w_var);
562 #endif // LIBMESH_ENABLE_DIRICHLET
577 equation_systems.
init();
589 #ifdef LIBMESH_HAVE_EXODUS_API
595 #endif // #ifdef LIBMESH_HAVE_EXODUS_API