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