OGS
ComputeIntersections.h File Reference
#include <vtkType.h>
#include <vtkUnstructuredGrid.h>
#include <Eigen/Core>
#include <range/v3/view/enumerate.hpp>
#include <vector>
#include "ComputeNaturalCoordsResult.h"
Include dependency graph for ComputeIntersections.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

struct  IntersectionResult
 Result of an intersection of a line with a cell. More...

Functions

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

◆ 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 136 of file ComputeIntersections.cpp.

141{
142 std::vector<std::vector<IntersectionResult>> ordered_intersections(
143 realcoords.rows() / 2);
144 assert(realcoords.rows() % 2 == 0);
145 assert(free_fraction.size() == realcoords.rows() / 2);
146 for (int i = 0; i < realcoords.rows(); i += 2)
147 {
148 ordered_intersections[i / 2] = findOrderedIntersections(
149 grid, realcoords.row(i), realcoords.row(i + 1),
150 free_fraction(i / 2), tol);
151 }
152 return ordered_intersections;
153}
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().

◆ 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 155 of file ComputeIntersections.cpp.

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

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().