34 using namespace Utilities;
38 using namespace FeatureUtilities;
39 namespace ContinentalPlateModels
45 min_depth(NaN::
DSNAN),
46 max_depth(NaN::
DSNAN),
47 top_temperature(NaN::
DSNAN),
48 bottom_temperature(NaN::
DSNAN),
52 this->
name =
"linear";
64 "Linear temperature model. Can be set to use an adiabatic temperature at the boundaries.");
68 "The depth in meters from which the composition of this feature is present.");
71 "The depth in meters to which the composition of this feature is present.");
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.");
78 "The temperature at the top in degree Kelvin of this feature. " 79 "If the value is below zero, an adiabatic temperature is used.");
90 WBAssert(max_depth >= min_depth,
"max depth needs to be larger or equal to min depth.");
101 const double gravity_norm,
103 const double feature_min_depth,
104 const double feature_max_depth)
const 106 if (depth <= max_depth && depth >=
min_depth)
110 if (depth <= max_depth_local && depth >= min_depth_local)
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);
116 if (top_temperature_local < 0)
120 this->world->specific_heat) * min_depth_local_local);
124 if (bottom_temperature_local < 0)
128 this->world->specific_heat) * max_depth_local_local);
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)));
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);
Objects::Surface max_depth_surface
void parse_entries(Parameters &prm, const std::vector< Point< 2 >> &coordinates) override final
double potential_mantle_temperature
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
Operations string_operations_to_enum(const std::string &operation)
#define WBAssert(condition, message)
static void declare_entries(Parameters &prm, const std::string &parent_name="")
Point< 2 > get_surface_point() const
Objects::Surface min_depth_surface
SurfaceValueInfo local_value(const Point< 2 > &check_point) const
double thermal_expansion_coefficient
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)
double bottom_temperature
T get(const std::string &name)
#define WB_REGISTER_FEATURE_CONTINENTAL_PLATE_TEMPERATURE_MODEL(classname, name)
WorldBuilder::World * world