honeycomb_kernels/grid_generation/
mod.rs

1//! grid generation implementation
2//!
3//! This module contains the routines and interface used to implement 2D and 3D grid generation.
4
5mod internals;
6
7#[cfg(test)]
8mod tests;
9
10use honeycomb_core::attributes::AttributeBind;
11use honeycomb_core::cmap::{CMap2, CMap3, CMapBuilder};
12use honeycomb_core::geometry::{CoordsFloat, Vertex2, Vertex3};
13
14/// # Builder-level error enum
15///
16/// This enum is used to describe all non-panic errors that can occur when using the [`GridBuilder`]
17/// structure.
18#[derive(thiserror::Error, Debug, PartialEq, Eq)]
19pub enum GridBuilderError {
20    // grid-related variants
21    /// One or multiple of the specified grid characteristics are invalid.
22    #[error("invalid grid parameters - {0}")]
23    InvalidGridParameters(&'static str),
24    /// The builder is missing one or multiple parameters to generate the grid.
25    #[error("insufficient parameters - please specify at least 2")]
26    MissingGridParameters,
27}
28
29/// # Grid builder structure
30///
31/// ## Generics
32///
33/// - `const D: usize` -- Dimension of the grid. Should be 2 or 3.
34/// - `T: CoordsFloat` -- Generic FP type that will be used by the map's vertices.
35///
36/// ## Example
37///
38/// ```rust
39/// # use honeycomb_kernels::grid_generation::{GridBuilderError};
40/// # fn main() -> Result<(), GridBuilderError> {
41/// use honeycomb_core::cmap::CMap3;
42/// use honeycomb_kernels::grid_generation::{GridBuilder};
43///
44/// let map: CMap3<f64> = GridBuilder::<3, f64>::default()
45///     .n_cells([2, 3, 2])
46///     .len_per_cell([1.0; 3])
47///     .build()?;
48///
49/// # Ok(())
50/// # }
51/// ```
52pub struct GridBuilder<const D: usize, T: CoordsFloat> {
53    pub(crate) map_builder: CMapBuilder<D>,
54    pub(crate) origin: [T; D],
55    pub(crate) n_cells: Option<[usize; D]>,
56    pub(crate) len_per_cell: Option<[T; D]>,
57    pub(crate) lens: Option<[T; D]>,
58    pub(crate) split_cells: bool,
59}
60
61impl<const D: usize, T: CoordsFloat> Default for GridBuilder<D, T> {
62    fn default() -> Self {
63        Self {
64            map_builder: CMapBuilder::<D>::from_n_darts(1),
65            origin: [T::zero(); D],
66            n_cells: None,
67            len_per_cell: None,
68            lens: None,
69            split_cells: false,
70        }
71    }
72}
73
74impl<const D: usize, T: CoordsFloat> GridBuilder<D, T> {
75    /// Set values for all dimensions
76    #[must_use = "unused builder object"]
77    pub fn n_cells(mut self, n_cells: [usize; D]) -> Self {
78        self.n_cells = Some(n_cells);
79        self
80    }
81
82    /// Set values for all dimensions
83    #[must_use = "unused builder object"]
84    pub fn len_per_cell(mut self, len_per_cell: [T; D]) -> Self {
85        self.len_per_cell = Some(len_per_cell);
86        self
87    }
88
89    /// Set values for all dimensions
90    #[must_use = "unused builder object"]
91    pub fn lens(mut self, lens: [T; D]) -> Self {
92        self.lens = Some(lens);
93        self
94    }
95
96    /// Set origin (most bottom-left vertex) of the grid
97    #[must_use = "unused builder object"]
98    pub fn origin(mut self, origin: [T; D]) -> Self {
99        self.origin = origin;
100        self
101    }
102
103    /// Indicate whether to split cells of the grid
104    ///
105    /// In 2D, this will result in triangular cells.
106    ///
107    /// In 3D, this will result in tetrahedral cells.
108    #[must_use = "unused builder object"]
109    pub fn split_cells(mut self, split: bool) -> Self {
110        self.split_cells = split;
111        self
112    }
113
114    /// Add the attribute `A` to the attributes the created map will contain.
115    ///
116    /// # Usage
117    ///
118    /// Each attribute must be uniquely typed, i.e. a single type or struct cannot be added twice
119    /// to the builder / map. This includes type aliases as these are not distinct from the
120    /// compiler's perspective.
121    ///
122    /// If you have multiple attributes that are represented using the same data type, you may want
123    /// to look into the **Newtype** pattern
124    /// [here](https://rust-unofficial.github.io/patterns/patterns/behavioural/newtype.html)
125    /// and [here](https://doc.rust-lang.org/rust-by-example/generics/new_types.html)
126    #[must_use = "unused builder object"]
127    pub fn add_attribute<A: AttributeBind + 'static>(mut self) -> Self {
128        let tmp = self.map_builder.add_attribute::<A>();
129        self.map_builder = tmp;
130        self
131    }
132}
133
134// --- parsing routine
135
136macro_rules! check_parameters {
137    ($id: ident, $msg: expr) => {
138        if $id.is_sign_negative() | $id.is_zero() {
139            return Err(GridBuilderError::InvalidGridParameters($msg));
140        }
141    };
142}
143
144impl<T: CoordsFloat> GridBuilder<2, T> {
145    /// Parse provided grid parameters to provide what's used to actually generate the grid.
146    #[allow(clippy::type_complexity)]
147    pub(crate) fn parse_2d(
148        self,
149    ) -> Result<(CMapBuilder<2>, Vertex2<T>, [usize; 2], [T; 2]), GridBuilderError> {
150        match (self.n_cells, self.len_per_cell, self.lens) {
151            // from # cells and lengths per cell
152            (Some([nx, ny]), Some([lpx, lpy]), lens) => {
153                if lens.is_some() {
154                    eprintln!(
155                        "W: All three grid parameters were specified, total lengths will be ignored"
156                    );
157                }
158                #[rustfmt::skip]
159                check_parameters!(lpx, "length per x cell is null or negative");
160                #[rustfmt::skip]
161                check_parameters!(lpy, "length per y cell is null or negative");
162                Ok((
163                    self.map_builder,
164                    Vertex2(self.origin[0], self.origin[1]),
165                    [nx, ny],
166                    [lpx, lpy],
167                ))
168            }
169            // from # cells and total lengths
170            (Some([nx, ny]), None, Some([lx, ly])) => {
171                #[rustfmt::skip]
172                check_parameters!(lx, "grid length along x is null or negative");
173                #[rustfmt::skip]
174                check_parameters!(ly, "grid length along y is null or negative");
175                Ok((
176                    self.map_builder,
177                    Vertex2(self.origin[0], self.origin[1]),
178                    [nx, ny],
179                    [lx / T::from(nx).unwrap(), ly / T::from(ny).unwrap()],
180                ))
181            }
182            // from lengths per cell and total lengths
183            (None, Some([lpx, lpy]), Some([lx, ly])) => {
184                #[rustfmt::skip]
185                check_parameters!(lpx, "length per x cell is null or negative");
186                #[rustfmt::skip]
187                check_parameters!(lpy, "length per y cell is null or negative");
188                #[rustfmt::skip]
189                check_parameters!(lx, "grid length along x is null or negative");
190                #[rustfmt::skip]
191                check_parameters!(ly, "grid length along y is null or negative");
192                Ok((
193                    self.map_builder,
194                    Vertex2(self.origin[0], self.origin[1]),
195                    [
196                        (lx / lpx).ceil().to_usize().unwrap(),
197                        (ly / lpy).ceil().to_usize().unwrap(),
198                    ],
199                    [lpx, lpy],
200                ))
201            }
202            (_, _, _) => Err(GridBuilderError::MissingGridParameters),
203        }
204    }
205}
206
207impl<T: CoordsFloat> GridBuilder<3, T> {
208    /// Parse provided grid parameters to provide what's used to actually generate the grid.
209    #[allow(clippy::type_complexity)]
210    pub(crate) fn parse_3d(
211        self,
212    ) -> Result<(CMapBuilder<3>, Vertex3<T>, [usize; 3], [T; 3]), GridBuilderError> {
213        match (self.n_cells, self.len_per_cell, self.lens) {
214            // from # cells and lengths per cell
215            (Some([nx, ny, nz]), Some([lpx, lpy, lpz]), lens) => {
216                if lens.is_some() {
217                    eprintln!(
218                        "W: All three grid parameters were specified, total lengths will be ignored"
219                    );
220                }
221                #[rustfmt::skip]
222                check_parameters!(lpx, "length per x cell is null or negative");
223                #[rustfmt::skip]
224                check_parameters!(lpy, "length per y cell is null or negative");
225                #[rustfmt::skip]
226                check_parameters!(lpz, "length per z cell is null or negative");
227                Ok((
228                    self.map_builder,
229                    Vertex3(self.origin[0], self.origin[1], self.origin[2]),
230                    [nx, ny, nz],
231                    [lpx, lpy, lpz],
232                ))
233            }
234            // from # cells and total lengths
235            (Some([nx, ny, nz]), None, Some([lx, ly, lz])) => {
236                #[rustfmt::skip]
237                check_parameters!(lx, "grid length along x is null or negative");
238                #[rustfmt::skip]
239                check_parameters!(ly, "grid length along y is null or negative");
240                #[rustfmt::skip]
241                check_parameters!(lz, "grid length along z is null or negative");
242                Ok((
243                    self.map_builder,
244                    Vertex3(self.origin[0], self.origin[1], self.origin[2]),
245                    [nx, ny, nz],
246                    [
247                        lx / T::from(nx).unwrap(),
248                        ly / T::from(ny).unwrap(),
249                        lz / T::from(nz).unwrap(),
250                    ],
251                ))
252            }
253            // from lengths per cell and total lengths
254            (None, Some([lpx, lpy, lpz]), Some([lx, ly, lz])) => {
255                #[rustfmt::skip]
256                check_parameters!(lpx, "length per x cell is null or negative");
257                #[rustfmt::skip]
258                check_parameters!(lpy, "length per y cell is null or negative");
259                #[rustfmt::skip]
260                check_parameters!(lpz, "length per z cell is null or negative");
261                #[rustfmt::skip]
262                check_parameters!(lx, "grid length along x is null or negative");
263                #[rustfmt::skip]
264                check_parameters!(ly, "grid length along y is null or negative");
265                #[rustfmt::skip]
266                check_parameters!(lz, "grid length along z is null or negative");
267                Ok((
268                    self.map_builder,
269                    Vertex3(self.origin[0], self.origin[1], self.origin[2]),
270                    [
271                        (lx / lpx).ceil().to_usize().unwrap(),
272                        (ly / lpy).ceil().to_usize().unwrap(),
273                        (lz / lpz).ceil().to_usize().unwrap(),
274                    ],
275                    [lpx, lpy, lpz],
276                ))
277            }
278            (_, _, _) => Err(GridBuilderError::MissingGridParameters),
279        }
280    }
281}
282
283impl<T: CoordsFloat> GridBuilder<2, T> {
284    #[allow(clippy::missing_panics_doc)]
285    /// Create a combinatorial map representing a 2D orthogonal grid.
286    ///
287    /// The map generated by this pre-definite value corresponds to an orthogonal mesh, with an
288    /// equal number of cells along each axis:
289    ///
290    /// ![`CMAP2_GRID`](https://lihpc-computational-geometry.github.io/honeycomb/user-guide/images/bg_grid.svg)
291    #[must_use = "unused builder object"]
292    pub fn unit_grid(n_cells_per_axis: usize) -> CMap2<T> {
293        GridBuilder::default()
294            .n_cells([n_cells_per_axis; 2])
295            .len_per_cell([T::one(); 2])
296            .build()
297            .expect("E: unreachable")
298    }
299
300    #[allow(clippy::missing_panics_doc)]
301    /// Create a combinatorial map representing a 2D orthogonal grid.
302    ///
303    /// The map generated by this pre-definite value corresponds to an orthogonal mesh, with an
304    /// equal number of cells along each axis. Each cell is split diagonally (top left to
305    /// bottom right) to form triangles:
306    ///
307    /// ![`CMAP2_GRID`](https://lihpc-computational-geometry.github.io/honeycomb/user-guide/images/bg_grid_tri.svg)
308    #[must_use = "unused builder object"]
309    pub fn unit_triangles(n_square: usize) -> CMap2<T> {
310        GridBuilder::default()
311            .n_cells([n_square; 2])
312            .len_per_cell([T::one(); 2])
313            .split_cells(true)
314            .build()
315            .expect("E: unreachable")
316    }
317
318    #[allow(clippy::missing_errors_doc)]
319    /// Consumes the builder and produce a combinatorial map object.
320    ///
321    /// This method is only available for `D == 2` or `D == 3`.
322    ///
323    /// # Return / Errors
324    ///
325    /// This method return a `Result` taking the following values:
326    /// - `Ok(map: CMap2<T>)` if generation was successful,
327    /// - `Err(GridBuilderError)` otherwise. See [`GridBuilderError`] for possible failures.
328    ///
329    /// # Panics
330    ///
331    /// This method may panic if type casting goes wrong during parameters parsing.
332    pub fn build(self) -> Result<CMap2<T>, GridBuilderError> {
333        let split = self.split_cells;
334        self.parse_2d().map(|(builder, origin, ns, lens)| {
335            if split {
336                internals::build_2d_splitgrid(builder, origin, ns, lens)
337            } else {
338                internals::build_2d_grid(builder, origin, ns, lens)
339            }
340        })
341    }
342}
343
344impl<T: CoordsFloat> GridBuilder<3, T> {
345    #[allow(clippy::missing_panics_doc)]
346    /// Create a combinatorial map representing a 3D orthogonal grid.
347    ///
348    /// The map generated by this pre-definite value corresponds to an orthogonal mesh, with an
349    /// equal number of cells along each axis:
350    ///
351    /// ![`CMAP2_GRID`](https://lihpc-computational-geometry.github.io/honeycomb/user-guide/images/hex.svg)
352    pub fn hex_grid(n_cells_per_axis: usize, cell_length: T) -> CMap3<T> {
353        GridBuilder::default()
354            .n_cells([n_cells_per_axis; 3])
355            .len_per_cell([cell_length; 3])
356            .build()
357            .expect("E: unreachable")
358    }
359
360    #[allow(clippy::missing_panics_doc)]
361    /// Create a combinatorial map representing a 3D orthogonal grid.
362    ///
363    /// The map generated by this pre-definite value corresponds to an orthogonal mesh, with an
364    /// equal number of cells along each axis. each hexahedral cell is cut into 5 tetrahedra;
365    /// this pattern is repeated with added symmetry to produce a conformal mesh:
366    ///
367    /// ![`CMAP2_GRID`](https://lihpc-computational-geometry.github.io/honeycomb/user-guide/images/tet.svg)
368    #[must_use = "unused builder object"]
369    pub fn tet_grid(n_cells_per_axis: usize, cell_length: T) -> CMap3<T> {
370        GridBuilder::default()
371            .n_cells([n_cells_per_axis; 3])
372            .len_per_cell([cell_length; 3])
373            .split_cells(true)
374            .build()
375            .expect("E: unreachable")
376    }
377
378    #[allow(clippy::missing_errors_doc)]
379    /// Consumes the builder and produce a combinatorial map object.
380    ///
381    /// This method is only available for `D == 2` or `D == 3`.
382    ///
383    /// # Return / Errors
384    ///
385    /// This method return a `Result` taking the following values:
386    /// - `Ok(map: CMap3<T>)` if generation was successful,
387    /// - `Err(GridBuilderError)` otherwise. See [`GridBuilderError`] for possible failures.
388    ///
389    /// # Panics
390    ///
391    /// This method may panic if type casting goes wrong during parameters parsing.
392    pub fn build(self) -> Result<CMap3<T>, GridBuilderError> {
393        let split = self.split_cells;
394        self.parse_3d().map(|(builder, origin, ns, lens)| {
395            if split {
396                internals::build_3d_tetgrid(builder, origin, ns, lens)
397            } else {
398                internals::build_3d_grid(builder, origin, ns, lens)
399            }
400        })
401    }
402}