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 51 of file CreateAnchors.cpp.

52{
53 std::vector<std::string> required_keys = {"anchor_stiffness",
54 "anchor_cross_sectional_area"};
55 std::vector<std::string> optional_keys = {
56 "maximum_anchor_stress", "initial_anchor_stress",
57 "residual_anchor_stress", "free_fraction"};
58 if (number_of_anchors == 0)
59 {
60 OGS_FATAL("No anchors found in the json.");
61 }
62 for (const auto& key : required_keys)
63 {
64 if (!data.contains(key))
65 {
66 OGS_FATAL("JSON file does not contain required key '{:s}'", key);
67 }
68 if (number_of_anchors != data[key].size())
69 {
71 "Number of anchor start points does not match the number of {} "
72 "entries.",
73 key);
74 }
75 for (size_t i = 0; i < number_of_anchors; ++i)
76 {
77 if (!data[key][i].is_number())
78 {
79 OGS_FATAL("Non-numeric element in JSON array for key {}!", key);
80 }
81 }
82 }
83 for (const auto& key : optional_keys)
84 {
85 if (data.contains(key))
86 {
87 if (number_of_anchors != data[key].size())
88 {
90 "Number of anchor start points does not match the number "
91 "of {} "
92 "entries.",
93 key);
94 }
95 for (size_t i = 0; i < number_of_anchors; ++i)
96 {
97 if (!data[key][i].is_number())
98 {
99 OGS_FATAL("Non-numeric element in JSON array for key {}!",
100 key);
101 }
102 if (key == "free_fraction")
103 {
104 double const free_fraction_number =
105 data["free_fraction"][i].get<double>();
106 if (free_fraction_number < 0.0 ||
107 free_fraction_number > 1.0)
108 {
109 OGS_FATAL(
110 "Free fraction must be in the range [0, 1] but is "
111 "{} for anchor #{}.",
112 free_fraction_number, i);
113 }
114 }
115 }
116 }
117 }
118}
#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 239 of file CreateAnchors.cpp.

240{
241 // parse cmdline -----------------------------------------------------------
242
243 TCLAP::CmdLine cmd(
244 "Computes natual coordinates from given real coordinates\n\n"
245
246 "OpenGeoSys-6 software, version " +
248 ".\n"
249 "Copyright (c) 2012-2025, OpenGeoSys Community "
250 "(http://www.opengeosys.org)",
252
253 auto log_level_arg = BaseLib::makeLogLevelArg();
254 cmd.add(log_level_arg);
255
256 TCLAP::ValueArg<std::string> input_filename_arg(
257 "i", "input", "Input (.vtu) bulk mesh file", true, "", "INPUT_FILE",
258 cmd);
259
260 TCLAP::ValueArg<std::string> output_filename_arg(
261 "o", "output", "Output (.vtu). Anchor mesh file", true, "",
262 "OUTPUT_FILE", cmd);
263
264 TCLAP::ValueArg<std::string> json_filename_arg(
265 "f", "json-file",
266 "Input (.json) JSON file containing anchor start and end points", true,
267 "", "INPUT_FILE", cmd);
268
269 TCLAP::ValueArg<double> tolerance_arg(
270 "", "tolerance", "Tolerance/Search length", false,
271 std::numeric_limits<double>::epsilon(), "float", cmd);
272
273 TCLAP::ValueArg<unsigned> max_iter_arg(
274 "", "max-iter",
275 "maximum number of iterations of the internal root-finding "
276 "algorithm, (min = 0)",
277 false, 5, "MAX_ITER", cmd);
278
279 cmd.parse(argc, argv);
280
281 BaseLib::MPI::Setup mpi_setup(argc, argv);
282 BaseLib::initOGSLogger(log_level_arg.getValue());
283
284 auto const& [json_data, free_fraction] = readJSON(json_filename_arg);
285 auto const bulk_mesh = readGrid(input_filename_arg.getValue());
286
287 auto split_anchor_coords =
288 getOrderedAnchorCoords(bulk_mesh, json_data.real_coords, free_fraction,
289 tolerance_arg.getValue());
290
291 AU::ComputeNaturalCoordsResult const anchor_data =
293 json_data);
295
296 // Compute natural coordinates
298 bulk_mesh, anchor_data, tolerance_arg.getValue(),
299 max_iter_arg.getValue());
300 // Write output
301 auto const output_mesh = AU::toVTKGrid(result);
302 writeGrid(output_mesh, output_filename_arg.getValue());
303
304 if (!result.success)
305 {
306 OGS_FATAL("Failed to compute natural coordinates.");
307 }
308
309 return EXIT_SUCCESS;
310}
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:64
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 27 of file CreateAnchors.cpp.

28{
29 vtkSmartPointer<vtkUnstructuredGrid> grid =
31 input_filename);
32 if (grid == nullptr)
33 {
34 OGS_FATAL("Could not open file '{}'", input_filename);
35 }
36 if (grid->GetNumberOfCells() == 0)
37 {
38 OGS_FATAL("Mesh file '{}' contains no cells.", input_filename);
39 }
40 return grid;
41}
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 123 of file CreateAnchors.cpp.

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

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

Referenced by main().