38 #include "libmesh/libmesh.h" 39 #include "libmesh/replicated_mesh.h" 40 #include "libmesh/mesh_refinement.h" 41 #include "libmesh/mesh_modification.h" 42 #include "libmesh/mesh_tools.h" 43 #include "libmesh/linear_implicit_system.h" 44 #include "libmesh/equation_systems.h" 45 #include "libmesh/fe.h" 46 #include "libmesh/quadrature.h" 47 #include "libmesh/node.h" 48 #include "libmesh/elem.h" 49 #include "libmesh/dof_map.h" 50 #include "libmesh/vector_value.h" 51 #include "libmesh/tensor_value.h" 52 #include "libmesh/dense_matrix.h" 53 #include "libmesh/dense_submatrix.h" 54 #include "libmesh/dense_vector.h" 55 #include "libmesh/dense_subvector.h" 56 #include "libmesh/sparse_matrix.h" 57 #include "libmesh/numeric_vector.h" 58 #include "libmesh/vtk_io.h" 59 #include "libmesh/exodusII_io.h" 60 #include "libmesh/enum_solver_package.h" 61 #include "libmesh/parallel.h" 62 #include "libmesh/getpot.h" 65 #include "libmesh/face_tri3_subdivision.h" 66 #include "libmesh/mesh_subdivision_support.h" 71 #if defined(LIBMESH_ENABLE_SECOND_DERIVATIVES) && LIBMESH_DIM > 2 76 const std::string & system_name);
80 int main (
int argc,
char ** argv)
87 "--enable-petsc, --enable-trilinos, or --enable-eigen");
90 libmesh_example_requires (3 == LIBMESH_DIM,
"3D support");
93 #ifndef LIBMESH_ENABLE_NODE_VALENCE 94 libmesh_example_requires (
false,
"--enable-node-valence");
98 #ifndef LIBMESH_ENABLE_AMR 99 libmesh_example_requires(
false,
"--enable-amr");
103 #ifndef LIBMESH_ENABLE_SECOND_DERIVATIVES 104 libmesh_example_requires(
false,
"--enable-second");
105 #elif LIBMESH_DIM > 2 119 const int n_refinements =
130 #if defined(LIBMESH_HAVE_VTK) 133 #if defined(LIBMESH_HAVE_EXODUS_API) 154 #if defined(LIBMESH_HAVE_VTK) 157 #if defined(LIBMESH_HAVE_EXODUS_API) 196 equation_systems.
init();
207 #if defined(LIBMESH_HAVE_VTK) 210 #if defined(LIBMESH_HAVE_EXODUS_API) 215 Node * center_node = 0;
217 for (
unsigned int nid=1; nid<
mesh.
n_nodes(); ++nid)
220 if (dist_sq < nearest_dist_sq)
222 nearest_dist_sq = dist_sq;
240 const Real D = E * h*h*h / (12*(1-nu*nu));
241 const Real w_analytic = 0.001265319 * L*L*L*L * q / D;
246 libMesh::out <<
"z-displacement of the center point: " << w << std::endl;
247 libMesh::out <<
"Analytic solution for pure bending: " << w_analytic << std::endl;
249 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES, LIBMESH_DIM > 2 251 #endif // #ifdef LIBMESH_ENABLE_AMR 257 #if defined(LIBMESH_ENABLE_SECOND_DERIVATIVES) && LIBMESH_DIM > 2 265 const std::string & libmesh_dbg_var(system_name))
269 libmesh_assert_equal_to (system_name,
"Shell");
285 const Real K = E * h / (1-nu*nu);
286 const Real D = E * h*h*h / (12*(1-nu*nu));
303 const int extraorder = 0;
307 fe->attach_quadrature_rule (qrule.get());
310 const std::vector<Real> & JxW = fe->get_JxW();
313 const std::vector<RealGradient> & dxyzdxi = fe->get_dxyzdxi();
314 const std::vector<RealGradient> & dxyzdeta = fe->get_dxyzdeta();
317 const std::vector<RealGradient> & d2xyzdxi2 = fe->get_d2xyzdxi2();
318 const std::vector<RealGradient> & d2xyzdeta2 = fe->get_d2xyzdeta2();
319 const std::vector<RealGradient> & d2xyzdxideta = fe->get_d2xyzdxideta();
323 const std::vector<std::vector<Real>> & phi = fe->get_phi();
324 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
325 const std::vector<std::vector<RealTensor>> & d2phi = fe->get_d2phi();
343 Kuu(Ke), Kuv(Ke), Kuw(Ke),
344 Kvu(Ke), Kvv(Ke), Kvw(Ke),
345 Kwu(Ke), Kwv(Ke), Kww(Ke);
355 std::vector<dof_id_type> dof_indices;
356 std::vector<dof_id_type> dof_indices_u;
357 std::vector<dof_id_type> dof_indices_v;
358 std::vector<dof_id_type> dof_indices_w;
362 for (
const auto & elem :
mesh.active_local_element_ptr_range())
369 if (sd_elem->is_ghost())
381 const std::size_t n_dofs = dof_indices.size();
382 const std::size_t n_u_dofs = dof_indices_u.size();
383 const std::size_t n_v_dofs = dof_indices_v.size();
384 const std::size_t n_w_dofs = dof_indices_w.size();
396 Ke.
resize (n_dofs, n_dofs);
412 Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
413 Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
414 Kuw.reposition (u_var*n_u_dofs, w_var*n_u_dofs, n_u_dofs, n_w_dofs);
416 Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
417 Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
418 Kvw.reposition (v_var*n_v_dofs, w_var*n_v_dofs, n_v_dofs, n_w_dofs);
420 Kwu.reposition (w_var*n_w_dofs, u_var*n_w_dofs, n_w_dofs, n_u_dofs);
421 Kwv.reposition (w_var*n_w_dofs, v_var*n_w_dofs, n_w_dofs, n_v_dofs);
422 Kww.
reposition (w_var*n_w_dofs, w_var*n_w_dofs, n_w_dofs, n_w_dofs);
424 Fu.reposition (u_var*n_u_dofs, n_u_dofs);
425 Fv.reposition (v_var*n_u_dofs, n_v_dofs);
429 for (
unsigned int qp=0; qp<qrule->n_points(); ++qp)
435 for (
unsigned int i=0; i<n_u_dofs; ++i)
436 Fw(i) += JxW[qp] * phi[i][qp] * q;
447 libmesh_assert_greater (jac, 0);
465 H(0,0) = a(1) * a(1);
466 H(0,1) = H(1,0) = nu * a(1) * a(0) + (1-nu) * a(2) * a(2);
467 H(0,2) = H(2,0) = -a(1) * a(2);
468 H(1,1) = a(0) * a(0);
469 H(1,2) = H(2,1) = -a(0) * a(2);
470 H(2,2) = 0.5 * ((1-nu) * a(1) * a(0) + (1+nu) * a(2) * a(2));
471 const Real det = a(0) * a(1) - a(2) * a(2);
472 libmesh_assert_not_equal_to (det * det, 0);
486 for (
unsigned int i=0; i<n_u_dofs; ++i)
488 for (
unsigned int j=0; j<n_u_dofs; ++j)
492 for (
unsigned int k=0; k<3; ++k)
494 MI(0,k) = dphi[i][qp](0) * a1(k);
495 MI(1,k) = dphi[i][qp](1) * a2(k);
496 MI(2,k) = dphi[i][qp](1) * a1(k)
497 + dphi[i][qp](0) * a2(k);
499 MJ(0,k) = dphi[j][qp](0) * a1(k);
500 MJ(1,k) = dphi[j][qp](1) * a2(k);
501 MJ(2,k) = dphi[j][qp](1) * a1(k)
502 + dphi[j][qp](0) * a2(k);
507 for (
unsigned int k=0; k<3; ++k)
509 const Real term_ik = dphi[i][qp](0) * a2xa3(k)
510 + dphi[i][qp](1) * a3xa1(k);
511 BI(0,k) = -d2phi[i][qp](0,0) * a3(k)
512 +(dphi[i][qp](0) * a11xa2(k)
513 + dphi[i][qp](1) * a1xa11(k)
514 + (a3*a11) * term_ik) / jac;
515 BI(1,k) = -d2phi[i][qp](1,1) * a3(k)
516 +(dphi[i][qp](0) * a22xa2(k)
517 + dphi[i][qp](1) * a1xa22(k)
518 + (a3*a22) * term_ik) / jac;
519 BI(2,k) = 2 * (-d2phi[i][qp](0,1) * a3(k)
520 +(dphi[i][qp](0) * a12xa2(k)
521 + dphi[i][qp](1) * a1xa12(k)
522 + (a3*a12) * term_ik) / jac);
524 const Real term_jk = dphi[j][qp](0) * a2xa3(k)
525 + dphi[j][qp](1) * a3xa1(k);
526 BJ(0,k) = -d2phi[j][qp](0,0) * a3(k)
527 +(dphi[j][qp](0) * a11xa2(k)
528 + dphi[j][qp](1) * a1xa11(k)
529 + (a3*a11) * term_jk) / jac;
530 BJ(1,k) = -d2phi[j][qp](1,1) * a3(k)
531 +(dphi[j][qp](0) * a22xa2(k)
532 + dphi[j][qp](1) * a1xa22(k)
533 + (a3*a22) * term_jk) / jac;
534 BJ(2,k) = 2 * (-d2phi[j][qp](0,1) * a3(k)
535 +(dphi[j][qp](0) * a12xa2(k)
536 + dphi[j][qp](1) * a1xa12(k)
537 + (a3*a12) * term_jk) / jac);
549 Kuu(i,j) += KIJ(0,0);
550 Kuv(i,j) += KIJ(0,1);
551 Kuw(i,j) += KIJ(0,2);
553 Kvu(i,j) += KIJ(1,0);
554 Kvv(i,j) += KIJ(1,1);
555 Kvw(i,j) += KIJ(1,2);
557 Kwu(i,j) += KIJ(2,0);
558 Kwv(i,j) += KIJ(2,1);
559 Kww(i,j) += KIJ(2,2);
580 for (
const auto & elem :
mesh.active_local_element_ptr_range())
586 if (!gh_elem->is_ghost())
591 for (
auto s : elem->side_index_range())
594 if (nb_elem ==
nullptr || nb_elem->
is_ghost())
612 const Node * nodes [4];
613 nodes[1] = gh_elem->node_ptr(s);
620 unsigned int n_int = 0;
622 while (nodes[0]->
id() == nodes[1]->
id() || nodes[0]->
id() == nodes[2]->
id())
623 nodes[0] = nb_elem->
node_ptr(++n_int);
626 const Real penalty = 1.e10;
632 for (
unsigned int n=0; n<4; ++n)
637 matrix.
add (u_dof, u_dof, penalty);
638 matrix.
add (v_dof, v_dof, penalty);
639 matrix.
add (w_dof, w_dof, penalty);
645 #endif // defined(LIBMESH_ENABLE_SECOND_DERIVATIVES) && LIBMESH_DIM > 2 class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
T command_line_next(std::string name, T default_value)
Use GetPot's search()/next() functions to get following arguments from the command line...
dof_id_type end_dof(const processor_id_type proc) const
This is the EquationSystems class.
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
A Node is like a Point, but with more information.
auto norm() const -> decltype(std::norm(T()))
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
TypeTensor< T > transpose() const
void resize(const unsigned int n)
Resize the vector.
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i]...
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
virtual void write(const std::string &) override
Output the mesh without solutions to a .pvtu file.
void assemble_shell(EquationSystems &es, const std::string &system_name)
NumericVector< Number > * rhs
The system matrix.
const Parallel::Communicator & comm() const
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
The libMesh namespace provides an interface to certain functionality in the library.
This class implements reading and writing meshes in the VTK format.
virtual void solve() override
Assembles & solves the linear system A*x=b.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Number current_solution(const dof_id_type global_dof_number) const
unsigned int variable_number(std::string_view var) const
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
Add value to the element (i,j).
Defines a dense subvector for use in finite element computations.
SolverPackage default_solver_package()
This class handles the numbering of degrees of freedom on a mesh.
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
void reposition(const unsigned int ioff, const unsigned int n)
Changes the location of the subvector in the parent vector.
unsigned int number() const
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
auto norm_sq() const -> decltype(std::norm(T()))
const T & get(std::string_view) const
Implements (adaptive) mesh refinement algorithms for a MeshBase.
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Defines a dense submatrix for use in Finite Element-type computations.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
The Tri3Subdivision element is a three-noded subdivision surface shell element used in mechanics calc...
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
const FEType & variable_type(const unsigned int i) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Node * node_ptr(const unsigned int i) const
virtual void write(const std::string &fname) override
This method implements writing a mesh to a specified file.
T & set(const std::string &)
const MeshBase & get_mesh() const
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Parameters parameters
Data structure holding arbitrary parameters.
void reposition(const unsigned int ioff, const unsigned int joff, const unsigned int new_m, const unsigned int new_n)
Changes the location of the submatrix in the parent matrix.
virtual const Point & point(const dof_id_type i) const =0
virtual void init()
Initialize all the systems.
dof_id_type first_dof(const processor_id_type proc) const
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
virtual const Node * node_ptr(const dof_id_type i) const =0
const DofMap & get_dof_map() const
int main(int argc, char **argv)
const SparseMatrix< Number > & get_system_matrix() const
virtual dof_id_type n_nodes() const =0
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.