Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_EB_STL_utils.H
Go to the documentation of this file.
1#ifndef AMREX_EB_STL_UTILS_H_
2#define AMREX_EB_STL_UTILS_H_
3
4#include <AMReX_Config.H>
5#include <AMReX_Geometry.H>
6#include <AMReX_MultiFab.H>
7#include <AMReX_Dim3.H>
8#include <AMReX_EB2_Graph.H>
9
10#include <algorithm>
11#include <cstdint>
12#include <limits>
13#include <utility>
14
21namespace amrex
22{
23
26{
27public:
28
30 struct Triangle {
31 XDim3 v1, v2, v3;
32
33 [[nodiscard]] Real cent (int d) const
34 {
35 static_assert(sizeof(XDim3) == sizeof(Real)*3);
36 return Real(1./3.)*((&v1.x)[d] + (&v2.x)[d] + (&v3.x)[d]); // NOLINT(clang-analyzer-security.ArrayBound)
37 }
38
39 [[nodiscard]] std::pair<Real,Real> minmax (int d) const
40 {
41 static_assert(sizeof(XDim3) == sizeof(Real)*3);
42 return std::minmax({(&v1.x)[d], (&v2.x)[d], (&v3.x)[d]}); // NOLINT(clang-analyzer-security.ArrayBound)
43 }
44 };
46
48 template <int M, int N>
49 struct BVHNodeT
50 {
51 RealBox boundingbox{AMREX_D_DECL(std::numeric_limits<Real>::max(),
52 std::numeric_limits<Real>::max(),
53 std::numeric_limits<Real>::max()),
54 AMREX_D_DECL(std::numeric_limits<Real>::lowest(),
55 std::numeric_limits<Real>::lowest(),
56 std::numeric_limits<Real>::lowest())};
57 STLtools::Triangle triangles[M];
58 XDim3 trinorm[M];
59 int children[N];
60 std::int8_t ntriangles = 0;
61 std::int8_t nchildren = 0;
62 };
64
65 static constexpr int m_bvh_max_size = 4; // max # of triangles in a leaf node
66 static constexpr int m_bvh_max_splits = 4; // max # of children
67 static constexpr int m_bvh_max_stack_size = 12; // max depth of the tree
68
69 using Node = BVHNodeT<m_bvh_max_size,m_bvh_max_splits>;
70
71 static constexpr int allregular = -1;
72 static constexpr int mixedcells = 0;
73 static constexpr int allcovered = 1;
74
76 void setBVHOptimization (bool flag) { m_bvh_optimization = flag; }
77
86 void read_stl_file (std::string const& fname, Real scale, Array<Real,3> const& center,
87 int reverse_normal);
88
92 void fill (MultiFab& mf, IntVect const& nghost, Geometry const& geom,
93 Real outside_value = -1._rt, Real inside_value = 1._rt) const;
94
96 [[nodiscard]] int getBoxType (Box const& box, Geometry const& geom, RunOn) const;
97
98 static constexpr bool isGPUable () noexcept { return true; }
99
101 void fillFab (BaseFab<Real>& levelset, const Geometry& geom, RunOn,
102 Box const& bounding_box) const;
103
105 void getIntercept (Array<Array4<Real>,AMREX_SPACEDIM> const& inter_arr,
106 Array<Array4<EB2::Type_t const>,AMREX_SPACEDIM> const& type_arr,
107 Array4<Real const> const& lst, Geometry const& geom, RunOn,
108 Box const& bounding_box) const;
109
111 static void updateIntercept (Array<Array4<Real>,AMREX_SPACEDIM> const& inter_arr,
112 Array<Array4<EB2::Type_t const>,AMREX_SPACEDIM> const& type_arr,
113 Array4<Real const> const& lst, Geometry const& geom) ;
114
116 void fillSignedDistance (MultiFab& mf, IntVect const& nghost, Geometry const& geom) const;
117
119 void prepare (Gpu::PinnedVector<Triangle> a_tri_pts); // public for cuda
120
121private:
122
123 bool m_bvh_optimization = true;
124
125 Gpu::DeviceVector<Triangle> m_tri_pts_d;
126 Gpu::DeviceVector<XDim3> m_tri_normals_d;
127 Gpu::DeviceVector<Node> m_bvh_nodes;
128
129 int m_num_tri=0;
130
131 XDim3 m_ptmin; // All triangles are inside the bounding box defined by
132 XDim3 m_ptmax; // m_ptmin and m_ptmax.
133 XDim3 m_ptref; // The reference point is slightly outside the bounding box.
134 bool m_boundry_is_outside; // Is the bounding box boundary outside or inside the object?
135
136 void read_ascii_stl_file (std::string const& fname, Real scale,
137 Array<Real,3> const& center, int reverse_normal,
138 Gpu::PinnedVector<Triangle>& a_tri_pts);
139 void read_binary_stl_file (std::string const& fname, Real scale,
140 Array<Real,3> const& center, int reverse_normal,
141 Gpu::PinnedVector<Triangle>& a_tri_pts);
142
143 static void build_bvh (Triangle* begin, Triangle * end, Gpu::PinnedVector<Node>& bvh_nodes);
144 static void bvh_size (int ntri, std::size_t& nnodes);
145};
146
147}
148#endif
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
A FortranArrayBox(FAB)-like object.
Definition AMReX_BaseFab.H:190
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
A Box with real dimensions.
Definition AMReX_RealBox.H:28
Utility class that loads STL meshes and samples them onto AMReX grids.
Definition AMReX_EB_STL_utils.H:26
static constexpr int allregular
Definition AMReX_EB_STL_utils.H:71
static constexpr int m_bvh_max_stack_size
Definition AMReX_EB_STL_utils.H:67
static void updateIntercept(Array< Array4< Real >, 3 > const &inter_arr, Array< Array4< EB2::Type_t const >, 3 > const &type_arr, Array4< Real const > const &lst, Geometry const &geom)
Update intercepts.
Definition AMReX_EB_STL_utils.cpp:1274
static constexpr int m_bvh_max_size
Definition AMReX_EB_STL_utils.H:65
static constexpr int m_bvh_max_splits
Definition AMReX_EB_STL_utils.H:66
void fillSignedDistance(MultiFab &mf, IntVect const &nghost, Geometry const &geom) const
Fill mf with signed-distance values sampled from the STL mesh.
Definition AMReX_EB_STL_utils.cpp:1331
BVHNodeT< m_bvh_max_size, m_bvh_max_splits > Node
Definition AMReX_EB_STL_utils.H:69
int getBoxType(Box const &box, Geometry const &geom, RunOn) const
Return whether the EB is regular/covered/mixed inside box.
Definition AMReX_EB_STL_utils.cpp:914
void prepare(Gpu::PinnedVector< Triangle > a_tri_pts)
Upload triangle data to device memory (exposed for CUDA workflows).
Definition AMReX_EB_STL_utils.cpp:559
void fill(MultiFab &mf, IntVect const &nghost, Geometry const &geom, Real outside_value=-1._rt, Real inside_value=1._rt) const
Fill a level-set MultiFab sampled from the STL mesh.
Definition AMReX_EB_STL_utils.cpp:838
void getIntercept(Array< Array4< Real >, 3 > const &inter_arr, Array< Array4< EB2::Type_t const >, 3 > const &type_arr, Array4< Real const > const &lst, Geometry const &geom, RunOn, Box const &bounding_box) const
Compute intercepts and face types for EB geometry based on the STL mesh.
Definition AMReX_EB_STL_utils.cpp:1090
void read_stl_file(std::string const &fname, Real scale, Array< Real, 3 > const &center, int reverse_normal)
Read an STL file, apply scaling/translation, and store triangles.
Definition AMReX_EB_STL_utils.cpp:406
void setBVHOptimization(bool flag)
Enable/disable BVH acceleration when sampling.
Definition AMReX_EB_STL_utils.H:76
static constexpr int mixedcells
Definition AMReX_EB_STL_utils.H:72
static constexpr bool isGPUable() noexcept
Definition AMReX_EB_STL_utils.H:98
static constexpr int allcovered
Definition AMReX_EB_STL_utils.H:73
void fillFab(BaseFab< Real > &levelset, const Geometry &geom, RunOn, Box const &bounding_box) const
Fill a BaseFab with level-set data (host/device aware).
Definition AMReX_EB_STL_utils.cpp:1022
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
std::array< T, N > Array
Definition AMReX_Array.H:26
Definition AMReX_Amr.cpp:49
RunOn
Definition AMReX_GpuControl.H:69
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2006
__host__ __device__ constexpr IntVectND< dim > scale(const IntVectND< dim > &p, int s) noexcept
Returns a IntVectND obtained by multiplying each of the components of this IntVectND by s.
Definition AMReX_IntVect.H:1105
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2015
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Definition AMReX_Dim3.H:13
Real x
Definition AMReX_Dim3.H:13