preCICE v3.2.0
Loading...
Searching...
No Matches
MeshContext.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <vector>
4#include "MappingContext.hpp"
5#include "SharedPointer.hpp"
7#include "mapping/Mapping.hpp"
11
12namespace precice::impl {
13
79
84
85inline void MeshContext::checkVerticesInsideAccessRegion(precice::span<const double> coordinates, const int meshDim, std::string_view functionName) const
86{
87
89 const auto nVertices = (coordinates.size() / meshDim);
90 Eigen::Map<const Eigen::MatrixXd> C(coordinates.data(), meshDim, nVertices);
91 Eigen::VectorXd minCoeffs = C.rowwise().minCoeff();
92 Eigen::VectorXd maxCoeffs = C.rowwise().maxCoeff();
93 bool minCheck = (minCoeffs.array() >= userDefinedAccessRegion->minCorner().array()).all();
94 bool maxCheck = (maxCoeffs.array() <= userDefinedAccessRegion->maxCorner().array()).all();
95 PRECICE_CHECK(minCheck && maxCheck, "The provided coordinates in \"{}\" are not within the access region defined with \"setMeshAccessRegion()\". "
96 "Minimum corner of the provided values is (x,y,z) = ({}), the minimum corner of the access region box is (x,y,z) = ({}). "
97 "Maximum corner of the provided values is (x,y,z) = ({}), the maximum corner of the access region box is (x,y,z) = ({}). ",
98 functionName, minCoeffs, userDefinedAccessRegion->minCorner(), maxCoeffs, userDefinedAccessRegion->maxCorner());
99 C.colwise().maxCoeff();
100 }
101}
102
104{
106 for (const auto &v : mesh->vertices()) {
107 // either the vertex lies within the region OR the user-defined region is not strictly necessary
109 // region is defined: only add if the vertex is inside the region
110 if (userDefinedAccessRegion->contains(v)) {
111 filteredVertices.push_back(std::cref(v));
112 }
113 } else if (!requiresBB) {
114 // region is not defined, so if filtering isn't required, add all vertices
115 filteredVertices.push_back(std::cref(v));
116 }
117 }
118 return filteredVertices;
119}
120
121} // namespace precice::impl
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
This class provides a lightweight logger.
Definition Logger.hpp:17
MeshRequirement
Specifies requirements for the input and output meshes of a mapping.
Definition Mapping.hpp:46
GeometricFilter
Defines the type of geometric filter used.
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
T max(T... args)
std::shared_ptr< Mesh > PtrMesh
std::shared_ptr< Partition > PtrPartition
T push_back(T... args)
T cref(T... args)
Stores a mesh and related objects and data.
mesh::PtrMesh mesh
Mesh holding the geometry data structure.
mapping::Mapping::MeshRequirement meshRequirement
Determines which mesh type has to be provided by the accessor.
void checkVerticesInsideAccessRegion(precice::span< const double > coordinates, const int meshDim, std::string_view functionName) const
partition::PtrPartition partition
Partition creating the parallel decomposition of the mesh.
partition::ReceivedPartition::GeometricFilter geoFilter
type of geometric filter
std::vector< MappingContext > toMappingContexts
Mappings used when mapping data to the mesh. Can be empty.
std::string receiveMeshFrom
Name of participant that creates the mesh.
bool provideMesh
True, if accessor does create the mesh.
std::vector< MappingContext > fromMappingContexts
Mappings used when mapping data from the mesh. Can be empty.
void require(mapping::Mapping::MeshRequirement requirement)
std::shared_ptr< mesh::BoundingBox > userDefinedAccessRegion
std::vector< std::reference_wrapper< const mesh::Vertex > > filterVerticesToLocalAccessRegion(bool requiresBB) const
double safetyFactor
bounding box to speed up decomposition of received mesh is increased by this safety factor