ASPECT
utilities.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 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_material_model_utilities_h
22 #define _aspect_material_model_utilities_h
23 
24 #include <aspect/global.h>
25 #include <deal.II/base/point.h>
26 #include <deal.II/base/symmetric_tensor.h>
27 #include <deal.II/fe/component_mask.h>
28 #include <deal.II/base/signaling_nan.h>
29 #include <deal.II/base/parameter_handler.h>
30 
31 namespace aspect
32 {
33  template <int dim> class SimulatorAccess;
34  namespace Utilities
35  {
50  using namespace ::Utilities;
51 
52  template <int dim>
53  class StructuredDataLookup;
54  }
55  namespace MaterialModel
56  {
57  template <int dim> class MaterialModelOutputs;
58  template <int dim> struct EquationOfStateOutputs;
59 
65  namespace MaterialUtilities
66  {
67  namespace Lookup
68  {
76  {
77  public:
78 
79  double
80  specific_heat(const double temperature,
81  const double pressure) const;
82 
83  double
84  density(const double temperature,
85  const double pressure) const;
86 
87  double
88  thermal_expansivity(const double temperature,
89  const double pressure) const;
90 
91  double
92  seismic_Vp(const double temperature,
93  const double pressure) const;
94 
95  double
96  seismic_Vs(const double temperature,
97  const double pressure) const;
98 
99  double
100  enthalpy(const double temperature,
101  const double pressure) const;
102 
108  double
109  dHdT (const double temperature,
110  const double pressure) const;
111 
117  double
118  dHdp (const double temperature,
119  const double pressure) const;
120 
133  std::array<std::pair<double, unsigned int>,2>
134  enthalpy_derivatives(const std::vector<double> &temperatures,
135  const std::vector<double> &pressures,
136  const unsigned int n_substeps = 1) const;
137 
138  double
139  dRhodp (const double temperature,
140  const double pressure) const;
141 
146  unsigned int
147  dominant_phase (const double temperature,
148  const double pressure) const;
149 
154  bool
155  has_dominant_phase() const;
156 
161  std::vector<std::string>
162  phase_volume_column_names() const;
163 
168  double
169  phase_volume_fraction(const int phase_id,
170  const double temperature,
171  const double pressure) const;
172 
177  std::array<double,2>
178  get_pT_steps() const;
179 
180 
187  const std::vector<std::string> &
188  get_dominant_phase_names() const;
189 
190  protected:
198  double
199  value (const double temperature,
200  const double pressure,
201  const Table<2, double> &values,
202  const bool interpol) const;
203 
209  unsigned int
210  value (const double temperature,
211  const double pressure,
212  const Table<2, unsigned int> &values) const;
213 
217  double get_nT(const double temperature) const;
218 
222  double get_np(const double pressure) const;
223 
224  ::Table<2,double> density_values;
225  ::Table<2,double> thermal_expansivity_values;
226  ::Table<2,double> specific_heat_values;
227  ::Table<2,double> vp_values;
228  ::Table<2,double> vs_values;
229  ::Table<2,double> enthalpy_values;
230  ::Table<2,unsigned int> dominant_phase_indices;
231 
238  std::vector<std::string> phase_column_names;
239  std::vector<::Table<2,double>> phase_volume_fractions;
240 
241  double delta_press;
242  double min_press;
243  double max_press;
244  double delta_temp;
245  double min_temp;
246  double max_temp;
247  unsigned int n_temperature;
248  unsigned int n_pressure;
249  unsigned int n_phases;
250  unsigned int n_columns;
253  std::vector<std::string> dominant_phase_names;
254  };
255 
261  {
262  public:
263  HeFESToReader(const std::string &material_filename,
264  const std::string &derivatives_filename,
265  const bool interpol,
266  const MPI_Comm comm);
267  };
268 
274  {
275  public:
276  PerplexReader(const std::string &filename,
277  const bool interpol,
278  const MPI_Comm comm);
279  };
280 
286  {
287  public:
288 
292  void
293  initialize(const MPI_Comm comm,
294  const std::string &data_directory,
295  const std::string &material_file_name);
296 
300  double
301  specific_heat(const double entropy,
302  const double pressure) const;
303 
307  double
308  density(const double entropy,
309  const double pressure) const;
310 
314  double
315  thermal_expansivity(const double entropy,
316  const double pressure) const;
317 
321  double
322  temperature(const double entropy,
323  const double pressure) const;
324 
328  double
329  seismic_vp(const double entropy,
330  const double pressure) const;
331 
335  double
336  seismic_vs(const double entropy,
337  const double pressure) const;
338 
342  Tensor<1, 2>
343  density_gradient(const double entropy,
344  const double pressure) const;
345 
346  private:
350  std::unique_ptr<Utilities::StructuredDataLookup<2>> material_lookup;
351  };
352  }
353 
369  std::vector<double>
370  compute_only_composition_fractions(const std::vector<double> &compositional_fields,
371  const std::vector<unsigned int> &indices_to_use);
372 
391  std::vector<double>
392  compute_composition_fractions(const std::vector<double> &compositional_fields,
393  const ComponentMask &field_mask = ComponentMask());
394 
402  std::vector<double>
403  compute_volumes_from_masses(const std::vector<double> &masses,
404  const std::vector<double> &densities,
405  const bool return_as_fraction);
406 
416  {
421  };
422 
423 
424 
431  parse_compositional_averaging_operation (const std::string &parameter_name,
432  const ParameterHandler &prm);
433 
434 
435 
455  double average_value (const std::vector<double> &volume_fractions,
456  const std::vector<double> &parameter_values,
457  const CompositionalAveragingOperation &average_type);
458 
459 
460 
475  template <int dim>
476  void
478  const std::vector<double> &mass_fractions,
479  const std::vector<double> &volume_fractions,
480  const unsigned int i,
482 
483 
484 
488  namespace PhaseUtilities
489  {
496  {
499  };
500  }
501 
518  double phase_average_value (const std::vector<double> &phase_function_values,
519  const std::vector<unsigned int> &n_phase_transitions_per_composition,
520  const std::vector<double> &parameter_values,
521  const unsigned int composition_index,
523 
524 
525 
531  template <int dim>
533  {
538  PhaseFunctionInputs(const double temperature,
539  const double pressure,
540  const double depth,
541  const double pressure_depth_derivative,
542  const unsigned int phase_transition_index);
543 
544  double temperature;
545  double pressure;
546  double depth;
548 
561  };
562 
571  template <int dim>
573  {
574  public:
575 
580  void initialize();
581 
588  double compute_value (const PhaseFunctionInputs<dim> &in) const;
589 
595  double compute_derivative () const;
596 
600  unsigned int n_phase_transitions () const;
601 
605  unsigned int n_phases () const;
606 
610  unsigned int n_phases_over_all_chemical_compositions () const;
611 
615  const std::vector<unsigned int> &
616  n_phase_transitions_for_each_chemical_composition () const;
617 
621  const std::vector<unsigned int> &
622  n_phases_for_each_chemical_composition () const;
623 
630  const std::vector<unsigned int> &
631  n_phase_transitions_for_each_composition () const;
632 
639  const std::vector<unsigned int> &
640  n_phases_for_each_composition () const;
641 
648  static
649  void
650  declare_parameters (ParameterHandler &prm);
651 
658  void
659  parse_parameters (ParameterHandler &prm);
660 
661 
662  private:
666  std::string data_directory;
667 
671  std::vector<std::string> material_file_names;
672 
676  std::vector<double> minimum_temperature;
677 
681  std::vector<double> maximum_temperature;
682 
686  std::vector<double> interval_temperature;
687 
691  std::vector<double> minimum_pressure;
692 
696  std::vector<double> maximum_pressure;
697 
701  std::vector<double> interval_pressure;
702 
707  std::vector<std::unique_ptr<Utilities::StructuredDataLookup<2>>> material_lookup;
708 
715  std::vector<unsigned int> transition_indicators;
716 
720  std::unique_ptr<std::vector<unsigned int>> n_phase_transitions_per_composition;
721 
725  std::vector<unsigned int> n_phases_per_composition;
726 
731 
735  std::vector<unsigned int> n_phases_per_chemical_composition;
736 
740  unsigned int n_phases_total;
741 
746  };
747 
755  template <int dim>
757  {
758  public:
764  double compute_value (const PhaseFunctionInputs<dim> &in) const;
765 
770  double compute_derivative (const PhaseFunctionInputs<dim> &in) const;
771 
775  unsigned int n_phase_transitions () const;
776 
780  unsigned int n_phases () const;
781 
785  unsigned int n_phases_over_all_chemical_compositions () const;
786 
791  double get_transition_slope (const unsigned int phase_transition_index) const;
792 
796  double get_transition_depth (const unsigned int phase_transition_index) const;
797 
801  const std::vector<unsigned int> &
802  n_phase_transitions_for_each_chemical_composition () const;
803 
807  const std::vector<unsigned int> &
808  n_phases_for_each_chemical_composition () const;
809 
816  const std::vector<unsigned int> &
817  n_phase_transitions_for_each_composition () const;
818 
825  const std::vector<unsigned int> &
826  n_phases_for_each_composition () const;
827 
834  static
835  void
836  declare_parameters (ParameterHandler &prm);
837 
844  void
845  parse_parameters (ParameterHandler &prm);
846 
847 
848  private:
853  std::vector<double> transition_depths;
854  std::vector<double> transition_pressures;
855  std::vector<double> transition_temperatures;
856  std::vector<double> transition_widths;
857  std::vector<double> transition_pressure_widths;
858  std::vector<double> transition_slopes;
861 
869 
873  std::unique_ptr<std::vector<unsigned int>> n_phase_transitions_per_composition;
874 
878  std::vector<unsigned int> n_phases_per_composition;
879 
884 
888  std::vector<unsigned int> n_phases_per_chemical_composition;
889 
893  unsigned int n_phases_total;
894 
899  };
900  }
901  }
902 }
903 
904 
905 #endif
void fill_averaged_equation_of_state_outputs(const EquationOfStateOutputs< dim > &eos_outputs, const std::vector< double > &mass_fractions, const std::vector< double > &volume_fractions, const unsigned int i, MaterialModelOutputs< dim > &out)
std::vector< unsigned int > n_phases_per_composition
Definition: utilities.h:878
std::unique_ptr< Utilities::StructuredDataLookup< 2 > > material_lookup
Definition: utilities.h:350
std::vector<::Table< 2, double > > phase_volume_fractions
Definition: utilities.h:239
void declare_parameters(ParameterHandler &prm)
std::unique_ptr< std::vector< unsigned int > > n_phase_transitions_per_composition
Definition: utilities.h:720
std::unique_ptr< std::vector< unsigned int > > n_phase_transitions_per_composition
Definition: utilities.h:873
std::vector< unsigned int > n_phase_transitions_per_chemical_composition
Definition: utilities.h:730
double average_value(const std::vector< double > &volume_fractions, const std::vector< double > &parameter_values, const CompositionalAveragingOperation &average_type)
std::vector< std::unique_ptr< Utilities::StructuredDataLookup< 2 > > > material_lookup
Definition: utilities.h:707
double phase_average_value(const std::vector< double > &phase_function_values, const std::vector< unsigned int > &n_phase_transitions_per_composition, const std::vector< double > &parameter_values, const unsigned int composition_index, const PhaseUtilities::PhaseAveragingOperation operation=PhaseUtilities::arithmetic)
std::vector< unsigned int > n_phase_transitions_per_chemical_composition
Definition: utilities.h:883
CompositionalAveragingOperation parse_compositional_averaging_operation(const std::string &parameter_name, const ParameterHandler &prm)
std::vector< unsigned int > n_phases_per_chemical_composition
Definition: utilities.h:888
Definition: compat.h:59
std::vector< double > compute_volumes_from_masses(const std::vector< double > &masses, const std::vector< double > &densities, const bool return_as_fraction)
std::vector< double > compute_only_composition_fractions(const std::vector< double > &compositional_fields, const std::vector< unsigned int > &indices_to_use)
std::vector< double > compute_composition_fractions(const std::vector< double > &compositional_fields, const ComponentMask &field_mask=ComponentMask())