ASPECT
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
crystal_preferred_orientation.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2022 - 2024 by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19  */
20 
21 #ifndef _aspect_particle_property_cpo_h
22 #define _aspect_particle_property_cpo_h
23 
26 #include <array>
27 
28 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
29 #include <boost/random.hpp>
30 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
31 
32 namespace aspect
33 {
34  namespace Particle
35  {
36  namespace Property
37  {
46  enum class DeformationType
47  {
49  };
50 
51 
66  {
68  };
69 
73  enum class AdvectionMethod
74  {
76  };
77 
85  {
87  };
88 
96  {
98  };
99 
129  template <int dim>
131  {
132  public:
136  CrystalPreferredOrientation() = default;
137 
142  void
143  initialize () override;
144 
156  void
157  initialize_one_particle_property (const Point<dim> &position,
158  std::vector<double> &particle_properties) const override;
159 
163  void
164  update_particle_properties (const ParticleUpdateInputs<dim> &inputs,
165  typename ParticleHandler<dim>::particle_iterator_range &particles) const override;
166 
172  need_update () const override;
173 
178  late_initialization_mode () const override;
179 
183  UpdateFlags
184  get_update_flags (const unsigned int component) const override;
185 
193  std::vector<std::pair<std::string, unsigned int>>
194  get_property_information() const override;
195 
214  std::pair<std::vector<double>, std::vector<Tensor<2,3>>>
215  compute_derivatives(const unsigned int cpo_index,
216  const ArrayView<double> &data,
217  const unsigned int mineral_i,
218  const SymmetricTensor<2,3> &strain_rate_3d,
219  const Tensor<2,3> &velocity_gradient_tensor,
220  const Point<dim> &position,
221  const typename DoFHandler<dim>::active_cell_iterator &cell,
222  const double temperature,
223  const double pressure,
224  const Tensor<1,dim> &velocity,
225  const std::vector<double> &compositions,
226  const SymmetricTensor<2,dim> &strain_rate,
227  const SymmetricTensor<2,dim> &deviatoric_strain_rate,
228  const double water_content) const;
229 
246  std::pair<std::vector<double>, std::vector<Tensor<2,3>>>
247  compute_derivatives_drex_2004(const DeformationType deformation_type,
248  const unsigned int cpo_index,
249  const ArrayView<double> &data,
250  const unsigned int mineral_i,
251  const SymmetricTensor<2,3> &strain_rate_3d,
252  const Tensor<2,3> &velocity_gradient_tensor,
253  const std::array<double,4> ref_resolved_shear_stress,
254  const bool prevent_nondimensionalization = false) const;
255 
256 
260  static
261  void
262  declare_parameters (ParameterHandler &prm);
263 
267  void
268  parse_parameters (ParameterHandler &prm) override;
269 
273  unsigned int
274  get_number_of_grains() const;
275 
279  unsigned int
280  get_number_of_minerals() const;
281 
292  determine_deformation_type(const DeformationTypeSelector deformation_type_selector,
293  const Point<dim> &position,
294  const typename DoFHandler<dim>::active_cell_iterator &cell,
295  const double temperature,
296  const double pressure,
297  const Tensor<1,dim> &velocity,
298  const std::vector<double> &compositions,
299  const SymmetricTensor<2,dim> &strain_rate,
300  const SymmetricTensor<2,dim> &deviatoric_strain_rate,
301  const double water_content) const;
302 
308  determine_deformation_type_karato_2008(const double stress,
309  const double water_content) const;
310 
320  std::array<double,4>
321  reference_resolved_shear_stress_from_deformation_type(DeformationType deformation_type,
322  double max_value = 1e60) const;
323 
331  inline
332  DeformationType get_deformation_type(const unsigned int cpo_data_position,
333  const ArrayView<double> &data,
334  const unsigned int mineral_i) const
335  {
336  return static_cast<DeformationType>(data[cpo_data_position + 0 + mineral_i * (n_grains * 10 + 2)]);
337  }
338 
347  inline
348  void set_deformation_type(const unsigned int cpo_data_position,
349  const ArrayView<double> &data,
350  const unsigned int mineral_i,
351  const DeformationType deformation_type) const
352  {
353  data[cpo_data_position + 0 + mineral_i * (n_grains * 10 + 2)] = static_cast<double>(deformation_type);
354  }
355 
363  inline
364  double get_volume_fraction_mineral(const unsigned int cpo_data_position,
365  const ArrayView<double> &data,
366  const unsigned int mineral_i) const
367  {
368  return data[cpo_data_position + 1 + mineral_i *(n_grains * 10 + 2)];
369  }
370 
379  inline
380  void set_volume_fraction_mineral(const unsigned int cpo_data_position,
381  const ArrayView<double> &data,
382  const unsigned int mineral_i,
383  const double volume_fraction_mineral) const
384  {
385  data[cpo_data_position + 1 + mineral_i *(n_grains * 10 + 2)] = volume_fraction_mineral;
386  }
387 
396  inline
397  double get_volume_fractions_grains(const unsigned int cpo_data_position,
398  const ArrayView<const double> &data,
399  const unsigned int mineral_i,
400  const unsigned int grain_i) const
401  {
402  return data[cpo_data_position + 2 + grain_i * 10 + mineral_i * (n_grains * 10 + 2)];
403  }
404 
414  inline
415  void set_volume_fractions_grains(const unsigned int cpo_data_position,
416  const ArrayView<double> &data,
417  const unsigned int mineral_i,
418  const unsigned int grain_i,
419  const double volume_fractions_grains) const
420  {
421  data[cpo_data_position + 2 + grain_i * 10 + mineral_i * (n_grains * 10 + 2)] = volume_fractions_grains;
422  }
423 
433  inline
434  Tensor<2,3> get_rotation_matrix_grains(const unsigned int cpo_data_position,
435  const ArrayView<const double> &data,
436  const unsigned int mineral_i,
437  const unsigned int grain_i) const
438  {
439  Tensor<2,3> rotation_matrix;
440  for (unsigned int i = 0; i < Tensor<2,3>::n_independent_components ; ++i)
441  {
442  const ::TableIndices<2> index = Tensor<2,3>::unrolled_to_component_indices(i);
443  rotation_matrix[index] = data[cpo_data_position + 3 + grain_i * 10 + mineral_i * (n_grains * 10 + 2) + i];
444  }
445  return rotation_matrix;
446  }
447 
457  inline
458  void set_rotation_matrix_grains(const unsigned int cpo_data_position,
459  const ArrayView<double> &data,
460  const unsigned int mineral_i,
461  const unsigned int grain_i,
462  const Tensor<2,3> &rotation_matrix) const
463  {
464  for (unsigned int i = 0; i < Tensor<2,3>::n_independent_components ; ++i)
465  {
466  const ::TableIndices<2> index = Tensor<2,3>::unrolled_to_component_indices(i);
467  data[cpo_data_position + 3 + grain_i * 10 + mineral_i * (n_grains * 10 + 2) + i] = rotation_matrix[index];
468  }
469  }
470 
471 
472  private:
476  void
477  compute_random_rotation_matrix(Tensor<2,3> &rotation_matrix) const;
478 
494  double
495  advect_forward_euler(const unsigned int cpo_data_position,
496  const ArrayView<double> &data,
497  const unsigned int mineral_i,
498  const double dt,
499  const std::pair<std::vector<double>, std::vector<Tensor<2,3>>> &derivatives) const;
500 
516  double
517  advect_backward_euler(const unsigned int cpo_data_position,
518  const ArrayView<double> &data,
519  const unsigned int mineral_i,
520  const double dt,
521  const std::pair<std::vector<double>, std::vector<Tensor<2,3>>> &derivatives) const;
522 
523 
530  std::pair<std::vector<double>, std::vector<Tensor<2,3>>>
531  compute_derivatives_spin_tensor(const Tensor<2,3> &velocity_gradient_tensor) const;
532 
536  mutable boost::mt19937 random_number_generator;
537  unsigned int random_number_seed;
538 
539  unsigned int n_grains;
540 
541  unsigned int n_minerals;
542 
546  unsigned int water_index;
547 
554  std::vector<DeformationTypeSelector> deformation_type_selector;
555 
559  std::vector<double> volume_fractions_minerals;
560 
565 
570 
576 
582 
591 
597 
601  double exponent_p;
602 
609 
614  double mobility;
615 
620 
632  std::vector<double> CPX_RRSS;
633 
643  std::vector<double> OlivineD_RRSS;
644 
647  };
648  }
649  }
650 }
651 
652 #endif
void set_volume_fraction_mineral(const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const double volume_fraction_mineral) const
Sets the value in the data array representing the volume fraction of a mineral.
DeformationTypeSelector
The type of deformation selector used by the CPO code.
CPODerivativeAlgorithm
The algorithm used to compute the derivatives of the grain size and rotation matrix used in the advec...
DeformationType
The type of deformation used by the CPO code.
double get_volume_fraction_mineral(const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i) const
Returns the value in the data array representing the volume fraction of a mineral.
void declare_parameters(ParameterHandler &prm)
AdvectionMethod
The type of Advection method used to advect the CPO properties.
CPOInitialGrainsModel
An enum used to determine how the initial grain sizes and orientations are set for all particles...
Tensor< 2, 3 > get_rotation_matrix_grains(const unsigned int cpo_data_position, const ArrayView< const double > &data, const unsigned int mineral_i, const unsigned int grain_i) const
Gets the rotation matrix for a grain in a mineral.
void set_deformation_type(const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const DeformationType deformation_type) const
Sets the value in the data array representing the deformation type.
DeformationType get_deformation_type(const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i) const
Returns the value in the data array representing the deformation type.
void set_rotation_matrix_grains(const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const unsigned int grain_i, const Tensor< 2, 3 > &rotation_matrix) const
Sets the rotation matrix for a grain in a mineral.
double get_volume_fractions_grains(const unsigned int cpo_data_position, const ArrayView< const double > &data, const unsigned int mineral_i, const unsigned int grain_i) const
Returns the value in the data array representing the volume fraction of a grain.
void set_volume_fractions_grains(const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const unsigned int grain_i, const double volume_fractions_grains) const
Sets the value in the data array representing the volume fraction of a grain.