ASPECT
elastic_tensor_decomposition.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2023 - 2024 by the authors of the ASPECT code.
3  This file is part of ASPECT.
4  ASPECT is free software; you can redistribute it and/or modify
5  it under the terms of the GNU General Public License as published by
6  the Free Software Foundation; either version 2, or (at your option)
7  any later version.
8  ASPECT is distributed in the hope that it will be useful,
9  but WITHOUT ANY WARRANTY; without even the implied warranty of
10  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11  GNU General Public License for more details.
12  You should have received a copy of the GNU General Public License
13  along with ASPECT; see the file LICENSE. If not see
14  <http://www.gnu.org/licenses/>.
15  */
16 
17 #ifndef _aspect_particle_property_elastic_tensor_decomposition_h
18 #define _aspect_particle_property_elastic_tensor_decomposition_h
19 
22 
23 namespace aspect
24 {
25  namespace Particle
26  {
27  namespace Property
28  {
29  namespace Utilities
30  {
35  std::array<unsigned int, 3>
36  indexed_even_permutation(const unsigned int index);
37 
44  SymmetricTensor<2,3>
45  compute_voigt_stiffness_tensor(const SymmetricTensor<2,6> &elastic_tensor);
46 
52  SymmetricTensor<2,3>
53  compute_dilatation_stiffness_tensor(const SymmetricTensor<2,6> &elastic_tensor);
54 
73  Tensor<2,3> compute_unpermutated_SCCS(const SymmetricTensor<2,3> &dilatation_stiffness_tensor,
74  const SymmetricTensor<2,3> &voigt_stiffness_tensor);
75 
87  std::array<std::array<double,3>,7>
89  const Tensor<2,3> &unpermutated_SCCS,
90  const SymmetricTensor<2,6> &elastic_matrix);
91 
92 
97  static const SymmetricTensor<2,21> projection_matrix_triclinic_to_monoclinic(
98  Tensor<2,21>(
99  {
100  {1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
101  {0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
102  {0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
103  {0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
104  {0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
105  {0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
106  {0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
107  {0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0},
108  {0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0},
109  {0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0},
110  {0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0},
111  {0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0},
112  {0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0},
113  {0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0},
114  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0},
115  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0},
116  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0},
117  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0},
118  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0},
119  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0},
120  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1},
121  })
122  );
123  static const SymmetricTensor<2,9> projection_matrix_monoclinic_to_orthorhombic(
124  Tensor<2,9>(
125  {
126  {1,0,0,0,0,0,0,0,0},
127  {0,1,0,0,0,0,0,0,0},
128  {0,0,1,0,0,0,0,0,0},
129  {0,0,0,1,0,0,0,0,0},
130  {0,0,0,0,1,0,0,0,0},
131  {0,0,0,0,0,1,0,0,0},
132  {0,0,0,0,0,0,1,0,0},
133  {0,0,0,0,0,0,0,1,0},
134  {0,0,0,0,0,0,0,0,1}
135  })
136  );
137  static const SymmetricTensor<2,9> projection_matrix_orthorhombic_to_tetragonal(
138  Tensor<2,9>(
139  {
140  {0.5,0.5,0.0,0.0,0.0,0.0,0.0,0.0,0.0},
141  {0.5,0.5,0.0,0.0,0.0,0.0,0.0,0.0,0.0},
142  {0.0,0.0,1.0,0.0,0.0,0.0,0.0,0.0,0.0},
143  {0.0,0.0,0.0,0.5,0.5,0.0,0.0,0.0,0.0},
144  {0.0,0.0,0.0,0.5,0.5,0.0,0.0,0.0,0.0},
145  {0.0,0.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0},
146  {0.0,0.0,0.0,0.0,0.0,0.0,0.5,0.5,0.0},
147  {0.0,0.0,0.0,0.0,0.0,0.0,0.5,0.5,0.0},
148  {0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,1.0}
149  })
150  );
151  static const SymmetricTensor<2,9> projection_matrix_tetragonal_to_hexagonal(
152  Tensor<2,9>(
153  {
154  {3./8. , 3./8. , 0.0, 0.0, 0.0, 1./(4.*sqrt(2.)) , 0.0, 0.0, 1./4. },
155  {3./8. , 3./8. , 0.0, 0.0, 0.0, 1./(4.*sqrt(2.)) , 0.0, 0.0, 1./4. },
156  {0.0 , 0.0 , 1.0, 0.0, 0.0, 0.0 , 0.0, 0.0, 0.0 },
157  {0.0 , 0.0 , 0.0, 0.5, 0.5, 0.0 , 0.0, 0.0, 0.0 },
158  {0.0 , 0.0 , 0.0, 0.5, 0.5, 0.0 , 0.0, 0.0, 0.0 },
159  {1./(4.*sqrt(2.)), 1./(4.*sqrt(2.)), 0.0, 0.0, 0.0, 3./4. , 0.0, 0.0, -1./(2.*sqrt(2.))},
160  {0.0 , 0.0 , 0.0, 0.0, 0.0, 0.0 , 0.5, 0.5, 0.0 },
161  {0.0 , 0.0 , 0.0, 0.0, 0.0, 0.0 , 0.5, 0.5, 0.0 },
162  {1./4. , 1./4. , 0.0, 0.0, 0.0, -1./(2.*sqrt(2.)) , 0.0, 0.0, 0.5 }
163  })
164  );
165  static const SymmetricTensor<2,9> projection_matrix_hexagonal_to_isotropic(
166  Tensor<2,9>(
167  {
168  {3./15. , 3./15. , 3./15. , sqrt(2.)/15. , sqrt(2.)/15. , sqrt(2.)/15. , 2./15. , 2./15. , 2./15. },
169  {3./15. , 3./15. , 3./15. , sqrt(2.)/15. , sqrt(2.)/15. , sqrt(2.)/15. , 2./15. , 2./15. , 2./15. },
170  {3./15. , 3./15. , 3./15. , sqrt(2.)/15. , sqrt(2.)/15. , sqrt(2.)/15. , 2./15. , 2./15. , 2./15. },
171  {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15. , 4./15. , 4./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. },
172  {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15. , 4./15. , 4./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. },
173  {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15. , 4./15. , 4./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. },
174  {2./15. , 2./15. , 2./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5. , 1./5. , 1./5. },
175  {2./15. , 2./15. , 2./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5. , 1./5. , 1./5. },
176  {2./15. , 2./15. , 2./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5. , 1./5. , 1./5. }
177  })
178  );
179  }
180 
187  template <int dim>
189  {
190  public:
194  ElasticTensorDecomposition() = default;
195 
200  void
201  initialize () override;
202 
214  void
215  initialize_one_particle_property (const Point<dim> &position,
216  std::vector<double> &particle_properties) const override;
217 
221  void
222  update_particle_properties (const ParticleUpdateInputs<dim> &inputs,
223  typename ParticleHandler<dim>::particle_iterator_range &particles) const override;
224 
230  need_update () const override;
231 
239  std::vector<std::pair<std::string, unsigned int>>
240  get_property_information() const override;
241 
242  private:
244  };
245  }
246  }
247 }
248 
249 #endif
std::array< unsigned int, 3 > indexed_even_permutation(const unsigned int index)
static const SymmetricTensor< 2, 9 > projection_matrix_tetragonal_to_hexagonal(Tensor< 2, 9 >({ {3./8., 3./8., 0.0, 0.0, 0.0, 1./(4.*sqrt(2.)), 0.0, 0.0, 1./4. }, {3./8., 3./8., 0.0, 0.0, 0.0, 1./(4.*sqrt(2.)), 0.0, 0.0, 1./4. }, {0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 }, {0.0, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0 }, {0.0, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0 }, {1./(4.*sqrt(2.)), 1./(4.*sqrt(2.)), 0.0, 0.0, 0.0, 3./4., 0.0, 0.0, -1./(2.*sqrt(2.))}, {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0 }, {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0 }, {1./4., 1./4., 0.0, 0.0, 0.0, -1./(2.*sqrt(2.)), 0.0, 0.0, 0.5 } }))
static const SymmetricTensor< 2, 9 > projection_matrix_monoclinic_to_orthorhombic(Tensor< 2, 9 >({ {1, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 1, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 1, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 1, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 1, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 1, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 1, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 1, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 1} }))
Tensor< 2, 3 > compute_unpermutated_SCCS(const SymmetricTensor< 2, 3 > &dilatation_stiffness_tensor, const SymmetricTensor< 2, 3 > &voigt_stiffness_tensor)
static const SymmetricTensor< 2, 21 > projection_matrix_triclinic_to_monoclinic(Tensor< 2, 21 >({ {1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0}, {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1}, }))
std::array< std::array< double, 3 >, 7 > compute_elastic_tensor_SCCS_decompositions(const Tensor< 2, 3 > &unpermutated_SCCS, const SymmetricTensor< 2, 6 > &elastic_matrix)
Definition: compat.h:59
SymmetricTensor< 2, 3 > compute_voigt_stiffness_tensor(const SymmetricTensor< 2, 6 > &elastic_tensor)
static const SymmetricTensor< 2, 9 > projection_matrix_hexagonal_to_isotropic(Tensor< 2, 9 >({ {3./15., 3./15., 3./15., sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 2./15., 2./15., 2./15. }, {3./15., 3./15., 3./15., sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 2./15., 2./15., 2./15. }, {3./15., 3./15., 3./15., sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 2./15., 2./15., 2./15. }, {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15., 4./15., 4./15., -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. }, {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15., 4./15., 4./15., -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. }, {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15., 4./15., 4./15., -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. }, {2./15., 2./15., 2./15., -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5., 1./5., 1./5. }, {2./15., 2./15., 2./15., -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5., 1./5., 1./5. }, {2./15., 2./15., 2./15., -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5., 1./5., 1./5. } }))
static const SymmetricTensor< 2, 9 > projection_matrix_orthorhombic_to_tetragonal(Tensor< 2, 9 >({ {0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, {0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, {0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, {0.0, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, {0.0, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, {0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0} }))
SymmetricTensor< 2, 3 > compute_dilatation_stiffness_tensor(const SymmetricTensor< 2, 6 > &elastic_tensor)