preCICE
Loading...
Searching...
No Matches
ReceivedMeshContext.cpp
Go to the documentation of this file.
2
3#include <Eigen/Core>
4
6#include "mesh/Mesh.hpp"
7#include "utils/assertion.hpp"
8
9namespace precice::impl {
10
11void ReceivedMeshContext::checkVerticesInsideAccessRegion(precice::span<const double> coordinates, int meshDim, std::string_view functionName) const
12{
14 return;
15 }
16 const auto nVertices = (coordinates.size() / meshDim);
17 Eigen::Map<const Eigen::MatrixXd> C(coordinates.data(), meshDim, nVertices);
18 Eigen::VectorXd minCoeffs = C.rowwise().minCoeff();
19 Eigen::VectorXd maxCoeffs = C.rowwise().maxCoeff();
20 bool minCheck = (minCoeffs.array() >= userDefinedAccessRegion->minCorner().array()).all();
21 bool maxCheck = (maxCoeffs.array() <= userDefinedAccessRegion->maxCorner().array()).all();
22 PRECICE_CHECK(minCheck && maxCheck, "The provided coordinates in \"{}\" are not within the access region defined with \"setMeshAccessRegion()\". "
23 "Minimum corner of the provided values is (x,y,z) = ({}), the minimum corner of the access region box is (x,y,z) = ({}). "
24 "Maximum corner of the provided values is (x,y,z) = ({}), the maximum corner of the access region box is (x,y,z) = ({}). ",
25 functionName, minCoeffs, userDefinedAccessRegion->minCorner(), maxCoeffs, userDefinedAccessRegion->maxCorner());
26 C.colwise().maxCoeff();
27}
28
29std::vector<std::reference_wrapper<const mesh::Vertex>> ReceivedMeshContext::filterVerticesToLocalAccessRegion(bool requiresBB) const
30{
31 std::vector<std::reference_wrapper<const mesh::Vertex>> filteredVertices;
32 for (const auto &v : mesh->vertices()) {
33 // either the vertex lies within the region OR the user-defined region is not strictly necessary
35 // region is defined: only add if the vertex is inside the region
36 if (userDefinedAccessRegion->contains(v)) {
37 filteredVertices.push_back(std::cref(v));
38 }
39 } else if (!requiresBB) {
40 // region is not defined, so if filtering isn't required, add all vertices
41 filteredVertices.push_back(std::cref(v));
42 }
43 }
44 return filteredVertices;
45}
46
47} // namespace precice::impl
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
constexpr pointer data() const noexcept
Definition span.hpp:500
constexpr size_type size() const noexcept
Definition span.hpp:469
mesh::PtrMesh mesh
Mesh holding the geometry data structure.
void checkVerticesInsideAccessRegion(precice::span< const double > coordinates, int meshDim, std::string_view functionName) const
std::shared_ptr< mesh::BoundingBox > userDefinedAccessRegion
std::vector< std::reference_wrapper< const mesh::Vertex > > filterVerticesToLocalAccessRegion(bool requiresBB) const
Filters vertices to those within the local access region.