OGS
CreateAnchors.cpp File Reference
#include <spdlog/spdlog.h>
#include <tclap/CmdLine.h>
#include <vtkXMLUnstructuredGridReader.h>
#include <vtkXMLUnstructuredGridWriter.h>
#include <fstream>
#include <nlohmann/json.hpp>
#include "BaseLib/FileTools.h"
#include "BaseLib/Logging.h"
#include "ComputeNaturalCoordsAlgorithm.h"
#include "InfoLib/GitInfo.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)
 
AU::ComputeNaturalCoordsResult 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 55 of file CreateAnchors.cpp.

56{
57 std::vector<std::string> required_keys = {"anchor_stiffness",
58 "anchor_cross_sectional_area"};
59 std::vector<std::string> optional_keys = {"maximum_anchor_stress",
60 "initial_anchor_stress",
61 "residual_anchor_stress"};
62 if (number_of_anchors == 0)
63 {
64 OGS_FATAL("No anchors found in the json.");
65 }
66 for (const auto& key : required_keys)
67 {
68 if (!data.contains(key))
69 {
70 OGS_FATAL("JSON file does not contain required key '{:s}'", key);
71 }
72 if (number_of_anchors != data[key].size())
73 {
75 "Number of anchor start points does not match the number of {} "
76 "entries.",
77 key);
78 }
79 for (size_t i = 0; i < number_of_anchors; ++i)
80 {
81 if (!data[key][i].is_number())
82 {
83 OGS_FATAL("Non-numeric element in JSON array for key {}!", key);
84 }
85 }
86 }
87 for (const auto& key : optional_keys)
88 {
89 if (data.contains(key))
90 {
91 if (number_of_anchors != data[key].size())
92 {
94 "Number of anchor start points does not match the number "
95 "of {} "
96 "entries.",
97 key);
98 }
99 for (size_t i = 0; i < number_of_anchors; ++i)
100 {
101 if (!data[key][i].is_number())
102 {
103 OGS_FATAL("Non-numeric element in JSON array for key {}!",
104 key);
105 }
106 }
107 }
108 }
109}
#define OGS_FATAL(...)
Definition Error.h:26
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 213 of file CreateAnchors.cpp.

214{
215 // parse cmdline -----------------------------------------------------------
216
217 TCLAP::CmdLine cmd(
218 "Computes natual coordinates from given real coordinates\n\n"
219
220 "OpenGeoSys-6 software, version " +
222 ".\n"
223 "Copyright (c) 2012-2025, OpenGeoSys Community "
224 "(http://www.opengeosys.org)",
226
227 TCLAP::ValueArg<std::string> log_level_arg(
228 "l", "log-level",
229 "the verbosity of logging messages: none, error, warn, info, debug, "
230 "all",
231 false,
232#ifdef NDEBUG
233 "info",
234#else
235 "all",
236#endif
237 "LOG_LEVEL", cmd);
238
239 TCLAP::ValueArg<std::string> input_filename_arg(
240 "i", "input", "Input bulk mesh", true, "", "VTU_FILE", cmd);
241
242 TCLAP::ValueArg<std::string> output_filename_arg(
243 "o", "output", "Anchor mesh", true, "", "VTU_FILE", cmd);
244
245 TCLAP::ValueArg<std::string> json_filename_arg(
246 "f", "json-file", "JSON file containing anchor start and end points",
247 true, "", "JSON_FILE", cmd);
248
249 TCLAP::ValueArg<double> tolerance_arg(
250 "", "tolerance", "Tolerance/Search length", false,
251 std::numeric_limits<double>::epsilon(), "float", cmd);
252
253 TCLAP::ValueArg<unsigned> max_iter_arg(
254 "", "max-iter",
255 "maximum number of iterations of the internal root-finding algorithm",
256 false, 5, "int", cmd);
257
258 cmd.parse(argc, argv);
259
260 BaseLib::setConsoleLogLevel(log_level_arg.getValue());
261 spdlog::set_pattern("%^%l:%$ %v");
262 spdlog::set_error_handler(
263 [](const std::string& msg)
264 {
265 std::cerr << "spdlog error: " << msg << std::endl;
266 OGS_FATAL("spdlog logger error occurred.");
267 });
268
269 AU::ComputeNaturalCoordsResult const json_data =
270 readJSON(json_filename_arg);
271 auto const bulk_mesh = readGrid(input_filename_arg.getValue());
272
273 // Compute natural coordinates
274 AU::ComputeNaturalCoordsResult const result =
275 AU::computeNaturalCoords(bulk_mesh, json_data, tolerance_arg.getValue(),
276 max_iter_arg.getValue());
277 // Write output
278 auto const output_mesh = AU::toVTKGrid(result);
279 writeGrid(output_mesh, output_filename_arg.getValue());
280
281 if (!result.success)
282 {
283 OGS_FATAL("Failed to compute natural coordinates.");
284 }
285
286 return EXIT_SUCCESS;
287}
AU::ComputeNaturalCoordsResult 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)
void setConsoleLogLevel(std::string const &level_string)
Definition Logging.cpp:37
GITINFOLIB_EXPORT const std::string ogs_version

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

◆ readGrid()

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

Definition at line 25 of file CreateAnchors.cpp.

26{
27 if (!BaseLib::IsFileExisting(input_filename))
28 {
29 OGS_FATAL("'{:s}' file does not exist.", input_filename);
30 }
31
32 vtkNew<vtkXMLUnstructuredGridReader> reader;
33 reader->SetFileName(input_filename.c_str());
34 reader->Update();
35
36 if (!reader->GetOutput())
37 {
38 OGS_FATAL("Could not open file '{}'", input_filename);
39 }
40 if (reader->GetOutput()->GetNumberOfCells() == 0)
41 {
42 OGS_FATAL("Mesh file '{}' contains no cells.", input_filename);
43 }
44 return reader->GetOutput();
45}
bool IsFileExisting(const std::string &strFilename)
Returns true if given file exists.
Definition FileTools.cpp:48

References BaseLib::IsFileExisting(), and OGS_FATAL.

Referenced by main().

◆ readJSON()

AU::ComputeNaturalCoordsResult readJSON ( TCLAP::ValueArg< std::string > const & input_filename)

Definition at line 111 of file CreateAnchors.cpp.

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

Referenced by main().

◆ writeGrid()

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

Definition at line 47 of file CreateAnchors.cpp.

48{
49 vtkNew<vtkXMLUnstructuredGridWriter> writer;
50 writer->SetFileName(output_filename.c_str());
51 writer->SetInputData(grid);
52 writer->Write();
53}

Referenced by main().