spade/delaunay_core/refinement.rs
1#[cfg(not(feature = "std"))]
2use hashbrown::{HashMap, HashSet};
3#[cfg(feature = "std")]
4use std::collections::{HashMap, HashSet};
5
6use alloc::collections::VecDeque;
7use alloc::vec::Vec;
8
9use num_traits::Float;
10
11use crate::{
12 delaunay_core::{dcel_operations::append_unconnected_vertex, math},
13 CdtEdge, ConstrainedDelaunayTriangulation, HasPosition, HintGenerator, Point2,
14 PositionInTriangulation, SpadeNum, Triangulation,
15};
16
17use super::{
18 DirectedEdgeHandle, FaceHandle, FixedFaceHandle, FixedUndirectedEdgeHandle, FixedVertexHandle,
19 InnerTag, TriangulationExt, UndirectedEdgeHandle,
20};
21
22/// Contains details about the outcome of a refinement procedure.
23///
24/// *See [ConstrainedDelaunayTriangulation::refine]*
25#[derive(Debug, Clone)]
26pub struct RefinementResult {
27 /// A `Vec` containing all outer faces that were excluded from refinement.
28 ///
29 /// This `Vec` will be empty unless [RefinementParameters::exclude_outer_faces] has been set.
30 /// In this case, the `Vec` contains all finite outer faces, including any additional outer faces
31 /// that were created during the refinement.
32 pub excluded_faces: Vec<FixedFaceHandle<InnerTag>>,
33
34 /// Set to `true` if the refinement could be completed regularly.
35 ///
36 /// This will be `false` if the refinement ran out of additional vertices
37 /// (see [RefinementParameters::with_max_additional_vertices]). Consider adapting the refinement parameters in this case,
38 /// either by using a higher additional vertex count or by e.g. lowering the [angle limit](RefinementParameters::with_angle_limit).
39 pub refinement_complete: bool,
40}
41
42/// Specifies the minimum allowed angle that should be kept after a refinement procedure.
43///
44/// The refinement algorithm will attempt to keep the *minimum angle in the triangulation* greater than
45/// an angle limit specified with this struct.
46///
47/// *See [ConstrainedDelaunayTriangulation::refine], [RefinementParameters::with_angle_limit]*
48#[derive(Copy, Clone, PartialEq, PartialOrd)]
49pub struct AngleLimit {
50 radius_to_shortest_edge_limit: f64,
51}
52
53impl AngleLimit {
54 /// Create a new angle limit from an angle given in degrees.
55 ///
56 /// Note that angles larger than 30 degrees will quickly lead to overrefinement as the algorithm
57 /// cannot necessarily guarantee termination (other than by limiting the number of additional inserted vertices).
58 ///
59 /// Defaults to 30°. An angle of 0 degrees will disable refining due to small angles.
60 ///
61 /// *See also [from_rad](AngleLimit::from_rad)*
62 pub fn from_deg(degree: f64) -> Self {
63 Self::from_rad(degree.to_radians())
64 }
65
66 /// Create a new angle limit from an angle given in radians.
67 ///
68 /// Note angles larger than 30 degrees (≈0.52rad = PI / 6) will quickly lead to poor refinement quality.
69 /// Passing in an angle of 0rad will disable refining due to small angles.
70 ///
71 /// *See also [from_deg](AngleLimit::from_deg)*
72 pub fn from_rad(rad: f64) -> Self {
73 let sin = rad.sin();
74 if sin == 0.0 {
75 Self::from_radius_to_shortest_edge_ratio(f64::INFINITY)
76 } else {
77 Self::from_radius_to_shortest_edge_ratio(0.5 / sin)
78 }
79 }
80
81 /// Returns the radius to shortest edge limit corresponding to this angle limit.
82 ///
83 /// See [from_radius_to_shortest_edge_ratio](AngleLimit::from_radius_to_shortest_edge_ratio) for more
84 /// information.
85 pub fn radius_to_shortest_edge_limit(&self) -> f64 {
86 self.radius_to_shortest_edge_limit
87 }
88
89 /// Creates a new angle limit by specifying the circumradius to shortest edge ratio that must be kept.
90 ///
91 /// For each face, this ratio is calculated by dividing the circumradius of the face by the length of its shortest
92 /// edge.
93 /// This ratio is related directly to the minimum allowed angle by the formula
94 /// `ratio = 1 / (2 sin * (min_angle))`.
95 /// The *larger* the allowed min angle is, the *smaller* ratio becomes.
96 ///
97 /// Larger ratio values will lead to a less refined triangulation. Passing in `f64::INFINITY` will disable
98 /// refining due to small angles.
99 ///
100 /// Defaults to 1.0 (30 degrees).
101 ///
102 /// # Example values
103 ///
104 /// | ratio | Bound on smallest angle (deg) | Bound on smallest angle (rad) |
105 /// |-------|-------------------------------|-------------------------------|
106 /// | 0.58 | 60.00° | 1.05 |
107 /// | 0.60 | 56.44° | 0.99 |
108 /// | 0.70 | 45.58° | 0.80 |
109 /// | 0.80 | 38.68° | 0.68 |
110 /// | 0.90 | 33.75° | 0.59 |
111 /// | 1.00 | 30.00° | 0.52 |
112 /// | 1.10 | 27.04° | 0.47 |
113 /// | 1.20 | 24.62° | 0.43 |
114 /// | 1.30 | 22.62° | 0.39 |
115 /// | +INF | 0° | 0 |
116 pub fn from_radius_to_shortest_edge_ratio(ratio: f64) -> Self {
117 Self {
118 radius_to_shortest_edge_limit: ratio,
119 }
120 }
121}
122
123impl alloc::fmt::Debug for AngleLimit {
124 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> alloc::fmt::Result {
125 f.debug_struct("AngleLimit")
126 .field(
127 "angle limit (deg)",
128 &(0.5 / self.radius_to_shortest_edge_limit)
129 .asin()
130 .to_degrees(),
131 )
132 .finish()
133 }
134}
135
136impl Default for AngleLimit {
137 fn default() -> Self {
138 Self::from_radius_to_shortest_edge_ratio(1.0)
139 }
140}
141
142#[derive(Debug, PartialEq, PartialOrd, Clone, Copy, Hash)]
143enum RefinementHint {
144 Ignore,
145 ShouldRefine,
146 MustRefine,
147}
148
149/// Controls how Delaunay refinement is performed.
150///
151/// Refer to [ConstrainedDelaunayTriangulation::refine] and methods implemented by this type for more details
152/// about which parameters are supported.
153///
154/// # Example
155///
156/// ```
157/// use spade::{AngleLimit, ConstrainedDelaunayTriangulation, Point2, RefinementParameters};
158///
159/// fn refine_cdt(cdt: &mut ConstrainedDelaunayTriangulation<Point2<f64>>) {
160/// let params = RefinementParameters::<f64>::new()
161/// .exclude_outer_faces(true)
162/// .keep_constraint_edges()
163/// .with_min_required_area(0.0001)
164/// .with_max_allowed_area(0.5)
165/// .with_angle_limit(AngleLimit::from_deg(25.0));
166///
167/// cdt.refine(params);
168/// }
169/// ```
170#[derive(Debug, PartialEq, Clone)]
171pub struct RefinementParameters<S: SpadeNum + Float> {
172 max_additional_vertices: Option<usize>,
173
174 angle_limit: AngleLimit,
175 min_area: Option<S>,
176 max_area: Option<S>,
177 keep_constraint_edges: bool,
178 exclude_outer_faces: bool,
179}
180
181impl<S: SpadeNum + Float> Default for RefinementParameters<S> {
182 fn default() -> Self {
183 Self {
184 max_additional_vertices: None,
185 angle_limit: AngleLimit::from_radius_to_shortest_edge_ratio(1.0),
186 min_area: None,
187 max_area: None,
188 exclude_outer_faces: false,
189 keep_constraint_edges: false,
190 }
191 }
192}
193
194impl<S: SpadeNum + Float> RefinementParameters<S> {
195 /// Creates a new set of `RefinementParameters`.
196 ///
197 /// The following values will be used by `new` and `Self::default`:
198 /// * `exclude_outer_faces`: disabled - all faces are used for refinement
199 /// * `keep_constraint_edges`: disabled
200 /// * `min_required_area`: disabled - no lower area limit is used
201 /// * `max_allowed_area`: disabled - no upper area limit is used
202 /// * `angle_limit`: 30 degrees by default.
203 /// * `num_additional_vertices`: 10 times the number of vertices in the triangulation
204 pub fn new() -> Self {
205 Self::default()
206 }
207
208 /// Specifies the smallest allowed inner angle in a refined triangulation.
209 ///
210 /// The refinement algorithm will attempt to insert additional points (called steiner points) until the
211 /// minimum angle is larger than the angle bound specified by the refinement parameters.
212 ///
213 /// Defaults to 30 degrees.
214 ///
215 /// Note that angle limits much larger than 30 degrees may not always terminate successfully - consider checking
216 /// [RefinementResult::refinement_complete] to make sure that the angle limit could actually be applied everywhere.
217 ///
218 /// # Examples of different angle limits
219 /// <table>
220 /// <tr><th>0° (no angle refinement)</th><th>20°</th><th>30°</th><th>34°</th></tr>
221 /// <tr><td>
222 #[doc = concat!(include_str!("../../images/angle_limit_00.svg"), "</td><td>")]
223 #[doc = concat!(include_str!("../../images/angle_limit_20.svg"), "</td><td>")]
224 #[doc = concat!(include_str!("../../images/angle_limit_30.svg"), "</td><td>")]
225 #[doc = concat!(include_str!("../../images/angle_limit_34.svg"), "</td></tr></table>")]
226 ///
227 /// *See also [ConstrainedDelaunayTriangulation::refine]*
228 pub fn with_angle_limit(mut self, angle_limit: AngleLimit) -> Self {
229 self.angle_limit = angle_limit;
230 self
231 }
232
233 /// Specifies a lower bound for a triangles area.
234 ///
235 /// The algorithm will attempt to ignore any triangle with an area below this limit. This can also prevent an
236 /// exhaustion of additionally available vertices (see [Self::with_max_additional_vertices]).
237 ///
238 /// Note that there is no guarantee that no face below this area bound will be kept intact - in some cases, a split
239 /// will still be required to restore the triangulation's Delaunay property. Also, this value does not specify a lower
240 /// bound for the smallest possible triangle in the triangulation.
241 ///
242 /// Should be set to something lower than [Self::with_max_allowed_area]. If this method is not called, no lower
243 /// bound check will be performed.
244 pub fn with_min_required_area(mut self, min_area: S) -> Self {
245 self.min_area = Some(min_area);
246 self
247 }
248
249 /// Specifies an upper bound for triangle areas in the triangulation.
250 ///
251 /// By default, the refinement tries to be conservative in how many vertices it adds. This will lead to an uneven
252 /// triangle size distribution - areas with larger input features will contain fewer, larger triangles whereas
253 /// regions with small features will contain more densely packed triangles.
254 /// By specifying an upper area bound for triangles, the resulting triangle sizes can be made more similar
255 /// as any large triangle above the bound will be split into smaller parts.
256 ///
257 /// # Examples of different maximum area values
258 #[doc = concat!(include_str!("../../images/refinement_maximum_area_no_limit.svg"))]
259 #[doc = concat!(include_str!("../../images/refinement_maximum_area_200.svg"))]
260 #[doc = concat!(include_str!("../../images/refinement_maximum_area_100.svg"))]
261 ///
262 /// Should be set to something larger than [Self::with_min_required_area]. If this method is not called, no upper area
263 /// bound check will be performed.
264 pub fn with_max_allowed_area(mut self, max_area: S) -> Self {
265 self.max_area = Some(max_area);
266 self
267 }
268
269 /// Specifies how many additional vertices may be inserted during Delaunay refinement.
270 ///
271 /// Refinement may, in some cases, fail to terminate if the angle limit is set too high
272 /// (see [with_angle_limit](Self::with_angle_limit)). Simply stopping the refinement after a certain number of vertices
273 /// has been inserted is an easy way to enforce termination. However, the resulting mesh may exhibit very poor quality
274 /// in this case - some areas may have become overly refined, others might be overlooked completely. Consider changing
275 /// the parameters (most notably the angle limit) if the refinement runs out of vertices.
276 ///
277 /// Use [RefinementResult::refinement_complete] to check if the number of additional vertices was sufficient.
278 pub fn with_max_additional_vertices(mut self, max_additional_vertices: usize) -> Self {
279 self.max_additional_vertices = Some(max_additional_vertices);
280 self
281 }
282
283 /// Prevents constraint edges from being split during refinement.
284 ///
285 /// By default, constraint edges may be split in order to restore the triangulation's Delaunay property.
286 /// The resulting two new edges will *become new constraint edges*, hence the original shape outlined by
287 /// constraint edges remains the same - no "gaps" or deviations are introduced.
288 ///
289 /// Enabling this option will, in general, reduce the quality of the resulting mesh - it is not necessarily
290 /// Delaunay anymore and faces adjacent to long constraint edges may violate the configured [AngleLimit].
291 pub fn keep_constraint_edges(mut self) -> Self {
292 self.keep_constraint_edges = true;
293 self
294 }
295
296 /// Allows to exclude outer faces from the refinement process.
297 ///
298 /// This is useful if the constraint edges form a *closed shape* with a clearly defined inner and outer part.
299 /// Spade will determine inner and outer faces by identifying which faces can be reached from the outer face
300 /// without "crossing" a constraint edge, similar to a flood fill algorithm.
301 ///
302 /// Any holes in the triangulation will also be excluded. More specifically, any point with an odd winding number
303 /// is considered to be inner (see e.g. [Wikipedia](https://en.wikipedia.org/wiki/Point_in_polygon#Winding_number_algorithm)).
304 ///
305 /// Note that excluded faces may still be subdivided if a neighboring edge needs to be split. However, they will never be the
306 /// *cause* for a subdivision - their angle and area is ignored.
307 ///
308 /// The resulting outer faces of the triangulation are returned by the call to [refine](ConstrainedDelaunayTriangulation::refine),
309 /// see [RefinementResult::excluded_faces].
310 ///
311 /// # Example
312 /// <table>
313 /// <tr><th>Unrefined</th><th>Refined</th></tr>
314 /// <tr><td>
315 #[doc = concat!(include_str!("../../images/exclude_outer_faces_unrefined.svg"), "</td><td>",include_str!("../../images/exclude_outer_faces_refined.svg"), " </td></tr></table>")]
316 ///
317 /// *A refinement operation configured to exclude outer faces. All colored faces are considered outer faces and are
318 /// ignored during refinement. Note that the inner part of the "A" shape forms a hole and is also excluded.*
319 pub fn exclude_outer_faces(mut self, exclude: bool) -> Self {
320 self.exclude_outer_faces = exclude;
321 self
322 }
323
324 fn get_refinement_hint<V, DE, UE, F>(
325 &self,
326 face: FaceHandle<InnerTag, V, DE, UE, F>,
327 ) -> RefinementHint
328 where
329 V: HasPosition<Scalar = S>,
330 {
331 if let Some(max_area) = self.max_area {
332 if face.area() > max_area {
333 return RefinementHint::MustRefine;
334 }
335 }
336
337 if let Some(min_area) = self.min_area {
338 if face.area() < min_area {
339 return RefinementHint::Ignore;
340 }
341 }
342
343 let (_, length2) = face.shortest_edge();
344 let (_, radius2) = face.circumcircle();
345
346 let ratio2 = radius2 / length2;
347 let angle_limit = self.angle_limit.radius_to_shortest_edge_limit;
348 if ratio2.into() > angle_limit * angle_limit {
349 RefinementHint::ShouldRefine
350 } else {
351 RefinementHint::Ignore
352 }
353 }
354}
355
356impl<V, DE, UE, F, L> ConstrainedDelaunayTriangulation<V, DE, UE, F, L>
357where
358 V: HasPosition + From<Point2<<V as HasPosition>::Scalar>>,
359 DE: Default,
360 UE: Default,
361 F: Default,
362 L: HintGenerator<<V as HasPosition>::Scalar>,
363 <V as HasPosition>::Scalar: Float,
364{
365 /// Refines a triangulation by inserting additional points to improve the quality of its mesh.
366 ///
367 /// *Mesh quality*, when applied to constrained delaunay triangulations (CDT), usually refers to how skewed its
368 /// triangles are. A skewed triangle is a triangle with very large or very small (acute) inner angles.
369 /// Some applications (e.g. interpolation and finite element methods) perform poorly in the presence of skewed triangles.
370 ///
371 /// Refining by inserting additional points (called "steiner points") may increase the minimum angle. The given
372 /// [RefinementParameters] should be used to modify the refinement behavior.
373 ///
374 /// # General usage
375 ///
376 /// The vertex type must implement `From<Point2<...>>` - otherwise, Spade cannot construct new steiner points at a
377 /// certain location. The refinement itself happens *in place* and will result in a valid CDT.
378 ///
379 /// ```
380 /// use spade::{ConstrainedDelaunayTriangulation, RefinementParameters, Point2, InsertionError, Triangulation};
381 ///
382 /// fn get_refined_triangulation(vertices: Vec<Point2<f64>>) ->
383 /// Result<ConstrainedDelaunayTriangulation<Point2<f64>>, InsertionError>
384 /// {
385 /// let mut cdt = ConstrainedDelaunayTriangulation::bulk_load(vertices)?;
386 /// let result = cdt.refine(RefinementParameters::default());
387 /// if !result.refinement_complete {
388 /// panic!("Refinement failed - I should consider using different parameters.")
389 /// }
390 ///
391 /// return Ok(cdt)
392 /// }
393 /// ```
394 ///
395 /// # Example image
396 ///
397 /// <table>
398 /// <tr><th>Unrefined</th><th>Refined</th></tr>
399 /// <tr><td>
400 #[doc = concat!(include_str!("../../images/unrefined.svg"), "</td><td>",include_str!("../../images/refined.svg"), " </td></tr></table>")]
401 ///
402 /// *A refinement example. The CDT on the left has some acute angles and skewed triangles.
403 /// The refined CDT on the right contains several additional points that prevents such triangles from appearing while keeping
404 /// all input vertices and constraint edges.*
405 ///
406 /// # Quality guarantees
407 ///
408 /// Refinement will ensure that the resulting triangulation fulfills a few properties:
409 /// - The triangulation's minimum angle will be larger than the angle specified by
410 /// [with_angle_limit](RefinementParameters::with_angle_limit).<br>
411 /// *Exception*: Acute input angles (small angles between initial constraint edges) cannot be maximized as the constraint edges
412 /// must be kept intact. The algorithm will, for the most part, leave those places unchanged.<br>
413 /// *Exception*: The refinement will often not be able to increase the minimal angle much above 30 degrees as every newly
414 /// inserted steiner point may create additional skewed triangles.
415 /// - The refinement will fullfil the *Delaunay Property*: Every triangle's circumcenter will not contain any other vertex.<br>
416 /// *Exception*: Refining with [keep_constraint_edges](RefinementParameters::keep_constraint_edges) cannot restore
417 /// the Delaunay property if doing so would require splitting a constraint edge.<br>
418 /// *Exception*: Refining with [exclude_outer_faces](RefinementParameters::exclude_outer_faces) will not
419 /// restore the Delaunay property of any outer face.
420 /// - Spade allows to specify a [maximum allowed triangle area](RefinementParameters::with_max_allowed_area).
421 /// The algorithm will attempt to subdivide any triangle with an area larger than this, independent of its smallest angle.
422 /// - Spade allows to specify a [minimum required triangle area](RefinementParameters::with_min_required_area).
423 /// The refinement will attempt to ignore any triangle with an area smaller than this parameter. This can prevent the
424 /// refinement algorithm from over-refining in some cases.
425 ///
426 /// # General limitations
427 ///
428 /// The algorithm may fail to terminate in some cases for a minimum angle limit larger than 30 degrees. Such a limit can
429 /// result in an endless loop: Every additionally inserted point creates more triangles that need to be refined.
430 ///
431 /// To prevent this, spade limits the number of additionally inserted steiner points
432 /// (see [RefinementParameters::with_max_additional_vertices]). However, this may leave the refinement in an incomplete state -
433 /// some areas of the input mesh may not have been triangulated at all, some will be overly refined.
434 /// Use [RefinementResult::refinement_complete] to identify if a refinement operation has succeeded without running out of
435 /// vertices.
436 ///
437 /// For mitigation, consider either lowering the minimum angle limit
438 /// (see [RefinementParameters::with_angle_limit]) or introduce a
439 /// [minimum required area](RefinementParameters::with_min_required_area).
440 ///
441 /// Meshes with very small input angles (angles between two constraint edges) may lead to poorly refined results.
442 /// Please consider providing a bug report if you encounter an input mesh which you think isn't refined well.
443 ///
444 /// # Stability guarantees
445 ///
446 /// While changing the interface of this method is considered to be a breaking change, changes to the specific
447 /// refinement process (e.g. which faces are split in which order) are not. Any patch release may change how
448 /// the same input mesh is being refined.
449 ///
450 /// # References
451 ///
452 /// This is an adaption of the classical refinement algorithms introduced by Jim Ruppert and Paul Chew.
453 ///
454 /// For a good introduction to the topic, refer to the slides from a short course at the thirteenth and fourteenth
455 /// International Meshing Roundtables (2005) by Jonathan Richard Shewchuk:
456 /// <https://people.eecs.berkeley.edu/~jrs/papers/imrtalk.pdf>
457 ///
458 /// Wikipedia: <https://en.wikipedia.org/wiki/Delaunay_refinement>
459 ///
460 ///
461 #[doc(alias = "Refinement")]
462 #[doc(alias = "Delaunay Refinement")]
463 pub fn refine(&mut self, parameters: RefinementParameters<V::Scalar>) -> RefinementResult {
464 use PositionInTriangulation::*;
465
466 let mut excluded_faces = if parameters.exclude_outer_faces {
467 calculate_outer_faces(self)
468 } else {
469 HashSet::new()
470 };
471
472 let mut legalize_edges_buffer = Vec::with_capacity(20);
473 let mut forcibly_split_segments_buffer = Vec::with_capacity(5);
474
475 // Maps each steiner point on an input edge onto the two vertices of that
476 // input edge. This helps in identifying when two steiner points share a common
477 // input angle
478 let mut constraint_edge_map = HashMap::new();
479
480 // Stores all edges that should be checked for encroachment
481 let mut encroached_segment_candidates =
482 VecDeque::with_capacity(self.num_constraints() + self.convex_hull_size());
483
484 encroached_segment_candidates.extend(
485 self.undirected_edges()
486 .filter(|edge| {
487 if parameters.keep_constraint_edges {
488 edge.is_part_of_convex_hull()
489 } else {
490 Self::is_fixed_edge(*edge)
491 }
492 })
493 .map(|edge| edge.fix()),
494 );
495
496 // Stores all faces that should be checked for their area and angles ("skinniness").
497 let mut skinny_triangle_candidates: VecDeque<_> = self.fixed_inner_faces().collect();
498
499 let num_initial_vertices: usize = self.num_vertices();
500 let num_additional_vertices = parameters
501 .max_additional_vertices
502 .unwrap_or(num_initial_vertices * 10);
503 let max_allowed_vertices =
504 usize::saturating_add(num_initial_vertices, num_additional_vertices);
505
506 let mut refinement_complete = true;
507
508 // Main loop of the algorithm
509 //
510 // Some terminology:
511 // - "Skinny triangle" refers to any triangle that has a minimum inner angle less than the allowed limit specified
512 // by the refinement parameters. The algorithm will attempt to insert steiner points to increase their minimal
513 // angle.
514 // - An edge is *encroached* by a point if that point lies in the diametral circle of the edge (the smallest circle
515 // fully containing the edge)
516 // - a "fixed" edge is a constraint edge or an edge of the convex hull. These are special as they may not be
517 // flipped - the input geometry must remain the same.
518 // - "input angle" is any angle between two fixed edges. Small input angles cannot be refined away as
519 // the input geometry must be kept intact.
520 // - "excluded faces" may exist if the triangulation's outer faces should not be refined. They are excluded from
521 // the third step in the main loop (see below). We don't simply delete these faces to keep the triangulation's
522 // convexity.
523 //
524 // Every iteration performs up to three checks:
525 // - First, check if any edges that must be split exists (`forcibly_split_segments_buffer`).
526 // - Second, check if any segment is encroached. If found, resolve the offending encroachment.
527 // Checking segments first makes sure that the algorithm
528 // restores the Delaunay property as quickly as possible.
529 // - Third, search for skinny triangles. Attempt to insert a new vertex at the triangle's circumcenter. If inserting
530 // such a vertex would encroach any fixed edge, add the encroached edge to the forcibly split segments buffer
531 // and revisit the face later.
532 //
533 // See method `resolve_encroachment` for more details on how step 1 and 2 manage to split edges in order to resolve
534 // an encroachment.
535 'main_loop: loop {
536 if self.num_vertices() >= max_allowed_vertices {
537 refinement_complete = false;
538 break;
539 }
540
541 // Step 1: Check for forcibly split segments.
542 if let Some(forcibly_split_segment) = forcibly_split_segments_buffer.pop() {
543 if !self.resolve_encroachment(
544 &mut encroached_segment_candidates,
545 &mut skinny_triangle_candidates,
546 &mut constraint_edge_map,
547 forcibly_split_segment,
548 &mut excluded_faces,
549 parameters.keep_constraint_edges,
550 ) {
551 // Re-processing the last skinny triangle only makes sense if the encroachment
552 // could actually be resolved. Otherwise, we'd retry the same failed split and
553 // potentially loop forever.
554 skinny_triangle_candidates.pop_back();
555 }
556 continue;
557 }
558
559 // Step 2: Check for encroached segments.
560 if let Some(segment_candidate) = encroached_segment_candidates.pop_front() {
561 // Check both adjacent faces of any candidate for encroachment.
562 for edge in segment_candidate.directed_edges() {
563 let edge = self.directed_edge(edge);
564
565 let is_excluded = edge
566 .face()
567 .as_inner()
568 .map(|face| excluded_faces.contains(&face.fix()))
569 .unwrap_or(true);
570
571 if is_excluded {
572 continue;
573 }
574
575 if let Some(opposite_position) = edge.opposite_position() {
576 if is_encroaching_edge(
577 edge.from().position(),
578 edge.to().position(),
579 opposite_position,
580 ) {
581 // The edge is encroaching
582 let _ = self.resolve_encroachment(
583 &mut encroached_segment_candidates,
584 &mut skinny_triangle_candidates,
585 &mut constraint_edge_map,
586 segment_candidate,
587 &mut excluded_faces,
588 parameters.keep_constraint_edges,
589 );
590 }
591 }
592 }
593
594 continue;
595 }
596
597 // Step 3: Take the next skinny triangle candidate
598 if let Some(face) = skinny_triangle_candidates.pop_front() {
599 if excluded_faces.contains(&face) {
600 continue;
601 }
602
603 let face = self.face(face);
604
605 let (shortest_edge, _) = face.shortest_edge();
606
607 let refinement_hint = parameters.get_refinement_hint(face);
608
609 if refinement_hint == RefinementHint::Ignore {
610 // Triangle is fine as is and can be skipped
611 continue;
612 }
613
614 if refinement_hint == RefinementHint::ShouldRefine
615 && !Self::is_fixed_edge(shortest_edge.as_undirected())
616 {
617 // Check if the shortest edge ends in two input edges that span a small
618 // input angle.
619 //
620 // Such an input angle cannot be maximized as that would require flipping at least one of its edges.
621 //
622 // See Miller, Gary; Pav, Steven; Walkington, Noel (2005). "When and why Delaunay refinement algorithms work".
623 // for more details on this idea.
624 let original_from = constraint_edge_map
625 .get(&shortest_edge.from().fix())
626 .copied();
627 let original_to = constraint_edge_map.get(&shortest_edge.to().fix()).copied();
628
629 for from_input_vertex in original_from.iter().flatten() {
630 for to_input_vertex in original_to.iter().flatten() {
631 if from_input_vertex == to_input_vertex {
632 // The two edges are input segments and join a common segment.
633 // Don't attempt to subdivide it any further, this is as good as we can get.
634 continue 'main_loop;
635 }
636 }
637 }
638 }
639
640 // Continue to resolve the skinny face
641 let circumcenter = face.circumcenter();
642
643 let locate_hint = face.vertices()[0].fix();
644
645 assert!(forcibly_split_segments_buffer.is_empty());
646 legalize_edges_buffer.clear();
647
648 // "Simulate" inserting the circumcenter by locating the insertion site and identifying which edges
649 // would need to be flipped (legalized) by the insertion. If any of these edges is fixed, an
650 // encroachment with this edge is found.
651 //
652 // First step: fill `legalize_edges_buffer` with the initial set of edges that would need to be legalized
653 // if the triangle's circumcenter would be inserted.
654 match self.locate_with_hint(circumcenter, locate_hint) {
655 OnEdge(edge) => {
656 let edge = self.directed_edge(edge);
657 if parameters.keep_constraint_edges && edge.is_constraint_edge() {
658 continue;
659 }
660
661 if edge.is_constraint_edge() {
662 // Splitting constraint edges may require updating the `excluded_faces` set.
663 // This is a little cumbersome, we'll re-use the existing implementation of edge
664 // splitting (see function resolve_encroachment).
665 forcibly_split_segments_buffer.push(edge.fix().as_undirected());
666 continue;
667 }
668
669 for edge in [edge, edge.rev()] {
670 if !edge.is_outer_edge() {
671 legalize_edges_buffer.extend([edge.next().fix(), edge.prev().fix()])
672 }
673 }
674 }
675 OnFace(face_under_circumcenter) => {
676 if excluded_faces.contains(&face_under_circumcenter) {
677 continue;
678 }
679 legalize_edges_buffer.extend(
680 self.face(face_under_circumcenter)
681 .adjacent_edges()
682 .map(|edge| edge.fix()),
683 );
684 }
685 OutsideOfConvexHull(_) => continue,
686 OnVertex(_) => continue,
687 NoTriangulation => unreachable!(),
688 };
689
690 let mut is_encroaching = false;
691
692 // Next step: Perform the regular legalization procedure by "simulating" edge flips
693 while let Some(edge) = legalize_edges_buffer.pop() {
694 let edge = self.directed_edge(edge);
695 let [from, to] = edge.as_undirected().positions();
696 if Self::is_fixed_edge(edge.as_undirected()) {
697 if is_encroaching_edge(from, to, circumcenter) {
698 // We found an encroaching edge! Makes sure that we won't attempt to
699 // insert the circumcenter.
700 is_encroaching = true;
701
702 if !parameters.keep_constraint_edges || !edge.is_constraint_edge() {
703 // New circumcenter would encroach a constraint edge. Don't insert the circumcenter
704 // but force splitting the segment
705 forcibly_split_segments_buffer.push(edge.as_undirected().fix());
706 }
707 }
708 continue; // Don't actually flip the edge as it's fixed - continue with any other edge instead.
709 }
710
711 // edge is not a fixed edge. Check if it needs to be legalized.
712 // We've already checked that this edge is not part of the convex hull - unwrap is safe
713 let opposite = edge.rev().opposite_position().unwrap();
714 let from = edge.from().position();
715 let to = edge.to().position();
716 let should_flip =
717 math::contained_in_circumference(opposite, to, from, circumcenter);
718
719 if should_flip {
720 let e1 = edge.rev().next().fix();
721 let e2 = edge.rev().prev().fix();
722
723 legalize_edges_buffer.push(e1);
724 legalize_edges_buffer.push(e2);
725 }
726 }
727
728 if !is_encroaching {
729 // The circumcenter doesn't encroach any segment. Continue really inserting it.
730 let new_vertex = self
731 .insert_with_hint(circumcenter.into(), locate_hint)
732 .expect("Failed to insert circumcenter, likely due to loss of precision. Consider refining with fewer additional vertices.");
733
734 // Add all new and changed faces to the skinny candidate list
735 skinny_triangle_candidates.extend(
736 self.vertex(new_vertex)
737 .out_edges()
738 .flat_map(|edge| edge.face().fix().as_inner()),
739 );
740 } else if !forcibly_split_segments_buffer.is_empty() {
741 // Revisit this face later. Since the encroached edge will have been split in the next iteration,
742 // inserting the circumcenter might succeed this time around.
743 skinny_triangle_candidates.push_back(face.fix());
744 }
745 } else {
746 // Done! This branch is reached if no skinny triangle could be identified anymore.
747 break;
748 }
749 }
750
751 RefinementResult {
752 excluded_faces: excluded_faces.iter().copied().collect(),
753 refinement_complete,
754 }
755 }
756
757 fn is_fixed_edge(edge: UndirectedEdgeHandle<V, DE, CdtEdge<UE>, F>) -> bool {
758 edge.is_constraint_edge() || edge.is_part_of_convex_hull()
759 }
760
761 fn resolve_encroachment(
762 &mut self,
763 encroached_segments_buffer: &mut VecDeque<FixedUndirectedEdgeHandle>,
764 encroached_faces_buffer: &mut VecDeque<FixedFaceHandle<InnerTag>>,
765 constraint_edge_map: &mut HashMap<FixedVertexHandle, [FixedVertexHandle; 2]>,
766 encroached_edge: FixedUndirectedEdgeHandle,
767 excluded_faces: &mut HashSet<FixedFaceHandle<InnerTag>>,
768 keep_constraint_edges: bool,
769 ) -> bool {
770 // Resolves an encroachment by splitting the encroached edge. Since this reduces the diametral circle, this will
771 // eventually get rid of the encroachment completely.
772 //
773 // There are a few details that make this more complicated:
774 //
775 // # Runaway encroachment
776 // Any input angle less than 45 degrees may lead to a "runaway encroachment". In such a situation, any of the
777 // angle's edges will encroach *the other* edge. This goes on forever, subdividing the edges infinitely.
778 //
779 // To work around this, spade will split edges at their center position only *once*.
780 // Any subsegment will not be split at its center position but *rounded towards the nearest power of 2*.
781 // With this behavior, neighboring edges will eventually share vertices equally far away from the offending angle's
782 // apex vertex. Points and edges in such a configuration cannot encroach each other. Refer to the original paper
783 // by Ruppert for more details.
784 //
785 // # Keeping track of which edges and faces have changed
786 // Since `resolve_encroachment` will create new edges and faces, we need to add those to the existing buffers as
787 // appropriate. This becomes a little convoluted when supporting all different refinement modes, e.g. excluded faces.
788
789 let segment = self.directed_edge(encroached_edge.as_directed());
790
791 let [v0, v1] = segment.vertices();
792
793 let half = Into::<V::Scalar>::into(0.5f32);
794
795 let v0_constraint_vertex = constraint_edge_map.get(&v0.fix()).copied();
796 let v1_constraint_vertex = constraint_edge_map.get(&v1.fix()).copied();
797
798 let (weight0, weight1) = match (v0_constraint_vertex, v1_constraint_vertex) {
799 (None, None) => {
800 // Split the segment exactly in the middle if it has not been split before.
801 (half, half)
802 }
803 _ => {
804 // One point is a steiner point, another point isn't. This will trigger rounding the distance to
805 // the nearest power of two to prevent runaway encroachment.
806
807 let half_length = segment.length_2().sqrt() * half;
808
809 let nearest_power_of_two = nearest_power_of_two(half_length);
810 let other_vertex_weight = half * nearest_power_of_two / half_length;
811 let original_vertex_weight = Into::<V::Scalar>::into(1.0) - other_vertex_weight;
812
813 if v0_constraint_vertex.is_none() {
814 // Orient the weight towards to original vertex. This makes sure that any edge participating in
815 // a runaway encroachment will end up with the same distance to the non-steiner (original) point.
816 (original_vertex_weight, other_vertex_weight)
817 } else {
818 (other_vertex_weight, original_vertex_weight)
819 }
820 }
821 };
822
823 let final_position = v0.position().mul(weight0).add(v1.position().mul(weight1));
824 if !validate_constructed_vertex(final_position, segment) {
825 return false;
826 }
827
828 let is_constraint_edge = segment.is_constraint_edge();
829 if keep_constraint_edges && is_constraint_edge {
830 return false;
831 }
832
833 let [is_left_side_excluded, is_right_side_excluded] =
834 [segment.face(), segment.rev().face()].map(|face| {
835 face.as_inner()
836 .is_some_and(|face| excluded_faces.contains(&face.fix()))
837 });
838
839 // Perform the actual split!
840 let segment = segment.fix();
841 let (v0, v1) = (v0.fix(), v1.fix());
842
843 let new_vertex = append_unconnected_vertex(self.s_mut(), final_position.into());
844 let [e1, e2] = self.insert_on_edge(segment, new_vertex);
845 self.hint_generator_mut()
846 .notify_vertex_inserted(new_vertex, final_position);
847
848 let original_vertices = v0_constraint_vertex
849 .or(v1_constraint_vertex)
850 .unwrap_or([v0, v1]);
851 constraint_edge_map.insert(new_vertex, original_vertices);
852
853 if is_constraint_edge {
854 // Make sure to update the constraint edges count as required.
855 self.handle_legal_edge_split([e1, e2]);
856 }
857
858 let (h1, h2) = (self.directed_edge(e1), self.directed_edge(e2));
859
860 if is_left_side_excluded {
861 // Any newly added face on the left becomes an excluded face
862 excluded_faces.insert(h1.face().fix().as_inner().unwrap());
863 excluded_faces.insert(h2.face().fix().as_inner().unwrap());
864 }
865
866 if is_right_side_excluded {
867 // Any newly added face on the right becomes an excluded face
868 excluded_faces.insert(h1.rev().face().fix().as_inner().unwrap());
869 excluded_faces.insert(h2.rev().face().fix().as_inner().unwrap());
870 }
871
872 self.legalize_vertex(new_vertex);
873
874 // Any of the faces that share an outgoing edge may be changed by the vertex insertion. Make sure that all of them
875 // will be revisited.
876 encroached_faces_buffer.extend(
877 self.vertex(new_vertex)
878 .out_edges()
879 .flat_map(|edge| edge.face().fix().as_inner()),
880 );
881
882 // Neighboring edges may have become encroached. Check if they need to be added to the encroached segment buffer.
883 encroached_segments_buffer.extend(
884 self.vertex(new_vertex)
885 .out_edges()
886 .filter(|edge| !edge.is_outer_edge())
887 .map(|edge| edge.next().as_undirected())
888 .filter(|edge| Self::is_fixed_edge(*edge))
889 .map(|edge| edge.fix()),
890 );
891
892 // Update encroachment candidates - any of the resulting edges may still be in an encroaching state.
893 encroached_segments_buffer.push_back(e1.as_undirected());
894 encroached_segments_buffer.push_back(e2.as_undirected());
895
896 true
897 }
898}
899
900/// Check if final_position would violate an ordering constraint. This is needed since final_position is constructed
901/// with imprecise calculations and may not even be representable in the underlying floating point type. In rare cases,
902/// this means that the newly formed triangles would not be ordered ccw.
903/// We'll simply skip these refinements steps as it should only happen for very bad input geometries.
904///
905/// Before (v0 = segment.from(), v1 = segment.to()):
906/// v2
907/// / \
908/// v0 --> v1
909/// \ /
910/// v3
911///
912/// After (before legalizing) - return if any face would be ordered cw
913/// v2
914/// / | \
915/// v0 -v-> v1
916/// \ | /
917/// v3
918fn validate_constructed_vertex<V, DE, UE, F>(
919 final_position: Point2<V::Scalar>,
920 segment: DirectedEdgeHandle<V, DE, UE, F>,
921) -> bool
922where
923 V: HasPosition,
924{
925 use math::is_ordered_ccw;
926 let [v0, v1] = segment.positions();
927
928 if math::validate_vertex(&final_position).is_err() {
929 return false;
930 }
931
932 if let Some(v2) = segment.opposite_position() {
933 if is_ordered_ccw(v0, v2, final_position) || is_ordered_ccw(v2, v1, final_position) {
934 return false;
935 }
936 }
937
938 if let Some(v3) = segment.rev().opposite_position() {
939 if is_ordered_ccw(v3, v0, final_position) || is_ordered_ccw(v1, v3, final_position) {
940 return false;
941 }
942 }
943
944 true
945}
946
947fn is_encroaching_edge<S: SpadeNum + Float>(
948 edge_from: Point2<S>,
949 edge_to: Point2<S>,
950 query_point: Point2<S>,
951) -> bool {
952 let edge_center = edge_from.add(edge_to).mul(0.5f32.into());
953 let radius_2 = edge_from.distance_2(edge_to) * 0.25.into();
954
955 query_point.distance_2(edge_center) < radius_2
956}
957
958fn nearest_power_of_two<S: Float + SpadeNum>(input: S) -> S {
959 input.log2().round().exp2()
960}
961
962fn calculate_outer_faces<V: HasPosition, DE: Default, UE: Default, F: Default, L>(
963 triangulation: &ConstrainedDelaunayTriangulation<V, DE, UE, F, L>,
964) -> HashSet<FixedFaceHandle<InnerTag>>
965where
966 L: HintGenerator<<V as HasPosition>::Scalar>,
967{
968 if triangulation.all_vertices_on_line() {
969 return HashSet::new();
970 }
971
972 // Determine excluded faces by "peeling of" outer layers and adding them to an outer layer set.
973 // This needs to be done repeatedly to also get inner "holes" within the triangulation
974 let mut inner_faces = HashSet::new();
975 let mut outer_faces = HashSet::new();
976
977 let mut current_todo_list: Vec<_> =
978 triangulation.convex_hull().map(|edge| edge.rev()).collect();
979 let mut next_todo_list = Vec::new();
980
981 let mut return_outer_faces = true;
982
983 loop {
984 // Every iteration of the outer while loop will peel of the outmost layer and pre-populate the
985 // next, inner layer.
986 while let Some(next_edge) = current_todo_list.pop() {
987 let (list, face_set) = if next_edge.is_constraint_edge() {
988 // We're crossing a constraint edge - add that face to the *next* todo list
989 (&mut next_todo_list, &mut inner_faces)
990 } else {
991 (&mut current_todo_list, &mut outer_faces)
992 };
993
994 if let Some(inner) = next_edge.face().as_inner() {
995 if face_set.insert(inner.fix()) {
996 list.push(next_edge.prev().rev());
997 list.push(next_edge.next().rev());
998 }
999 }
1000 }
1001
1002 if next_todo_list.is_empty() {
1003 break;
1004 }
1005 core::mem::swap(&mut inner_faces, &mut outer_faces);
1006 core::mem::swap(&mut next_todo_list, &mut current_todo_list);
1007
1008 return_outer_faces = !return_outer_faces;
1009 }
1010
1011 if return_outer_faces {
1012 outer_faces
1013 } else {
1014 inner_faces
1015 }
1016}
1017
1018#[cfg(test)]
1019mod test {
1020 use super::HashSet;
1021 use alloc::{vec, vec::Vec};
1022
1023 use crate::{
1024 handles::FixedVertexHandle,
1025 test_utilities::{random_points_with_seed, SEED},
1026 AngleLimit, ConstrainedDelaunayTriangulation, InsertionError, Point2, RefinementParameters,
1027 Triangulation as _,
1028 };
1029
1030 pub type Cdt = ConstrainedDelaunayTriangulation<Point2<f64>>;
1031
1032 #[test]
1033 fn test_zero_angle_limit_dbg() {
1034 let limit = AngleLimit::from_deg(0.0);
1035 let debug_string = alloc::format!("{:?}", limit);
1036 assert_eq!(debug_string, "AngleLimit { angle limit (deg): 0.0 }");
1037 }
1038
1039 #[test]
1040 fn test_zero_angle_limit() -> Result<(), InsertionError> {
1041 let limit = AngleLimit::from_deg(0.0);
1042
1043 assert_eq!(limit.radius_to_shortest_edge_limit(), f64::INFINITY);
1044
1045 let mut vertices = random_points_with_seed(20, SEED);
1046
1047 // Insert an artificial outer boundary that will prevent the convex hull from being encroached.
1048 // This should prevent any refinement.
1049 vertices.push(Point2::new(100.0, 100.0));
1050 vertices.push(Point2::new(100.0, 0.0));
1051 vertices.push(Point2::new(100.0, -100.0));
1052 vertices.push(Point2::new(0.0, -100.0));
1053 vertices.push(Point2::new(-100.0, -100.0));
1054 vertices.push(Point2::new(-100.0, 0.0));
1055 vertices.push(Point2::new(-100.0, 100.0));
1056 vertices.push(Point2::new(0.0, 100.0));
1057
1058 let mut cdt = Cdt::bulk_load(vertices)?;
1059
1060 let initial_num_vertices = cdt.num_vertices();
1061 cdt.refine(RefinementParameters::new().with_angle_limit(limit));
1062
1063 assert_eq!(initial_num_vertices, cdt.num_vertices());
1064
1065 cdt.refine(RefinementParameters::new());
1066 assert!(initial_num_vertices < cdt.num_vertices());
1067
1068 Ok(())
1069 }
1070
1071 #[test]
1072 fn test_nearest_power_of_two() {
1073 use super::nearest_power_of_two;
1074
1075 for i in 1..400u32 {
1076 let log = (i as f64).log2() as u32;
1077 let nearest = nearest_power_of_two(i as f64);
1078
1079 assert!((1 << log) as f64 == nearest || (1 << (log + 1)) as f64 == nearest);
1080 }
1081
1082 assert_eq!(0.25, nearest_power_of_two(0.25));
1083 assert_eq!(0.25, nearest_power_of_two(0.27));
1084 assert_eq!(0.5, nearest_power_of_two(0.5));
1085 assert_eq!(1.0, nearest_power_of_two(0.75));
1086 assert_eq!(2.0, nearest_power_of_two(1.5));
1087 assert_eq!(2.0, nearest_power_of_two(2.5));
1088 assert_eq!(4.0, nearest_power_of_two(3.231));
1089 assert_eq!(4.0, nearest_power_of_two(4.0));
1090 }
1091
1092 #[test]
1093 fn test_small_refinement() -> Result<(), InsertionError> {
1094 let vertices = random_points_with_seed(20, SEED);
1095 let mut cdt = Cdt::bulk_load(vertices)?;
1096
1097 let mut peekable = cdt.fixed_vertices().peekable();
1098 while let (Some(p0), Some(p1)) = (peekable.next(), peekable.peek()) {
1099 if !cdt.can_add_constraint(p0, *p1) {
1100 cdt.add_constraint(p0, *p1);
1101 }
1102 }
1103
1104 cdt.refine(Default::default());
1105
1106 cdt.cdt_sanity_check_with_params(false);
1107
1108 Ok(())
1109 }
1110
1111 #[test]
1112 fn test_sharp_angle_refinement() -> Result<(), InsertionError> {
1113 let mut cdt = Cdt::new();
1114
1115 // This sharp angle should only be subdivided once rather than infinitely often.
1116 cdt.add_constraint_edge(Point2::new(-40.0, 80.0), Point2::new(0.0, 0.0))?;
1117 cdt.add_constraint_edge(Point2::new(0.0, 0.0), Point2::new(4.0, 90.0))?;
1118
1119 let refinement_parameters = RefinementParameters::new()
1120 .with_max_additional_vertices(10)
1121 .with_angle_limit(AngleLimit::from_radius_to_shortest_edge_ratio(1.0));
1122
1123 let result = cdt.refine(refinement_parameters);
1124
1125 assert!(result.refinement_complete);
1126 cdt.cdt_sanity_check();
1127 Ok(())
1128 }
1129
1130 #[test]
1131 fn test_simple_exclude_outer_faces() -> Result<(), InsertionError> {
1132 let mut cdt = Cdt::new();
1133 let vertices = [
1134 Point2::new(-1.0, 1.0),
1135 Point2::new(0.0, 0.5),
1136 Point2::new(1.0, 1.0),
1137 Point2::new(1.0, -1.0),
1138 Point2::new(-1.0, -1.0),
1139 Point2::new(-1.0, 1.0),
1140 ];
1141
1142 for points in vertices.windows(2) {
1143 let p1 = points[0];
1144 let p2 = points[1];
1145 cdt.add_constraint_edge(p1, p2)?;
1146 }
1147
1148 let excluded_faces = super::calculate_outer_faces(&cdt);
1149
1150 let v2 = cdt.locate_vertex(Point2::new(1.0, 1.0)).unwrap().fix();
1151 let v0 = cdt.locate_vertex(Point2::new(-1.0, 1.0)).unwrap().fix();
1152
1153 let excluded_face = cdt
1154 .get_edge_from_neighbors(v2, v0)
1155 .and_then(|edge| edge.face().as_inner())
1156 .unwrap()
1157 .fix();
1158
1159 assert_eq!(
1160 excluded_faces,
1161 HashSet::from_iter(core::iter::once(excluded_face))
1162 );
1163
1164 Ok(())
1165 }
1166
1167 fn test_shape() -> [Point2<f64>; 6] {
1168 [
1169 Point2::new(-1.0, 1.0),
1170 Point2::new(0.0, 0.7),
1171 Point2::new(1.0, 1.0),
1172 Point2::new(2.0, 0.0),
1173 Point2::new(0.5, -1.0),
1174 Point2::new(-0.5, -2.0),
1175 ]
1176 }
1177
1178 #[test]
1179 fn test_exclude_outer_faces() -> Result<(), InsertionError> {
1180 let mut cdt = Cdt::new();
1181
1182 let scale_factors = [1.0, 2.0, 3.0, 4.0];
1183
1184 let mut current_excluded_faces = super::calculate_outer_faces(&cdt);
1185
1186 assert!(current_excluded_faces.is_empty());
1187
1188 for factor in scale_factors {
1189 cdt.add_constraint_edges(test_shape().iter().map(|p| p.mul(factor)), true)?;
1190
1191 let next_excluded_faces = super::calculate_outer_faces(&cdt);
1192
1193 assert!(current_excluded_faces.len() < next_excluded_faces.len());
1194
1195 current_excluded_faces = next_excluded_faces;
1196 assert!(current_excluded_faces.len() < cdt.num_inner_faces());
1197 }
1198
1199 Ok(())
1200 }
1201
1202 #[test]
1203 fn test_exclude_outer_faces_with_non_closed_mesh() -> Result<(), InsertionError> {
1204 let mut cdt = Cdt::new();
1205 cdt.add_constraint_edges(test_shape(), false)?;
1206
1207 let refinement_result = cdt.refine(
1208 RefinementParameters::new()
1209 .exclude_outer_faces(true)
1210 .with_angle_limit(AngleLimit::from_radius_to_shortest_edge_ratio(
1211 f64::INFINITY,
1212 )),
1213 );
1214
1215 assert_eq!(
1216 refinement_result.excluded_faces.len(),
1217 cdt.num_inner_faces()
1218 );
1219
1220 cdt.refine(RefinementParameters::new().exclude_outer_faces(true));
1221
1222 Ok(())
1223 }
1224
1225 #[test]
1226 fn test_failing_refinement() -> Result<(), InsertionError> {
1227 // f32 is important - only then, rounding errors will lead to violating the ccw property when
1228 // the refinement splits an edge. See issue #96
1229 let mut cdt = ConstrainedDelaunayTriangulation::<Point2<f32>>::new();
1230
1231 #[rustfmt::skip]
1232 let vert = [
1233 [Point2 { x: -50.023544f32, y: -25.29227 }, Point2 { x: 754.23883, y: -25.29227 }],
1234 [Point2 { x: 754.23883, y: 508.68994 }, Point2 { x: -50.023544, y: 508.68994 }],
1235 [Point2 { x: -44.20742, y: 316.5185 }, Point2 { x: -50.023544, y: 318.19534 }],
1236 [Point2 { x: 11.666269, y: 339.4947 }, Point2 { x: 15.110367, y: 335.44138 }],
1237 [Point2 { x: 335.06403, y: 122.91455 }, Point2 { x: 340.15436, y: 132.04283 }],
1238 [Point2 { x: 446.92468, y: -6.7025666 }, Point2 { x: 458.70944, y: 14.341333 }],
1239 [Point2 { x: 458.70944, y: 14.341333 }, Point2 { x: 471.58313, y: 7.1453195 }],
1240 [Point2 { x: 467.80966, y: 0.40460825 }, Point2 { x: 468.6454, y: -0.061800003 }],
1241 [Point2 { x: 464.55957, y: -7.3688636 }, Point2 { x: 465.48816, y: -7.890797 }],
1242 [Point2 { x: 465.48816, y: -7.890797 }, Point2 { x: 461.57117, y: -14.898027 }],
1243 [Point2 { x: 465.42877, y: 10.587858 }, Point2 { x: 453.93112, y: 17.01763 }],
1244 ];
1245
1246 for [start, end] in vert {
1247 cdt.add_constraint_edge(start, end)?;
1248 }
1249
1250 cdt.refine(Default::default());
1251 cdt.cdt_sanity_check();
1252
1253 Ok(())
1254 }
1255
1256 #[test]
1257 fn repro_refine_hang_regression() -> Result<(), InsertionError> {
1258 let points = repro_refine_hang_points();
1259 let constraints = repro_refine_hang_constraints();
1260
1261 let mut cdt = Cdt::bulk_load_cdt(points, vec![])?;
1262
1263 for edge in &constraints {
1264 cdt.add_constraint_and_split(
1265 FixedVertexHandle::from_index(edge[0]),
1266 FixedVertexHandle::from_index(edge[1]),
1267 |v| v,
1268 );
1269 }
1270
1271 let initial_num_constraints = cdt.num_constraints();
1272
1273 let result = cdt.refine(
1274 RefinementParameters::<f64>::new()
1275 .keep_constraint_edges()
1276 .with_max_allowed_area(0.4 * 0.32 * 0.32)
1277 .with_max_additional_vertices(1_000_000),
1278 );
1279
1280 assert!(
1281 result.refinement_complete,
1282 "Refinement ran out of additional vertices for regression input"
1283 );
1284 cdt.cdt_sanity_check();
1285
1286 assert_eq!(initial_num_constraints, cdt.num_constraints());
1287
1288 Ok(())
1289 }
1290
1291 fn repro_refine_hang_points() -> Vec<Point2<f64>> {
1292 vec![
1293 Point2::new(0.0, 0.0),
1294 Point2::new(1.0, 0.0),
1295 Point2::new(1.0, 1.0),
1296 Point2::new(0.0, 1.0),
1297 Point2::new(0.5, 0.5),
1298 Point2::new(0.5040435515666191, 0.48484848484848486),
1299 Point2::new(0.5080871031332381, 0.4696969696969697),
1300 Point2::new(0.5121306546998572, 0.45454545454545453),
1301 Point2::new(0.5161742062664761, 0.4393939393939394),
1302 Point2::new(0.5202177578330952, 0.42424242424242425),
1303 Point2::new(0.5242613093997143, 0.40909090909090906),
1304 Point2::new(0.5283048609663333, 0.3939393939393939),
1305 Point2::new(0.5323484125329524, 0.3787878787878788),
1306 Point2::new(0.5363919640995715, 0.36363636363636365),
1307 Point2::new(0.5404355156661904, 0.3484848484848485),
1308 Point2::new(0.5444790672328095, 0.33333333333333337),
1309 Point2::new(0.5485226187994285, 0.3181818181818182),
1310 Point2::new(0.5525661703660476, 0.30303030303030304),
1311 Point2::new(0.5566097219326667, 0.2878787878787879),
1312 Point2::new(0.5606532734992857, 0.2727272727272727),
1313 Point2::new(0.5646968250659048, 0.25757575757575757),
1314 Point2::new(0.5687403766325237, 0.24242424242424243),
1315 Point2::new(0.5727839281991428, 0.2272727272727273),
1316 Point2::new(0.5768274797657619, 0.2121212121212121),
1317 Point2::new(0.5808710313323809, 0.19696969696969702),
1318 Point2::new(0.584914582899, 0.18181818181818188),
1319 Point2::new(0.5889581344656191, 0.16666666666666669),
1320 Point2::new(0.593001686032238, 0.15151515151515155),
1321 Point2::new(0.5970452375988571, 0.13636363636363635),
1322 Point2::new(0.6010887891654761, 0.12121212121212122),
1323 Point2::new(0.6051323407320952, 0.10606060606060608),
1324 Point2::new(0.6091758922987143, 0.09090909090909088),
1325 Point2::new(0.6132194438653333, 0.0757575757575758),
1326 Point2::new(0.6172629954319524, 0.06060606060606066),
1327 Point2::new(0.6213065469985714, 0.04545454545454547),
1328 Point2::new(0.6253500985651904, 0.030303030303030387),
1329 Point2::new(0.6289250162663348, 0.01690752419363295),
1330 Point2::new(0.6334372016984285, 5.551115123125783e-17),
1331 Point2::new(0.5156817958365265, 0.5),
1332 Point2::new(0.531363591673053, 0.5),
1333 Point2::new(0.5470453875095795, 0.5),
1334 Point2::new(0.562727183346106, 0.5),
1335 Point2::new(0.5784089791826326, 0.5),
1336 Point2::new(0.5940907750191591, 0.5),
1337 Point2::new(0.6097725708556856, 0.5),
1338 Point2::new(0.6254543666922121, 0.5),
1339 Point2::new(0.6411361625287386, 0.5),
1340 Point2::new(0.6568179583652651, 0.5),
1341 Point2::new(0.6724997542017918, 0.5),
1342 Point2::new(0.6881815500383183, 0.5),
1343 Point2::new(0.7038633458748448, 0.5),
1344 Point2::new(0.7195451417113713, 0.5),
1345 Point2::new(0.7352269375478978, 0.5),
1346 Point2::new(0.7509087333844243, 0.5),
1347 Point2::new(0.7665905292209508, 0.5),
1348 Point2::new(0.7822723250574773, 0.5),
1349 Point2::new(0.797954120894004, 0.5),
1350 Point2::new(0.8136359167305304, 0.5),
1351 Point2::new(0.8293177125670569, 0.5),
1352 Point2::new(0.8449995084035835, 0.5),
1353 Point2::new(0.86068130424011, 0.5),
1354 Point2::new(0.8763631000766365, 0.5),
1355 Point2::new(0.892044895913163, 0.5),
1356 Point2::new(0.9077266917496896, 0.5),
1357 Point2::new(0.9234084875862161, 0.5),
1358 Point2::new(0.9390902834227426, 0.5),
1359 Point2::new(0.9547720792592691, 0.5),
1360 Point2::new(0.9704538750957956, 0.5),
1361 Point2::new(0.9861356709323221, 0.5),
1362 Point2::new(1.0, 0.5),
1363 ]
1364 }
1365
1366 fn repro_refine_hang_constraints() -> Vec<[usize; 2]> {
1367 vec![
1368 [4, 5],
1369 [5, 6],
1370 [6, 7],
1371 [7, 8],
1372 [8, 9],
1373 [9, 10],
1374 [10, 11],
1375 [11, 12],
1376 [12, 13],
1377 [13, 14],
1378 [14, 15],
1379 [15, 16],
1380 [16, 17],
1381 [17, 18],
1382 [18, 19],
1383 [19, 20],
1384 [20, 21],
1385 [21, 22],
1386 [22, 23],
1387 [23, 24],
1388 [24, 25],
1389 [25, 26],
1390 [26, 27],
1391 [27, 28],
1392 [28, 29],
1393 [29, 30],
1394 [30, 31],
1395 [31, 32],
1396 [32, 33],
1397 [33, 34],
1398 [34, 35],
1399 [35, 36],
1400 [36, 37],
1401 [4, 38],
1402 [38, 39],
1403 [39, 40],
1404 [40, 41],
1405 [41, 42],
1406 [42, 43],
1407 [43, 44],
1408 [44, 45],
1409 [45, 46],
1410 [46, 47],
1411 [47, 48],
1412 [48, 49],
1413 [49, 50],
1414 [50, 51],
1415 [51, 52],
1416 [52, 53],
1417 [53, 54],
1418 [54, 55],
1419 [55, 56],
1420 [56, 57],
1421 [57, 58],
1422 [58, 59],
1423 [59, 60],
1424 [60, 61],
1425 [61, 62],
1426 [62, 63],
1427 [63, 64],
1428 [64, 65],
1429 [65, 66],
1430 [66, 67],
1431 [67, 68],
1432 [68, 69],
1433 ]
1434 }
1435}