PRISMS-PF Manual
Loading...
Searching...
No Matches
mechanics.h
Go to the documentation of this file.
1#pragma once
2
3#include <deal.II/base/tensor.h>
4
5#include <prismspf/config.h>
6
8
9namespace Mechanics
10{
11
15 template <unsigned int dim>
16 constexpr unsigned int voigt_tensor_size = (2 * dim) - 1 + (dim / 3);
17
22 template <unsigned int dim, typename T>
23 inline DEAL_II_ALWAYS_INLINE void
24 compute_stress(const dealii::Tensor<2, voigt_tensor_size<dim>, T> &elasticity_tensor,
25 const dealii::Tensor<1, voigt_tensor_size<dim>, T> &strain,
26 dealii::Tensor<1, voigt_tensor_size<dim>, T> &stress)
27 {
28 stress = elasticity_tensor * strain;
29 }
30
36 template <unsigned int dim, typename T>
37 inline DEAL_II_ALWAYS_INLINE void
38 compute_stress(const dealii::Tensor<2, voigt_tensor_size<dim>, T> &elasticity_tensor,
39 const dealii::Tensor<2, dim, T> &strain,
40 dealii::Tensor<2, dim, T> &stress)
41 {
42 dealii::Tensor<1, voigt_tensor_size<dim>, T> sigma;
43 dealii::Tensor<1, voigt_tensor_size<dim>, T> epsilon;
44
45 if constexpr (dim == 3)
46 {
47 const int xx_dir = 0;
48 const int yy_dir = 1;
49 const int zz_dir = 2;
50 const int yz_dir = 3;
51 const int xz_dir = 4;
52 const int xy_dir = 5;
53
54 epsilon[xx_dir] = strain[xx_dir][xx_dir];
55 epsilon[yy_dir] = strain[yy_dir][yy_dir];
56 epsilon[zz_dir] = strain[zz_dir][zz_dir];
57
58 // In Voigt notation: epsilon are engineering shear strains
59 epsilon[yz_dir] = strain[yy_dir][zz_dir] + strain[zz_dir][yy_dir];
60 epsilon[xz_dir] = strain[xx_dir][zz_dir] + strain[zz_dir][xx_dir];
61 epsilon[xy_dir] = strain[xx_dir][yy_dir] + strain[yy_dir][xx_dir];
62
63 // Multiply elasticity_tensor and epsilon to get sigma
64 sigma = elasticity_tensor * epsilon;
65
66 stress[xx_dir][xx_dir] = sigma[xx_dir];
67 stress[yy_dir][yy_dir] = sigma[yy_dir];
68 stress[zz_dir][zz_dir] = sigma[zz_dir];
69
70 stress[yy_dir][zz_dir] = sigma[yz_dir];
71 stress[zz_dir][yy_dir] = sigma[yz_dir];
72
73 stress[xx_dir][zz_dir] = sigma[xz_dir];
74 stress[zz_dir][xx_dir] = sigma[xz_dir];
75
76 stress[xx_dir][yy_dir] = sigma[xy_dir];
77 stress[yy_dir][xx_dir] = sigma[xy_dir];
78 }
79 else if constexpr (dim == 2)
80 {
81 const int xx_dir = 0;
82 const int yy_dir = 1;
83 const int xy_dir = 2;
84
85 epsilon[xx_dir] = strain[xx_dir][xx_dir];
86 epsilon[yy_dir] = strain[yy_dir][yy_dir];
87
88 // In Voigt notation: epsilon are engineering shear strains
89 epsilon[xy_dir] = strain[xx_dir][yy_dir] + strain[yy_dir][xx_dir];
90
91 // Multiply elasticity_tensor and epsilon to get sigma
92 sigma = elasticity_tensor * epsilon;
93
94 stress[xx_dir][xx_dir] = sigma[xx_dir];
95 stress[yy_dir][yy_dir] = sigma[yy_dir];
96 stress[xx_dir][yy_dir] = sigma[xy_dir];
97 stress[yy_dir][xx_dir] = sigma[xy_dir];
98 }
99 else
100 {
101 const int xx_dir = 0;
102
103 stress[xx_dir][xx_dir] =
104 elasticity_tensor[xx_dir][xx_dir] * strain[xx_dir][xx_dir];
105 }
106 }
107
108} // namespace Mechanics
109
110PRISMS_PF_END_NAMESPACE
Definition mechanics.h:10
constexpr unsigned int voigt_tensor_size
Voigt notation index range.
Definition mechanics.h:16
DEAL_II_ALWAYS_INLINE void compute_stress(const dealii::Tensor< 2, voigt_tensor_size< dim >, T > &elasticity_tensor, const dealii::Tensor< 1, voigt_tensor_size< dim >, T > &strain, dealii::Tensor< 1, voigt_tensor_size< dim >, T > &stress)
Compute the stress with a given displacement and elasticity tensor. This assumes that the provided pa...
Definition mechanics.h:24
Definition conditional_ostreams.cc:20