ASPECT
crystal_preferred_orientation.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2022 - 2023 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 
118  template <int dim>
120  {
121  public:
126 
131  void
132  initialize () override;
133 
145  void
146  initialize_one_particle_property (const Point<dim> &position,
147  std::vector<double> &particle_properties) const override;
148 
170  void
171  update_one_particle_property (const unsigned int data_position,
172  const Point<dim> &position,
173  const Vector<double> &solution,
174  const std::vector<Tensor<1,dim>> &gradients,
175  const ArrayView<double> &particle_properties) const override;
176 
182  need_update () const override;
183 
188  late_initialization_mode () const override;
189 
194  UpdateFlags
195  get_needed_update_flags () const override;
196 
204  std::vector<std::pair<std::string, unsigned int>>
205  get_property_information() const override;
206 
224  std::pair<std::vector<double>, std::vector<Tensor<2,3>>>
225  compute_derivatives(const unsigned int cpo_index,
226  const ArrayView<double> &data,
227  const unsigned int mineral_i,
228  const SymmetricTensor<2,3> &strain_rate_3d,
229  const Tensor<2,3> &velocity_gradient_tensor,
230  const Point<dim> &position,
231  const double temperature,
232  const double pressure,
233  const Tensor<1,dim> &velocity,
234  const std::vector<double> &compositions,
235  const SymmetricTensor<2,dim> &strain_rate,
236  const SymmetricTensor<2,dim> &deviatoric_strain_rate,
237  const double water_content) const;
238 
254  std::pair<std::vector<double>, std::vector<Tensor<2,3>>>
255  compute_derivatives_drex_2004(const unsigned int cpo_index,
256  const ArrayView<double> &data,
257  const unsigned int mineral_i,
258  const SymmetricTensor<2,3> &strain_rate_3d,
259  const Tensor<2,3> &velocity_gradient_tensor,
260  const std::array<double,4> ref_resolved_shear_stress,
261  const bool prevent_nondimensionalization = false) const;
262 
263 
267  static
268  void
269  declare_parameters (ParameterHandler &prm);
270 
274  void
275  parse_parameters (ParameterHandler &prm) override;
276 
280  unsigned int
281  get_number_of_grains() const;
282 
286  unsigned int
287  get_number_of_minerals() const;
288 
299  determine_deformation_type(const DeformationTypeSelector deformation_type_selector,
300  const Point<dim> &position,
301  const double temperature,
302  const double pressure,
303  const Tensor<1,dim> &velocity,
304  const std::vector<double> &compositions,
305  const SymmetricTensor<2,dim> &strain_rate,
306  const SymmetricTensor<2,dim> &deviatoric_strain_rate,
307  const double water_content) const;
308 
314  determine_deformation_type_karato_2008(const double stress,
315  const double water_content) const;
316 
326  std::array<double,4>
327  reference_resolved_shear_stress_from_deformation_type(DeformationType deformation_type,
328  double max_value = 1e60) const;
329 
337  inline
338  double get_deformation_type(const unsigned int cpo_data_position,
339  const ArrayView<double> &data,
340  const unsigned int mineral_i) const
341  {
342  return data[cpo_data_position + 0 + mineral_i * (n_grains * 10 + 2)];
343  }
344 
353  inline
354  void set_deformation_type(const unsigned int cpo_data_position,
355  const ArrayView<double> &data,
356  const unsigned int mineral_i,
357  const double deformation_type) const
358  {
359  data[cpo_data_position + 0 + mineral_i * (n_grains * 10 + 2)] = deformation_type;
360  }
361 
369  inline
370  double get_volume_fraction_mineral(const unsigned int cpo_data_position,
371  const ArrayView<double> &data,
372  const unsigned int mineral_i) const
373  {
374  return data[cpo_data_position + 1 + mineral_i *(n_grains * 10 + 2)];
375  }
376 
385  inline
386  void set_volume_fraction_mineral(const unsigned int cpo_data_position,
387  const ArrayView<double> &data,
388  const unsigned int mineral_i,
389  const double volume_fraction_mineral) const
390  {
391  data[cpo_data_position + 1 + mineral_i *(n_grains * 10 + 2)] = volume_fraction_mineral;
392  }
393 
402  inline
403  double get_volume_fractions_grains(const unsigned int cpo_data_position,
404  const ArrayView<const double> &data,
405  const unsigned int mineral_i,
406  const unsigned int grain_i) const
407  {
408  return data[cpo_data_position + 2 + grain_i * 10 + mineral_i * (n_grains * 10 + 2)];
409  }
410 
420  inline
421  void set_volume_fractions_grains(const unsigned int cpo_data_position,
422  const ArrayView<double> &data,
423  const unsigned int mineral_i,
424  const unsigned int grain_i,
425  const double volume_fractions_grains) const
426  {
427  data[cpo_data_position + 2 + grain_i * 10 + mineral_i * (n_grains * 10 + 2)] = volume_fractions_grains;
428  }
429 
439  inline
440  Tensor<2,3> get_rotation_matrix_grains(const unsigned int cpo_data_position,
441  const ArrayView<const double> &data,
442  const unsigned int mineral_i,
443  const unsigned int grain_i) const
444  {
445  Tensor<2,3> rotation_matrix;
446  for (unsigned int i = 0; i < Tensor<2,3>::n_independent_components ; ++i)
447  {
448  const ::TableIndices<2> index = Tensor<2,3>::unrolled_to_component_indices(i);
449  rotation_matrix[index] = data[cpo_data_position + 3 + grain_i * 10 + mineral_i * (n_grains * 10 + 2) + i];
450  }
451  return rotation_matrix;
452  }
453 
463  inline
464  void set_rotation_matrix_grains(const unsigned int cpo_data_position,
465  const ArrayView<double> &data,
466  const unsigned int mineral_i,
467  const unsigned int grain_i,
468  const Tensor<2,3> &rotation_matrix) const
469  {
470  for (unsigned int i = 0; i < Tensor<2,3>::n_independent_components ; ++i)
471  {
472  const ::TableIndices<2> index = Tensor<2,3>::unrolled_to_component_indices(i);
473  data[cpo_data_position + 3 + grain_i * 10 + mineral_i * (n_grains * 10 + 2) + i] = rotation_matrix[index];
474  }
475  }
476 
477 
478  private:
482  void
483  compute_random_rotation_matrix(Tensor<2,3> &rotation_matrix) const;
484 
500  double
501  advect_forward_euler(const unsigned int cpo_data_position,
502  const ArrayView<double> &data,
503  const unsigned int mineral_i,
504  const double dt,
505  const std::pair<std::vector<double>, std::vector<Tensor<2,3>>> &derivatives) const;
506 
522  double
523  advect_backward_euler(const unsigned int cpo_data_position,
524  const ArrayView<double> &data,
525  const unsigned int mineral_i,
526  const double dt,
527  const std::pair<std::vector<double>, std::vector<Tensor<2,3>>> &derivatives) const;
528 
529 
536  std::pair<std::vector<double>, std::vector<Tensor<2,3>>>
537  compute_derivatives_spin_tensor(const Tensor<2,3> &velocity_gradient_tensor) const;
538 
542  mutable boost::mt19937 random_number_generator;
543  unsigned int random_number_seed;
544 
545  unsigned int n_grains;
546 
547  unsigned int n_minerals;
548 
552  unsigned int water_index;
553 
560  std::vector<DeformationTypeSelector> deformation_type_selector;
561 
565  std::vector<double> volume_fractions_minerals;
566 
571 
576 
582 
588 
597 
603 
607  double exponent_p;
608 
615 
620  double mobility;
621 
624  };
625  }
626  }
627 }
628 
629 #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)
double 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.
AdvectionMethod
The type of Advection method used to advect the CPO properties.
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 double deformation_type) const
Sets 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.