ASPECT
latent_heat_melt.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2013 - 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_material_model_latent_heat_melt_h
22 #define _aspect_material_model_latent_heat_melt_h
23 
26 #include <aspect/melt.h>
27 
28 namespace aspect
29 {
30  namespace MaterialModel
31  {
39  template <int dim>
42  public ::aspect::SimulatorAccess<dim>
43  {
44  public:
50  MaterialModel::MaterialModelOutputs<dim> &out) const override;
51 
52 
70  bool is_compressible () const override;
76  std::vector<double> &melt_fractions) const override;
77 
78 
86  static
87  void
88  declare_parameters (ParameterHandler &prm);
89 
93  void
94  parse_parameters (ParameterHandler &prm) override;
99  private:
101  double reference_T;
102  double eta;
109 
113  double k_value;
114 
116 
121  // for the solidus temperature
122  double A1; // °C
123  double A2; // °C/Pa
124  double A3; // °C/(Pa^2)
125 
126  // for the lherzolite liquidus temperature
127  double B1; // °C
128  double B2; // °C/Pa
129  double B3; // °C/(Pa^2)
130 
131  // for the liquidus temperature
132  double C1; // °C
133  double C2; // °C/Pa
134  double C3; // °C/(Pa^2)
135 
136  // for the reaction coefficient of pyroxene
137  double r1; // cpx/melt
138  double r2; // cpx/melt/GPa
139  double M_cpx; // mass fraction of pyroxene
140 
141  // melt fraction exponent
142  double beta;
143 
144  // entropy change upon melting
146 
151  // for the melting temperature
152  double D1; // °C
153  double D2; // °C/Pa
154  double D3; // °C/(Pa^2)
155  // for the melt-fraction dependence of productivity
156  double E1;
157  double E2;
158 
159  // for the maximum melt fraction of pyroxenite
160  double F_px_max;
161 
162  // the relative density of molten material (compared to solid)
164 
166 
171  virtual
172  double
173  melt_fraction (const double temperature,
174  const double pressure,
175  const std::vector<double> &compositional_fields,
176  const Point<dim> &position) const;
177 
178  virtual
179  double
181  const double pressure,
182  const std::vector<double> &compositional_fields,
183  const Point<dim> &position) const;
184 
185  virtual
186  double
188  const double pressure,
189  const std::vector<double> &compositional_fields,
190  const Point<dim> &position) const;
191 
192  double
193  entropy_derivative ( const double temperature,
194  const double pressure,
195  const std::vector<double> &compositional_fields,
196  const Point<dim> &position,
197  const NonlinearDependence::Dependence dependence) const;
198  };
199 
200  }
201 }
202 
203 #endif
double entropy_derivative(const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position, const NonlinearDependence::Dependence dependence) const
static void declare_parameters(ParameterHandler &prm)
bool is_compressible() const override
virtual double pyroxenite_melt_fraction(const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const
void parse_parameters(ParameterHandler &prm) override
Definition: compat.h:59
void evaluate(const MaterialModel::MaterialModelInputs< dim > &in, MaterialModel::MaterialModelOutputs< dim > &out) const override
void melt_fractions(const MaterialModel::MaterialModelInputs< dim > &in, std::vector< double > &melt_fractions) const override
virtual double melt_fraction(const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const
virtual double peridotite_melt_fraction(const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const