libMesh
Public Member Functions | Private Attributes | List of all members
LinearElasticity Class Reference
Inheritance diagram for LinearElasticity:
[legend]

Public Member Functions

 LinearElasticity (EquationSystems &es_in)
 
Real kronecker_delta (unsigned int i, unsigned int j)
 Kronecker delta function. More...
 
Real elasticity_tensor (unsigned int i, unsigned int j, unsigned int k, unsigned int l)
 Evaluate the fourth order tensor (C_ijkl) that relates stress to strain. More...
 
void assemble ()
 Assemble the system matrix and right-hand side vector. More...
 
void compute_stresses ()
 
 LinearElasticity (EquationSystems &es_in)
 
Real kronecker_delta (unsigned int i, unsigned int j)
 Kronecker delta function. More...
 
Real elasticity_tensor (unsigned int i, unsigned int j, unsigned int k, unsigned int l)
 Evaluate the fourth order tensor (C_ijkl) that relates stress to strain. More...
 
void assemble ()
 Assemble the system matrix and right-hand side vector. More...
 

Private Attributes

EquationSystemses
 

Detailed Description

Definition at line 114 of file systems_of_equations_ex6.C.

Constructor & Destructor Documentation

◆ LinearElasticity() [1/2]

LinearElasticity::LinearElasticity ( EquationSystems es_in)
inline

Definition at line 121 of file systems_of_equations_ex6.C.

121  :
122  es(es_in)
123  {}

◆ LinearElasticity() [2/2]

LinearElasticity::LinearElasticity ( EquationSystems es_in)
inline

Definition at line 206 of file systems_of_equations_ex9.C.

206  :
207  es(es_in)
208  {}

Member Function Documentation

◆ assemble() [1/2]

void LinearElasticity::assemble ( )
inlinevirtual

Assemble the system matrix and right-hand side vector.

Implements libMesh::System::Assembly.

Definition at line 157 of file systems_of_equations_ex6.C.

158  {
159  const MeshBase & mesh = es.get_mesh();
160 
161  const unsigned int dim = mesh.mesh_dimension();
162 
163  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Elasticity");
164 
165  const unsigned int u_var = system.variable_number ("u");
166 
167  const DofMap & dof_map = system.get_dof_map();
168  FEType fe_type = dof_map.variable_type(u_var);
169  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
170  QGauss qrule (dim, fe_type.default_quadrature_order());
171  fe->attach_quadrature_rule (&qrule);
172 
173  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
174  QGauss qface(dim-1, fe_type.default_quadrature_order());
175  fe_face->attach_quadrature_rule (&qface);
176 
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();
180 
182  DenseSubMatrix<Number> Ke_var[3][3] =
183  {
187  };
188 
190 
191  DenseSubVector<Number> Fe_var[3] =
195 
196  std::vector<dof_id_type> dof_indices;
197  std::vector<std::vector<dof_id_type>> dof_indices_var(3);
198 
199  for (const auto & elem : mesh.active_local_element_ptr_range())
200  {
201  dof_map.dof_indices (elem, dof_indices);
202  for (unsigned int var=0; var<3; var++)
203  dof_map.dof_indices (elem, dof_indices_var[var], var);
204 
205  const unsigned int n_dofs = dof_indices.size();
206  const unsigned int n_var_dofs = dof_indices_var[0].size();
207 
208  fe->reinit (elem);
209 
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);
214 
215  Fe.resize (n_dofs);
216  for (unsigned int var=0; var<3; var++)
217  Fe_var[var].reposition (var*n_var_dofs, n_var_dofs);
218 
219  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
220  {
221  // assemble \int_Omega C_ijkl u_k,l v_i,j \dx
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) +=
229  JxW[qp] * elasticity_tensor(i,j,k,l) * dphi[dof_j][qp](l) * dphi[dof_i][qp](j);
230 
231  // assemble \int_Omega f_i v_i \dx
232  DenseVector<Number> f_vec(3);
233  f_vec(0) = 0.;
234  f_vec(1) = 0.;
235  f_vec(2) = -1.;
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]);
239  }
240 
241  // assemble \int_\Gamma g_i v_i \ds
242  DenseVector<Number> g_vec(3);
243  g_vec(0) = 0.;
244  g_vec(1) = 0.;
245  g_vec(2) = -1.;
246  {
247  for (auto side : elem->side_index_range())
248  if (elem->neighbor_ptr(side) == nullptr)
249  {
250  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
251  const std::vector<Real> & JxW_face = fe_face->get_JxW();
252 
253  fe_face->reinit(elem, side);
254 
255  // Apply a traction
256  for (unsigned int qp=0; qp<qface.n_points(); qp++)
257  if (mesh.get_boundary_info().has_boundary_id(elem, side, BOUNDARY_ID_MAX_X))
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]);
261  }
262  }
263 
264  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
265 
266  system.matrix->add_matrix (Ke, dof_indices);
267  system.rhs->add_vector (Fe, dof_indices);
268  }
269  }

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), libMesh::FEType::default_quadrature_order(), dim, libMesh::DofMap::dof_indices(), elasticity_tensor(), libMesh::MeshBase::get_boundary_info(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::BoundaryInfo::has_boundary_id(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, libMesh::System::variable_number(), and libMesh::DofMap::variable_type().

◆ assemble() [2/2]

void LinearElasticity::assemble ( )
inlinevirtual

Assemble the system matrix and right-hand side vector.

Implements libMesh::System::Assembly.

Definition at line 242 of file systems_of_equations_ex9.C.

243  {
244  const MeshBase & mesh = es.get_mesh();
245 
246  const unsigned int dim = mesh.mesh_dimension();
247 
248  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Elasticity");
249 
250  const unsigned int u_var = system.variable_number ("u");
251 
252  const DofMap & dof_map = system.get_dof_map();
253  FEType fe_type = dof_map.variable_type(u_var);
254  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
255  QGauss qrule (dim, fe_type.default_quadrature_order());
256  fe->attach_quadrature_rule (&qrule);
257 
258  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
259  QGauss qface(dim-1, fe_type.default_quadrature_order());
260  fe_face->attach_quadrature_rule (&qface);
261 
262  const std::vector<Real> & JxW = fe->get_JxW();
263  const std::vector<std::vector<Real>> & phi = fe->get_phi();
264  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
265 
267  DenseSubMatrix<Number> Ke_var[3][3] =
268  {
272  };
273 
275 
276  DenseSubVector<Number> Fe_var[3] =
280 
281  std::vector<dof_id_type> dof_indices;
282  std::vector<std::vector<dof_id_type>> dof_indices_var(3);
283 
284  for (const auto & elem : mesh.active_local_element_ptr_range())
285  {
286  dof_map.dof_indices (elem, dof_indices);
287  for (unsigned int var=0; var<3; var++)
288  dof_map.dof_indices (elem, dof_indices_var[var], var);
289 
290  const unsigned int n_dofs = dof_indices.size();
291  const unsigned int n_var_dofs = dof_indices_var[0].size();
292 
293  fe->reinit (elem);
294 
295  Ke.resize (n_dofs, n_dofs);
296  for (unsigned int var_i=0; var_i<3; var_i++)
297  for (unsigned int var_j=0; var_j<3; var_j++)
298  Ke_var[var_i][var_j].reposition (var_i*n_var_dofs, var_j*n_var_dofs, n_var_dofs, n_var_dofs);
299 
300  Fe.resize (n_dofs);
301  for (unsigned int var=0; var<3; var++)
302  Fe_var[var].reposition (var*n_var_dofs, n_var_dofs);
303 
304  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
305  {
306  // assemble \int_Omega C_ijkl u_k,l v_i,j \dx
307  for (unsigned int dof_i=0; dof_i<n_var_dofs; dof_i++)
308  for (unsigned int dof_j=0; dof_j<n_var_dofs; dof_j++)
309  for (unsigned int i=0; i<3; i++)
310  for (unsigned int j=0; j<3; j++)
311  for (unsigned int k=0; k<3; k++)
312  for (unsigned int l=0; l<3; l++)
313  Ke_var[i][k](dof_i,dof_j) +=
314  JxW[qp] * elasticity_tensor(i,j,k,l) * dphi[dof_j][qp](l) * dphi[dof_i][qp](j);
315 
316  // assemble \int_Omega f_i v_i \dx
317  DenseVector<Number> f_vec(3);
318  if(elem->subdomain_id() == 101)
319  {
320  f_vec(0) = 1.;
321  f_vec(1) = 1.;
322  f_vec(2) = 0.;
323  }
324  else if(elem->subdomain_id() == 1)
325  {
326  f_vec(0) = 0.36603;
327  f_vec(1) = 1.36603;
328  f_vec(2) = 0.;
329  }
330  for (unsigned int dof_i=0; dof_i<n_var_dofs; dof_i++)
331  for (unsigned int i=0; i<3; i++)
332  Fe_var[i](dof_i) += JxW[qp] * (f_vec(i) * phi[dof_i][qp]);
333  }
334 
335  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
336 
337  system.matrix->add_matrix (Ke, dof_indices);
338  system.rhs->add_vector (Fe, dof_indices);
339  }
340  }

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), libMesh::FEType::default_quadrature_order(), dim, libMesh::DofMap::dof_indices(), elasticity_tensor(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, libMesh::System::variable_number(), and libMesh::DofMap::variable_type().

◆ compute_stresses()

void LinearElasticity::compute_stresses ( )
inline

Definition at line 272 of file systems_of_equations_ex6.C.

273  {
274  const MeshBase & mesh = es.get_mesh();
275  const unsigned int dim = mesh.mesh_dimension();
276 
277  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Elasticity");
278 
279  unsigned int displacement_vars[3];
280  displacement_vars[0] = system.variable_number ("u");
281  displacement_vars[1] = system.variable_number ("v");
282  displacement_vars[2] = system.variable_number ("w");
283  const unsigned int u_var = system.variable_number ("u");
284 
285  const DofMap & dof_map = system.get_dof_map();
286  FEType fe_type = dof_map.variable_type(u_var);
287  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
288  QGauss qrule (dim, fe_type.default_quadrature_order());
289  fe->attach_quadrature_rule (&qrule);
290 
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();
294 
295  // Also, get a reference to the ExplicitSystem
296  ExplicitSystem & stress_system = es.get_system<ExplicitSystem>("StressSystem");
297  const DofMap & stress_dof_map = stress_system.get_dof_map();
298  unsigned int sigma_vars[6];
299  sigma_vars[0] = stress_system.variable_number ("sigma_00");
300  sigma_vars[1] = stress_system.variable_number ("sigma_01");
301  sigma_vars[2] = stress_system.variable_number ("sigma_02");
302  sigma_vars[3] = stress_system.variable_number ("sigma_11");
303  sigma_vars[4] = stress_system.variable_number ("sigma_12");
304  sigma_vars[5] = stress_system.variable_number ("sigma_22");
305  unsigned int vonMises_var = stress_system.variable_number ("vonMises");
306 
307  // Storage for the stress dof indices on each element
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;
311 
312  for (const auto & elem : mesh.active_local_element_ptr_range())
313  {
314  for (unsigned int var=0; var<3; var++)
315  dof_map.dof_indices (elem, dof_indices_var[var], displacement_vars[var]);
316 
317  const unsigned int n_var_dofs = dof_indices_var[0].size();
318 
319  fe->reinit (elem);
320 
321  std::vector<DenseMatrix<Number>> stress_tensor_qp(qrule.n_points());
322  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
323  {
324  stress_tensor_qp[qp].resize(3,3);
325 
326  // Row is variable u1, u2, or u3, column is x, y, or z
327  DenseMatrix<Number> grad_u(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]);
332 
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);
338  }
339 
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);
344 
345  // Below we project each component of the stress tensor onto a L2_LAGRANGE discretization.
346  // Note that this gives a discontinuous stress plot on element boundaries, which is
347  // appropriate. We then also get the von Mises stress from the projected stress tensor.
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++)
351  {
352  stress_dof_map.dof_indices (elem, stress_dof_indices_var, sigma_vars[stress_var_index]);
353 
354  const unsigned int n_proj_dofs = stress_dof_indices_var.size();
355 
356  DenseMatrix<Real> Me(n_proj_dofs, n_proj_dofs);
357  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
358  {
359  for(unsigned int i=0; i<n_proj_dofs; i++)
360  for(unsigned int j=0; j<n_proj_dofs; j++)
361  {
362  Me(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp]);
363  }
364  }
365 
366  DenseVector<Number> Fe(n_proj_dofs);
367  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
368  for(unsigned int i=0; i<n_proj_dofs; i++)
369  {
370  Fe(i) += JxW[qp] * stress_tensor_qp[qp](var_i,var_j) * phi[i][qp];
371  }
372 
373  DenseVector<Number> projected_data;
374  Me.cholesky_solve(Fe, projected_data);
375 
376  for(unsigned int index=0; index<n_proj_dofs; index++)
377  {
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));
382 
383  elem_sigma_vec[index](var_i,var_j) = projected_data(index);
384  }
385 
386  stress_var_index++;
387  }
388 
389  for (std::size_t index=0; index<elem_sigma_vec.size(); index++)
390  {
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);
394 
395  // Get the von Mises stress from the projected stress tensor
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.))));
402 
403  dof_id_type dof_index = vonmises_dof_indices_var[index];
404 
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);
408  }
409  }
410 
411  // Should call close and update when we set vector entries directly
412  stress_system.solution->close();
413  stress_system.update();
414  }

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DenseMatrix< T >::cholesky_solve(), libMesh::System::current_solution(), libMesh::FEType::default_quadrature_order(), dim, libMesh::DofMap::dof_indices(), elasticity_tensor(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::System::n_vars(), std::pow(), libMesh::System::solution, std::sqrt(), libMesh::System::update(), libMesh::System::variable_number(), and libMesh::DofMap::variable_type().

Referenced by main().

◆ elasticity_tensor() [1/2]

Real LinearElasticity::elasticity_tensor ( unsigned int  i,
unsigned int  j,
unsigned int  k,
unsigned int  l 
)
inline

Evaluate the fourth order tensor (C_ijkl) that relates stress to strain.

Definition at line 137 of file systems_of_equations_ex6.C.

141  {
142  // Hard code material parameters for the sake of simplicity
143  const Real poisson_ratio = 0.3;
144  const Real young_modulus = 1.;
145 
146  // Define the Lame constants
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));
149 
150  return lambda_1 * kronecker_delta(i, j) * kronecker_delta(k, l) +
151  lambda_2 * (kronecker_delta(i, k) * kronecker_delta(j, l) + kronecker_delta(i, l) * kronecker_delta(j, k));
152  }

References kronecker_delta(), and libMesh::Real.

◆ elasticity_tensor() [2/2]

Real LinearElasticity::elasticity_tensor ( unsigned int  i,
unsigned int  j,
unsigned int  k,
unsigned int  l 
)
inline

Evaluate the fourth order tensor (C_ijkl) that relates stress to strain.

Definition at line 222 of file systems_of_equations_ex9.C.

226  {
227  // Hard code material parameters for the sake of simplicity
228  const Real poisson_ratio = 0.3;
229  const Real young_modulus = 1.;
230 
231  // Define the Lame constants
232  const Real lambda_1 = (young_modulus*poisson_ratio)/((1.+poisson_ratio)*(1.-2.*poisson_ratio));
233  const Real lambda_2 = young_modulus/(2.*(1.+poisson_ratio));
234 
235  return lambda_1 * kronecker_delta(i, j) * kronecker_delta(k, l) +
236  lambda_2 * (kronecker_delta(i, k) * kronecker_delta(j, l) + kronecker_delta(i, l) * kronecker_delta(j, k));
237  }

References kronecker_delta(), and libMesh::Real.

◆ kronecker_delta() [1/2]

Real LinearElasticity::kronecker_delta ( unsigned int  i,
unsigned int  j 
)
inline

Kronecker delta function.

Definition at line 128 of file systems_of_equations_ex6.C.

130  {
131  return i == j ? 1. : 0.;
132  }

◆ kronecker_delta() [2/2]

Real LinearElasticity::kronecker_delta ( unsigned int  i,
unsigned int  j 
)
inline

Kronecker delta function.

Definition at line 213 of file systems_of_equations_ex9.C.

215  {
216  return i == j ? 1. : 0.;
217  }

Member Data Documentation

◆ es

EquationSystems & LinearElasticity::es
private

Definition at line 117 of file systems_of_equations_ex6.C.


The documentation for this class was generated from the following files:
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::System::n_vars
unsigned int n_vars() const
Definition: system.h:2155
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::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::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
std::sqrt
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::DenseMatrix< Number >
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
LinearElasticity::elasticity_tensor
Real elasticity_tensor(unsigned int i, unsigned int j, unsigned int k, unsigned int l)
Evaluate the fourth order tensor (C_ijkl) that relates stress to strain.
Definition: systems_of_equations_ex6.C:137
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::System::current_solution
Number current_solution(const dof_id_type global_dof_number) const
Definition: system.C:194
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::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
std::pow
double pow(double a, int b)
Definition: libmesh_augment_std_namespace.h:58
libMesh::DofMap::variable_type
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:1924
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::ExplicitSystem
Manages consistently variables, degrees of freedom, and coefficient vectors for explicit systems.
Definition: explicit_system.h:48
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
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
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::System::variable_number
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1232
LinearElasticity::es
EquationSystems & es
Definition: systems_of_equations_ex6.C:117
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
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
libMesh::BoundaryInfo::has_boundary_id
bool has_boundary_id(const Node *const node, const boundary_id_type id) const
Definition: boundary_info.C:972
libMesh::System::update
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:408
LinearElasticity::kronecker_delta
Real kronecker_delta(unsigned int i, unsigned int j)
Kronecker delta function.
Definition: systems_of_equations_ex6.C:128
libMesh::DenseVector< Number >