Skip to content

Commit

Permalink
prevent repeated edge splits
Browse files Browse the repository at this point in the history
  • Loading branch information
Ahdhn committed Nov 5, 2023
1 parent c2cac54 commit d69c19c
Show file tree
Hide file tree
Showing 3 changed files with 86 additions and 34 deletions.
57 changes: 42 additions & 15 deletions apps/Remesh/remesh_kernels.cuh
Original file line number Diff line number Diff line change
Expand Up @@ -36,6 +36,7 @@ __global__ static void compute_average_edge_length(
template <typename T, uint32_t blockThreads>
__global__ static void edge_split(rxmesh::Context context,
const rxmesh::VertexAttribute<T> coords,
rxmesh::EdgeAttribute<int8_t> updated,
const T high_edge_len_sq)
{
// EV for calc edge len
Expand All @@ -54,25 +55,30 @@ __global__ static void edge_split(rxmesh::Context context,
return;
}

Bitmask is_updated(cavity.patch_info().edges_capacity[0], shrd_alloc);

auto should_split = [&](const EdgeHandle& eh, const VertexIterator& iter) {
const Vec3<T> v0(
coords(iter[0], 0), coords(iter[0], 1), coords(iter[0], 2));
const Vec3<T> v1(
coords(iter[1], 0), coords(iter[1], 1), coords(iter[1], 2));
if (updated(eh) == 0) {
const Vec3<T> v0(
coords(iter[0], 0), coords(iter[0], 1), coords(iter[0], 2));
const Vec3<T> v1(
coords(iter[1], 0), coords(iter[1], 1), coords(iter[1], 2));

const T edge_len = glm::distance2(v0, v1);
const T edge_len = glm::distance2(v0, v1);

if (edge_len > high_edge_len_sq) {
cavity.create(eh);
if (edge_len > high_edge_len_sq) {
cavity.create(eh);
}
}
};

Query<blockThreads> query(context, cavity.patch_id());
query.dispatch<Op::EV>(block, shrd_alloc, should_split);
block.sync();

if (cavity.prologue(block, shrd_alloc, coords)) {
if (cavity.prologue(block, shrd_alloc, coords, updated)) {

is_updated.reset(block);

cavity.for_each_cavity(block, [&](uint16_t c, uint16_t size) {
assert(size == 4);
Expand All @@ -92,14 +98,22 @@ __global__ static void edge_split(rxmesh::Context context,
cavity.add_edge(new_v, cavity.get_cavity_vertex(c, 0));
const DEdgeHandle e_init = e0;

is_updated.set(e0.local_id(), true);

if (e0.is_valid()) {
for (uint16_t i = 0; i < size; ++i) {
const DEdgeHandle e = cavity.get_cavity_edge(c, i);

is_updated.set(e.local_id(), true);

const DEdgeHandle e1 =
(i == size - 1) ?
e_init.get_flip_dedge() :
cavity.add_edge(
cavity.get_cavity_vertex(c, i + 1), new_v);

is_updated.set(e1.local_id(), true);

if (!e1.is_valid()) {
break;
}
Expand All @@ -115,6 +129,14 @@ __global__ static void edge_split(rxmesh::Context context,
}

cavity.epilogue(block);

if (cavity.is_successful()) {
detail::for_each_edge(cavity.patch_info(), [&](EdgeHandle eh) {
if (is_updated(eh.local_id())) {
updated(eh) = 1;
}
});
}
}

template <typename T, uint32_t blockThreads>
Expand Down Expand Up @@ -158,13 +180,18 @@ __global__ static void edge_collapse(rxmesh::Context context,
// check if the edge is short and should be collapsed
auto should_collapse = [&](const EdgeHandle& eh,
const VertexIterator& iter) {
const VertexHandle v0(iter[0]), v1(iter[1]);

const Vec3<T> p0(coords(v0, 0), coords(v0, 1), coords(v0, 2));
const Vec3<T> p1(coords(v1, 0), coords(v1, 1), coords(v1, 2));
const T edge_len_sq = glm::distance2(p0, p1);
if (edge_len_sq < low_edge_len_sq) {
cavity.create(eh);
// only when both the two end of the edge are not tri-valent, we may do
// the collapse. Otherwise, we will run into inconsistent topology
// problem i.e., foldover
if (!is_tri_valent(iter.local(0)) && !is_tri_valent(iter.local(1))) {
const VertexHandle v0(iter[0]), v1(iter[1]);

const Vec3<T> p0(coords(v0, 0), coords(v0, 1), coords(v0, 2));
const Vec3<T> p1(coords(v1, 0), coords(v1, 1), coords(v1, 2));
const T edge_len_sq = glm::distance2(p0, p1);
if (edge_len_sq < low_edge_len_sq) {
cavity.create(eh);
}
}
};
query.dispatch<Op::EV>(block, shrd_alloc, should_collapse);
Expand Down
54 changes: 35 additions & 19 deletions apps/Remesh/remesh_rxmesh.cuh
Original file line number Diff line number Diff line change
Expand Up @@ -5,9 +5,10 @@
#include "remesh_kernels.cuh"

template <typename T>
inline void split_long_edges(rxmesh::RXMeshDynamic& rx,
rxmesh::VertexAttribute<T>* coords,
const T high_edge_len_sq)
inline void split_long_edges(rxmesh::RXMeshDynamic& rx,
rxmesh::VertexAttribute<T>* coords,
rxmesh::EdgeAttribute<int8_t>* updated,
const T high_edge_len_sq)
{
using namespace rxmesh;

Expand All @@ -20,16 +21,24 @@ inline void split_long_edges(rxmesh::RXMeshDynamic& rx,
// rx.get_context().m_patch_scheduler.print_list();

LaunchBox<blockThreads> launch_box;
rx.prepare_launch_box(
{Op::EV}, launch_box, (void*)edge_split<T, blockThreads>);
rx.prepare_launch_box({Op::EV},
launch_box,
(void*)edge_split<T, blockThreads>,
true,
false,
false,
[&](uint32_t v, uint32_t e, uint32_t f) {
return detail::mask_num_bytes(e) +
ShmemAllocator::default_alignment;
});

edge_split<T, blockThreads><<<DIVIDE_UP(launch_box.blocks, 2),
launch_box.num_threads,
launch_box.smem_bytes_dyn>>>(
rx.get_context(), *coords, high_edge_len_sq);
rx.get_context(), *coords, *updated, high_edge_len_sq);

rx.cleanup();
rx.slice_patches(*coords);
rx.slice_patches(*coords, *updated);
rx.cleanup();

// RXMESH_TRACE(" ");
Expand All @@ -43,11 +52,14 @@ inline void split_long_edges(rxmesh::RXMeshDynamic& rx,
// bool show = true;
// if (show) {
// coords->move(DEVICE, HOST);
// updated->move(DEVICE, HOST);
// rx.update_polyscope();
// auto ps_mesh = rx.get_polyscope_mesh();
// ps_mesh->updateVertexPositions(*coords);
// ps_mesh->setEnabled(false);
//
// ps_mesh->addEdgeScalarQuantity("updated", *updated);
//
// rx.render_vertex_patch();
// rx.render_edge_patch();
// rx.render_face_patch();
Expand All @@ -57,10 +69,11 @@ inline void split_long_edges(rxmesh::RXMeshDynamic& rx,
}

template <typename T>
inline void collapse_short_edges(rxmesh::RXMeshDynamic& rx,
rxmesh::VertexAttribute<T>* coords,
const T low_edge_len_sq,
const T high_edge_len_sq)
inline void collapse_short_edges(rxmesh::RXMeshDynamic& rx,
rxmesh::VertexAttribute<T>* coords,
rxmesh::EdgeAttribute<int8_t>* updated,
const T low_edge_len_sq,
const T high_edge_len_sq)
{
using namespace rxmesh;

Expand All @@ -76,7 +89,7 @@ inline void collapse_short_edges(rxmesh::RXMeshDynamic& rx,
false,
true,
[=](uint32_t v, uint32_t e, uint32_t f) {
return detail::mask_num_bytes(e) +
return detail::mask_num_bytes(v) +
ShmemAllocator::default_alignment;
});

Expand All @@ -86,14 +99,15 @@ inline void collapse_short_edges(rxmesh::RXMeshDynamic& rx,
rx.get_context(), *coords, low_edge_len_sq, high_edge_len_sq);

rx.cleanup();
rx.slice_patches(*coords);
rx.slice_patches(*coords, *updated);
rx.cleanup();
}
}

template <typename T>
inline void equalize_valences(rxmesh::RXMeshDynamic& rx,
rxmesh::VertexAttribute<T>* coords)
inline void equalize_valences(rxmesh::RXMeshDynamic& rx,
rxmesh::EdgeAttribute<int8_t>* updated,
rxmesh::VertexAttribute<T>* coords)
{
using namespace rxmesh;

Expand All @@ -115,7 +129,7 @@ inline void equalize_valences(rxmesh::RXMeshDynamic& rx,
launch_box.smem_bytes_dyn>>>(rx.get_context(), *coords);

rx.cleanup();
rx.slice_patches(*coords);
rx.slice_patches(*coords, *updated);
rx.cleanup();
}
}
Expand Down Expand Up @@ -156,6 +170,7 @@ inline void remesh_rxmesh(rxmesh::RXMeshDynamic& rx)

auto coords = rx.get_input_vertex_coordinates();
auto new_coords = rx.add_vertex_attribute<float>("newCoords", 3);
auto updated = rx.add_edge_attribute<int8_t>("Updated", 1);

// compute average edge length
float* average_edge_len;
Expand Down Expand Up @@ -189,10 +204,11 @@ inline void remesh_rxmesh(rxmesh::RXMeshDynamic& rx)
timer.start();

for (uint32_t iter = 0; iter < Arg.num_iter; ++iter) {
split_long_edges(rx, coords.get(), high_edge_len_sq);
updated->reset(0, DEVICE);
split_long_edges(rx, coords.get(), updated.get(), high_edge_len_sq);
collapse_short_edges(
rx, coords.get(), low_edge_len_sq, high_edge_len_sq);
equalize_valences(rx, coords.get());
rx, coords.get(), updated.get(), low_edge_len_sq, high_edge_len_sq);
equalize_valences(rx, updated.get(), coords.get());
tangential_relaxation(rx, coords.get(), new_coords.get());
std::swap(new_coords, coords);
}
Expand Down
9 changes: 9 additions & 0 deletions include/rxmesh/cavity_manager.cuh
Original file line number Diff line number Diff line change
Expand Up @@ -86,6 +86,15 @@ struct CavityManager
__device__ __inline__ bool is_successful(HandleT seed);


/**
* @brief check if all cavities created in this CavityManager will be
* written to global memory
*/
__device__ __inline__ bool is_successful()
{
return m_write_to_gmem && !m_s_remove_fill_in[0];
}

/**
* @brief returns a handle to the mesh element that has created a give
* cavity
Expand Down

0 comments on commit d69c19c

Please sign in to comment.