tritet 2.0.0

Triangle and tetrahedron mesh generators
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
# Triangle and tetrahedron mesh generators <!-- omit from toc --> 

[![Test](https://github.com/cpmech/tritet/actions/workflows/test_and_coverage.yml/badge.svg)](https://github.com/cpmech/tritet/actions/workflows/test_and_coverage.yml)
[![Windows & macOS](https://github.com/cpmech/tritet/actions/workflows/windows_and_macos.yml/badge.svg)](https://github.com/cpmech/tritet/actions/workflows/windows_and_macos.yml)
[![Test on Arch Linux](https://github.com/cpmech/tritet/actions/workflows/test_on_arch_linux.yml/badge.svg)](https://github.com/cpmech/tritet/actions/workflows/test_on_arch_linux.yml)

## Contents <!-- omit from toc --> 

- [Introduction]#introduction
- [Installation]#installation
- [Setting Cargo.toml]#setting-cargotoml
- [Examples]#examples
  - [Mesh generation using JSON input files]#mesh-generation-using-json-input-files
  - [2D Delaunay triangulation]#2d-delaunay-triangulation
  - [2D Voronoi tessellation]#2d-voronoi-tessellation
  - [2D mesh generation using input data structure]#2d-mesh-generation-using-input-data-structure
  - [2D mesh generation using setup functions]#2d-mesh-generation-using-setup-functions
  - [3D Delaunay triangulation]#3d-delaunay-triangulation
  - [3D mesh generation using the input data structure]#3d-mesh-generation-using-the-input-data-structure
- [Definitions for Triangle (By J. R. Shewchuk)]#definitions-for-triangle-by-j-r-shewchuk
- [For developers]#for-developers

## Introduction

This crate implements Triangle and Tetrahedron mesh generators by wrapping the best tools around, namely, [Triangle](https://www.cs.cmu.edu/~quake/triangle.html) and [Tetgen](http://tetgen.org/).

Here, all the data structures accessed by the C/C++ codes are allocated on the "C-side" by (carefully) using "malloc/new." 😅 We then make use of [Valgrind](https://valgrind.org/) and tests to make sure that there are no leaks. In this way, there is no performance loss of the C-code while enabling the convenience of Rust.

The resulting Rust interface to Triangle and Tetgen is somewhat low-level. However, other projects could use this interface to make higher-level functions.

The code works in multithreaded applications---not exhaustively verified but tested. See, for example, the comprehensive tests in [mem_check_triangle_build.rs](https://github.com/cpmech/tritet/blob/main/src/bin/mem_check_triangle_build.rs) and [mem_check_tetgen_build.rs](https://github.com/cpmech/tritet/blob/main/src/bin/mem_check_tetgen_build.rs)

A higher-level crate is available for mesh generation (and more): [Gemlab: Geometry, meshes, and numerical integration for finite element analyses](https://github.com/cpmech/gemlab).

See the documentation for further information:

- [Tritet documentation]https://docs.rs/tritet - Contains the API reference and examples

## Installation

Install some libraries:

```bash
sudo apt install build-essential
```

## Setting Cargo.toml

[![Crates.io](https://img.shields.io/crates/v/tritet.svg)](https://crates.io/crates/tritet)

👆 Check the crate version and update your Cargo.toml accordingly:

```toml
[dependencies]
tritet = "*"
```

## Examples

Note: set `SAVE_FIGURE` to true to generate the figures.

### Mesh generation using JSON input files

Tritet contains two executables to generate meshes:

1. `trigen_mesh` to generate quality triangle meshes with boundary markers and volume and angle constraints
2. `tetgen_mesh` to generate quality tetrahedral meshes with boundary markers and volume and angle constraints

Below is a JSON input for `trigen_mesh`:

```json
{
  "points": [
    [0, 0.0, 0.0],
    [0, 1.0, 0.0],
    [0, 1.0, 1.0],
    [0, 0.0, 1.0],
    [0, 0.2, 0.2],
    [0, 0.8, 0.2],
    [0, 0.8, 0.8],
    [0, 0.2, 0.8],
    [0, 0.0, 0.5],
    [0, 0.2, 0.5],
    [0, 0.8, 0.5],
    [0, 1.0, 0.5]
  ],
  "segments": [
    [-1,  0,  1],
    [-1,  1,  2],
    [-1,  2,  3],
    [-1,  3,  0],
    [-1,  4,  5],
    [-1,  5,  6],
    [-1,  6,  7],
    [-1,  7,  4],
    [-1,  8,  9],
    [-1, 10, 11]
  ],
  "holes": [
    [0.5, 0.5]
  ],
  "regions": [
    [1, 0.1, 0.1, null],
    [2, 0.1, 0.9, null]
  ]
}
```

Which can be used as follows:

```bash
cargo run --bin trigen_mesh -- data/input/example_tri_input.json /tmp/tritet -s -v0.01
```

Where `-s` indicates SVG file generation (if Python/Matplotlib is available; otherwise an error arises),
and `v0.01` indicates an area (volume) constraint of 0.01. The output is shown below:

![example_tri_input.svg](https://raw.githubusercontent.com/cpmech/tritet/main/data/figures/example_tri_input.svg)

Below is a JSON input for `tetgen_mesh`:

```json
{
  "points": [
    [0, 0.0, 1.0, 0.0],
    [0, 0.0, 0.0, 0.0],
    [0, 1.0, 1.0, 0.0],
    [0, 0.0, 1.0, 1.0]
  ],
  "facets": [
    [0, [0, 2, 1]],
    [0, [0, 1, 3]],
    [0, [0, 3, 2]],
    [0, [1, 2, 3]]
  ],
  "holes": [],
  "regions": [
    [1, 0.1, 0.9, 0.1, null]
  ]
}
```

Which can be used as follows:

```bash
cargo run --bin tetgen_mesh -- data/input/example_tet_input.json /tmp/tritet -s -v0.1
```

Where `-s` indicates SVG file generation (if Python/Matplotlib is available; otherwise an error arises),
and `v0.1` indicates an volume constraint of 0.1. The output is shown below:

![example_tet_input.svg](https://raw.githubusercontent.com/cpmech/tritet/main/data/figures/example_tet_input.svg)

### 2D Delaunay triangulation

```rust
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](https://raw.githubusercontent.com/cpmech/tritet/main/data/figures/doc_triangle_delaunay_1.svg)

### 2D Voronoi tessellation

```rust
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](https://raw.githubusercontent.com/cpmech/tritet/main/data/figures/doc_triangle_voronoi_1.svg)

### 2D mesh generation using input data structure

```rust
use plotpy::Plot;
use tritet::{InputDataTriMesh, StrError, Trigen};

const SAVE_FIGURE: bool = false;

fn main() -> Result<(), StrError> {
    // set input data
    let input_data = InputDataTriMesh {
        points: vec![
            (0, 0.0, 0.0), // boundary marker, x, y
            (0, 1.0, 0.0),
            (0, 1.0, 1.0),
            (0, 0.0, 1.0),
            (0, 0.2, 0.2),
            (0, 0.8, 0.2),
            (0, 0.8, 0.8),
            (0, 0.2, 0.8),
            (0, 0.0, 0.5),
            (0, 0.2, 0.5),
            (0, 0.8, 0.5),
            (0, 1.0, 0.5),
        ],
        segments: vec![
            (-1, 0, 1), // boundary marker, point indices
            (-1, 1, 2),
            (-1, 2, 3),
            (-1, 3, 0),
            (-1, 4, 5),
            (-1, 5, 6),
            (-1, 6, 7),
            (-1, 7, 4),
            (-1, 8, 9),
            (-1, 10, 11),
        ],
        holes: vec![
            (0.5, 0.5), // x, y
        ],
        regions: vec![
            (1, 0.1, 0.1, None), // marker, x, y, max area
            (2, 0.1, 0.9, None),
        ],
    };

    // allocate generator from input data
    let trigen = Trigen::from_input_data(&input_data)?;

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

    // 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(())
}
```

### 2D mesh generation using setup functions

```rust
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](https://raw.githubusercontent.com/cpmech/tritet/main/data/figures/doc_triangle_mesh_1.svg)

### 3D Delaunay triangulation

```rust
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](https://raw.githubusercontent.com/cpmech/tritet/main/data/figures/example_tetgen_delaunay_1.svg)

### 3D mesh generation using the input data structure

```rust
use plotpy::Plot;
use tritet::{InputDataTetMesh, StrError, Tetgen};

const SAVE_FIGURE: bool = false;

fn main() -> Result<(), StrError> {
    // set input data
    let input_data = InputDataTetMesh {
        points: vec![
            (0, 0.0, 1.0, 0.0), // marker, x, y, z
            (0, 0.0, 0.0, 0.0),
            (0, 1.0, 1.0, 0.0),
            (0, 0.0, 1.0, 1.0),
        ],
        facets: vec![
            (0, vec![0, 2, 1]), // marker, point indices
            (0, vec![0, 1, 3]),
            (0, vec![0, 3, 2]),
            (0, vec![1, 2, 3]),
        ],
        holes: vec![],                           // no holes
        regions: vec![(1, 0.1, 0.9, 0.1, None)], // region marker, x, y, z, max volume
    };

    // allocate generator from input data
    let tetgen = Tetgen::from_input_data(&input_data)?;

    // generate mesh
    let global_max_volume = Some(0.5);
    tetgen.generate_mesh(false, false, global_max_volume, None)?;

    // 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/doc_tetgen_mesh_2.svg")?;
    }

    assert_eq!(tetgen.out_ncell(), 7);
    assert_eq!(tetgen.out_npoint(), 10);
    Ok(())
}
```

#g## 3D mesh generation using setup functions

Note: set `SAVE_VTU_FILE` to true to generate Paraview file.

```rust
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](https://raw.githubusercontent.com/cpmech/tritet/main/data/figures/example_tetgen_mesh_1.svg)

## Definitions for Triangle (By J. R. Shewchuk)

This section is part of the text written by J. R. Shewchuk in [Triangle](https://www.cs.cmu.edu/~quake/triangle.html). The figures in this section are also Copyright by J. R. Shewchuk.

A Delaunay triangulation (Figure 1) of a vertex set (Figure 2) is a triangulation of the vertex set with the property that no vertex in the vertex set falls in the interior of the circumcircle (circle that passes through all three vertices) of any triangle in the triangulation.

**Figure 1**. Delaunay triangulation:

![Delaunay triangulation](https://raw.githubusercontent.com/cpmech/tritet/main/data/shewchuk/dots.1.ele.gif)

**Figure 2**. Vertex set (cloud of points):

![Vertex set (cloud of points).](https://raw.githubusercontent.com/cpmech/tritet/main/data/shewchuk/dots.node.gif)

A Voronoi diagram (Figure 3) of a vertex set is a subdivision of the plane into polygonal regions (some of which may be infinite), where each region is the set of points in the plane that are closer to some input vertex than to any other input vertex. (The Voronoi diagram is the geometric dual of the Delaunay triangulation.)

**Figure 3**. Voronoi diagram:

![Voronoi diagram.](https://raw.githubusercontent.com/cpmech/tritet/main/data/shewchuk/dots.1.v.edge.gif)

A Planar Straight Line Graph (**PSLG**, Figure 4) is a collection of vertices and segments. Segments are edges whose endpoints are vertices in the PSLG, and whose presence in any mesh generated from the PSLG is enforced.

**Figure 4**. Planar Straight Line Graph (PSLG):

![Planar Straight Line Graph (PSLG).](https://raw.githubusercontent.com/cpmech/tritet/main/data/shewchuk/A.poly.gif)

A constrained Delaunay triangulation of a PSLG (Figure 5) is similar to a Delaunay triangulation, but each PSLG segment is present as a single edge in the triangulation. A constrained Delaunay triangulation is not truly a Delaunay triangulation. Some of its triangles might not be Delaunay, but they are all constrained Delaunay.

**Figure 5**. Constrained Delaunay triangulation of a PSLG:

![Constrained Delaunay triangulation of a PSLG.](https://raw.githubusercontent.com/cpmech/tritet/main/data/shewchuk/A.constrain.gif)

A conforming Delaunay triangulation (CDT, Figure 6) of a PSLG is a true Delaunay triangulation in which each PSLG segment may have been subdivided into several edges by the insertion of additional vertices, called Steiner points. Steiner points are necessary to allow the segments to exist in the mesh while maintaining the Delaunay property. Steiner points are also inserted to meet constraints on the minimum angle and maximum triangle area.

**Figure 6**. Conforming Delaunay triangulation (CDT):

![Conforming Delaunay triangulation (CDT)](https://raw.githubusercontent.com/cpmech/tritet/main/data/shewchuk/A.conform.gif)

A constrained conforming Delaunay triangulation (CCDT, Figure 7) of a PSLG is a constrained Delaunay triangulation that includes Steiner points. It usually takes fewer vertices to make a good-quality CCDT than a good-quality CDT, because the triangles do not need to be Delaunay (although they still must be constrained Delaunay).

**Figure 7**. Constrained conforming Delaunay triangulation (CCDT):

![Constrained conforming Delaunay triangulation (CCDT)](https://raw.githubusercontent.com/cpmech/tritet/main/data/shewchuk/A.ccdt.gif)

## For developers

Install cargo-valgrind:

```bash
cargo install cargo-valgrind
```

Then check for memory leaks (none ;-):

```bash
bash memcheck.bash
```