ASPECT
S40RTS_perturbation.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 
22 #ifndef _aspect_initial_temperature_S40RTS_perturbation_h
23 #define _aspect_initial_temperature_S40RTS_perturbation_h
24 
26 #include <aspect/utilities.h>
27 
28 namespace aspect
29 {
30  namespace InitialTemperature
31  {
32  namespace internal
33  {
34  namespace S40RTS
35  {
37  {
38  public:
39  SphericalHarmonicsLookup(const std::string &filename,
40  const MPI_Comm comm);
41 
43  const std::vector<double> &
44  cos_coeffs() const;
45 
47  const std::vector<double> &
48  sin_coeffs() const;
49 
50  unsigned int maxdegree() const;
51 
52  private:
53  unsigned int order;
54  std::vector<double> a_lm;
55  std::vector<double> b_lm;
56  };
57 
59  {
60  public:
61  SplineDepthsLookup(const std::string &filename,
62  const MPI_Comm comm);
63 
64  const std::vector<double> &
65  spline_depths() const;
66 
67  private:
68  std::vector<double> depths;
69  };
70  }
71  }
72 
73  template <int dim>
75 
76 
86  template <int dim>
87  class S40RTSPerturbation : public Interface<dim>, public ::aspect::SimulatorAccess<dim>
88  {
89  public:
94 
99  void
100  initialize () override;
101 
105  virtual
106  double get_Vs (const Point<dim> &position) const;
107 
111  double initial_temperature (const Point<dim> &position) const override;
112 
116  static
117  void
118  declare_parameters (ParameterHandler &prm);
119 
123  void
124  parse_parameters (ParameterHandler &prm) override;
125 
126  private:
127 
132  {
135  };
136 
141 
145  std::string data_directory;
147 
156 
171 
178 
184 
189  unsigned int specified_max_degree;
190 
197 
202  std::unique_ptr<internal::S40RTS::SphericalHarmonicsLookup> spherical_harmonics_lookup;
203 
208  std::unique_ptr<internal::S40RTS::SplineDepthsLookup> spline_depths_lookup;
209 
214 
218  unsigned int vs_to_density_index;
219 
224 
225  template <int dim2> friend class PatchOnS40RTS;
226  };
227 
228  }
229 }
230 
231 #endif
SphericalHarmonicsLookup(const std::string &filename, const MPI_Comm comm)
const std::vector< double > & sin_coeffs() const
Declare a function that returns the sine coefficients.
void declare_parameters(ParameterHandler &prm)
aspect::Utilities::AsciiDataProfile< dim > profile
std::unique_ptr< internal::S40RTS::SphericalHarmonicsLookup > spherical_harmonics_lookup
std::unique_ptr< internal::S40RTS::SplineDepthsLookup > spline_depths_lookup
Definition: compat.h:59
const std::vector< double > & cos_coeffs() const
Declare a function that returns the cosine coefficients.