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