World Builder  1.1.0-pre
A geodynamic initial conditions generator
subducting_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 
23 #include "glm/glm.h"
30 #include "world_builder/world.h"
31 #include <algorithm>
32 #include <array>
33 
34 using namespace std;
35 
36 namespace WorldBuilder
37 {
38  using namespace Utilities;
39 
40  namespace Features
41  {
42  SubductingPlate::SubductingPlate(WorldBuilder::World *world_)
43  :
44  reference_point(0,0,cartesian)
45  {
46  this->world = world_;
47  this->name = "subducting plate";
48  }
49 
51  = default;
52 
53 
54 
56  {
57  using namespace rapidjson;
58  Document &declarations = prm.declarations;
59 
60  const std::string path = prm.get_full_json_path();
61 
62  /*
63  ideally:
64  {
65  "model": "fault",
66  "name": "${1:default name}",
67  "dip point":[0.0,0.0],
68  "coordinates": [[0.0,0.0]],
69  "segments": [],
70  "sections": [],
71  "temperature models":[{"model":"uniform", "temperature":600.0}],
72  "composition models":[{"model":"uniform", "compositions": [0], "fractions":[1.0]}]
73  }
74  */
75 
76  Pointer((path + "/body").c_str()).Set(declarations,"object");
77  Pointer((path + "/body/model").c_str()).Set(declarations,"subducting plate");
78  Pointer((path + "/body/name").c_str()).Set(declarations,"${1:My Subducting Plate}");
79  Pointer((path + "/body/coordinates").c_str()).Create(declarations).SetArray();
80  Pointer((path + "/body/temperature models").c_str()).Create(declarations).SetArray();
81  Pointer((path + "/body/composition models").c_str()).Create(declarations).SetArray();
82  }
83 
84 
85 
86  void
88  const std::string &parent_name,
89  const std::vector<std::string> &required_entries)
90  {
91  // This statement is needed because of the recursion associated with
92  // the sections entry.
93  if (parent_name == "items")
94  {
95  prm.enter_subsection("properties");
96  }
97  else
98  {
99  prm.declare_entry("", Types::Object(required_entries), "Subducting slab object. Requires properties `model` and `coordinates`.");
100  }
101 
102 
103  prm.declare_entry("min depth", Types::Double(0),
104  "The depth to which this feature is present");
105  prm.declare_entry("max depth", Types::Double(std::numeric_limits<double>::max()),
106  "The depth to which this feature is present");
107  prm.declare_entry("dip point", Types::Point<2>(),
108  "The depth to which this feature is present");
109 
110  /*prm.declare_entry("segments", Types::Array(Types::Segment(0,Point<2>(0,0,invalid),Point<2>(0,0,invalid),Point<2>(0,0,invalid),
111  Types::PluginSystem("", Features::SubductingPlateModels::Temperature::Interface::declare_entries, {"model"}),
112  Types::PluginSystem("", Features::SubductingPlateModels::Composition::Interface::declare_entries, {"model"}),
113  Types::PluginSystem("", Features::SubductingPlateModels::Grains::Interface::declare_entries, {"model"}))),
114  "The depth to which this feature is present");*/
120  "The depth to which this feature is present");
121 
122  prm.declare_entry("temperature models",
124  "A list of temperature models.");
125  prm.declare_entry("composition models",
127  "A list of composition models.");
128  prm.declare_entry("grains models",
130  "A list of grains models.");
131  prm.declare_entry("velocity models",
133  "A list of velocity models.");
134 
135  if (parent_name != "items")
136  {
137  // This only happens if we are not in sections
138  prm.declare_entry("sections", Types::Array(Types::PluginSystem("",Features::SubductingPlate::declare_entries, {"coordinate"}, false)),"A list of feature properties for a coordinate.");
139  }
140  else
141  {
142 
143  // this only happens in sections
144  prm.declare_entry("coordinate", Types::UnsignedInt(0),
145  "The coordinate which should be overwritten");
146 
147  prm.leave_subsection();
148  }
149  }
150 
151  void
153  {
154  const CoordinateSystem coordinate_system = prm.coordinate_system->natural_coordinate_system();
155 
156  this->name = prm.get<std::string>("name");
157 
158  std::string tag = prm.get<std::string>("tag");
159  if (tag == "")
160  {
161  tag = "subducting plate";
162  }
164 
165  this->get_coordinates("coordinates", prm, coordinate_system);
166 
167 
168  starting_depth = prm.get<double>("min depth");
169  maximum_depth = prm.get<double>("max depth");
170 
171  const size_t n_sections = this->original_number_of_coordinates;
172 
173  reference_point = prm.get<Point<2> >("dip point");
174 
175  if (coordinate_system == spherical)
176  {
177  // When spherical, input is in degrees, so change to radians for internal use.
178  reference_point *= (Consts::PI/180.);
179  }
180 
181  default_temperature_models.resize(0);
182  default_composition_models.resize(0);
183  default_grains_models.resize(0);
184  default_velocity_models.resize(0);
189 
190  // get the default segments.
195 
196 
197  // This vector stores segments to this coordinate/section.
198  //First used (raw) pointers to the segment relevant to this coordinate/section,
199  // but I do not trust it won't fail when memory is moved. So storing the all the data now.
200  segment_vector.resize(0);
201  segment_vector.resize(n_sections, default_segment_vector);
202 
203 
204  // now search whether a section is present, if so, replace the default segments.
205  std::vector<std::unique_ptr<Features::SubductingPlate> > sections_vector;
206  prm.get_unique_pointers("sections", sections_vector);
207 
208  prm.enter_subsection("sections");
209  for (unsigned int i_section = 0; i_section < n_sections; ++i_section)
210  {
211  // first check whether this section/coordinate has a a special overwrite
212  for (unsigned int i_sector = 0; i_sector < sections_vector.size(); ++i_sector)
213  {
214  prm.enter_subsection(std::to_string(i_sector));
215  {
216  const unsigned int change_coord_number = prm.get<unsigned int>("coordinate");
217 
218  WBAssertThrow(segment_vector.size() > change_coord_number, "Error: for subducting plate with name: '" << this->name
219  << "', trying to change the section of coordinate " << change_coord_number
220  << " while only " << segment_vector.size() << " coordinates are defined.");
221 
222  std::vector<std::shared_ptr<Features::SubductingPlateModels::Temperature::Interface> > local_default_temperature_models;
223  std::vector<std::shared_ptr<Features::SubductingPlateModels::Composition::Interface> > local_default_composition_models;
224  std::vector<std::shared_ptr<Features::SubductingPlateModels::Grains::Interface> > local_default_grains_models;
225  std::vector<std::shared_ptr<Features::SubductingPlateModels::Velocity::Interface> > local_default_velocity_models;
226 
227  if (!prm.get_shared_pointers<Features::SubductingPlateModels::Temperature::Interface>("temperature models", local_default_temperature_models))
228  {
229  // no local temperature model, use global default
230  local_default_temperature_models = default_temperature_models;
231  }
232 
233  if (!prm.get_shared_pointers<Features::SubductingPlateModels::Composition::Interface>("composition models", local_default_composition_models))
234  {
235  // no local composition model, use global default
236  local_default_composition_models = default_composition_models;
237  }
238 
239  if (!prm.get_shared_pointers<Features::SubductingPlateModels::Grains::Interface>("grains models", local_default_grains_models))
240  {
241  // no local composition model, use global default
242  local_default_grains_models = default_grains_models;
243  }
244 
245  if (!prm.get_shared_pointers<Features::SubductingPlateModels::Velocity::Interface>("velocity models", local_default_velocity_models))
246  {
247  // no local composition model, use global default
248  local_default_velocity_models = default_velocity_models;
249  }
250 
251  segment_vector[change_coord_number] = prm.get_vector<Objects::Segment<Features::SubductingPlateModels::Temperature::Interface,
252  Features::SubductingPlateModels::Composition::Interface,
253  Features::SubductingPlateModels::Grains::Interface,
254  Features::SubductingPlateModels::Velocity::Interface> >("segments", local_default_temperature_models, local_default_composition_models, local_default_grains_models, local_default_velocity_models);
255 
256 
257  WBAssertThrow(segment_vector[change_coord_number].size() == default_segment_vector.size(),
258  "Error: There are not the same amount of segments in section with coordinate " << change_coord_number
259  << " (" << segment_vector[change_coord_number].size() << " segments) as in the default segment ("
260  << default_segment_vector.size() << " segments). This is not allowed.");
261 
262  prm.enter_subsection("segments");
263  {
264  for (unsigned int i = 0; i < segment_vector[change_coord_number].size(); ++i)
265  {
266  prm.enter_subsection(std::to_string(i));
267  {
268  prm.enter_subsection("temperature models");
269  {
270  for (unsigned int j = 0; j < segment_vector[change_coord_number][i].temperature_systems.size(); ++j)
271  {
272  prm.enter_subsection(std::to_string(j));
273  {
274  segment_vector[change_coord_number][i].temperature_systems[j]->parse_entries(prm);
275  }
276  prm.leave_subsection();
277  }
278  }
279  prm.leave_subsection();
280 
281 
282  prm.enter_subsection("composition models");
283  {
284  for (unsigned int j = 0; j < segment_vector[change_coord_number][i].composition_systems.size(); ++j)
285  {
286  prm.enter_subsection(std::to_string(j));
287  {
288  segment_vector[change_coord_number][i].composition_systems[j]->parse_entries(prm);
289  }
290  prm.leave_subsection();
291  }
292  }
293  prm.leave_subsection();
294 
295  prm.enter_subsection("grains models");
296  {
297  for (unsigned int j = 0; j < segment_vector[change_coord_number][i].grains_systems.size(); ++j)
298  {
299  prm.enter_subsection(std::to_string(j));
300  {
301  segment_vector[change_coord_number][i].grains_systems[j]->parse_entries(prm);
302  }
303  prm.leave_subsection();
304  }
305  }
306  prm.leave_subsection();
307 
308  prm.enter_subsection("velocity models");
309  {
310  for (unsigned int j = 0; j < segment_vector[change_coord_number][i].velocity_systems.size(); ++j)
311  {
312  prm.enter_subsection(std::to_string(j));
313  {
314  segment_vector[change_coord_number][i].velocity_systems[j]->parse_entries(prm);
315  }
316  prm.leave_subsection();
317  }
318  }
319  prm.leave_subsection();
320  }
321  prm.leave_subsection();
322  }
323  }
324  prm.leave_subsection();
325 
326  }
327  prm.leave_subsection();
328  }
329 
330  }
331  prm.leave_subsection();
332 
333 
334  prm.enter_subsection("segments");
335  {
336  for (unsigned int i = 0; i < default_segment_vector.size(); ++i)
337  {
338  prm.enter_subsection(std::to_string(i));
339  {
340  prm.enter_subsection("temperature models");
341  {
342  for (unsigned int j = 0; j < default_segment_vector[i].temperature_systems.size(); ++j)
343  {
344  prm.enter_subsection(std::to_string(j));
345  {
346  default_segment_vector[i].temperature_systems[j]->parse_entries(prm);
347  }
348  prm.leave_subsection();
349  }
350  }
351  prm.leave_subsection();
352 
353 
354  prm.enter_subsection("composition models");
355  {
356  for (unsigned int j = 0; j < default_segment_vector[i].composition_systems.size(); ++j)
357  {
358  prm.enter_subsection(std::to_string(j));
359  {
360  default_segment_vector[i].composition_systems[j]->parse_entries(prm);
361  }
362  prm.leave_subsection();
363  }
364  }
365  prm.leave_subsection();
366 
367 
368  prm.enter_subsection("grains models");
369  {
370  for (unsigned int j = 0; j < default_segment_vector[i].grains_systems.size(); ++j)
371  {
372  prm.enter_subsection(std::to_string(j));
373  {
374  default_segment_vector[i].grains_systems[j]->parse_entries(prm);
375  }
376  prm.leave_subsection();
377  }
378  }
379  prm.leave_subsection();
380 
381  prm.enter_subsection("velocity models");
382  {
383  for (unsigned int j = 0; j < default_segment_vector[i].velocity_systems.size(); ++j)
384  {
385  prm.enter_subsection(std::to_string(j));
386  {
387  default_segment_vector[i].velocity_systems[j]->parse_entries(prm);
388  }
389  prm.leave_subsection();
390  }
391  }
392  prm.leave_subsection();
393  }
394  prm.leave_subsection();
395  }
396  }
397  prm.leave_subsection();
398 
406 
407  for (unsigned int i = 0; i < segment_vector.size(); ++i)
408  {
409  double local_total_slab_length = 0;
410  slab_segment_lengths[i].resize(segment_vector[i].size());
411  slab_segment_thickness[i].resize(segment_vector[i].size(), Point<2>(invalid));
413  slab_segment_angles[i].resize(segment_vector[i].size(), Point<2>(invalid));
414  for (unsigned int j = 0; j < segment_vector[i].size(); ++j)
415  {
416  slab_segment_lengths[i][j] = segment_vector[i][j].value_length;
417  local_total_slab_length += segment_vector[i][j].value_length;
418 
419  slab_segment_thickness[i][j] = segment_vector[i][j].value_thickness;
422  slab_segment_top_truncation[i][j] = segment_vector[i][j].value_top_truncation;
423 
424  slab_segment_angles[i][j] = segment_vector[i][j].value_angle * (Consts::PI/180);
425  }
426  total_slab_length[i] = local_total_slab_length;
427  maximum_total_slab_length = std::max(maximum_total_slab_length, local_total_slab_length);
428  }
429 
430  // Here, we compute the spherical bounding box using the two extreme points of the box containing all the surface
431  // coordinates and an additional buffer zone that accounts for the slab thickness and length. The first and second
432  // points correspond to the lower left and the upper right corners of the bounding box, respectively (see the
433  // documentation in include/bounding_box.h).
434  // For the spherical system, the buffer zone along the longitudal direction is calculated using the
435  // corresponding latitude points.
436 
437  // Find minimal and maximal coordinates. Do this by finding the
438  // leftmost/rightmost point with regard to either the [0] or [1]
439  // coordinate, and then takes its [0] or [1] element.
440  auto compare_x_coordinate = [](auto p1, auto p2)
441  {
442  return p1[0]<p2[0];
443  };
444 
445  min_along_x = (*std::min_element(coordinates.begin(), coordinates.end(), compare_x_coordinate)) [0];
446  max_along_x = (*std::max_element(coordinates.begin(), coordinates.end(), compare_x_coordinate)) [0];
447 
448 
449  auto compare_y_coordinate = [](auto p1, auto p2)
450  {
451  return p1[1]<p2[1];
452  };
453 
454  min_along_y = (*std::min_element(coordinates.begin(), coordinates.end(), compare_y_coordinate)) [1];
455  max_along_y = (*std::max_element(coordinates.begin(), coordinates.end(), compare_y_coordinate)) [1];
456 
459 
461 
462  // Compute the surface bounding box
464  if (world->parameters.coordinate_system->natural_coordinate_system() == CoordinateSystem::spherical)
465  {
466  const double starting_radius_inv = 1 / (world->parameters.coordinate_system->max_model_depth());
467  std::pair<Point<2>, Point<2> > &spherical_bounding_box = surface_bounding_box.get_boundary_points();
468 
469  const double buffer_around_slab_spherical = 2 * Consts::PI * buffer_around_slab_cartesian * starting_radius_inv;
470 
471  spherical_bounding_box.first = {(min_along_x - buffer_around_slab_spherical * min_lat_cos_inv) ,
472  (min_along_y - buffer_around_slab_spherical), spherical
473  } ;
474 
475  spherical_bounding_box.second = {(max_along_x + buffer_around_slab_spherical * max_lat_cos_inv) ,
476  (max_along_y + buffer_around_slab_spherical), spherical
477  };
478  }
479  else if (world->parameters.coordinate_system->natural_coordinate_system() == CoordinateSystem::cartesian)
480  {
481  std::pair<Point<2>, Point<2> > &bounding_box = surface_bounding_box.get_boundary_points();
482  bounding_box.first = {min_along_x, min_along_y, cartesian};
483  bounding_box.second = {max_along_x, max_along_y, cartesian};
485  }
486  }
487 
488 
489 
490  const BoundingBox<2> &
492  {
493  return surface_bounding_box;
494  }
495 
496 
497  void
498  SubductingPlate::properties(const Point<3> &position_in_cartesian_coordinates,
499  const Objects::NaturalCoordinate &position_in_natural_coordinates,
500  const double depth,
501  const std::vector<std::array<unsigned int,3>> &properties,
502  const double gravity_norm,
503  const std::vector<size_t> &entry_in_output,
504  std::vector<double> &output) const
505  {
506  // The depth variable is the distance from the surface to the position, the depth
507  // coordinate is the distance from the bottom of the model to the position and
508  // the starting radius is the distance from the bottom of the model to the surface.
509  const double starting_radius = position_in_natural_coordinates.get_depth_coordinate() + depth - starting_depth;
510 
511  WBAssert(std::abs(starting_radius) > std::numeric_limits<double>::epsilon(), "World Builder error: starting_radius can not be zero. "
512  << "Position = " << position_in_cartesian_coordinates[0] << ':' << position_in_cartesian_coordinates[1] << ':' << position_in_cartesian_coordinates[2]
513  << ", position_in_natural_coordinates.get_depth_coordinate() = " << position_in_natural_coordinates.get_depth_coordinate()
514  << ", depth = " << depth
515  << ", starting_depth " << starting_depth
516  );
517 
518  // todo: explain and check -starting_depth
519  if (depth <= maximum_depth && depth >= starting_depth && depth <= maximum_total_slab_length + maximum_slab_thickness &&
520  get_surface_bounding_box().point_inside(Point<2>(position_in_natural_coordinates.get_surface_coordinates(),
521  world->parameters.coordinate_system->natural_coordinate_system())))
522  {
523  /*WBAssert(coordinates.size() == slab_segment_lengths.size(),
524  "Internal error: The size of coordinates (" << coordinates.size()
525  << ") and slab_segment_lengths (" << slab_segment_lengths.size() << ") are different.");
526  WBAssert(coordinates.size() == slab_segment_angles.size(),
527  "Internal error: The size of coordinates (" << coordinates.size()
528  << ") and slab_segment_angles (" << slab_segment_angles.size() << ") are different.");
529  WBAssert(coordinates.size() == slab_segment_angles.size(),
530  "Internal error: The size of coordinates (" << coordinates.size()
531  << ") and one_dimensional_coordinates (" << one_dimensional_coordinates.size() << ") are different.");*/
532  // todo: explain
533  const WorldBuilder::Utilities::PointDistanceFromCurvedPlanes distance_from_planes =
534  WorldBuilder::Utilities::distance_point_from_curved_planes(position_in_cartesian_coordinates,
535  position_in_natural_coordinates,
537  coordinates,
540  starting_radius,
542  false,
543  this->bezier_curve);
544 
545  const double distance_from_plane = distance_from_planes.distance_from_plane;
546  const double distance_along_plane = distance_from_planes.distance_along_plane;
547  const double section_fraction = distance_from_planes.fraction_of_section;
548  const size_t current_section = distance_from_planes.section;
549  const size_t next_section = current_section + 1;
550  const size_t current_segment = distance_from_planes.segment; // the original value was a unsigned in, converting it back.
551  //const size_t next_segment = current_segment + 1;
552  const double segment_fraction = distance_from_planes.fraction_of_segment;
553 
554  if (abs(distance_from_plane) < std::numeric_limits<double>::infinity() || (distance_along_plane) < std::numeric_limits<double>::infinity())
555  {
556  // We want to do both section (horizontal) and segment (vertical) interpolation.
557  // first for thickness
558  const double thickness_up = slab_segment_thickness[current_section][current_segment][0]
559  + section_fraction
560  * (slab_segment_thickness[next_section][current_segment][0]
561  - slab_segment_thickness[current_section][current_segment][0]);
562  const double thickness_down = slab_segment_thickness[current_section][current_segment][1]
563  + section_fraction
564  * (slab_segment_thickness[next_section][current_segment][1]
565  - slab_segment_thickness[current_section][current_segment][1]);
566  const double thickness_local = thickness_up + segment_fraction * (thickness_down - thickness_up);
567 
568  // if the thickness is zero, we don't need to compute anything, so return.
569  if (std::fabs(thickness_local) < 2.0 * std::numeric_limits<double>::epsilon())
570  return;
571 
572  // secondly for top truncation
573  const double top_truncation_up = slab_segment_top_truncation[current_section][current_segment][0]
574  + section_fraction
575  * (slab_segment_top_truncation[next_section][current_segment][0]
576  - slab_segment_top_truncation[current_section][current_segment][0]);
577  const double top_truncation_down = slab_segment_top_truncation[current_section][current_segment][1]
578  + section_fraction
579  * (slab_segment_top_truncation[next_section][current_segment][1]
580  - slab_segment_top_truncation[current_section][current_segment][1]);
581  const double top_truncation_local = top_truncation_up + segment_fraction * (top_truncation_down - top_truncation_up);
582 
583  // if the thickness is smaller than what is truncated off at the top, we don't need to compute anything, so return.
584  if (thickness_local < top_truncation_local)
585  return;
586 
587  const double max_slab_length = total_slab_length[current_section] +
588  section_fraction *
589  (total_slab_length[next_section] - total_slab_length[current_section]);
590 
591  if (distance_from_plane >= top_truncation_local &&
592  distance_from_plane <= thickness_local &&
593  distance_along_plane >= 0 &&
594  distance_along_plane <= max_slab_length)
595  {
596  // Inside the slab!
597  const Features::AdditionalParameters additional_parameters = {max_slab_length,thickness_local};
598  for (unsigned int i_property = 0; i_property < properties.size(); ++i_property)
599  {
600  switch (properties[i_property][0])
601  {
602  case 1: // temperature
603  {
604  double temperature_current_section = output[entry_in_output[i_property]];
605  double temperature_next_section = output[entry_in_output[i_property]];
606 
607  for (const auto &temperature_model: segment_vector[current_section][current_segment].temperature_systems)
608  {
609  temperature_current_section = temperature_model->get_temperature(position_in_cartesian_coordinates,
610  depth,
611  gravity_norm,
612  temperature_current_section,
613  starting_depth,
615  distance_from_planes,
616  additional_parameters);
617 
618  WBAssert(!std::isnan(temperature_current_section), "Temperature is not a number: " << temperature_current_section
619  << ", based on a temperature model with the name " << temperature_model->get_name() << ", in feature " << this->name);
620  WBAssert(std::isfinite(temperature_current_section), "Temperature is not a finite: " << temperature_current_section
621  << ", based on a temperature model with the name " << temperature_model->get_name() << ", in feature " << this->name);
622 
623  }
624 
625  for (const auto &temperature_model: segment_vector[next_section][current_segment].temperature_systems)
626  {
627  temperature_next_section = temperature_model->get_temperature(position_in_cartesian_coordinates,
628  depth,
629  gravity_norm,
630  temperature_next_section,
631  starting_depth,
633  distance_from_planes,
634  additional_parameters);
635 
636  WBAssert(!std::isnan(temperature_next_section), "Temperature is not a number: " << temperature_next_section
637  << ", based on a temperature model with the name " << temperature_model->get_name() << ", in feature " << this->name);
638  WBAssert(std::isfinite(temperature_next_section), "Temperature is not a finite: " << temperature_next_section
639  << ", based on a temperature model with the name " << temperature_model->get_name() << ", in feature " << this->name);
640 
641  }
642 
643  // linear interpolation between current and next section temperatures
644  output[entry_in_output[i_property]] = temperature_current_section + section_fraction * (temperature_next_section - temperature_current_section);
645  break;
646  }
647  case 2: // composition
648  {
649  double composition_current_section = output[entry_in_output[i_property]];
650  double composition_next_section = output[entry_in_output[i_property]];
651 
652  for (const auto &composition_model: segment_vector[current_section][current_segment].composition_systems)
653  {
654  composition_current_section = composition_model->get_composition(position_in_cartesian_coordinates,
655  depth,
656  properties[i_property][1],
657  composition_current_section,
658  starting_depth,
660  distance_from_planes,
661  additional_parameters);
662 
663  WBAssert(!std::isnan(composition_current_section), "Composition_current_section is not a number: " << composition_current_section
664  << ", based on a composition model with the name " << composition_model->get_name() << ", in feature " << this->name);
665  WBAssert(std::isfinite(composition_current_section), "Composition_current_section is not a finite: " << composition_current_section
666  << ", based on a composition model with the name " << composition_model->get_name() << ", in feature " << this->name);
667 
668  }
669 
670  for (const auto &composition_model: segment_vector[next_section][current_segment].composition_systems)
671  {
672  composition_next_section = composition_model->get_composition(position_in_cartesian_coordinates,
673  depth,
674  properties[i_property][1],
675  composition_next_section,
676  starting_depth,
678  distance_from_planes,
679  additional_parameters);
680 
681  WBAssert(!std::isnan(composition_next_section), "Composition_next_section is not a number: " << composition_next_section
682  << ", based on a composition model with the name " << composition_model->get_name() << ", in feature " << this->name);
683  WBAssert(std::isfinite(composition_next_section), "Composition_next_section is not a finite: " << composition_next_section
684  << ", based on a composition model with the name " << composition_model->get_name() << ", in feature " << this->name);
685 
686  }
687 
688  // linear interpolation between current and next section temperatures
689  output[entry_in_output[i_property]] = composition_current_section + section_fraction * (composition_next_section - composition_current_section);
690  break;
691  }
692  case 3: // grains
693  {
694  WorldBuilder::grains grains(output,properties[i_property][2],entry_in_output[i_property]);
695  WorldBuilder::grains grains_current_section = grains;
696  WorldBuilder::grains grains_next_section = grains;
697 
698  for (const auto &grains_model: segment_vector[current_section][current_segment].grains_systems)
699  {
700  grains_current_section = grains_model->get_grains(position_in_cartesian_coordinates,
701  depth,
702  properties[i_property][1],
703  grains_current_section,
704  starting_depth,
706  distance_from_planes,
707  additional_parameters);
708 
709  }
710 
711  for (const auto &grains_model: segment_vector[next_section][current_segment].grains_systems)
712  {
713  grains_next_section = grains_model->get_grains(position_in_cartesian_coordinates,
714  depth,
715  properties[i_property][1],
716  grains_next_section,
717  starting_depth,
719  distance_from_planes,
720  additional_parameters);
721 
722  }
723 
724  // linear interpolation between current and next section temperatures
725  for (size_t i = 0; i < grains.sizes.size(); i++)
726  {
727  grains.sizes[i] = grains_current_section.sizes[i] + section_fraction * (grains_next_section.sizes[i] - grains_current_section.sizes[i]);
728  }
729 
730  // average two rotations matrices through quaternions.
731  for (size_t i = 0; i < grains_current_section.rotation_matrices.size(); i++)
732  {
733  const glm::quaternion::quat quat_current = glm::quaternion::quat_cast(grains_current_section.rotation_matrices[i]);
734  const glm::quaternion::quat quat_next = glm::quaternion::quat_cast(grains_next_section.rotation_matrices[i]);
735 
736  const glm::quaternion::quat quat_average = glm::quaternion::slerp(quat_current,quat_next,section_fraction);
737 
738  grains.rotation_matrices[i] = glm::quaternion::mat3_cast(quat_average);
739  }
740 
741  grains.unroll_into(output,entry_in_output[i_property]);
742  break;
743  }
744  case 4:
745  {
746  output[entry_in_output[i_property]] = static_cast<double>(tag_index);
747  break;
748  }
749  case 5: // velocity
750  {
751  std::array<double,3> velocity_current_section;// = Point<3>(cartesian);
752  velocity_current_section[0] = output[entry_in_output[i_property]];
753  velocity_current_section[1] = output[entry_in_output[i_property]+1];
754  velocity_current_section[2] = output[entry_in_output[i_property]]+2;
755  std::array<double,3> velocity_next_section = velocity_current_section;// output[entry_in_output[i_property]];
756 
757  for (const auto &velocity_model: segment_vector[current_section][current_segment].velocity_systems)
758  {
759  velocity_current_section = velocity_model->get_velocity(position_in_cartesian_coordinates,
760  depth,
761  gravity_norm,
762  velocity_current_section,
763  starting_depth,
765  distance_from_planes,
766  additional_parameters);
767 
768  //WBAssert(!std::isnan(velocity_current_section), "Velocity is not a number: " << velocity_current_section
769  // << ", based on a velocity model with the name " << velocity_model->get_name() << ", in feature " << this->name);
770  //WBAssert(std::isfinite(velocity_current_section), "Velocity is not a finite: " << velocity_current_section
771  // << ", based on a velocity model with the name " << velocity_model->get_name() << ", in feature " << this->name);
772 
773  }
774 
775  for (const auto &velocity_model: segment_vector[next_section][current_segment].velocity_systems)
776  {
777  velocity_next_section = velocity_model->get_velocity(position_in_cartesian_coordinates,
778  depth,
779  gravity_norm,
780  velocity_next_section,
781  starting_depth,
783  distance_from_planes,
784  additional_parameters);
785 
786  //WBAssert(!std::isnan(velocity_next_section), "Velocity is not a number: " << velocity_next_section
787  // << ", based on a velocity model with the name " << velocity_model->get_name() << ", in feature " << this->name);
788  //WBAssert(std::isfinite(velocity_next_section), "Velocity is not a finite: " << velocity_next_section
789  // << ", based on a velocity model with the name " << velocity_model->get_name() << ", in feature " << this->name);
790 
791  }
792 
793  // linear interpolation between current and next section velocitys
794  output[entry_in_output[i_property]] = velocity_current_section[0] + section_fraction * (velocity_next_section[0] - velocity_current_section[0]);
795  output[entry_in_output[i_property]+1] = velocity_current_section[1] + section_fraction * (velocity_next_section[1] - velocity_current_section[1]);
796  output[entry_in_output[i_property]+2] = velocity_current_section[2] + section_fraction * (velocity_next_section[2] - velocity_current_section[2]);
797  break;
798  }
799  default:
800  {
801  WBAssertThrow(false,
802  "Internal error: Unimplemented property provided. " <<
803  "Only temperature (1), composition (2), grains (3), tag (4) or velocity (5) are allowed. "
804  "Provided property number was: " << properties[i_property][0]);
805  }
806  }
807  }
808 
809  }
810  }
811  }
812 
813  }
814 
816  SubductingPlate::distance_to_feature_plane(const Point<3> &position_in_cartesian_coordinates,
817  const Objects::NaturalCoordinate &position_in_natural_coordinates,
818  const double depth) const
819  {
820  // The depth variable is the distance from the surface to the position, the depth
821  // coordinate is the distance from the bottom of the model to the position and
822  // the starting radius is the distance from the bottom of the model to the surface.
823  const double starting_radius = position_in_natural_coordinates.get_depth_coordinate() + depth - starting_depth;
824 
825  WBAssert(std::abs(starting_radius) > std::numeric_limits<double>::epsilon(), "World Builder error: starting_radius can not be zero. "
826  << "Position = " << position_in_cartesian_coordinates[0] << ':' << position_in_cartesian_coordinates[1] << ':' << position_in_cartesian_coordinates[2]
827  << ", natural_coordinate.get_depth_coordinate() = " << position_in_natural_coordinates.get_depth_coordinate()
828  << ", depth = " << depth
829  << ", starting_depth " << starting_depth
830  );
831 
832  const WorldBuilder::Utilities::PointDistanceFromCurvedPlanes distance_from_planes =
833  WorldBuilder::Utilities::distance_point_from_curved_planes(position_in_cartesian_coordinates,
834  position_in_natural_coordinates,
836  coordinates,
839  starting_radius,
841  false,
842  this->bezier_curve);
843 
844  Objects::PlaneDistances plane_distances(distance_from_planes.distance_from_plane, distance_from_planes.distance_along_plane);
845  return plane_distances;
846  }
847 
851  WB_REGISTER_FEATURE(SubductingPlate, subducting plate)
852  } // namespace Features
853 } // namespace WorldBuilder
854 
void parse_entries(Parameters &prm) override final
std::array< double, 2 > get_surface_coordinates() const
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)
Definition: interface.cc:41
std::vector< std::shared_ptr< Features::SubductingPlateModels::Grains::Interface > > default_grains_models
std::vector< std::vector< Point< 2 > > > slab_segment_angles
bool get_shared_pointers(const std::string &name, std::vector< std::shared_ptr< T > > &)
Definition: parameters.cc:1845
std::vector< double > sizes
Definition: grains.h:47
static void make_snippet(Parameters &prm)
std::vector< std::shared_ptr< Features::SubductingPlateModels::Composition::Interface > > default_composition_models
std::vector< std::vector< Point< 2 > > > slab_segment_top_truncation
std::pair< Point< spacedim >, Point< spacedim > > & get_boundary_points()
Definition: bounding_box.h:318
std::vector< std::vector< Objects::Segment< Features::SubductingPlateModels::Temperature::Interface, Features::SubductingPlateModels::Composition::Interface, Features::SubductingPlateModels::Grains::Interface, Features::SubductingPlateModels::Velocity::Interface > > > segment_vector
void extend(const double amount)
Definition: bounding_box.h:355
#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::shared_ptr< Features::SubductingPlateModels::Velocity::Interface > > default_velocity_models
double cos(const double angle)
Definition: point.h:97
static void declare_entries(Parameters &prm, const std::string &parent_name, const std::vector< std::string > &required_entries)
Definition: interface.cc:41
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
PointDistanceFromCurvedPlanes distance_point_from_curved_planes(const Point< 3 > &check_point, const Objects::NaturalCoordinate &natural_coordinate, const Point< 2 > &reference_point, const std::vector< Point< 2 > > &point_list, const std::vector< std::vector< double > > &plane_segment_lengths, const std::vector< std::vector< Point< 2 > > > &plane_segment_angles, const double start_radius, const std::unique_ptr< CoordinateSystems::Interface > &coordinate_system, const bool only_positive, const Objects::BezierCurve &bezier_curve)
Definition: utilities.cc:380
std::vector< Objects::Segment< Features::SubductingPlateModels::Temperature::Interface, Features::SubductingPlateModels::Composition::Interface, Features::SubductingPlateModels::Grains::Interface, Features::SubductingPlateModels::Velocity::Interface > > default_segment_vector
std::vector< std::string > feature_tags
Definition: world.h:308
WorldBuilder::World * world
Definition: interface.h:127
Objects::PlaneDistances distance_to_feature_plane(const Point< 3 > &position_in_cartesian_coordinates, const Objects::NaturalCoordinate &position_in_natural_coordinates, const double depth) const override
#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::unique_ptr< WorldBuilder::CoordinateSystems::Interface > coordinate_system
Definition: parameters.h:268
std::vector< std::vector< Point< 2 > > > slab_segment_thickness
std::vector< std::vector< double > > slab_segment_lengths
static void declare_entries(Parameters &prm, const std::string &parent_name, const std::vector< std::string > &required_entries)
Definition: interface.cc:42
const BoundingBox< 2 > & get_surface_bounding_box() const
static void declare_entries(Parameters &prm, const std::string &parent_name="", const std::vector< std::string > &required_entries={})
std::size_t original_number_of_coordinates
Definition: interface.h:149
std::vector< std::array< std::array< double, 3 >, 3 > > rotation_matrices
Definition: grains.h:51
constexpr double PI
Definition: consts.h:30
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
std::vector< T > get_vector(const std::string &name)
std::vector< std::shared_ptr< Features::SubductingPlateModels::Temperature::Interface > > default_temperature_models
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< Point< 2 > > coordinates
Definition: interface.h:154
bool point_inside(const Point< spacedim > &p, const double tolerance=std::numeric_limits< double >::epsilon()) const
Definition: bounding_box.h:371
static void declare_entries(Parameters &prm, const std::string &parent_name, const std::vector< std::string > &required_entries)
Definition: interface.cc:41