libMesh
systems_of_equations_ex5.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // <h1> Systems Example 5 - Linear Elastic Cantilever with Constraint </h1>
21 // \author David Knezevic
22 // \date 2012
23 //
24 // In this example we extend systems_of_equations_ex4 to enforce a constraint.
25 // We apply a uniform load on the top surface of the cantilever, and we
26 // determine the traction on the right boundary in order to obtain zero
27 // average vertical displacement on the right boundary of the domain.
28 //
29 // This constraint is enforced via a Lagrange multiplier (SCALAR variable).
30 // The system we solve, therefore, is of the form:
31 // a(u,v) + \lambda g(v) = f(v)
32 // g(u) = 0
33 // Here \lambda tells us the traction required to satisfy the constraint.
34 
35 // C++ include files that we need
36 #include <iostream>
37 #include <algorithm>
38 #include <math.h>
39 
40 // libMesh includes
41 #include "libmesh/libmesh.h"
42 #include "libmesh/mesh.h"
43 #include "libmesh/mesh_generation.h"
44 #include "libmesh/exodusII_io.h"
45 #include "libmesh/linear_implicit_system.h"
46 #include "libmesh/equation_systems.h"
47 #include "libmesh/fe.h"
48 #include "libmesh/quadrature_gauss.h"
49 #include "libmesh/dof_map.h"
50 #include "libmesh/sparse_matrix.h"
51 #include "libmesh/numeric_vector.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/perf_log.h"
57 #include "libmesh/elem.h"
58 #include "libmesh/boundary_info.h"
59 #include "libmesh/zero_function.h"
60 #include "libmesh/dirichlet_boundaries.h"
61 #include "libmesh/string_to_enum.h"
62 #include "libmesh/getpot.h"
63 #include "libmesh/enum_solver_package.h"
64 
65 // Bring in everything from the libMesh namespace
66 using namespace libMesh;
67 
68 // Matrix and right-hand side assemble
70  const std::string & system_name);
71 
72 // Define the elasticity tensor, which is a fourth-order tensor
73 // i.e. it has four indices i, j, k, l
74 Real eval_elasticity_tensor(unsigned int i,
75  unsigned int j,
76  unsigned int k,
77  unsigned int l);
78 
79 // Begin the main program.
80 int main (int argc, char ** argv)
81 {
82  // Initialize libMesh and any dependent libraries
83  LibMeshInit init (argc, argv);
84 
85  // This example requires a linear solver package.
86  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
87  "--enable-petsc, --enable-trilinos, or --enable-eigen");
88 
89  // This example NaNs with the Eigen sparse linear solvers
90  libmesh_example_requires(libMesh::default_solver_package() != EIGEN_SOLVERS, "--enable-petsc or --enable-laspack");
91 
92  // Initialize the cantilever mesh
93  const unsigned int dim = 2;
94 
95  // Skip this 2D example if libMesh was compiled as 1D-only.
96  libmesh_example_requires(dim <= LIBMESH_DIM, "2D support");
97 
98  // We use Dirichlet boundary conditions here
99 #ifndef LIBMESH_ENABLE_DIRICHLET
100  libmesh_example_requires(false, "--enable-dirichlet");
101 #endif
102 
103  // Create a 2D mesh distributed across the default MPI communicator.
104  Mesh mesh(init.comm(), dim);
106  50, 10,
107  0., 1.,
108  0., 0.2,
109  QUAD9);
110 
111  // Print information about the mesh to the screen.
112  mesh.print_info();
113 
114  // Create an equation systems object.
115  EquationSystems equation_systems (mesh);
116 
117  // Declare the system and its variables.
118  // Create a system named "Elasticity"
119  LinearImplicitSystem & system =
120  equation_systems.add_system<LinearImplicitSystem> ("Elasticity");
121 
122 #ifdef LIBMESH_ENABLE_DIRICHLET
123  // Add two displacement variables, u and v, to the system
124  unsigned int u_var = system.add_variable("u", SECOND, LAGRANGE);
125  unsigned int v_var = system.add_variable("v", SECOND, LAGRANGE);
126 
127  // Add a SCALAR variable for the Lagrange multiplier to enforce our constraint
128  system.add_variable("lambda", FIRST, SCALAR);
129 
131 
132  // Construct a Dirichlet boundary condition object
133  // We impose a "clamped" boundary condition on the
134  // "left" boundary, i.e. bc_id = 3
135  std::set<boundary_id_type> boundary_ids;
136  boundary_ids.insert(3);
137 
138  // Create a vector storing the variable numbers which the BC applies to
139  std::vector<unsigned int> variables(2);
140  variables[0] = u_var;
141  variables[1] = v_var;
142 
143  // Create a ZeroFunction to initialize dirichlet_bc
144  ZeroFunction<> zf;
145 
146  // Most DirichletBoundary users will want to supply a "locally
147  // indexed" functor
148  DirichletBoundary dirichlet_bc(boundary_ids, variables, zf,
150 
151  // We must add the Dirichlet boundary condition _before_
152  // we call equation_systems.init()
153  system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
154 #endif // LIBMESH_ENABLE_DIRICHLET
155 
156  // Initialize the data structures for the equation system.
157  equation_systems.init();
158 
159  // Print information about the system to the screen.
160  equation_systems.print_info();
161 
162  // Solve the system
163  system.solve();
164 
165  // Plot the solution
166 #ifdef LIBMESH_HAVE_EXODUS_API
167  ExodusII_IO (mesh).write_equation_systems("displacement.e", equation_systems);
168 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
169 
170  // All done.
171  return 0;
172 }
173 
174 
176  const std::string & libmesh_dbg_var(system_name))
177 {
178  libmesh_assert_equal_to (system_name, "Elasticity");
179 
180  const MeshBase & mesh = es.get_mesh();
181 
182  const unsigned int dim = mesh.mesh_dimension();
183 
184  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Elasticity");
185 
186  const unsigned int u_var = system.variable_number ("u");
187  const unsigned int v_var = system.variable_number ("v");
188  const unsigned int lambda_var = system.variable_number ("lambda");
189 
190  const DofMap & dof_map = system.get_dof_map();
191  FEType fe_type = dof_map.variable_type(0);
192  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
193  QGauss qrule (dim, fe_type.default_quadrature_order());
194  fe->attach_quadrature_rule (&qrule);
195 
196  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
197  QGauss qface(dim-1, fe_type.default_quadrature_order());
198  fe_face->attach_quadrature_rule (&qface);
199 
200  const std::vector<Real> & JxW = fe->get_JxW();
201  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
202 
205 
207  Kuu(Ke), Kuv(Ke),
208  Kvu(Ke), Kvv(Ke);
209  DenseSubMatrix<Number> Klambda_v(Ke), Kv_lambda(Ke);
210 
212  Fu(Fe),
213  Fv(Fe);
214 
215  std::vector<dof_id_type> dof_indices;
216  std::vector<dof_id_type> dof_indices_u;
217  std::vector<dof_id_type> dof_indices_v;
218  std::vector<dof_id_type> dof_indices_lambda;
219 
220  for (const auto & elem : mesh.active_local_element_ptr_range())
221  {
222  dof_map.dof_indices (elem, dof_indices);
223  dof_map.dof_indices (elem, dof_indices_u, u_var);
224  dof_map.dof_indices (elem, dof_indices_v, v_var);
225  dof_map.dof_indices (elem, dof_indices_lambda, lambda_var);
226 
227  const unsigned int n_dofs = dof_indices.size();
228  const unsigned int n_u_dofs = dof_indices_u.size();
229  const unsigned int n_v_dofs = dof_indices_v.size();
230  const unsigned int n_lambda_dofs = dof_indices_lambda.size();
231 
232  fe->reinit (elem);
233 
234  Ke.resize (n_dofs, n_dofs);
235  Fe.resize (n_dofs);
236 
237  Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
238  Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
239 
240  Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
241  Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
242 
243  // Also, add a row and a column to enforce the constraint
244  Kv_lambda.reposition (v_var*n_u_dofs, v_var*n_u_dofs+n_v_dofs, n_v_dofs, 1);
245  Klambda_v.reposition (v_var*n_v_dofs+n_v_dofs, v_var*n_v_dofs, 1, n_v_dofs);
246 
247  Fu.reposition (u_var*n_u_dofs, n_u_dofs);
248  Fv.reposition (v_var*n_u_dofs, n_v_dofs);
249 
250  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
251  {
252  for (unsigned int i=0; i<n_u_dofs; i++)
253  for (unsigned int j=0; j<n_u_dofs; j++)
254  {
255  // Tensor indices
256  unsigned int C_i, C_j, C_k, C_l;
257  C_i=0, C_k=0;
258 
259  C_j=0, C_l=0;
260  Kuu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
261 
262  C_j=1, C_l=0;
263  Kuu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
264 
265  C_j=0, C_l=1;
266  Kuu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
267 
268  C_j=1, C_l=1;
269  Kuu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
270  }
271 
272  for (unsigned int i=0; i<n_u_dofs; i++)
273  for (unsigned int j=0; j<n_v_dofs; j++)
274  {
275  // Tensor indices
276  unsigned int C_i, C_j, C_k, C_l;
277  C_i=0, C_k=1;
278 
279  C_j=0, C_l=0;
280  Kuv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
281 
282  C_j=1, C_l=0;
283  Kuv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
284 
285  C_j=0, C_l=1;
286  Kuv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
287 
288  C_j=1, C_l=1;
289  Kuv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
290  }
291 
292  for (unsigned int i=0; i<n_v_dofs; i++)
293  for (unsigned int j=0; j<n_u_dofs; j++)
294  {
295  // Tensor indices
296  unsigned int C_i, C_j, C_k, C_l;
297  C_i=1, C_k=0;
298 
299  C_j=0, C_l=0;
300  Kvu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
301 
302  C_j=1, C_l=0;
303  Kvu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
304 
305  C_j=0, C_l=1;
306  Kvu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
307 
308  C_j=1, C_l=1;
309  Kvu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
310  }
311 
312  for (unsigned int i=0; i<n_v_dofs; i++)
313  for (unsigned int j=0; j<n_v_dofs; j++)
314  {
315  // Tensor indices
316  unsigned int C_i, C_j, C_k, C_l;
317  C_i=1, C_k=1;
318 
319  C_j=0, C_l=0;
320  Kvv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
321 
322  C_j=1, C_l=0;
323  Kvv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
324 
325  C_j=0, C_l=1;
326  Kvv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
327 
328  C_j=1, C_l=1;
329  Kvv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
330  }
331  }
332 
333  {
334  std::vector<boundary_id_type> bc_ids;
335  for (auto side : elem->side_index_range())
336  if (elem->neighbor_ptr(side) == nullptr)
337  {
338  mesh.get_boundary_info().boundary_ids (elem, side, bc_ids);
339 
340  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
341  const std::vector<Real> & JxW_face = fe_face->get_JxW();
342 
343  fe_face->reinit(elem, side);
344 
345  for (std::vector<boundary_id_type>::const_iterator b =
346  bc_ids.begin(); b != bc_ids.end(); ++b)
347  {
348  const boundary_id_type bc_id = *b;
349  for (unsigned int qp=0; qp<qface.n_points(); qp++)
350  {
351  // Add the loading
352  if (bc_id == 2)
353  for (unsigned int i=0; i<n_v_dofs; i++)
354  Fv(i) += JxW_face[qp] * (-1.) * phi_face[i][qp];
355 
356  // Add the constraint contributions
357  if (bc_id == 1)
358  {
359  for (unsigned int i=0; i<n_v_dofs; i++)
360  for (unsigned int j=0; j<n_lambda_dofs; j++)
361  Kv_lambda(i,j) += JxW_face[qp] * (-1.) * phi_face[i][qp];
362 
363  for (unsigned int i=0; i<n_lambda_dofs; i++)
364  for (unsigned int j=0; j<n_v_dofs; j++)
365  Klambda_v(i,j) += JxW_face[qp] * (-1.) * phi_face[j][qp];
366  }
367  }
368  }
369  }
370  }
371 
372  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
373 
374  system.matrix->add_matrix (Ke, dof_indices);
375  system.rhs->add_vector (Fe, dof_indices);
376  }
377 }
378 
380  unsigned int j,
381  unsigned int k,
382  unsigned int l)
383 {
384  // Define the Poisson ratio
385  const Real nu = 0.3;
386 
387  // Define the Lame constants (lambda_1 and lambda_2) based on Poisson ratio
388  const Real lambda_1 = nu / ((1. + nu) * (1. - 2.*nu));
389  const Real lambda_2 = 0.5 / (1 + nu);
390 
391  // Define the Kronecker delta functions that we need here
392  Real delta_ij = (i == j) ? 1. : 0.;
393  Real delta_il = (i == l) ? 1. : 0.;
394  Real delta_ik = (i == k) ? 1. : 0.;
395  Real delta_jl = (j == l) ? 1. : 0.;
396  Real delta_jk = (j == k) ? 1. : 0.;
397  Real delta_kl = (k == l) ? 1. : 0.;
398 
399  return lambda_1 * delta_ij * delta_kl + lambda_2 * (delta_ik * delta_jl + delta_il * delta_jk);
400 }
libMesh::BoundaryInfo::boundary_ids
std::vector< boundary_id_type > boundary_ids(const Node *node) const
Definition: boundary_info.C:985
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
libMesh::EquationSystems::add_system
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
Definition: equation_systems.C:345
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::DofMap::dof_indices
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:1967
libMesh::DofMap::add_dirichlet_boundary
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
Definition: dof_map_constraints.C:4390
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::DenseSubMatrix
Defines a dense submatrix for use in Finite Element-type computations.
Definition: dense_submatrix.h:45
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
assemble_elasticity
void assemble_elasticity(EquationSystems &es, const std::string &system_name)
libMesh::DenseMatrix< Number >
libMesh::LinearImplicitSystem::solve
virtual void solve() override
Assembles & solves the linear system A*x=b.
Definition: linear_implicit_system.C:108
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
libMesh::SECOND
Definition: enum_order.h:43
libMesh::boundary_id_type
int8_t boundary_id_type
Definition: id_types.h:51
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::LOCAL_VARIABLE_ORDER
Definition: dirichlet_boundaries.h:62
libMesh::MeshTools::Generation::build_square
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
Definition: mesh_generation.C:1501
libMesh::System::attach_assemble_function
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.
Definition: system.C:1755
libMesh::NumericVector::add_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].
Definition: numeric_vector.C:363
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libMesh::System::add_variable
unsigned int add_variable(const std::string &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.
Definition: system.C:1069
libMesh::EquationSystems::init
virtual void init()
Initialize all the systems.
Definition: equation_systems.C:96
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::DofMap::variable_type
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:1924
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::ZeroFunction
ConstFunction that simply returns 0.
Definition: zero_function.h:36
libMesh::EquationSystems::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Definition: equation_systems.C:1314
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::DofMap::constrain_element_matrix_and_vector
void constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix and vector.
Definition: dof_map.h:2034
libMesh::DenseSubVector< Number >
libMesh::MeshOutput::write_equation_systems
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 ...
Definition: mesh_output.C:31
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::SparseMatrix::add_matrix
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.
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
main
int main(int argc, char **argv)
Definition: systems_of_equations_ex5.C:80
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
libMesh::DirichletBoundary
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
Definition: dirichlet_boundaries.h:88
libMesh::System::variable_number
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1232
libMesh::DenseSubVector::reposition
void reposition(const unsigned int ioff, const unsigned int n)
Changes the location of the subvector in the parent vector.
Definition: dense_subvector.h:174
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEType::default_quadrature_order
Order default_quadrature_order() const
Definition: fe_type.h:353
libMesh::LinearImplicitSystem
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
Definition: linear_implicit_system.h:55
eval_elasticity_tensor
Real eval_elasticity_tensor(unsigned int i, unsigned int j, unsigned int k, unsigned int l)
Definition: systems_of_equations_ex5.C:379
libMesh::FIRST
Definition: enum_order.h:42
libMesh::EIGEN_SOLVERS
Definition: enum_solver_package.h:40
libMesh::DenseVector< Number >
libMesh::DenseSubMatrix::reposition
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.
Definition: dense_submatrix.h:173