ASPECT
dynamic_core.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 doc/COPYING. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 
22 
23 #ifndef _aspect_boundary_temperature_dynamic_core_h
24 #define _aspect_boundary_temperature_dynamic_core_h
25 
28 
29 #include <tuple>
30 
31 namespace aspect
32 {
33  namespace BoundaryTemperature
34  {
35 
39  namespace internal
40  {
41  struct CoreData
42  {
47  double Qs,Qr,Qg,Qk,Ql;
48 
56  double Es,Er,Eg,Ek,El,Eh;
57 
64  double Ri,Ti,Xi;
65 
69  double Q,H;
70 
74  double dt;
75 
80  double dR_dt,dT_dt,dX_dt;
81 
85  double Q_OES;
86 
88  };
89 
94  using SolveTimeStepResult = std::tuple<double, double, double>;
95  }
96 
97 
105  template <int dim>
106  class DynamicCore : public Interface<dim>, public aspect::SimulatorAccess<dim>
107  {
108  public:
112  DynamicCore();
113 
118  void
119  update() override;
120 
124  const internal::CoreData &
125  get_core_data() const;
126 
132  bool
133  is_OES_used() const;
134 
145  double
146  boundary_temperature (const types::boundary_id boundary_indicator,
147  const Point<dim> &location) const override;
148 
156  double
157  minimal_temperature (const std::set<types::boundary_id> &fixed_boundary_ids) const override;
158 
166  double
167  maximal_temperature (const std::set<types::boundary_id> &fixed_boundary_ids) const override;
168 
173  static
174  void
175  declare_parameters (ParameterHandler &prm);
176 
181  void
182  parse_parameters (ParameterHandler &prm) override;
183 
184  private:
185 
191 
196 
201 
205  types::boundary_id inner_boundary_id;
206  types::boundary_id outer_boundary_id;
207 
211  double init_dT_dt;
212 
216  double init_dR_dt;
217 
221  double init_dX_dt;
222 
227 
231  double Rc;
232 
236  double X_init;
237 
241  double Delta;
242 
246  double P_CMB;
247 
255  double Tm0;
256  double Tm1;
257  double Tm2;
258  double Theta;
260 
264  bool use_bw11;
265 
266  //Variables for formulation of \cite NPB+04
270  double K0;
271 
275  double Alpha;
276 
280  double Rho_0;
281 
285  double Rho_cen;
286 
290  double Lh;
291 
295  double Rh;
296 
300  double Beta_c;
301 
305  double k_c;
306 
310  double Cp;
311 
316 
320  std::vector<double> heating_rate;
321 
325  std::vector<double> half_life;
326 
330  std::vector<double> initial_concentration;
331 
335  double L;
336  double D;
337 
341  double Mc;
342 
346  unsigned int max_steps;
347 
351  double dTa;
352 
356  std::string name_OES;
358  {
359  double t;
360  double w;
361  };
362  std::vector<struct str_data_OES> data_OES;
363  void read_data_OES();
364  double compute_OES(double t) const;
365 
392  internal::SolveTimeStepResult solve_time_step() const;
393 
398  double compute_dT(const double r) const;
399 
404  double compute_Tc(const double r) const;
405 
410  double compute_Ts(const double r) const;
411 
416  double compute_solidus(const double X, const double pressure) const;
417 
422  double compute_initial_Ri(const double T) const;
423 
428  double compute_X(const double r) const;
429 
433  double compute_mass(const double r) const;
434 
438  double fun_Sn(const double B, const double R, const unsigned int n) const;
439 
443  double compute_rho(const double r) const;
444 
449  double compute_T(const double Tc, const double r) const;
450 
454  double compute_pressure(const double r) const;
455 
459  double compute_gravity_potential(const double r) const;
460 
469  std::pair<double,double>
470  compute_specific_heating(const double Tc) const;
471 
479  std::pair<double,double>
480  compute_radio_heating(const double Tc) const;
481 
491  std::pair<double,double>
492  compute_gravity_heating(const double Tc, const double r, const double X) const;
493 
501  std::pair<double,double>
502  compute_adiabatic_heating(const double Tc) const;
503 
512  std::pair<double,double>
513  compute_latent_heating(const double Tc, const double r) const;
514 
519  double
520  compute_heat_solution(const double Tc, const double r, const double X) const;
521 
525  double compute_radioheating_rate() const;
526 
531  void update_core_data();
532 
533  };
534  }
535 }
536 
537 
538 #endif
std::tuple< double, double, double > SolveTimeStepResult
Definition: dynamic_core.h:94
std::vector< struct str_data_OES > data_OES
Definition: dynamic_core.h:362
void declare_parameters(ParameterHandler &prm)
std::vector< double > initial_concentration
Definition: dynamic_core.h:330