6 个版本 (3 个破坏性更新)

0.4.0 2023 年 10 月 15 日
0.3.1 2023 年 10 月 10 日
0.2.1 2023 年 7 月 3 日
0.2.0 2022 年 6 月 28 日
0.1.0 2022 年 6 月 16 日

#789 in 数学


2 个 Crates 中使用 (通过 gemlab)

MIT 许可证

2MB
42K SLoC

C++ 27K SLoC // 0.2% comments C 12K SLoC // 0.2% comments Rust 3.5K SLoC // 0.0% comments BASH 10 SLoC // 0.2% comments

三角形和四面体网格生成器

codecov Test & Coverage Windows & macOS

内容

简介

这个包通过封装最好的工具,即 TriangleTetgen,实现了三角形和四面体网格生成器。

在这里,所有由 C/C++ 代码访问的数据结构都是通过小心地使用 "malloc/new" 在 "C-side" 上分配的。😅 然后我们使用 Valgrind 和测试来确保没有泄漏。这样,C 代码没有性能损失,同时使得 Rust 方便易用。

Triangle 和 Tetgen 的 Rust 接口相对较低级。然而,其他项目可以使用此接口来创建更高级别的函数。

该代码适用于多线程应用程序——尽管没有彻底验证,但已测试。例如,请参阅 mem_check_triangle_build.rsmem_check_tetgen_build.rs 中的综合测试。

还有一个更高层次的 Crates 用于网格生成(等等):Gemlab:有限元分析中的几何、网格和数值积分

请参阅文档以获取更多信息

安装

安装一些库

sudo apt install build-essential

设置 Cargo.toml

Crates.io

👆 检查 Crate 版本并相应更新您的 Cargo.toml

[dependencies]
tritet = "*"

示例

注意:将 SAVE_FIGURE 设置为 true 以生成图形。

2D Delaunay 三角剖分

use plotpy::Plot;
use tritet::{StrError, Trigen};

const SAVE_FIGURE: bool = false;

fn main() -> Result<(), StrError> {
    // allocate data for 10 points
    let mut trigen = Trigen::new(10, None, None, None)?;

    // set points
    trigen
        .set_point(0, 0, 0.478554, 0.00869692)?
        .set_point(1, 0, 0.13928, 0.180603)?
        .set_point(2, 0, 0.578587, 0.760349)?
        .set_point(3, 0, 0.903726, 0.975904)?
        .set_point(4, 0, 0.0980015, 0.981755)?
        .set_point(5, 0, 0.133721, 0.348832)?
        .set_point(6, 0, 0.648071, 0.369534)?
        .set_point(7, 0, 0.230951, 0.558482)?
        .set_point(8, 0, 0.0307942, 0.459123)?
        .set_point(9, 0, 0.540745, 0.331184)?;

    // generate Delaunay triangulation
    trigen.generate_delaunay(false)?;

    // draw triangles
    if SAVE_FIGURE {
        let mut plot = Plot::new();
        trigen.draw_triangles(&mut plot, true, true, true, true, None, None, None);
        plot.set_equal_axes(true)
            .set_figure_size_points(600.0, 600.0)
            .save("/tmp/tritet/doc_triangle_delaunay_1.svg")?;
    }
    Ok(())
}

doc_triangle_delaunay_1.svg

2D Voronoi 剖分

use plotpy::Plot;
use tritet::{StrError, Trigen};

const SAVE_FIGURE: bool = false;

fn main() -> Result<(), StrError> {
    // allocate data for 10 points
    let mut trigen = Trigen::new(10, None, None, None)?;

    // set points
    trigen
        .set_point(0, 0, 0.478554, 0.00869692)?
        .set_point(1, 0, 0.13928, 0.180603)?
        .set_point(2, 0, 0.578587, 0.760349)?
        .set_point(3, 0, 0.903726, 0.975904)?
        .set_point(4, 0, 0.0980015, 0.981755)?
        .set_point(5, 0, 0.133721, 0.348832)?
        .set_point(6, 0, 0.648071, 0.369534)?
        .set_point(7, 0, 0.230951, 0.558482)?
        .set_point(8, 0, 0.0307942, 0.459123)?
        .set_point(9, 0, 0.540745, 0.331184)?;

    // generate Voronoi tessellation
    trigen.generate_voronoi(false)?;

    // draw Voronoi diagram
    if SAVE_FIGURE {
        let mut plot = Plot::new();
        trigen.draw_voronoi(&mut plot);
        plot.set_equal_axes(true)
            .set_figure_size_points(600.0, 600.0)
            .save("/tmp/tritet/doc_triangle_voronoi_1.svg")?;
    }
    Ok(())
}

doc_triangle_voronoi_1.svg

2D 网格生成

use plotpy::Plot;
use tritet::{StrError, Trigen};

const SAVE_FIGURE: bool = false;

fn main() -> Result<(), StrError> {
    // allocate data for 12 points, 10 segments, 2 regions, and 1 hole
    let mut trigen = Trigen::new(12, Some(10), Some(2), Some(1))?;

    // set points
    trigen
        .set_point(0, 0, 0.0, 0.0)?
        .set_point(1, 0, 1.0, 0.0)?
        .set_point(2, 0, 1.0, 1.0)?
        .set_point(3, 0, 0.0, 1.0)?
        .set_point(4, 0, 0.2, 0.2)?
        .set_point(5, 0, 0.8, 0.2)?
        .set_point(6, 0, 0.8, 0.8)?
        .set_point(7, 0, 0.2, 0.8)?
        .set_point(8, 0, 0.0, 0.5)?
        .set_point(9, 0, 0.2, 0.5)?
        .set_point(10, 0, 0.8, 0.5)?
        .set_point(11, 0, 1.0, 0.5)?;

    // set segments
    trigen
        .set_segment(0, -1, 0, 1)?
        .set_segment(1, -1, 1, 2)?
        .set_segment(2, -1, 2, 3)?
        .set_segment(3, -1, 3, 0)?
        .set_segment(4, -1, 4, 5)?
        .set_segment(5, -1, 5, 6)?
        .set_segment(6, -1, 6, 7)?
        .set_segment(7, -1, 7, 4)?
        .set_segment(8, -1, 8, 9)?
        .set_segment(9, -1, 10, 11)?;

    // set regions
    trigen
        .set_region(0, 1, 0.1, 0.1, None)?
        .set_region(1, 2, 0.1, 0.9, None)?;

    // set holes
    trigen.set_hole(0, 0.5, 0.5)?;

    // generate o2 mesh without constraints
    trigen.generate_mesh(false, true, false, None, None)?;

    // draw mesh
    if SAVE_FIGURE {
        let mut plot = Plot::new();
        trigen.draw_triangles(&mut plot, true, true, true, true, None, None, None);
        plot.set_equal_axes(true)
            .set_figure_size_points(600.0, 600.0)
            .save("/tmp/tritet/doc_triangle_mesh_1.svg")?;
    }
    Ok(())
}

doc_triangle_mesh_1.svg

3D Delaunay 三角剖分

use plotpy::Plot;
use tritet::{StrError, Tetgen};

const SAVE_FIGURE: bool = false;

fn main() -> Result<(), StrError> {
    // allocate data for 8 points
    let mut tetgen = Tetgen::new(8, None, None, None)?;

    // set points
    tetgen
        .set_point(0, 0, 0.0, 0.0, 0.0)?
        .set_point(1, 0, 1.0, 0.0, 0.0)?
        .set_point(2, 0, 1.0, 1.0, 0.0)?
        .set_point(3, 0, 0.0, 1.0, 0.0)?
        .set_point(4, 0, 0.0, 0.0, 1.0)?
        .set_point(5, 0, 1.0, 0.0, 1.0)?
        .set_point(6, 0, 1.0, 1.0, 1.0)?
        .set_point(7, 0, 0.0, 1.0, 1.0)?;

    // generate Delaunay triangulation
    tetgen.generate_delaunay(false)?;

    // draw edges of tetrahedra
    if SAVE_FIGURE {
        let mut plot = Plot::new();
        tetgen.draw_wireframe(&mut plot, true, true, true, true, None, None, None);
        plot.set_equal_axes(true)
            .set_figure_size_points(600.0, 600.0)
            .save("/tmp/tritet/example_tetgen_delaunay_1.svg")?;
    }
    Ok(())
}

example_tetgen_delaunay_1.svg

3D 网格生成

注意:将 SAVE_VTU_FILE 设置为 true 以生成 Paraview 文件。

use plotpy::Plot;
use tritet::{StrError, Tetgen};

const SAVE_VTU_FILE: bool = false;
const SAVE_FIGURE: bool = false;

fn main() -> Result<(), StrError> {
    // allocate data for 16 points and 12 facets
    // (one cube/hole inside another cube)
    let mut tetgen = Tetgen::new(
        16,
        Some(vec![
            4, 4, 4, 4, 4, 4, // inner cube
            4, 4, 4, 4, 4, 4, // outer cube
        ]),
        Some(1),
        Some(1),
    )?;

    // inner cube
    tetgen
        .set_point(0, 0, 0.0, 0.0, 0.0)?
        .set_point(1, 0, 1.0, 0.0, 0.0)?
        .set_point(2, 0, 1.0, 1.0, 0.0)?
        .set_point(3, 0, 0.0, 1.0, 0.0)?
        .set_point(4, 0, 0.0, 0.0, 1.0)?
        .set_point(5, 0, 1.0, 0.0, 1.0)?
        .set_point(6, 0, 1.0, 1.0, 1.0)?
        .set_point(7, 0, 0.0, 1.0, 1.0)?;

    // outer cube
    tetgen
        .set_point(8,  0, -1.0, -1.0, -1.0)?
        .set_point(9,  0, 2.0, -1.0, -1.0)?
        .set_point(10, 0, 2.0, 2.0, -1.0)?
        .set_point(11, 0, -1.0, 2.0, -1.0)?
        .set_point(12, 0, -1.0, -1.0, 2.0)?
        .set_point(13, 0, 2.0, -1.0, 2.0)?
        .set_point(14, 0, 2.0, 2.0, 2.0)?
        .set_point(15, 0, -1.0, 2.0, 2.0)?;

    // inner cube
    tetgen
        .set_facet_point(0, 0, 0)?
        .set_facet_point(0, 1, 4)?
        .set_facet_point(0, 2, 7)?
        .set_facet_point(0, 3, 3)?;
    tetgen
        .set_facet_point(1, 0, 1)?
        .set_facet_point(1, 1, 2)?
        .set_facet_point(1, 2, 6)?
        .set_facet_point(1, 3, 5)?;
    tetgen
        .set_facet_point(2, 0, 0)?
        .set_facet_point(2, 1, 1)?
        .set_facet_point(2, 2, 5)?
        .set_facet_point(2, 3, 4)?;
    tetgen
        .set_facet_point(3, 0, 2)?
        .set_facet_point(3, 1, 3)?
        .set_facet_point(3, 2, 7)?
        .set_facet_point(3, 3, 6)?;
    tetgen
        .set_facet_point(4, 0, 0)?
        .set_facet_point(4, 1, 3)?
        .set_facet_point(4, 2, 2)?
        .set_facet_point(4, 3, 1)?;
    tetgen
        .set_facet_point(5, 0, 4)?
        .set_facet_point(5, 1, 5)?
        .set_facet_point(5, 2, 6)?
        .set_facet_point(5, 3, 7)?;

    // outer cube
    tetgen
        .set_facet_point(6, 0, 8 + 0)?
        .set_facet_point(6, 1, 8 + 4)?
        .set_facet_point(6, 2, 8 + 7)?
        .set_facet_point(6, 3, 8 + 3)?;
    tetgen
        .set_facet_point(7, 0, 8 + 1)?
        .set_facet_point(7, 1, 8 + 2)?
        .set_facet_point(7, 2, 8 + 6)?
        .set_facet_point(7, 3, 8 + 5)?;
    tetgen
        .set_facet_point(8, 0, 8 + 0)?
        .set_facet_point(8, 1, 8 + 1)?
        .set_facet_point(8, 2, 8 + 5)?
        .set_facet_point(8, 3, 8 + 4)?;
    tetgen
        .set_facet_point(9, 0, 8 + 2)?
        .set_facet_point(9, 1, 8 + 3)?
        .set_facet_point(9, 2, 8 + 7)?
        .set_facet_point(9, 3, 8 + 6)?;
    tetgen
        .set_facet_point(10, 0, 8 + 0)?
        .set_facet_point(10, 1, 8 + 3)?
        .set_facet_point(10, 2, 8 + 2)?
        .set_facet_point(10, 3, 8 + 1)?;
    tetgen
        .set_facet_point(11, 0, 8 + 4)?
        .set_facet_point(11, 1, 8 + 5)?
        .set_facet_point(11, 2, 8 + 6)?
        .set_facet_point(11, 3, 8 + 7)?;

    // set region and hole
    tetgen.set_region(0, 1, -0.9, -0.9, -0.9, None)?;
    tetgen.set_hole(0, 0.5, 0.5, 0.5)?;

    // generate mesh
    tetgen.generate_mesh(false, false, None, None)?;

    // generate file for Paraview
    if SAVE_VTU_FILE {
        tetgen.write_vtu("/tmp/tritet/example_tetgen_mesh_1.vtu")?;
    }

    // draw edges of tetrahedra
    if SAVE_FIGURE {
        let mut plot = Plot::new();
        tetgen.draw_wireframe(&mut plot, true, true, true, true, None, None, None);
        plot.set_equal_axes(true)
            .set_figure_size_points(600.0, 600.0)
            .save("/tmp/tritet/example_tetgen_mesh_1.svg")?;
    }
    Ok(())
}

example_tetgen_mesh_1.svg

针对开发者

安装 cargo-valgrind

cargo install cargo-valgrind

然后检查内存泄漏(没有问题;-)

bash memcheck.bash

依赖项

~0.6–0.8MB