World Builder  1.1.0-pre
A geodynamic initial conditions generator
linear.cc
Go to the documentation of this file.
1 /*
2  Copyright (C) 2018-2024 by the authors of the World Builder code.
3 
4  This file is part of the World Builder.
5 
6  This program is free software: you can redistribute it and/or modify
7  it under the terms of the GNU Lesser General Public License as published
8  by the Free Software Foundation, either version 2 of the License, or
9  (at your option) any later version.
10 
11  This program 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 Lesser General Public License for more details.
15 
16  You should have received a copy of the GNU Lesser General Public License
17  along with this program. If not, see <https://www.gnu.org/licenses/>.
18 */
19 
21 
22 
23 #include "world_builder/nan.h"
27 #include "world_builder/world.h"
28 
29 
30 namespace WorldBuilder
31 {
32  using namespace Utilities;
33 
34  namespace Features
35  {
36  namespace SubductingPlateModels
37  {
38  namespace Temperature
39  {
41  :
42  min_depth(NaN::DSNAN),
43  max_depth(NaN::DSNAN),
44  top_temperature(NaN::DSNAN),
45  bottom_temperature(NaN::DSNAN),
46  operation(Operations::REPLACE)
47  {
48  this->world = world_;
49  this->name = "linear";
50  }
51 
53  = default;
54 
55  void
56  Linear::declare_entries(Parameters &prm, const std::string & /*unused*/)
57  {
58  // Document plugin and require entries if needed.
59  // Add `max distance slab top` center to the required parameters.
60  prm.declare_entry("", Types::Object({"max distance slab top"}),
61  "Linear temperature model. Can be set to use an adiabatic temperature at the boundaries.");
62 
63  // Declare entries of this plugin
64  prm.declare_entry("min distance slab top", Types::Double(0),
65  "todo The depth in meters from which the composition of this feature is present.");
66 
67  prm.declare_entry("max distance slab top", Types::Double(std::numeric_limits<double>::max()),
68  "todo The depth in meters to which the composition of this feature is present.");
69 
70  prm.declare_entry("top temperature", Types::Double(293.15),
71  "The temperature at the top in degree Kelvin of this feature."
72  "If the value is below zero, the an adiabatic temperature is used.");
73 
74  prm.declare_entry("bottom temperature", Types::Double(-1),
75  "The temperature at the bottom in degree Kelvin of this feature. "
76  "If the value is below zero, an adiabatic temperature is used.");
77  }
78 
79  void
81  {
82  min_depth = prm.get<double>("min distance slab top");
83  max_depth = prm.get<double>("max distance slab top");
84  WBAssert(max_depth >= min_depth, "max depth needs to be larger or equal to min depth.");
85  operation = string_operations_to_enum(prm.get<std::string>("operation"));
86  top_temperature = prm.get<double>("top temperature");
87  bottom_temperature = prm.get<double>("bottom temperature");
88  }
89 
90 
91  double
92  Linear::get_temperature(const Point<3> & /*position_in_cartesian_coordinates*/,
93  const double /*depth*/,
94  const double gravity_norm,
95  double temperature_,
96  const double /*feature_min_depth*/,
97  const double /*feature_max_depth*/,
99  const AdditionalParameters & /*additional_parameters*/) const
100  {
101  if (distance_from_plane.distance_from_plane <= max_depth && distance_from_plane.distance_from_plane >= min_depth)
102  {
103  const double min_depth_local = min_depth;
104  const double max_depth_local = max_depth;
105 
106 
107  double top_temperature_local = top_temperature;
108  if (top_temperature_local < 0)
109  {
110  top_temperature_local = this->world->potential_mantle_temperature *
111  std::exp(((this->world->thermal_expansion_coefficient * gravity_norm) /
112  this->world->specific_heat) * min_depth_local);
113  }
114 
115  double bottom_temperature_local = bottom_temperature;
116  if (bottom_temperature_local < 0)
117  {
118  bottom_temperature_local = this->world->potential_mantle_temperature *
119  std::exp(((this->world->thermal_expansion_coefficient * gravity_norm) /
120  this->world->specific_heat) * max_depth_local);
121  }
122 
123  const double new_temperature = top_temperature_local +
124  (distance_from_plane.distance_from_plane - min_depth_local) *
125  ((bottom_temperature_local - top_temperature_local) / (max_depth_local - min_depth_local));
126  return apply_operation(operation,temperature_,new_temperature);
127 
128 
129  }
130  return temperature_;
131  }
132 
134  } // namespace Temperature
135  } // namespace SubductingPlateModels
136  } // namespace Features
137 } // namespace WorldBuilder
138 
const double DSNAN
Definition: nan.h:41
double potential_mantle_temperature
Definition: world.h:268
double get_temperature(const Point< 3 > &position, const double depth, const double gravity, double temperature, const double feature_min_depth, const double feature_max_depth, const WorldBuilder::Utilities::PointDistanceFromCurvedPlanes &distance_from_planes, const AdditionalParameters &additional_parameters) const override final
Definition: linear.cc:92
Operations string_operations_to_enum(const std::string &operation)
static void declare_entries(Parameters &prm, const std::string &parent_name="")
Definition: linear.cc:56
#define WBAssert(condition, message)
Definition: assert.h:27
#define WB_REGISTER_FEATURE_SUBDUCTING_PLATE_TEMPERATURE_MODEL(classname, name)
Definition: interface.h:156
double thermal_expansion_coefficient
Definition: world.h:283
void parse_entries(Parameters &prm) override final
Definition: linear.cc:80
double apply_operation(const Operations operation, const double old_value, const double new_value)
void declare_entry(const std::string &name, const Types::Interface &type, const std::string &documentation)
Definition: parameters.cc:197
T get(const std::string &name)