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 "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)
 
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 50 of file CreateAnchors.cpp.

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

209{
210 // parse cmdline -----------------------------------------------------------
211
212 TCLAP::CmdLine cmd(
213 "Computes natual coordinates from given real coordinates\n\n"
214
215 "OpenGeoSys-6 software, version " +
217 ".\n"
218 "Copyright (c) 2012-2025, OpenGeoSys Community "
219 "(http://www.opengeosys.org)",
221
222 auto log_level_arg = BaseLib::makeLogLevelArg();
223 cmd.add(log_level_arg);
224
225 TCLAP::ValueArg<std::string> input_filename_arg(
226 "i", "input", "Input (.vtu) bulk mesh file", true, "", "INPUT_FILE",
227 cmd);
228
229 TCLAP::ValueArg<std::string> output_filename_arg(
230 "o", "output", "Output (.vtu). Anchor mesh file", true, "",
231 "OUTPUT_FILE", cmd);
232
233 TCLAP::ValueArg<std::string> json_filename_arg(
234 "f", "json-file",
235 "Input (.json) JSON file containing anchor start and end points", true,
236 "", "INPUT_FILE", cmd);
237
238 TCLAP::ValueArg<double> tolerance_arg(
239 "", "tolerance", "Tolerance/Search length", false,
240 std::numeric_limits<double>::epsilon(), "float", cmd);
241
242 TCLAP::ValueArg<unsigned> max_iter_arg(
243 "", "max-iter",
244 "maximum number of iterations of the internal root-finding "
245 "algorithm, (min = 0)",
246 false, 5, "MAX_ITER", cmd);
247
248 cmd.parse(argc, argv);
249
250 BaseLib::MPI::Setup mpi_setup(argc, argv);
251 BaseLib::initOGSLogger(log_level_arg.getValue());
252
253 AU::ComputeNaturalCoordsResult const json_data =
254 readJSON(json_filename_arg);
255 auto const bulk_mesh = readGrid(input_filename_arg.getValue());
256
257 // Compute natural coordinates
258 AU::ComputeNaturalCoordsResult const result =
259 AU::computeNaturalCoords(bulk_mesh, json_data, tolerance_arg.getValue(),
260 max_iter_arg.getValue());
261 // Write output
262 auto const output_mesh = AU::toVTKGrid(result);
263 writeGrid(output_mesh, output_filename_arg.getValue());
264
265 if (!result.success)
266 {
267 OGS_FATAL("Failed to compute natural coordinates.");
268 }
269
270 return EXIT_SUCCESS;
271}
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)
TCLAP::ValueArg< std::string > makeLogLevelArg()
void initOGSLogger(std::string const &log_level)
Definition Logging.cpp:64
GITINFOLIB_EXPORT const std::string ogs_version

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

◆ readGrid()

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

Definition at line 26 of file CreateAnchors.cpp.

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

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

Referenced by main().

◆ readJSON()

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

Definition at line 106 of file CreateAnchors.cpp.

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

43{
44 vtkNew<vtkXMLUnstructuredGridWriter> writer;
45 writer->SetFileName(output_filename.c_str());
46 writer->SetInputData(grid);
47 writer->Write();
48}

Referenced by main().