Skip to content

ai-sci-computing/quads

Folders and files

NameName
Last commit message
Last commit date

Latest commit

 

History

1 Commit
 
 
 
 
 
 
 
 
 
 
 
 
 
 

Repository files navigation

quads

Quad mesh of the Fertility model

Quad meshing of closed triangle surfaces by stripe-pattern parametrization. Cross field → 4-cover stripe pattern → integer-grid extraction → boundary cleanup → Steiner-fill holes.

Algorithm

input triangle mesh (closed, manifold)
        │
        ▼
  cross field   (n=4 direction field, smoothed via Globally Optimal
        │       Direction Fields; Knöppel et al. SIGGRAPH 2013)
        ▼
  stripe pattern on the 4-fold branched cover
        │       (smallest eigenvector of a stripe-energy / mass system;
        │       Knöppel et al. "Stripe Patterns on Surfaces", SIGGRAPH 2015)
        ▼
  per-face (u, v) texture coords + planar placement triangles
        │
        ▼
  integer-grid extraction
        │       per placement triangle, enumerate integer (u, v) lattice
        │       points + boundary crossings, lift to 3D, dedup, emit
        │       quads where all four corners exist
        ▼
  raw quad mesh (mostly quads, some n-gons / branch-incident garbage)
        │
        ▼
  purify + boundary cleanup
        │       drop non-quads and branch-incident quads (creates holes);
        │       iteratively drop ear/floater quads (≥3 boundary edges)
        │       and reflex-corner quads (any vertex with k=1 incidence)
        │       until the hole boundaries are smooth
        ▼
  Steiner hole fill
        │       per hole cycle: decompose at pinch vertices into simple
        │       sub-cycles; choose path by hole size:
        │         small  (n_rings=1): recursive shortest-diagonal split,
        │                              no Steiner vertices
        │         medium (n_rings=2): single centroid + fan
        │         large  (n_rings≥3): concentric rings interpolated to
        │                              centroid + recursive split of innermost
        │       all fill quads use reversed winding for manifold consistency
        ▼
  tangential smoothing (10 iterations, face-barycenter Laplacian)
        │       smooth Steiner vertices + cycle boundary; project back
        │       to input surface via inline triangle BVH
        ▼
final quad mesh

The quad-extraction code is in src/quad_extract.{h,cpp}. The cross-field solve is in src/cross_field.cpp, the stripe-pattern solve in src/stripe_pattern.cpp. Detailed math notes are in docs/design.md. A Doxyfile is included for API doc generation:

doxygen Doxyfile      # writes HTML to docs/api/html/

Build

Requires CMake ≥ 3.10, a C++17 compiler, and Eigen. Optional but recommended: SuiteSparse (CHOLMOD), used for the sparse Cholesky factorization of the stripe energy matrix — without it the build falls back to Eigen's slower SimplicialLDLT. On macOS:

brew install suite-sparse
cmake -S . -B build -DCMAKE_BUILD_TYPE=Release
cmake --build build -j

Polyscope, GLFW, ImGui, GLM, and Spectra are bundled via FetchContent.

Usage

GUI (Polyscope viewer with cross field, stripe pattern, layout, quad mesh quantities and a frequency slider):

./build/quads path/to/mesh.obj

Editing singularities (GUI)

The tool selector in the side panel offers three interactive cross-field edits. After each edit the stripe pattern and quad mesh are recomputed.

  • Drag singularity — press-and-hold the left mouse on a singularity sphere and release on a target face to relocate it. While dragging, an orange streamline overlay shows how the cross field would respond to the proposed location. The streamlines are re-traced on every hover from a trivial-connection update plus one face-Poisson solve (both O(|V|)); the full stripe eigensolve is too slow to run per frame and is deferred until release.
  • Place +1/−1 pair — click a face for the +1, then a face for the −1.
  • Cancel opposite pair — click two opposite-sign singularities to cancel them.

Stripe solve toggles

Two checkboxes in the side panel change what the stripe solver sees:

  • Hodge project omega to closed — replace the cross-field 1-form ω with its Hodge projection onto closed forms on the 4-fold cover before the eigensolve. Removes any curl on non-branch faces (topologically required curl stays in the branch hexagons), which cleans up the integer-grid behavior in irregular regions. Adds one face-Laplacian solve per stripe pattern compute.
  • Kill non-branch stripe defects — after the eigenpair (ψ, φ) is solved, some non-branch faces can still carry per-face stripe winding Σσ = 2π·n_ijk with n_ijk ≠ 0 — the "rounding tension" that breaks isoline matching across shared edges. This step solves the same face-Laplacian to find a coexact 1-form η that cancels the defects exactly, then transports it through (ψ, φ) so the stripe phases follow. Branch faces keep their topologically required winding.

Headless test

Runs the full pipeline and dumps an OBJ + per-stage timing breakdown + watertight check:

./build/quads_extract_test path/to/mesh.obj <frequency> <out.obj>

frequency is the target stripe density (≈ output cells per integration length); higher → more, smaller quads.

Timings

All numbers wall-clock, single-threaded, average of 3 runs on an Apple Silicon laptop. Instant Meshes ≈ Jakob et al., SIGGRAPH Asia 2015, built from github.com/wjakob/instant-meshes. IM (-D, 1 thread) is single-threaded dominant-mode (no 1→4 subdivision); IM parallel is default pure-quad mode using all cores.

input vertices output quads quads (1 thread) IM (1 thread) IM (parallel)
1.4k 40 20 ms 10 ms 10 ms
6.5k 800 120 ms 130 ms 70 ms
14k 2.4k 530 ms 290 ms 140 ms
~370k (~750k F) 600k 7.7 s 147 s 10.7 s

Improvement on large meshes — both vs single-threaded IM (~19×) and vs parallel IM (~1.4×). On small meshes IM is 2-4× faster.

The cross-over is around 20-50k output quads. Below it, the stripe pattern's sparse Cholesky factorization dominates our pipeline (60-73% of small-mesh time); IM has no equivalent global eigensolve. Above it, factorize is amortized out and our O(F) extraction + Steiner fill beats IM's per-output-cell cost.

Strengths

  • Linear scaling on large meshes. The dominant cost above ~50k quads is extract_quads and fill_holes_steiner, both O(F_output) with small constants. The eigensolve runs once on the input mesh, so its cost stops mattering as output density grows.
  • Watertight + manifold output on most inputs. Hole boundaries are smoothed (ears and reflex corners dropped) and fill quads are emitted with reversed winding for half-edge consistency.
  • No external mesh-processing dependencies. Inline AABB tree for closest-point projection, no libigl.

Weaknesses

  • Slow on small meshes. The stripe-pattern eigensolve is structural; CHOLMOD supernodal LLT is already the fastest sparse SPD solver available, so this is an algorithmic floor rather than an optimization gap. Algorithms without a global eigensolve (Instant Meshes, MIQ) are faster below ~50k output quads.
  • Triangles can appear at odd-cycle holes. Per-cycle parity is topologically locked (E_bdy = 4F − 2E_int is always even but individual cycle lengths can be odd). The Steiner fan for an odd cycle emits exactly one triangle. We've moved that triangle from the cycle centroid to the most acute corner of the inner-ring polygon so it doesn't dominate the eye, but it's still there. Ways to eliminate remaining triangles — widening odd cycles to merge with adjacent ones — are deferred future work.
  • Some non-manifold edges remain on tricky inputs. Bumpy or near-degenerate surfaces with many small holes can leave on the order of tens of non-manifold edges out of ~1M total — typically from interactions between the recursive inner-ring split and the pinch-cycle decomposition. The watertight check in tests/test_quad_extract.cpp reports them.
  • Closed surfaces only. Boundary handling is not implemented.
  • Single-threaded. Most steps would parallelize naturally (extract_quads per placement triangle, BVH queries per Steiner vertex, smoothing iterations) but currently aren't.

References

Algorithmic basis:

  • Knöppel, Crane, Pinkall, Schröder. Globally Optimal Direction Fields. ACM Transactions on Graphics (SIGGRAPH 2013), 32(4):59:1–59:10. — the smooth cross-field stage.
  • Knöppel, Crane, Pinkall, Schröder. Stripe Patterns on Surfaces. ACM Transactions on Graphics (SIGGRAPH 2015), 34(4):39:1–39:11. — the n=2 line-field stripe pattern that we extend to the 4-fold cover for n=4.

Comparison baseline:

  • Jakob, Tarini, Panozzo, Sorkine-Hornung. Instant Field-Aligned Meshes. ACM Transactions on Graphics (SIGGRAPH Asia 2015), 34(6):189:1–189:15. — code: github.com/wjakob/instant-meshes.

Tools and libraries:

  • Sharp, Soliman, Crane. Polyscope: A C++ & Python viewer for 3D data. polyscope.run — used for the GUI.
  • Davis. SuiteSparse / CHOLMOD. Sparse direct solvers; supernodal LLT used for the stripe energy factorization.
  • Guennebaud, Jacob, et al. Eigen v3. eigen.tuxfamily.org — sparse and dense linear algebra throughout.

License

Free for personal, academic, and research use. No warranty. Commercial use requires written permission from the author.

About

Stripe pattern approach to quad meshing

Resources

Stars

Watchers

Forks

Releases

No releases published

Packages

 
 
 

Contributors