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 #include "world_builder/nan.h"
28 #include "world_builder/world.h"
29 
30 
31 namespace WorldBuilder
32 {
33  using namespace Utilities;
34 
35  namespace Features
36  {
37  namespace OceanicPlateModels
38  {
39  namespace Temperature
40  {
42  :
43  min_depth(NaN::DSNAN),
44  max_depth(NaN::DSNAN),
45  top_temperature(NaN::DSNAN),
46  bottom_temperature(NaN::DSNAN),
47  operation(Operations::REPLACE)
48  {
49  this->world = world_;
50  this->name = "linear";
51  }
52 
54  = default;
55 
56  void
57  Linear::declare_entries(Parameters &prm, const std::string & /*unused*/)
58  {
59  // Document plugin and require entries if needed.
60  // Add max depth to the required parameters.
61  prm.declare_entry("", Types::Object({"max depth"}),
62  "Linear temperature model. Can be set to use an adiabatic temperature at the boundaries.");
63 
64  // Declare entries of this plugin
66  "The depth in meters from which the temperature of this feature is present.");
67 
68  prm.declare_entry("max depth", Types::OneOf(Types::Double(std::numeric_limits<double>::max()),Types::Array(Types::ValueAtPoints(std::numeric_limits<double>::max(), 2.))),
69  "The depth in meters to which the temperature of this feature is present.");
70 
71  prm.declare_entry("top temperature", Types::Double(293.15),
72  "The temperature at the top in degree Kelvin of this feature."
73  "If the value is below zero, the an adiabatic temperature is used.");
74 
75  prm.declare_entry("bottom temperature", Types::Double(-1),
76  "The temperature at the top in degree Kelvin of this feature. "
77  "If the value is below zero, an adiabatic temperature is used.");
78  }
79 
80  void
81  Linear::parse_entries(Parameters &prm, const std::vector<Point<2>> &coordinates)
82  {
83  min_depth_surface = Objects::Surface(prm.get("min depth",coordinates));
85  max_depth_surface = Objects::Surface(prm.get("max depth",coordinates));
87  WBAssert(max_depth >= min_depth, "max depth needs to be larger or equal to min depth.");
88  operation = string_operations_to_enum(prm.get<std::string>("operation"));
89  top_temperature = prm.get<double>("top temperature");
90  bottom_temperature = prm.get<double>("bottom temperature");
91  }
92 
93 
94  double
95  Linear::get_temperature(const Point<3> & /*position_in_cartesian_coordinates*/,
96  const Objects::NaturalCoordinate &position_in_natural_coordinates,
97  const double depth,
98  const double gravity_norm,
99  double temperature_,
100  const double feature_min_depth,
101  const double feature_max_depth) const
102  {
103  if (depth <= max_depth && depth >= min_depth)
104  {
105  const double min_depth_local = min_depth_surface.constant_value ? min_depth : min_depth_surface.local_value(position_in_natural_coordinates.get_surface_point()).interpolated_value;
106  const double max_depth_local = max_depth_surface.constant_value ? max_depth : max_depth_surface.local_value(position_in_natural_coordinates.get_surface_point()).interpolated_value;
107  if (depth <= max_depth_local && depth >= min_depth_local)
108  {
109  const double min_depth_local_local = std::max(feature_min_depth, min_depth_local);
110  const double max_depth_local_local = std::min(feature_max_depth, max_depth_local);
111 
112  double top_temperature_local = top_temperature;
113  if (top_temperature_local < 0)
114  {
115  top_temperature_local = this->world->potential_mantle_temperature *
116  std::exp(((this->world->thermal_expansion_coefficient * gravity_norm) /
117  this->world->specific_heat) * min_depth_local_local);
118  }
119 
120  double bottom_temperature_local = bottom_temperature;
121  if (bottom_temperature_local < 0)
122  {
123  bottom_temperature_local = this->world->potential_mantle_temperature *
124  std::exp(((this->world->thermal_expansion_coefficient * gravity_norm) /
125  this->world->specific_heat) * max_depth_local_local);
126  }
127 
128  const double new_temperature = top_temperature_local + (max_depth_local_local - min_depth_local_local < 10.0*std::numeric_limits<double>::epsilon() ? 0.0 :
129  (depth - min_depth_local_local) *
130  ((bottom_temperature_local - top_temperature_local) / (max_depth_local_local - min_depth_local_local)));
131 
132  return apply_operation(operation,temperature_,new_temperature);
133  }
134  }
135 
136 
137  return temperature_;
138  }
139 
141  } // namespace Temperature
142  } // namespace OceanicPlateModels
143  } // namespace Features
144 } // namespace WorldBuilder
145 
const double DSNAN
Definition: nan.h:41
double potential_mantle_temperature
Definition: world.h:268
static void declare_entries(Parameters &prm, const std::string &parent_name="")
Definition: linear.cc:57
Operations string_operations_to_enum(const std::string &operation)
#define WBAssert(condition, message)
Definition: assert.h:27
#define WB_REGISTER_FEATURE_OCEANIC_PLATE_TEMPERATURE_MODEL(classname, name)
Definition: interface.h:151
void parse_entries(Parameters &prm, const std::vector< Point< 2 >> &coordinates) override final
Definition: linear.cc:81
SurfaceValueInfo local_value(const Point< 2 > &check_point) const
Definition: surface.cc:149
double get_temperature(const Point< 3 > &position, const Objects::NaturalCoordinate &position_in_natural_coordinates, const double depth, const double gravity, double temperature, const double feature_min_depth, const double feature_max_depth) const override final
Definition: linear.cc:95
double thermal_expansion_coefficient
Definition: world.h:283
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)