OGS
ComputeIntersections.cpp File Reference
#include "ComputeIntersections.h"
#include <spdlog/spdlog.h>
#include <vtkCellData.h>
#include <vtkCellLocator.h>
#include <vtkExtractEdges.h>
#include <vtkIdTypeArray.h>
#include <vtkLine.h>
#include "BaseLib/Logging.h"
Include dependency graph for ComputeIntersections.cpp:

Go to the source code of this file.

Functions

int GetGridDimension (vtkUnstructuredGrid *grid)
bool isPointClose (const Eigen::Vector3d &a, const Eigen::Vector3d &b, double tol)
bool compareByT (const IntersectionResult &a, const IntersectionResult &b)
std::vector< IntersectionResultfindOrderedIntersections (vtkUnstructuredGrid *grid, Eigen::Vector3d const &p0, Eigen::Vector3d const &p1, const double free_fraction, double const tol)
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. The line segment is defined by consecutive entries in realcoords.
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.

Function Documentation

◆ compareByT()

bool compareByT ( const IntersectionResult & a,
const IntersectionResult & b )

Definition at line 42 of file ComputeIntersections.cpp.

43{
44 return a.t < b.t;
45}

References IntersectionResult::t.

Referenced by findOrderedIntersections().

◆ findOrderedIntersections()

std::vector< IntersectionResult > findOrderedIntersections ( vtkUnstructuredGrid * grid,
Eigen::Vector3d const & p0,
Eigen::Vector3d const & p1,
const double free_fraction,
double const tol )

Definition at line 47 of file ComputeIntersections.cpp.

53{
54 std::vector<IntersectionResult> intersections;
55
56 intersections.push_back(IntersectionResult(p0, 0.0));
57 intersections.push_back(IntersectionResult(p1, 1.0));
58 Eigen::Vector3d transition_point = p0 + (p1 - p0) * free_fraction;
59 if (!isPointClose(p0, transition_point, tol) &&
60 !isPointClose(p1, transition_point, tol))
61 {
62 INFO("Adding transition point at t = {:.2f} for free fraction {:.2f}",
63 free_fraction, free_fraction);
64 intersections.push_back(
65 IntersectionResult(transition_point, free_fraction));
66 }
67 // For cell intersection
68 vtkSmartPointer<vtkGenericCell> cell =
69 vtkSmartPointer<vtkGenericCell>::New();
70
71 vtkNew<vtkIdList> filteredCellIds;
72 vtkNew<vtkCellLocator> locator;
73 vtkSmartPointer<vtkDataSet> cells_to_be_intersected;
74 if (GetGridDimension(grid) < 3)
75 {
76 vtkSmartPointer<vtkExtractEdges> extractEdges =
77 vtkSmartPointer<vtkExtractEdges>::New();
78 extractEdges->SetInputData(grid);
79 extractEdges->Update();
80
81 cells_to_be_intersected = extractEdges->GetOutput();
82 }
83 else
84 {
85 cells_to_be_intersected = grid;
86 }
87 locator->SetDataSet(cells_to_be_intersected);
88 locator->BuildLocator();
89 locator->FindCellsAlongLine(p0.data(), p1.data(), tol, filteredCellIds);
90 for (vtkIdType i = 0; i < filteredCellIds->GetNumberOfIds(); ++i)
91 {
92 vtkIdType cellId = filteredCellIds->GetId(i);
93 cells_to_be_intersected->GetCell(cellId, cell);
94
95 double t;
96 double x[3];
97 double pcoords[3];
98 int subId;
99
100 int intersected = cell->IntersectWithLine(p0.data(), p1.data(), tol, t,
101 x, pcoords, subId);
102
103 Eigen::Vector3d x_vec(x[0], x[1], x[2]);
104
105 if (!intersected)
106 {
107 continue;
108 }
109 bool isDuplicate = false;
110 for (const auto& existing : intersections)
111 {
112 if (isPointClose(existing.point, x_vec, tol))
113 {
114 isDuplicate = true;
115 break;
116 }
117 }
118 if (!isDuplicate && (t < tol || t > free_fraction))
119 {
120 IntersectionResult result;
121 result.t = t;
122 result.point = x_vec;
123 intersections.push_back(result);
124 }
125 }
126 std::sort(intersections.begin(), intersections.end(), compareByT);
127 return intersections;
128}
bool compareByT(const IntersectionResult &a, const IntersectionResult &b)
bool isPointClose(const Eigen::Vector3d &a, const Eigen::Vector3d &b, double tol)
int GetGridDimension(vtkUnstructuredGrid *grid)
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
Result of an intersection of a line with a cell.

References compareByT(), GetGridDimension(), INFO(), isPointClose(), IntersectionResult::point, and IntersectionResult::t.

Referenced by getOrderedAnchorCoords().

◆ GetGridDimension()

int GetGridDimension ( vtkUnstructuredGrid * grid)

Definition at line 15 of file ComputeIntersections.cpp.

16{
17 vtkIdType numCells = grid->GetNumberOfCells();
18 int maxDim = 0;
19
20 vtkSmartPointer<vtkGenericCell> cell =
21 vtkSmartPointer<vtkGenericCell>::New();
22
23 for (vtkIdType i = 0; i < numCells; ++i)
24 {
25 grid->GetCell(i, cell);
26 int dim = cell->GetCellDimension();
27 if (dim > maxDim)
28 {
29 maxDim = dim;
30 }
31 }
32 return maxDim;
33}

Referenced by findOrderedIntersections().

◆ getOrderedAnchorCoords()

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. The line segment is defined by consecutive entries in realcoords.

Parameters
gridThe vtkUnstructuredGrid to intersect with.
realcoordsoriginal coordinates of the anchor start and end points.
free_fractionThe fraction of the line segment to consider for intersections.
tolTolerance for considering points as identical.
Returns
a vector of vectors containing the intersection results which need to be converted using setPhysicalPropertiesForIntersectionPoints() to fill in physical properties and to obtain the final anchor format

Definition at line 130 of file ComputeIntersections.cpp.

135{
136 std::vector<std::vector<IntersectionResult>> ordered_intersections(
137 realcoords.rows() / 2);
138 assert(realcoords.rows() % 2 == 0);
139 assert(free_fraction.size() == realcoords.rows() / 2);
140 for (int i = 0; i < realcoords.rows(); i += 2)
141 {
142 ordered_intersections[i / 2] = findOrderedIntersections(
143 grid, realcoords.row(i), realcoords.row(i + 1),
144 free_fraction(i / 2), tol);
145 }
146 return ordered_intersections;
147}
std::vector< IntersectionResult > findOrderedIntersections(vtkUnstructuredGrid *grid, Eigen::Vector3d const &p0, Eigen::Vector3d const &p1, const double free_fraction, double const tol)

References findOrderedIntersections().

Referenced by main().

◆ isPointClose()

bool isPointClose ( const Eigen::Vector3d & a,
const Eigen::Vector3d & b,
double tol )

Definition at line 35 of file ComputeIntersections.cpp.

37{
38 return (a - b).squaredNorm() <= tol * tol;
39}

Referenced by findOrderedIntersections().

◆ setPhysicalPropertiesForIntersectionPoints()

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.

Parameters
anchor_coordsThe intersection points of the anchors with the bulk mesh.
original_anchor_dataThe original anchor data read from the JSON file.
Returns
struct containing the anchor data which will be stored in the unstructured grid

Definition at line 149 of file ComputeIntersections.cpp.

152{
153 std::vector<IntersectionResult> all_intersections;
154 std::vector<int> anchorids;
155 for (auto&& [anchor_id, intersections] :
156 ranges::views::enumerate(anchor_coords))
157 {
158 for (auto&& [index, result] : ranges::views::enumerate(intersections))
159 {
160 if (index > 0 && index < intersections.size() - 1)
161 {
162 anchorids.push_back(anchor_id);
163 all_intersections.push_back(result);
164 }
165 anchorids.push_back(anchor_id);
166 all_intersections.push_back(result);
167 }
168 }
169 assert(all_intersections.size() == anchorids.size());
170 assert(all_intersections.size() % 2 == 0);
171 Eigen::MatrixXd realcoords(anchorids.size(), 3);
172 auto const number_of_anchors = anchorids.size() / 2;
173 Eigen::VectorXd initial_stress(number_of_anchors);
174 Eigen::VectorXd maximum_stress(number_of_anchors);
175 Eigen::VectorXd residual_stress(number_of_anchors);
176 Eigen::VectorXd cross_sectional_area(number_of_anchors);
177 Eigen::VectorXd stiffness(number_of_anchors);
178 for (std::size_t i = 0; i < all_intersections.size(); ++i)
179 {
180 realcoords.row(i).noalias() = all_intersections[i].point;
181 if (i % 2 == 0)
182 {
183 initial_stress(i / 2) =
184 original_anchor_data.initial_anchor_stress(anchorids[i]);
185 maximum_stress(i / 2) =
186 original_anchor_data.maximum_anchor_stress(anchorids[i]);
187 residual_stress(i / 2) =
188 original_anchor_data.residual_anchor_stress(anchorids[i]);
189 cross_sectional_area(i / 2) =
190 original_anchor_data.anchor_cross_sectional_area(anchorids[i]);
191 stiffness(i / 2) =
192 original_anchor_data.anchor_stiffness(anchorids[i]);
193 }
194 }
196 anchor_data.real_coords = realcoords;
197 anchor_data.initial_anchor_stress = initial_stress;
198 anchor_data.maximum_anchor_stress = maximum_stress;
199 anchor_data.residual_anchor_stress = residual_stress;
200 anchor_data.anchor_cross_sectional_area = cross_sectional_area;
201 anchor_data.anchor_stiffness = stiffness;
202 return anchor_data;
203}

References ApplicationUtils::ComputeNaturalCoordsResult::anchor_cross_sectional_area, ApplicationUtils::ComputeNaturalCoordsResult::anchor_stiffness, ApplicationUtils::ComputeNaturalCoordsResult::initial_anchor_stress, ApplicationUtils::ComputeNaturalCoordsResult::maximum_anchor_stress, ApplicationUtils::ComputeNaturalCoordsResult::real_coords, and ApplicationUtils::ComputeNaturalCoordsResult::residual_anchor_stress.

Referenced by main().