World Builder  1.1.0-pre
A geodynamic initial conditions generator
smooth.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 #include <world_builder/assert.h>
23 #include <world_builder/nan.h>
26 
32 
33 
34 namespace WorldBuilder
35 {
36  using namespace Utilities;
37 
38  namespace Features
39  {
40  namespace FaultModels
41  {
42  namespace Composition
43  {
44  Smooth::Smooth(WorldBuilder::World *world_)
45  :
46  min_distance(NaN::DSNAN),
47  side_distance(NaN::DSNAN),
48  operation(Operations::REPLACE)
49  {
50  this->world = world_;
51  this->name = "smooth";
52  }
53 
55  = default;
56 
57  void
58  Smooth::declare_entries(Parameters &prm, const std::string & /*unused*/)
59  {
60  // Add compositions to the required parameters.
61  prm.declare_entry("", Types::Object({"compositions"}), "Compositional model object");
62 
63  prm.declare_entry("min distance fault center", Types::Double(0),
64  "The distance in meters from which the composition of this feature is present.");
65  prm.declare_entry("side distance fault center", Types::Double(std::numeric_limits<double>::max()),
66  "The distance over which the composition is reduced from 1 to 0.");
67  prm.declare_entry("center fractions", Types::Array(Types::Double(1.0),1),
68  "The composition fraction at the center of the fault.");
69  prm.declare_entry("side fractions", Types::Array(Types::Double(0.0),1),
70  "The composition fraction at the sides of this feature.");
71  prm.declare_entry("compositions", Types::Array(Types::UnsignedInt(),0),
72  "A list with the labels of the composition which are present there.");
73  prm.declare_entry("operation", Types::String("replace", std::vector<std::string> {"replace", "replace defined only", "add", "subtract"}),
74  "Whether the value should replace any value previously defined at this location (replace) or "
75  "add the value to the previously define value. Replacing implies that all compositions not "
76  "explicitly defined are set to zero. To only replace the defined compositions use the replace only defined option.");
77  }
78 
79  void
81  {
82  min_distance = prm.get<double>("min distance fault center");
83  side_distance = prm.get<double>("side distance fault center");
84  WBAssert(side_distance >= min_distance, "distance at the side needs to be larger or equal than the min distance.");
85  operation = string_operations_to_enum(prm.get<std::string>("operation"));
86  center_fraction = prm.get_vector<double>("center fractions");
87  side_fraction = prm.get_vector<double>("side fractions");
88  compositions = prm.get_vector<unsigned int>("compositions");
89  }
90 
91 
92  double
93  Smooth::get_composition( const Point<3> & /*position*/,
94  const double /*depth*/,
95  const unsigned int composition_number,
96  double composition_,
97  const double /*feature_min_depth*/,
98  const double /*feature_max_depth*/,
100  const AdditionalParameters & /*additional_parameters*/) const
101  {
102  double composition = composition_;
103 
104  for (unsigned int i =0; i < compositions.size(); ++i)
105  {
106  if (compositions[i] == composition_number)
107  {
108 
109  // Hyperbolic tangent goes from 0 to 1 over approximately x=(0, 2) without any arguments. The function is written
110  // so that the composition returned 1 to 0 over the side_distance on either sides.
111  composition = (center_fraction[i] - side_fraction[i]) * ( 1 - std::tanh(10 * (distance_from_planes.distance_from_plane
112  - side_distance/2)/side_distance ) )/2 ;
113 
114  return apply_operation(operation,composition_,composition);
115  }
116  }
117 
118  if (operation == Operations::REPLACE)
119  return 0.0;
120 
121  return composition;
122  }
123 
125  } // namespace Composition
126  } // namespace FaultModels
127  } // namespace Features
128 } // namespace WorldBuilder
129 
void parse_entries(Parameters &prm) override final
Definition: smooth.cc:80
static void declare_entries(Parameters &prm, const std::string &parent_name="")
Definition: smooth.cc:58
double get_composition(const Point< 3 > &position, const double depth, const unsigned int composition_number, double composition, 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: smooth.cc:93
const double DSNAN
Definition: nan.h:41
Operations string_operations_to_enum(const std::string &operation)
#define WBAssert(condition, message)
Definition: assert.h:27
#define WB_REGISTER_FEATURE_FAULT_COMPOSITION_MODEL(classname, name)
Definition: interface.h:154
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)
Definition: parameters.cc:197
std::vector< T > get_vector(const std::string &name)
T get(const std::string &name)