World Builder  1.1.0-pre
A geodynamic initial conditions generator
oceanic_plate.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 
27 #include "world_builder/nan.h"
34 #include "world_builder/world.h"
35 
36 
37 namespace WorldBuilder
38 {
39  using namespace Utilities;
40 
41  namespace Features
42  {
43  OceanicPlate::OceanicPlate(WorldBuilder::World *world_)
44  :
45  min_depth(NaN::DSNAN),
46  max_depth(NaN::DSNAN)
47  {
48  this->world = world_;
49  this->name = "oceanic plate";
50  }
51 
53  = default;
54 
55 
56 
58  {
59  using namespace rapidjson;
60  Document &declarations = prm.declarations;
61 
62  const std::string path = prm.get_full_json_path();
63 
64  /*
65  ideally:
66  {
67  "model": "fault",
68  "name": "${1:default name}",
69  "dip point":[0.0,0.0],
70  "coordinates": [[0.0,0.0]],
71  "segments": [],
72  "sections": [],
73  "temperature models":[{"model":"uniform", "temperature":600.0}],
74  "composition models":[{"model":"uniform", "compositions": [0], "fractions":[1.0]}]
75  }
76  */
77 
78  Pointer((path + "/body").c_str()).Set(declarations,"object");
79  Pointer((path + "/body/model").c_str()).Set(declarations,"oceanic plate");
80  Pointer((path + "/body/name").c_str()).Set(declarations,"${1:My Oceanic Plate}");
81  Pointer((path + "/body/coordinates").c_str()).Create(declarations).SetArray();
82  Pointer((path + "/body/temperature models").c_str()).Create(declarations).SetArray();
83  Pointer((path + "/body/composition models").c_str()).Create(declarations).SetArray();
84  }
85 
86 
87 
88  void
90  const std::string & /*unused*/,
91  const std::vector<std::string> &required_entries)
92  {
93  prm.declare_entry("", Types::Object(required_entries), "Oceanic plate object. Requires properties `model` and `coordinates`.");
94 
96  "The depth from which this feature is present");
97  prm.declare_entry("max depth", Types::OneOf(Types::Double(std::numeric_limits<double>::max()),Types::Array(Types::ValueAtPoints(std::numeric_limits<double>::max(), 2.))),
98  "The depth to which this feature is present");
99  prm.declare_entry("temperature models",
101  "A list of temperature models.");
102  prm.declare_entry("composition models",
104  "A list of composition models.");
105  prm.declare_entry("grains models",
107  "A list of grains models.");
108  prm.declare_entry("velocity models",
110  "A list of velocity models.");
111  }
112 
113  void
115  {
116 
117  const CoordinateSystem coordinate_system = prm.coordinate_system->natural_coordinate_system();
118 
119  this->name = prm.get<std::string>("name");
120 
121  std::string tag = prm.get<std::string>("tag");
122  if (tag == "")
123  {
124  tag = "oceanic plate";
125  }
127 
128  this->get_coordinates("coordinates", prm, coordinate_system);
129 
132 
135 
136 
138 
139  prm.enter_subsection("temperature models");
140  {
141  for (unsigned int i = 0; i < temperature_models.size(); ++i)
142  {
143  prm.enter_subsection(std::to_string(i));
144  {
145  temperature_models[i]->parse_entries(prm,coordinates);
146  }
147  prm.leave_subsection();
148  }
149  }
150  prm.leave_subsection();
151 
152 
154 
155  prm.enter_subsection("composition models");
156  {
157  for (unsigned int i = 0; i < composition_models.size(); ++i)
158  {
159  prm.enter_subsection(std::to_string(i));
160  {
161  composition_models[i]->parse_entries(prm,coordinates);
162  }
163  prm.leave_subsection();
164  }
165  }
166  prm.leave_subsection();
167 
168 
170 
171  prm.enter_subsection("grains models");
172  {
173  for (unsigned int i = 0; i < grains_models.size(); ++i)
174  {
175  prm.enter_subsection(std::to_string(i));
176  {
177  grains_models[i]->parse_entries(prm,coordinates);
178  }
179  prm.leave_subsection();
180  }
181  }
182  prm.leave_subsection();
183 
185 
186  prm.enter_subsection("velocity models");
187  {
188  for (unsigned int i = 0; i < velocity_models.size(); ++i)
189  {
190  prm.enter_subsection(std::to_string(i));
191  {
192  velocity_models[i]->parse_entries(prm,coordinates);
193  }
194  prm.leave_subsection();
195  }
196  }
197  prm.leave_subsection();
198  }
199 
200 
201  void
202  OceanicPlate::properties(const Point<3> &position_in_cartesian_coordinates,
203  const Objects::NaturalCoordinate &position_in_natural_coordinates,
204  const double depth,
205  const std::vector<std::array<unsigned int,3>> &properties,
206  const double gravity_norm,
207  const std::vector<size_t> &entry_in_output,
208  std::vector<double> &output) const
209  {
210  if (depth <= max_depth && depth >= min_depth &&
212  world->parameters.coordinate_system->natural_coordinate_system())))
213  {
214  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;
215  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;
216  if (depth <= max_depth_local && depth >= min_depth_local)
217  {
218  for (unsigned int i_property = 0; i_property < properties.size(); ++i_property)
219  {
220  switch (properties[i_property][0])
221  {
222  case 1: // temperature
223  {
224  for (const auto &temperature_model: temperature_models)
225  {
226  output[entry_in_output[i_property]] = temperature_model->get_temperature(position_in_cartesian_coordinates,
227  position_in_natural_coordinates,
228  depth,
229  gravity_norm,
230  output[entry_in_output[i_property]],
231  min_depth_local,
232  max_depth_local);
233 
234  WBAssert(!std::isnan(output[entry_in_output[i_property]]), "Temperature is not a number: " << output[entry_in_output[i_property]]
235  << ", based on a temperature model with the name " << temperature_model->get_name() << ", in feature " << this->name);
236  WBAssert(std::isfinite(output[entry_in_output[i_property]]), "Temperature is not a finite: " << output[entry_in_output[i_property]]
237  << ", based on a temperature model with the name " << temperature_model->get_name() << ", in feature " << this->name);
238 
239  }
240  break;
241  }
242  case 2: // composition
243  {
244  for (const auto &composition_model: composition_models)
245  {
246  output[entry_in_output[i_property]] = composition_model->get_composition(position_in_cartesian_coordinates,
247  position_in_natural_coordinates,
248  depth,
249  properties[i_property][1],
250  output[entry_in_output[i_property]],
251  min_depth_local,
252  max_depth_local);
253 
254  WBAssert(!std::isnan(output[entry_in_output[i_property]]), "Composition is not a number: " << output[entry_in_output[i_property]]
255  << ", based on a composition model with the name " << composition_model->get_name() << ", in feature " << this->name);
256  WBAssert(std::isfinite(output[entry_in_output[i_property]]), "Composition is not a finite: " << output[entry_in_output[i_property]]
257  << ", based on a composition model with the name " << composition_model->get_name() << ", in feature " << this->name);
258 
259  }
260 
261  break;
262  }
263  case 3: // grains
264  {
265  WorldBuilder::grains grains(output,properties[i_property][2],entry_in_output[i_property]);
266  for (const auto &grains_model: grains_models)
267  {
268  grains = grains_model->get_grains(position_in_cartesian_coordinates,
269  position_in_natural_coordinates,
270  depth,
271  properties[i_property][1],
272  grains,
273  min_depth_local,
274  max_depth_local);
275 
276  }
277  grains.unroll_into(output,entry_in_output[i_property]);
278  break;
279  }
280  case 4:
281  {
282  output[entry_in_output[i_property]] = static_cast<double>(tag_index);
283  break;
284  }
285  case 5: // velocity
286  {
287  std::array<double, 3> velocity = {{0,0,0}};
288  for (const auto &velocity_model: velocity_models)
289  {
290  velocity = velocity_model->get_velocity(position_in_cartesian_coordinates,
291  position_in_natural_coordinates,
292  depth,
293  gravity_norm,
294  velocity,
295  min_depth_local,
296  max_depth_local);
297 
298  //WBAssert(!std::isnan(output[entry_in_output[i_property]]), "Velocity is not a number: " << output[entry_in_output[i_property]]
299  // << ", based on a velocity model with the name " << velocity_model->get_name() << ", in feature " << this->name);
300  //WBAssert(std::isfinite(output[entry_in_output[i_property]]), "Velocity is not a finite: " << output[entry_in_output[i_property]]
301  // << ", based on a velocity model with the name " << velocity_model->get_name() << ", in feature " << this->name);
302  }
303 
304  output[entry_in_output[i_property]] = velocity[0];
305  output[entry_in_output[i_property]+1] = velocity[1];
306  output[entry_in_output[i_property]+2] = velocity[2];
307  break;
308  }
309  break;
310  default:
311  {
312  WBAssertThrow(false,
313  "Internal error: Unimplemented property provided. " <<
314  "Only temperature (1), composition (2), grains (3), tag (4) or velocity (5) are allowed. "
315  "Provided property number was: " << properties[i_property][0]);
316  }
317  }
318  }
319  }
320  }
321  }
322 
326  WB_REGISTER_FEATURE(OceanicPlate, oceanic plate)
327  } // namespace Features
328 } // namespace WorldBuilder
329 
static void declare_entries(Parameters &prm, const std::string &parent_name, const std::vector< std::string > &required_entries)
Definition: interface.cc:42
std::array< double, 2 > get_surface_coordinates() const
const double DSNAN
Definition: nan.h:41
std::vector< std::unique_ptr< Features::OceanicPlateModels::Velocity::Interface > > velocity_models
size_t add_vector_unique(std::vector< std::string > &vector, const std::string &add_string)
#define WB_REGISTER_FEATURE(classname, name)
Definition: interface.h:212
void enter_subsection(const std::string &name)
Definition: parameters.cc:1871
void properties(const Point< 3 > &position_in_cartesian_coordinates, const Objects::NaturalCoordinate &position_in_natural_coordinates, const double depth, const std::vector< std::array< unsigned int, 3 >> &properties, const double gravity, const std::vector< size_t > &entry_in_output, std::vector< double > &output) const override final
std::vector< std::unique_ptr< Features::OceanicPlateModels::Composition::Interface > > composition_models
#define WBAssert(condition, message)
Definition: assert.h:27
std::string get_full_json_path(size_t max_size=std::numeric_limits< size_t >::max()) const
Definition: parameters.cc:1933
void get_coordinates(const std::string &name, Parameters &prm, const CoordinateSystem coordinate_system)
Definition: interface.cc:140
std::vector< std::unique_ptr< Features::OceanicPlateModels::Grains::Interface > > grains_models
static void make_snippet(Parameters &prm)
static void declare_entries(Parameters &prm, const std::string &parent_name="", const std::vector< std::string > &required_entries={})
void parse_entries(Parameters &prm) override final
std::vector< std::string > feature_tags
Definition: world.h:308
WorldBuilder::World * world
Definition: interface.h:127
#define WBAssertThrow(condition, message)
Definition: assert.h:40
bool get_unique_pointers(const std::string &name, std::vector< std::unique_ptr< T > > &vector)
Definition: parameters.cc:1770
static void declare_entries(Parameters &prm, const std::string &parent_name, const std::vector< std::string > &required_entries)
Definition: interface.cc:40
SurfaceValueInfo local_value(const Point< 2 > &check_point) const
Definition: surface.cc:149
std::unique_ptr< WorldBuilder::CoordinateSystems::Interface > coordinate_system
Definition: parameters.h:268
static void declare_entries(Parameters &prm, const std::string &parent_name, const std::vector< std::string > &required_entries)
Definition: interface.cc:41
bool polygon_contains_point(const std::vector< Point< 2 > > &point_list, const Point< 2 > &point)
Definition: utilities.cc:46
void declare_entry(const std::string &name, const Types::Interface &type, const std::string &documentation)
Definition: parameters.cc:197
Parameters parameters
Definition: world.h:253
T get(const std::string &name)
rapidjson::Document declarations
Definition: parameters.h:248
void unroll_into(std::vector< double > &vector, const size_t start_entry=0) const
Definition: grains.cc:58
std::vector< std::unique_ptr< Features::OceanicPlateModels::Temperature::Interface > > temperature_models
static void declare_entries(Parameters &prm, const std::string &parent_name, const std::vector< std::string > &required_entries)
Definition: interface.cc:41
std::vector< Point< 2 > > coordinates
Definition: interface.h:154