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  {
44  using namespace aspect::Utilities;
45 
50  std::array<unsigned int, 3>
51  indexed_even_permutation(const unsigned int index);
52 
59  SymmetricTensor<2,3>
60  compute_voigt_stiffness_tensor(const SymmetricTensor<2,6> &elastic_tensor);
61 
67  SymmetricTensor<2,3>
68  compute_dilatation_stiffness_tensor(const SymmetricTensor<2,6> &elastic_tensor);
69 
88  Tensor<2,3> compute_unpermutated_SCCS(const SymmetricTensor<2,3> &dilatation_stiffness_tensor,
89  const SymmetricTensor<2,3> &voigt_stiffness_tensor);
90 
102  std::array<std::array<double,3>,7>
104  const Tensor<2,3> &unpermutated_SCCS,
105  const SymmetricTensor<2,6> &elastic_matrix);
106 
107 
112  static const SymmetricTensor<2,21> projection_matrix_triclinic_to_monoclinic(
113  Tensor<2,21>(
114  {
115  {1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
116  {0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
117  {0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
118  {0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
119  {0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
120  {0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
121  {0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0,0},
122  {0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0,0},
123  {0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0,0},
124  {0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0,0},
125  {0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0,0},
126  {0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0,0},
127  {0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0,0},
128  {0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0,0},
129  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0,0},
130  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0,0},
131  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0,0},
132  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0,0},
133  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0,0},
134  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1,0},
135  {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,1},
136  })
137  );
138  static const SymmetricTensor<2,9> projection_matrix_monoclinic_to_orthorhombic(
139  Tensor<2,9>(
140  {
141  {1,0,0,0,0,0,0,0,0},
142  {0,1,0,0,0,0,0,0,0},
143  {0,0,1,0,0,0,0,0,0},
144  {0,0,0,1,0,0,0,0,0},
145  {0,0,0,0,1,0,0,0,0},
146  {0,0,0,0,0,1,0,0,0},
147  {0,0,0,0,0,0,1,0,0},
148  {0,0,0,0,0,0,0,1,0},
149  {0,0,0,0,0,0,0,0,1}
150  })
151  );
152  static const SymmetricTensor<2,9> projection_matrix_orthorhombic_to_tetragonal(
153  Tensor<2,9>(
154  {
155  {0.5,0.5,0.0,0.0,0.0,0.0,0.0,0.0,0.0},
156  {0.5,0.5,0.0,0.0,0.0,0.0,0.0,0.0,0.0},
157  {0.0,0.0,1.0,0.0,0.0,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  {0.0,0.0,0.0,0.5,0.5,0.0,0.0,0.0,0.0},
160  {0.0,0.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0},
161  {0.0,0.0,0.0,0.0,0.0,0.0,0.5,0.5,0.0},
162  {0.0,0.0,0.0,0.0,0.0,0.0,0.5,0.5,0.0},
163  {0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,1.0}
164  })
165  );
166  static const SymmetricTensor<2,9> projection_matrix_tetragonal_to_hexagonal(
167  Tensor<2,9>(
168  {
169  {3./8. , 3./8. , 0.0, 0.0, 0.0, 1./(4.*sqrt(2.)) , 0.0, 0.0, 1./4. },
170  {3./8. , 3./8. , 0.0, 0.0, 0.0, 1./(4.*sqrt(2.)) , 0.0, 0.0, 1./4. },
171  {0.0 , 0.0 , 1.0, 0.0, 0.0, 0.0 , 0.0, 0.0, 0.0 },
172  {0.0 , 0.0 , 0.0, 0.5, 0.5, 0.0 , 0.0, 0.0, 0.0 },
173  {0.0 , 0.0 , 0.0, 0.5, 0.5, 0.0 , 0.0, 0.0, 0.0 },
174  {1./(4.*sqrt(2.)), 1./(4.*sqrt(2.)), 0.0, 0.0, 0.0, 3./4. , 0.0, 0.0, -1./(2.*sqrt(2.))},
175  {0.0 , 0.0 , 0.0, 0.0, 0.0, 0.0 , 0.5, 0.5, 0.0 },
176  {0.0 , 0.0 , 0.0, 0.0, 0.0, 0.0 , 0.5, 0.5, 0.0 },
177  {1./4. , 1./4. , 0.0, 0.0, 0.0, -1./(2.*sqrt(2.)) , 0.0, 0.0, 0.5 }
178  })
179  );
180  static const SymmetricTensor<2,9> projection_matrix_hexagonal_to_isotropic(
181  Tensor<2,9>(
182  {
183  {3./15. , 3./15. , 3./15. , sqrt(2.)/15. , sqrt(2.)/15. , sqrt(2.)/15. , 2./15. , 2./15. , 2./15. },
184  {3./15. , 3./15. , 3./15. , sqrt(2.)/15. , sqrt(2.)/15. , sqrt(2.)/15. , 2./15. , 2./15. , 2./15. },
185  {3./15. , 3./15. , 3./15. , sqrt(2.)/15. , sqrt(2.)/15. , sqrt(2.)/15. , 2./15. , 2./15. , 2./15. },
186  {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15. , 4./15. , 4./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. },
187  {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15. , 4./15. , 4./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. },
188  {sqrt(2.)/15., sqrt(2.)/15., sqrt(2.)/15., 4./15. , 4./15. , 4./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15. },
189  {2./15. , 2./15. , 2./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5. , 1./5. , 1./5. },
190  {2./15. , 2./15. , 2./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5. , 1./5. , 1./5. },
191  {2./15. , 2./15. , 2./15. , -sqrt(2.)/15., -sqrt(2.)/15., -sqrt(2.)/15., 1./5. , 1./5. , 1./5. }
192  })
193  );
194  }
195 
202  template <int dim>
204  {
205  public:
209  ElasticTensorDecomposition() = default;
210 
215  void
216  initialize () override;
217 
229  void
230  initialize_one_particle_property (const Point<dim> &position,
231  std::vector<double> &particle_properties) const override;
232 
236  void
237  update_particle_properties (const ParticleUpdateInputs<dim> &inputs,
238  typename ParticleHandler<dim>::particle_iterator_range &particles) const override;
239 
245  need_update () const override;
246 
254  std::vector<std::pair<std::string, unsigned int>>
255  get_property_information() const override;
256 
257  private:
259  };
260  }
261  }
262 }
263 
264 #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)