OGS
ComputeNaturalCoordsAlgorithm.h
Go to the documentation of this file.
1
10#pragma once
11
12#include <vtkCell.h>
13#include <vtkLine.h>
14#include <vtkPointData.h>
15
16#include <algorithm>
17
20#include "FindCellsForPoints.h"
24
26{
27
29{
30public:
31 void append(Eigen::Vector3d const& natural_coords,
32 Eigen::Vector3d const& real_coords,
33 vtkIdType bulk_element_id,
34 vtkIdType point_cloud_node_id)
35 {
36 natural_coordss.push_back(natural_coords);
37 real_coordss.push_back(real_coords);
38 bulk_element_ids.push_back(bulk_element_id);
39 point_cloud_node_ids.push_back(point_cloud_node_id);
40 }
41
42 void fail() { success = false; }
43
45 Eigen::VectorXd const& initial_anchor_stress,
46 Eigen::VectorXd const& maximum_anchor_stress,
47 Eigen::VectorXd const& residual_anchor_stress,
48 Eigen::VectorXd const& anchor_cross_sectional_area,
49 Eigen::VectorXd const& anchor_stiffness)
50 {
51 auto bei = Eigen::Map<Eigen::VectorX<vtkIdType>>(
52 bulk_element_ids.data(), bulk_element_ids.size());
53 auto pcni = Eigen::Map<Eigen::VectorX<vtkIdType>>(
55 return {conv(natural_coordss),
57 initial_anchor_stress,
58 maximum_anchor_stress,
59 residual_anchor_stress,
60 anchor_cross_sectional_area,
61 anchor_stiffness,
62 bei,
63 pcni,
64 success};
65 }
66
67private:
68 static Eigen::MatrixXd conv(std::vector<Eigen::Vector3d> const& in)
69 {
70 Eigen::MatrixXd out(in.size(), 3);
71 for (std::size_t i = 0; i < in.size(); ++i)
72 {
73 out.row(i).noalias() = in[i].transpose();
74 }
75 return out;
76 }
77
78 std::vector<Eigen::Vector3d> natural_coordss;
79 std::vector<Eigen::Vector3d> real_coordss;
80 std::vector<vtkIdType> bulk_element_ids;
81 std::vector<vtkIdType> point_cloud_node_ids;
82 bool success = true;
83};
84
86 vtkUnstructuredGrid* const bulk_mesh,
87 ComputeNaturalCoordsResult const& json_data, double const tolerance,
88 int const max_iter)
89{
90 // Check inputs ------------------------------------------------------------
91
92 if (json_data.real_coords.cols() != 3)
93 {
94 OGS_FATAL("Wrong number of input coordinates");
95 }
96
97 // General definitions -----------------------------------------------------
98 using SolverRegistry =
101
102 double const real_coords_tolerance = tolerance;
103
104 // Initialization ----------------------------------------------------------
105 FindCellsForPoints findCellsForPoints;
106 findCellsForPoints.initialize(bulk_mesh);
107
109
110 // Do computation for each point ------------------------------------------
111 for (Eigen::Index rc_idx = 0; rc_idx < json_data.real_coords.rows();
112 ++rc_idx)
113 {
114 Eigen::RowVector3d const real_coords_single_point =
115 json_data.real_coords.row(rc_idx);
116
117 // Find cells for point
118 auto const filtered_bulk_mesh =
119 findCellsForPoints.find(real_coords_single_point, tolerance);
120
121 // Check multiplicity
122 int const actual_multiplicity = filtered_bulk_mesh->GetNumberOfCells();
123 if (actual_multiplicity == 0)
124 {
125 OGS_FATAL(
126 "Did not find any cell for the point {} (coordinates: {})",
127 rc_idx, real_coords_single_point);
128 }
129 assert(actual_multiplicity > 0);
130 if (actual_multiplicity != 1)
131 {
132 INFO(
133 "Found more then one cell (namely {}) for point #{} "
134 "(coordinates: {})",
135 actual_multiplicity, rc_idx, real_coords_single_point);
136 }
137
138 // Convert to OGS
139 std::unique_ptr<MeshLib::Mesh> ogs_mesh(
141 filtered_bulk_mesh, "filtered_bulk_mesh"));
142
143 auto const* bulk_element_ids = MeshLib::bulkElementIDs(*ogs_mesh);
144
145 auto const& elements = ogs_mesh->getElements();
146
147 // Compute natural coordinates
148 // Found a solution, no need to test other elements ignoring
149 // multiplicity.
150 constexpr std::size_t elt_idx = 0;
151 auto const& element = *elements[elt_idx];
152 auto const bulk_element_id = bulk_element_ids->getComponent(elt_idx, 0);
153
154 auto const& solver = SolverRegistry::getFor(element);
155 auto const opt_natural_coords = solver.solve(
156 element, real_coords_single_point, max_iter, real_coords_tolerance);
157
158 result.append(opt_natural_coords, real_coords_single_point,
159 bulk_element_id, rc_idx);
160 }
161
162 return result.finished(
163 json_data.initial_anchor_stress, json_data.maximum_anchor_stress,
165 json_data.anchor_stiffness);
166}
167
168template <typename T>
169void addPointData(vtkUnstructuredGrid* grid, std::string const& name,
170 Eigen::MatrixX<T> const& data)
171{
172 INFO("converting {}", name);
173 vtkNew<vtkAOSDataArrayTemplate<T>> array;
174 array->SetName(name.c_str());
175 array->SetNumberOfComponents(data.cols());
176 array->SetNumberOfTuples(grid->GetNumberOfPoints());
177
178 if (grid->GetNumberOfPoints() != data.rows())
179 {
180 OGS_FATAL(
181 "Got {} rows in the table but expected {} rows, same as number of "
182 "points in the grid.",
183 data.rows(), grid->GetNumberOfPoints());
184 }
185 for (Eigen::Index i = 0; i < data.rows(); ++i)
186 {
187 // copy to contiguous storage
188 Eigen::RowVectorX<T> const row = data.row(i);
189
190 array->SetTypedTuple(i, row.data());
191 DBUG("> {}", row);
192 }
193
194 grid->GetPointData()->AddArray(array);
195}
196
197template <typename T>
198void addCellData(vtkUnstructuredGrid* grid, std::string const& name,
199 Eigen::MatrixX<T> const& data)
200{
201 INFO("converting {}", name);
202 vtkNew<vtkAOSDataArrayTemplate<T>> array;
203 array->SetName(name.c_str());
204 array->SetNumberOfComponents(data.cols());
205 array->SetNumberOfTuples(grid->GetNumberOfCells());
206
207 if (grid->GetNumberOfCells() != data.rows())
208 {
209 OGS_FATAL(
210 "Got {} rows in the table but expected {} rows, same as number of "
211 "cells in the grid.",
212 data.rows(), grid->GetNumberOfCells());
213 }
214 for (Eigen::Index i = 0; i < data.rows(); ++i)
215 {
216 // copy to contiguous storage
217 Eigen::RowVectorX<T> const row = data.row(i);
218
219 array->SetTypedTuple(i, row.data());
220 DBUG("> {}", row);
221 }
222
223 grid->GetCellData()->AddArray(array);
224}
225
232vtkSmartPointer<vtkUnstructuredGrid> toVTKGrid(
233 ComputeNaturalCoordsResult const& result)
234{
235 INFO("converting points");
236 auto const n_pts = result.natural_coords.rows();
237 if (n_pts == 0)
238 {
239 OGS_FATAL(
240 "No point was found in the mesh. Please check whether all anchors "
241 "are within the model region.");
242 }
243 else if (n_pts == 1)
244 {
245 OGS_FATAL(
246 "Only one point was found in the mesh. Please check whether all "
247 "anchors are within the model region.");
248 }
249 else if (n_pts % 2 != 0)
250 {
251 // number should be even if multiplicity is ignored
252 OGS_FATAL(
253 "Number of points is not even. Anchors are believed to consist of "
254 "a start and an end point.");
255 }
256
257 vtkNew<vtkUnstructuredGrid> grid;
258
259 // Points ------------------------------------------------------------------
260 vtkNew<vtkPoints> points;
261 points->SetDataTypeToDouble();
262 points->SetNumberOfPoints(n_pts);
263
264 auto const& css = result.real_coords;
265 for (Eigen::Index i = 0; i < css.rows(); ++i)
266 {
267 points->SetPoint(i, css(i, 0), css(i, 1), css(i, 2));
268
269 if ((i % 2) == 1)
270 {
271 vtkNew<vtkLine> line;
272 line->GetPointIds()->SetId(0, i - 1);
273 line->GetPointIds()->SetId(1, i);
274 grid->InsertNextCell(line->GetCellType(), line->GetPointIds());
275 }
276 }
277
278 grid->SetPoints(points);
279
280 // Point data --------------------------------------------------------------
281 addPointData<double>(grid, "natural_coordinates", result.natural_coords);
282 addPointData<std::size_t>(grid, "bulk_element_ids",
283 result.bulk_element_ids.cast<std::size_t>());
284 addPointData<std::size_t>(grid, "point_cloud_node_ids",
285 result.point_cloud_node_ids.cast<std::size_t>());
286 addCellData<double>(grid, "initial_anchor_stress",
287 result.initial_anchor_stress);
288 addCellData<double>(grid, "maximum_anchor_stress",
289 result.maximum_anchor_stress);
290 addCellData<double>(grid, "residual_anchor_stress",
292 addCellData<double>(grid, "anchor_cross_sectional_area",
294 addCellData<double>(grid, "anchor_stiffness", result.anchor_stiffness);
295 return grid;
296}
297} // namespace ApplicationUtils
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:36
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Definition of the VtkMeshConverter class.
ComputeNaturalCoordsResult finished(Eigen::VectorXd const &initial_anchor_stress, Eigen::VectorXd const &maximum_anchor_stress, Eigen::VectorXd const &residual_anchor_stress, Eigen::VectorXd const &anchor_cross_sectional_area, Eigen::VectorXd const &anchor_stiffness)
static Eigen::MatrixXd conv(std::vector< Eigen::Vector3d > const &in)
void append(Eigen::Vector3d const &natural_coords, Eigen::Vector3d const &real_coords, vtkIdType bulk_element_id, vtkIdType point_cloud_node_id)
vtkSmartPointer< vtkUnstructuredGrid > find(Eigen::Vector3d const &coords, double const tolerance)
void initialize(vtkUnstructuredGrid *const bulk_mesh)
static MeshLib::Mesh * convertUnstructuredGrid(vtkUnstructuredGrid *grid, bool const compute_element_neighbors=false, std::string const &mesh_name="vtkUnstructuredGrid")
Converts a vtkUnstructuredGrid object to a Mesh.
void addPointData(vtkUnstructuredGrid *grid, std::string const &name, Eigen::MatrixX< T > const &data)
void addCellData(vtkUnstructuredGrid *grid, std::string const &name, Eigen::MatrixX< T > const &data)
vtkSmartPointer< vtkUnstructuredGrid > toVTKGrid(ComputeNaturalCoordsResult const &result)
ComputeNaturalCoordsResult computeNaturalCoords(vtkUnstructuredGrid *const bulk_mesh, ComputeNaturalCoordsResult const &json_data, double const tolerance, int const max_iter)
PropertyVector< std::size_t > const * bulkElementIDs(Mesh const &mesh)
Definition Mesh.cpp:301