World Builder  1.1.0-pre
A geodynamic initial conditions generator
uniform.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 
30 
31 namespace WorldBuilder
32 {
33 
34  using namespace Utilities;
35 
36  namespace Features
37  {
38  namespace FaultModels
39  {
40  namespace Grains
41  {
42  Uniform::Uniform(WorldBuilder::World *world_)
43  :
44  min_depth(NaN::DSNAN),
45  max_depth(NaN::DSNAN)
46 
47  {
48  this->world = world_;
49  this->name = "uniform";
50  }
51 
53  = default;
54 
55  void
56  Uniform::declare_entries(Parameters &prm, const std::string & /*unused*/)
57  {
58  // Document plugin and require entries if needed.
59  // Add compositions the required parameters.
60  prm.declare_entry("", Types::Object({"compositions"}),
61  "Uniform grains model. All grains start exactly the same.");
62 
63  // Declare entries of this plugin
64  prm.declare_entry("min distance fault center", Types::Double(0),
65  "The distance from the fault center in meters from which the composition of this feature is present.");
66  prm.declare_entry("max distance fault center", Types::Double(std::numeric_limits<double>::max()),
67  "The distance from the fault in meters to which the composition of this feature is present.");
68 
69  prm.declare_entry("compositions", Types::Array(Types::UnsignedInt(),0),
70  "A list with the integer labels of the composition which are present there.");
71 
72  prm.declare_entry("rotation matrices", Types::Array(Types::Array(Types::Array(Types::Double(0),3,3),3,3),0),
73  "A list with the labels of the grains which are present there for each compositions.");
74 
75  prm.declare_entry("Euler angles z-x-z", Types::Array(Types::Array(Types::Double(0),3,3),0),
76  "A list with the z-x-z Euler angles of the grains which are present there for each compositions.");
77 
78  prm.declare_entry("orientation operation", Types::String("replace", std::vector<std::string> {"replace"}),
79  "Whether the value should replace any value previously defined at this location (replace) or "
80  "add the value to the previously define value (add, not implemented). Replacing implies that all values not "
81  "explicitly defined are set to zero.");
82 
83  prm.declare_entry("grain sizes",
85  "A list of the size of all of the grains in each composition. If set to <0, the size will be set so that the total is equal to 1.");
86 
87 
88  }
89 
90  void
92  {
93  min_depth = prm.get<double>("min distance fault center");
94  max_depth = prm.get<double>("max distance fault center");
95  compositions = prm.get_vector<unsigned int>("compositions");
96 
97  const bool set_euler_angles = prm.check_entry("Euler angles z-x-z");
98  const bool set_rotation_matrices = prm.check_entry("rotation matrices");
99 
100  WBAssertThrow(!(set_euler_angles == true && set_rotation_matrices == true),
101  "Only Euler angles or Rotation matrices may be set, but both are set for " << prm.get_full_json_path());
102 
103 
104  WBAssertThrow(!(set_euler_angles == false && set_rotation_matrices == false),
105  "Euler angles or Rotation matrices have to be set, but neither are set for " << prm.get_full_json_path());
106 
107  if (set_euler_angles)
108  {
109  std::vector<std::array<double,3> > euler_angles_vector = prm.get_vector<std::array<double,3> >("Euler angles z-x-z");
110  rotation_matrices.resize(euler_angles_vector.size());
111  for (size_t i = 0; i<euler_angles_vector.size(); ++i)
112  {
113  rotation_matrices[i] = WorldBuilder::Utilities::euler_angles_to_rotation_matrix(euler_angles_vector[i][0],euler_angles_vector[i][1],euler_angles_vector[i][2]);
114  }
115 
116  }
117  else
118  {
119  rotation_matrices = prm.get_vector<std::array<std::array<double,3>,3> >("rotation matrices");
120  }
121 
122  operation = prm.get<std::string>("orientation operation");
123  grain_sizes = prm.get_vector<double>("grain sizes");
124 
125 
126  WBAssertThrow(compositions.size() == rotation_matrices.size(),
127  "There are not the same amount of compositions (" << compositions.size()
128  << ") and rotation_matrices (" << rotation_matrices.size() << ").");
129  WBAssertThrow(compositions.size() == grain_sizes.size(),
130  "There are not the same amount of compositions (" << compositions.size()
131  << ") and grain_sizes (" << grain_sizes.size() << ").");
132  }
133 
134 
136  Uniform::get_grains(const Point<3> & /*position_in_cartesian_coordinates*/,
137  const double /*depth*/,
138  const unsigned int composition_number,
139  WorldBuilder::grains grains_,
140  const double /*feature_min_depth*/,
141  const double /*feature_max_depth*/,
142  const WorldBuilder::Utilities::PointDistanceFromCurvedPlanes &distance_from_planes,
143  const AdditionalParameters & /*additional_parameters*/) const
144  {
145  WorldBuilder::grains grains_local = grains_;
146  if (std::fabs(distance_from_planes.distance_from_plane) <= max_depth && std::fabs(distance_from_planes.distance_from_plane) >= min_depth)
147  {
148  for (unsigned int i =0; i < compositions.size(); ++i)
149  {
150  if (compositions[i] == composition_number)
151  {
152  std::fill(grains_local.rotation_matrices.begin(),grains_local.rotation_matrices.end(),rotation_matrices[i]);
153 
154  const double size = grain_sizes[i] < 0 ? 1.0/static_cast<double>(grains_local.sizes.size()) : grain_sizes[i];
155  std::fill(grains_local.sizes.begin(),grains_local.sizes.end(),size);
156 
157  return grains_local;
158  }
159  }
160  }
161  return grains_local;
162  }
164  } // namespace Grains
165  } // namespace FaultModels
166  } // namespace Features
167 } // namespace WorldBuilder
168 
169 
bool check_entry(const std::string &name) const
Definition: parameters.cc:205
std::vector< std::array< std::array< double, 3 >, 3 > > rotation_matrices
Definition: uniform.h:104
const double DSNAN
Definition: nan.h:41
std::vector< double > sizes
Definition: grains.h:47
std::string get_full_json_path(size_t max_size=std::numeric_limits< size_t >::max()) const
Definition: parameters.cc:1933
#define WB_REGISTER_FEATURE_FAULT_GRAINS_MODEL(classname, name)
Definition: interface.h:155
static void declare_entries(Parameters &prm, const std::string &parent_name="")
Definition: uniform.cc:56
#define WBAssertThrow(condition, message)
Definition: assert.h:40
void parse_entries(Parameters &prm) override final
Definition: uniform.cc:91
WorldBuilder::grains get_grains(const Point< 3 > &position, const double depth, const unsigned int composition_number, WorldBuilder::grains grains, const double feature_min_depth, const double feature_max_depth, const WorldBuilder::Utilities::PointDistanceFromCurvedPlanes &distance_from_planes, const AdditionalParameters &additional_parameters) const override final
Definition: uniform.cc:136
std::array< std::array< double, 3 >, 3 > euler_angles_to_rotation_matrix(double phi1_d, double theta_d, double phi2_d)
Definition: utilities.cc:1147
std::vector< std::array< std::array< double, 3 >, 3 > > rotation_matrices
Definition: grains.h:51
void declare_entry(const std::string &name, const Types::Interface &type, const std::string &documentation)
Definition: parameters.cc:197
std::vector< T > get_vector(const std::string &name)
T get(const std::string &name)