OGS
CreateAnchors.cpp File Reference
#include <tclap/CmdLine.h>
#include <vtkXMLUnstructuredGridWriter.h>
#include <fstream>
#include <nlohmann/json.hpp>
#include "BaseLib/FileTools.h"
#include "BaseLib/Logging.h"
#include "BaseLib/MPI.h"
#include "BaseLib/TCLAPArguments.h"
#include "ComputeIntersections.h"
#include "ComputeNaturalCoordsAlgorithm.h"
#include "InfoLib/GitInfo.h"
#include "MeshLib/IO/VtkIO/VtuInterface.h"
Include dependency graph for CreateAnchors.cpp:

Go to the source code of this file.

Functions

vtkSmartPointer< vtkUnstructuredGrid > readGrid (std::string const &input_filename)
void writeGrid (vtkUnstructuredGrid *grid, std::string const &output_filename)
void checkJSONEntries (nlohmann::json const &data, size_t number_of_anchors)
std::tuple< AU::ComputeNaturalCoordsResult, Eigen::VectorXd > readJSON (TCLAP::ValueArg< std::string > const &input_filename)
int main (int argc, char **argv)

Function Documentation

◆ checkJSONEntries()

void checkJSONEntries ( nlohmann::json const & data,
size_t number_of_anchors )

Definition at line 45 of file CreateAnchors.cpp.

46{
47 std::vector<std::string> required_keys = {"anchor_stiffness",
48 "anchor_cross_sectional_area"};
49 std::vector<std::string> optional_keys = {
50 "maximum_anchor_stress", "initial_anchor_stress",
51 "residual_anchor_stress", "free_fraction"};
52 if (number_of_anchors == 0)
53 {
54 OGS_FATAL("No anchors found in the json.");
55 }
56 for (const auto& key : required_keys)
57 {
58 if (!data.contains(key))
59 {
60 OGS_FATAL("JSON file does not contain required key '{:s}'", key);
61 }
62 if (number_of_anchors != data[key].size())
63 {
65 "Number of anchor start points does not match the number of {} "
66 "entries.",
67 key);
68 }
69 for (size_t i = 0; i < number_of_anchors; ++i)
70 {
71 if (!data[key][i].is_number())
72 {
73 OGS_FATAL("Non-numeric element in JSON array for key {}!", key);
74 }
75 }
76 }
77 for (const auto& key : optional_keys)
78 {
79 if (data.contains(key))
80 {
81 if (number_of_anchors != data[key].size())
82 {
84 "Number of anchor start points does not match the number "
85 "of {} "
86 "entries.",
87 key);
88 }
89 for (size_t i = 0; i < number_of_anchors; ++i)
90 {
91 if (!data[key][i].is_number())
92 {
93 OGS_FATAL("Non-numeric element in JSON array for key {}!",
94 key);
95 }
96 if (key == "free_fraction")
97 {
98 double const free_fraction_number =
99 data["free_fraction"][i].get<double>();
100 if (free_fraction_number < 0.0 ||
101 free_fraction_number > 1.0)
102 {
103 OGS_FATAL(
104 "Free fraction must be in the range [0, 1] but is "
105 "{} for anchor #{}.",
106 free_fraction_number, i);
107 }
108 }
109 }
110 }
111 }
112}
#define OGS_FATAL(...)
Definition Error.h:19
constexpr int size(int const displacement_dim)
Vectorized tensor size for given displacement dimension.

References OGS_FATAL.

Referenced by readJSON().

◆ main()

int main ( int argc,
char ** argv )

Definition at line 233 of file CreateAnchors.cpp.

234{
235 // parse cmdline -----------------------------------------------------------
236
237 TCLAP::CmdLine cmd(
238 "Computes natual coordinates from given real coordinates\n\n"
239
240 "OpenGeoSys-6 software, version " +
242 ".\n"
243 "Copyright (c) 2012-2026, OpenGeoSys Community "
244 "(http://www.opengeosys.org)",
246
247 auto log_level_arg = BaseLib::makeLogLevelArg();
248 cmd.add(log_level_arg);
249
250 TCLAP::ValueArg<std::string> input_filename_arg(
251 "i", "input", "Input (.vtu) bulk mesh file", true, "", "INPUT_FILE",
252 cmd);
253
254 TCLAP::ValueArg<std::string> output_filename_arg(
255 "o", "output", "Output (.vtu). Anchor mesh file", true, "",
256 "OUTPUT_FILE", cmd);
257
258 TCLAP::ValueArg<std::string> json_filename_arg(
259 "f", "json-file",
260 "Input (.json) JSON file containing anchor start and end points", true,
261 "", "INPUT_FILE", cmd);
262
263 TCLAP::ValueArg<double> tolerance_arg(
264 "", "tolerance", "Tolerance/Search length", false,
265 std::numeric_limits<double>::epsilon(), "float", cmd);
266
267 TCLAP::ValueArg<unsigned> max_iter_arg(
268 "", "max-iter",
269 "maximum number of iterations of the internal root-finding "
270 "algorithm, (min = 0)",
271 false, 5, "MAX_ITER", cmd);
272
273 cmd.parse(argc, argv);
274
275 BaseLib::MPI::Setup mpi_setup(argc, argv);
276 BaseLib::initOGSLogger(log_level_arg.getValue());
277
278 auto const& [json_data, free_fraction] = readJSON(json_filename_arg);
279 auto const bulk_mesh = readGrid(input_filename_arg.getValue());
280
281 auto split_anchor_coords =
282 getOrderedAnchorCoords(bulk_mesh, json_data.real_coords, free_fraction,
283 tolerance_arg.getValue());
284
285 AU::ComputeNaturalCoordsResult const anchor_data =
287 json_data);
289
290 // Compute natural coordinates
292 bulk_mesh, anchor_data, tolerance_arg.getValue(),
293 max_iter_arg.getValue());
294 // Write output
295 auto const output_mesh = AU::toVTKGrid(result);
296 writeGrid(output_mesh, output_filename_arg.getValue());
297
298 if (!result.success)
299 {
300 OGS_FATAL("Failed to compute natural coordinates.");
301 }
302
303 return EXIT_SUCCESS;
304}
std::vector< std::vector< IntersectionResult > > getOrderedAnchorCoords(vtkUnstructuredGrid *grid, Eigen::MatrixX3d const &realcoords, Eigen::VectorXd const &free_fraction, double const tol)
Finds intersection points of a line segment with the cells of a vtkUnstructuredGrid....
AU::ComputeNaturalCoordsResult setPhysicalPropertiesForIntersectionPoints(std::vector< std::vector< IntersectionResult > > const &anchor_coords, AU::ComputeNaturalCoordsResult const &original_anchor_data)
fills the physical properties of the intersection points based on the original anchor data.
std::tuple< AU::ComputeNaturalCoordsResult, Eigen::VectorXd > readJSON(TCLAP::ValueArg< std::string > const &input_filename)
vtkSmartPointer< vtkUnstructuredGrid > readGrid(std::string const &input_filename)
void writeGrid(vtkUnstructuredGrid *grid, std::string const &output_filename)
vtkSmartPointer< vtkUnstructuredGrid > toVTKGrid(ComputeNaturalCoordsResult const &result)
ComputeNaturalCoordsResult computeNaturalCoords(vtkUnstructuredGrid *const bulk_mesh, ComputeNaturalCoordsResult const &json_data, double const tolerance, int const max_iter)
TCLAP::ValueArg< std::string > makeLogLevelArg()
void initOGSLogger(std::string const &log_level)
Definition Logging.cpp:56
GITINFOLIB_EXPORT const std::string ogs_version

References ApplicationUtils::computeNaturalCoords(), getOrderedAnchorCoords(), BaseLib::initOGSLogger(), BaseLib::makeLogLevelArg(), OGS_FATAL, GitInfoLib::GitInfo::ogs_version, readGrid(), readJSON(), setPhysicalPropertiesForIntersectionPoints(), ApplicationUtils::ComputeNaturalCoordsResult::success, ApplicationUtils::toVTKGrid(), and writeGrid().

◆ readGrid()

vtkSmartPointer< vtkUnstructuredGrid > readGrid ( std::string const & input_filename)

Definition at line 21 of file CreateAnchors.cpp.

22{
23 vtkSmartPointer<vtkUnstructuredGrid> grid =
25 input_filename);
26 if (grid == nullptr)
27 {
28 OGS_FATAL("Could not open file '{}'", input_filename);
29 }
30 if (grid->GetNumberOfCells() == 0)
31 {
32 OGS_FATAL("Mesh file '{}' contains no cells.", input_filename);
33 }
34 return grid;
35}
static vtkSmartPointer< vtkUnstructuredGrid > readVtuFileToVtkUnstructuredGrid(std::string const &file_name)

References OGS_FATAL, and MeshLib::IO::VtuInterface::readVtuFileToVtkUnstructuredGrid().

Referenced by main().

◆ readJSON()

std::tuple< AU::ComputeNaturalCoordsResult, Eigen::VectorXd > readJSON ( TCLAP::ValueArg< std::string > const & input_filename)

Reads a JSON file containing anchor start and end points and physical properties of the anchors. Also returns the free fraction of each anchor. If the free fraction is not given, it is set to 1.0 for all anchors.

Definition at line 117 of file CreateAnchors.cpp.

119{
120 using json = nlohmann::json;
121 std::ifstream f(input_filename.getValue());
122 if (!f)
123 {
124 OGS_FATAL("Could not open file '{:s}'", input_filename.getValue());
125 }
126 json const data = json::parse(f);
127
128 auto const number_of_anchors = data["anchor_start_points"].size();
129 if (number_of_anchors != data["anchor_end_points"].size())
130 {
131 OGS_FATAL(
132 "Number of anchor start points does not match the number of "
133 "anchor end points.");
134 }
135 checkJSONEntries(data, number_of_anchors);
136
137 Eigen::MatrixX3d realcoords(2 * number_of_anchors, 3);
138
139 auto get_coordinates = [&data](const char* key, std::size_t const i)
140 {
141 const auto& v = data[key][i].get_ref<json::array_t const&>();
142 if (v.size() != 3)
143 {
144 OGS_FATAL(
145 "Expected a vector of length three for {:s} {}. Got vector of "
146 "size {}.",
147 key, i, v.size());
148 }
149 return Eigen::RowVector3d{v[0].get<double>(), v[1].get<double>(),
150 v[2].get<double>()};
151 };
152
153 for (std::size_t i = 0; i < number_of_anchors; i++)
154 {
155 realcoords.row(2 * i).noalias() =
156 get_coordinates("anchor_start_points", i);
157 realcoords.row(2 * i + 1).noalias() =
158 get_coordinates("anchor_end_points", i);
159 }
160 Eigen::VectorXd initial_anchor_stress(number_of_anchors);
161 Eigen::VectorXd maximum_anchor_stress(number_of_anchors);
162 Eigen::VectorXd residual_anchor_stress(number_of_anchors);
163 Eigen::VectorXd anchor_cross_sectional_area(number_of_anchors);
164 Eigen::VectorXd anchor_stiffness(number_of_anchors);
165 Eigen::VectorXd free_fraction(number_of_anchors);
166 if (!data.contains("initial_anchor_stress"))
167 {
168 initial_anchor_stress.setZero();
169 }
170 else
171 {
172 for (size_t i = 0; i < number_of_anchors; ++i)
173 {
174 initial_anchor_stress(i) =
175 data["initial_anchor_stress"][i].get<double>();
176 }
177 }
178 if (!data.contains("maximum_anchor_stress"))
179 {
180 maximum_anchor_stress = Eigen::VectorXd::Constant(
181 number_of_anchors, std::numeric_limits<double>::max());
182 }
183 else
184 {
185 for (size_t i = 0; i < number_of_anchors; ++i)
186 {
187 maximum_anchor_stress(i) =
188 data["maximum_anchor_stress"][i].get<double>();
189 }
190 }
191 if (!data.contains("residual_anchor_stress"))
192 {
193 residual_anchor_stress = Eigen::VectorXd::Constant(
194 number_of_anchors, std::numeric_limits<double>::max());
195 }
196 else
197 {
198 for (size_t i = 0; i < number_of_anchors; ++i)
199 {
200 residual_anchor_stress(i) =
201 data["residual_anchor_stress"][i].get<double>();
202 }
203 }
204 if (!data.contains("free_fraction"))
205 {
206 free_fraction = Eigen::VectorXd::Constant(number_of_anchors, 1.0);
207 }
208 else
209 {
210 for (size_t i = 0; i < number_of_anchors; ++i)
211 {
212 free_fraction(i) = data["free_fraction"][i].get<double>();
213 }
214 }
215
216 for (size_t i = 0; i < number_of_anchors; ++i)
217 {
218 anchor_cross_sectional_area(i) =
219 data["anchor_cross_sectional_area"][i].get<double>();
220 anchor_stiffness(i) = data["anchor_stiffness"][i].get<double>();
221 }
223 json_data.real_coords = realcoords;
224 json_data.initial_anchor_stress = initial_anchor_stress;
225 json_data.maximum_anchor_stress = maximum_anchor_stress;
226 json_data.residual_anchor_stress = residual_anchor_stress;
227 json_data.anchor_cross_sectional_area = anchor_cross_sectional_area;
228 json_data.anchor_stiffness = anchor_stiffness;
229 json_data.success = true;
230 return {json_data, free_fraction};
231}
void checkJSONEntries(nlohmann::json const &data, size_t number_of_anchors)

References ApplicationUtils::ComputeNaturalCoordsResult::anchor_cross_sectional_area, ApplicationUtils::ComputeNaturalCoordsResult::anchor_stiffness, checkJSONEntries(), ApplicationUtils::ComputeNaturalCoordsResult::initial_anchor_stress, ApplicationUtils::ComputeNaturalCoordsResult::maximum_anchor_stress, OGS_FATAL, ApplicationUtils::ComputeNaturalCoordsResult::real_coords, ApplicationUtils::ComputeNaturalCoordsResult::residual_anchor_stress, and ApplicationUtils::ComputeNaturalCoordsResult::success.

Referenced by main().

◆ writeGrid()

void writeGrid ( vtkUnstructuredGrid * grid,
std::string const & output_filename )

Definition at line 37 of file CreateAnchors.cpp.

38{
39 vtkNew<vtkXMLUnstructuredGridWriter> writer;
40 writer->SetFileName(output_filename.c_str());
41 writer->SetInputData(grid);
42 writer->Write();
43}

Referenced by main().