Generating meshes with Wildmeshing

Generating meshes with Wildmeshing#

In this demo, we will show how to generate a mesh from STL/ply files based on MRI images. We use the files from the paper In-silico molecular enrichment and clearance of the human intracranial space

We start by importing all the packages required for this task.

from mpi4py import MPI
import wildmeshing
import dolfinx
from pathlib import Path
import json
import numpy as np
import pyvista
import ufl
import basix.ufl
import os

Next, we use the already downloaded stl files from this chapter. We will use WildMeshing to combine these. We create unique volume markers for the lateral ventricles (LV.ply) and the third and fourth ventricle (V34.ply). Next we embed the ventricular system within the brain parenchyma (parenchyma_incl_ventr.ply), before we add the dura interface (skull.ply), which is the interfacebetween the CSF and the skull.

We extract these files from the "WILDFENICS_DATA_PATH". If you have not set this variable, please change it to the appropriate location on your system.

folder = Path(os.environ["WILDFENICS_DATA_PATH"])
assert folder.exists(), "Could not find surface files"
tree = {
    "operation": "union",
    "left": str((folder / "skull.ply").absolute().as_posix()),
    "right": {
        "operation": "union",
        "left": str((folder / "parenchyma_incl_ventr.ply").absolute().as_posix()),
        "right": {
            "operation": "union",
            "left": str((folder / "LV.ply").absolute().as_posix()),
            "right": str((folder / "V34.ply").absolute().as_posix()),
        },
    },
}

Next, we can set up a 3D meshing instance and load in the instruction-set

tetra = wildmeshing.Tetrahedralizer(
    stop_quality=10, max_its=30, edge_length_r=0.03, epsilon=0.00225
)
tetra.load_csg_tree(json.dumps(tree))
TBB threads 4
o-[I/O         ] Loading file /root/data/mesh/standard/surfaces/skull.ply...
                 (FP64) nb_v:138154 nb_e:0 nb_f:276304 nb_b:0 tri:1 dim:3
                 Attributes on vertices: point[3]
                 Loading file /root/data/mesh/standard/surfaces/parenchyma_incl_ventr.ply...
                 (FP64) nb_v:317019 nb_e:0 nb_f:635594 nb_b:0 tri:1 dim:3
                 Attributes on vertices: point[3]
                 Loading file /root/data/mesh/standard/surfaces/LV.ply...
                 (FP64) nb_v:17572 nb_e:0 nb_f:35144 nb_b:0 tri:1 dim:3
                 Attributes on vertices: point[3]
                 Loading file /root/data/mesh/standard/surfaces/V34.ply...
                 (FP64) nb_v:4608 nb_e:0 nb_f:9212 nb_b:0 tri:1 dim:3
                 Attributes on vertices: point[3]
bbox_diag_length = 0.266221
ideal_edge_length = 0.00798662
stage = 2
eps_input = 0.000598997
eps = 0.000331599
eps_simplification = 0.000265279
eps_coplanar = 2.66221e-07
dd = 0.000399331
dd_simplification = 0.000319465

Then, we can generate the mesh

# Create mesh
tetra.tetrahedralize()
[2025-04-23 06:49:43.508] [float-tetwild] [info] remove duplicates: 
[2025-04-23 06:49:43.508] [float-tetwild] [info] #v: 477353 -> 477353
[2025-04-23 06:49:43.508] [float-tetwild] [info] #f: 956254 -> 956254
collapsing 20.9093
swapping 0.245065
[2025-04-23 06:50:04.894] [float-tetwild] [info] remove duplicates: 
[2025-04-23 06:50:04.894] [float-tetwild] [info] #v: 49442 -> 49442
[2025-04-23 06:50:04.894] [float-tetwild] [info] #f: 100369 -> 100369
[2025-04-23 06:50:04.895] [float-tetwild] [info] #v = 49442
[2025-04-23 06:50:04.895] [float-tetwild] [info] #f = 100256
#boundary_e1 = 309
#boundary_e2 = 120
[2025-04-23 06:50:05.259] [float-tetwild] [info] preprocessing 21.9628s
[2025-04-23 06:50:05.259] [float-tetwild] [info] 
[2025-04-23 06:50:05.489] [float-tetwild] [info] #v = 50230
[2025-04-23 06:50:05.490] [float-tetwild] [info] #t = 333917
[2025-04-23 06:50:05.490] [float-tetwild] [info] tetrahedralizing 0.230555s
[2025-04-23 06:50:05.490] [float-tetwild] [info] 
[2025-04-23 06:50:05.491] [float-tetwild] [info] triangle insertion start, #f = 100256, #v = 50230, #t = 333917
[2025-04-23 06:50:05.627] [float-tetwild] [info] matched #f = 74494, uninserted #f = 25762
[2025-04-23 06:50:14.508] [float-tetwild] [info] insert_one_triangle * n done, #v = 153689, #t = 867204
[2025-04-23 06:50:14.508] [float-tetwild] [info] uninserted #f = 0/25762
[2025-04-23 06:50:14.508] [float-tetwild] [info] total timing: 0s
[2025-04-23 06:50:14.942] [float-tetwild] [info] pair_track_surface_fs done
#boundary_e1 = 309
#boundary_e2 = 120
[2025-04-23 06:50:14.991] [float-tetwild] [info] find_boundary_edges done
[2025-04-23 06:50:15.040] [float-tetwild] [info] time1 = 0.048998
[2025-04-23 06:50:15.052] [float-tetwild] [info] uninsert boundary #e = 0/429
[2025-04-23 06:50:15.052] [float-tetwild] [info] time2 = 3.8e-05
[2025-04-23 06:50:15.052] [float-tetwild] [info] time3 = 0.011282
[2025-04-23 06:50:15.052] [float-tetwild] [info] time4 = 0.000279
[2025-04-23 06:50:15.052] [float-tetwild] [info] time5 = 7.1e-05
[2025-04-23 06:50:15.052] [float-tetwild] [info] time6 = 3e-06
[2025-04-23 06:50:15.058] [float-tetwild] [info] uninserted #f = 0/25762
known_surface_fs.size = 0
known_not_surface_fs.size = 0
[2025-04-23 06:50:17.727] [float-tetwild] [info] mark_surface_fs done
[2025-04-23 06:50:17.748] [float-tetwild] [info] #b_edge1 = 429, #b_edges2 = 0
[2025-04-23 06:50:17.779] [float-tetwild] [info] cutting 12.2882s
[2025-04-23 06:50:17.779] [float-tetwild] [info] 
initializing...
edge collapsing...
fixed 0 tangled element
success(env) = 85433
success = 100369(816835)
success(env) = 3244
success = 3457(533025)
success(env) = 265
success = 276(225638)
success(env) = 38
success = 42(29639)
success(env) = 3
success = 3(4127)
success(env) = 0
success = 0(602)
edge collapsing done!
time = 13.6636s
#v = 44083
#t = 288305
max_energy = 3666.9
avg_energy = 7.73081
//////////////// pass 0 ////////////////
edge splitting...
fixed 0 tangled element
success = 43023(43023)
edge splitting done!
time = 0.321668s
#v = 87106
#t = 503543
max_energy = 3666.9
avg_energy = 7.63438
edge collapsing...
fixed 0 tangled element
success(env) = 189
success = 29304(622753)
success(env) = 76
success = 1438(158340)
success(env) = 16
success = 197(48181)
success(env) = 4
success = 31(8477)
success(env) = 1
success = 10(1504)
success(env) = 0
success = 3(545)
success(env) = 0
success = 0(29)
edge collapsing done!
time = 5.11634s
#v = 56123
#t = 348998
max_energy = 3666.9
avg_energy = 6.65935
edge swapping...
fixed 0 tangled element
success3 = 33780
success4 = 39739
success5 = 3923
success = 77442(344669)
edge swapping done!
time = 2.43238s
#v = 56123
#t = 319141
max_energy = 1442.78
avg_energy = 5.56402
vertex smoothing...
success = 30244(52230)
vertex smoothing done!
time = 1.99705s
#v = 56123
#t = 319141
max_energy = 974.658
avg_energy = 5.21819
//////////////// pass 1 ////////////////
edge splitting...
fixed 0 tangled element
success = 5168(5168)
edge splitting done!
time = 0.197898s
#v = 61291
#t = 344429
max_energy = 974.658
avg_energy = 5.16486
edge collapsing...
fixed 0 tangled element
success(env) = 2686
success = 6960(520830)
success(env) = 176
success = 372(264695)
success(env) = 22
success = 47(26316)
success(env) = 1
success = 5(3545)
success(env) = 0
success = 1(299)
success(env) = 0
success = 0(25)
edge collapsing done!
time = 4.60197s
#v = 53906
#t = 306339
max_energy = 974.658
avg_energy = 5.12321
edge swapping...
fixed 0 tangled element
success3 = 3992
success4 = 14119
success5 = 1264
success = 19375(200759)
edge swapping done!
time = 1.35485s
#v = 53906
#t = 303611
max_energy = 974.658
avg_energy = 5.01678
vertex smoothing...
success = 25483(50012)
vertex smoothing done!
time = 1.77253s
#v = 53906
#t = 303611
max_energy = 974.658
avg_energy = 4.90812
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 0
enlarge envelope, eps = 0.000368443
//////////////// pass 2 ////////////////
edge splitting...
fixed 0 tangled element
success = 59431(59431)
edge splitting done!
time = 0.438149s
#v = 113337
#t = 620301
max_energy = 974.658
avg_energy = 5.03669
edge collapsing...
fixed 0 tangled element
success(env) = 2460
success = 42557(543439)
success(env) = 284
success = 2908(305370)
success(env) = 43
success = 319(84253)
success(env) = 10
success = 61(12622)
success(env) = 4
success = 19(2063)
success(env) = 0
success = 1(644)
success(env) = 0
success = 2(18)
success(env) = 0
success = 0(26)
edge collapsing done!
time = 6.0712s
#v = 67470
#t = 385764
max_energy = 974.658
avg_energy = 4.49476
edge swapping...
fixed 0 tangled element
success3 = 5515
success4 = 22806
success5 = 2129
success = 30450(285269)
edge swapping done!
time = 2.37577s
#v = 67470
#t = 382378
max_energy = 974.658
avg_energy = 4.32659
vertex smoothing...
success = 41956(63569)
vertex smoothing done!
time = 2.03473s
#v = 67470
#t = 382378
max_energy = 974.658
avg_energy = 4.13635
//////////////// pass 3 ////////////////
edge splitting...
fixed 0 tangled element
success = 12772(12772)
edge splitting done!
time = 0.267526s
#v = 80242
#t = 444465
max_energy = 974.658
avg_energy = 4.16244
edge collapsing...
fixed 0 tangled element
success(env) = 1227
success = 12557(440586)
success(env) = 91
success = 535(182942)
success(env) = 11
success = 62(18062)
success(env) = 2
success = 13(2352)
success(env) = 0
success = 1(480)
success(env) = 0
success = 1(17)
success(env) = 0
success = 1(45)
success(env) = 0
success = 0(60)
edge collapsing done!
time = 3.8073s
#v = 67072
#t = 380126
max_energy = 974.658
avg_energy = 4.12651
edge swapping...
fixed 0 tangled element
success3 = 1736
success4 = 11633
success5 = 856
success = 14225(249764)
edge swapping done!
time = 1.66611s
#v = 67072
#t = 379246
max_energy = 974.658
avg_energy = 4.08461
vertex smoothing...
success = 39205(63168)
vertex smoothing done!
time = 1.89019s
#v = 67072
#t = 379246
max_energy = 974.658
avg_energy = 4.00851
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 0
//////////////// pass 4 ////////////////
edge splitting...
fixed 0 tangled element
success = 7677(7677)
edge splitting done!
time = 0.244518s
#v = 74749
#t = 418733
max_energy = 974.658
avg_energy = 4.06409
edge collapsing...
fixed 0 tangled element
success(env) = 995
success = 6358(649398)
success(env) = 79
success = 412(137544)
success(env) = 7
success = 44(12273)
success(env) = 1
success = 3(1823)
success(env) = 0
success = 0(85)
edge collapsing done!
time = 4.14441s
#v = 67932
#t = 384316
max_energy = 974.658
avg_energy = 4.00984
edge swapping...
fixed 0 tangled element
success3 = 833
success4 = 6718
success5 = 382
success = 7933(240471)
edge swapping done!
time = 1.51871s
#v = 67932
#t = 383865
max_energy = 974.658
avg_energy = 3.9907
vertex smoothing...
success = 37809(64029)
vertex smoothing done!
time = 1.81357s
#v = 67932
#t = 383865
max_energy = 974.658
avg_energy = 3.94688
//////////////// pass 5 ////////////////
edge splitting...
fixed 0 tangled element
success = 3155(3155)
edge splitting done!
time = 0.226985s
#v = 71087
#t = 399119
max_energy = 974.658
avg_energy = 3.96332
edge collapsing...
fixed 0 tangled element
success(env) = 547
success = 3478(632621)
success(env) = 23
success = 130(76101)
success(env) = 2
success = 13(5014)
success(env) = 0
success = 1(359)
success(env) = 0
success = 0(50)
edge collapsing done!
time = 3.51978s
#v = 67465
#t = 381345
max_energy = 974.658
avg_energy = 3.95017
edge swapping...
fixed 0 tangled element
success3 = 349
success4 = 3633
success5 = 209
success = 4191(231341)
edge swapping done!
time = 1.32499s
#v = 67465
#t = 381205
max_energy = 974.658
avg_energy = 3.94393
vertex smoothing...
success = 35415(63556)
vertex smoothing done!
time = 1.74964s
#v = 67465
#t = 381205
max_energy = 974.658
avg_energy = 3.9229
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 0
//////////////// pass 6 ////////////////
edge splitting...
fixed 0 tangled element
success = 757(757)
edge splitting done!
time = 0.214329s
#v = 68222
#t = 384815
max_energy = 974.658
avg_energy = 3.92833
edge collapsing...
fixed 0 tangled element
success(env) = 211
success = 1026(722852)
success(env) = 9
success = 40(43557)
success(env) = 0
success = 2(2755)
success(env) = 0
success = 0(116)
edge collapsing done!
time = 3.62429s
#v = 67154
#t = 379485
max_energy = 974.658
avg_energy = 3.92441
edge swapping...
fixed 0 tangled element
success3 = 171
success4 = 1890
success5 = 87
success = 2148(226163)
edge swapping done!
time = 1.26286s
#v = 67154
#t = 379401
max_energy = 974.658
avg_energy = 3.92146
vertex smoothing...
success = 33546(63250)
vertex smoothing done!
time = 1.67114s
#v = 67154
#t = 379401
max_energy = 974.658
avg_energy = 3.91177
//////////////// pass 7 ////////////////
edge splitting...
fixed 0 tangled element
success = 577(577)
edge splitting done!
time = 0.212084s
#v = 67731
#t = 382029
max_energy = 974.658
avg_energy = 3.91399
edge collapsing...
fixed 0 tangled element
success(env) = 123
success = 739(719574)
success(env) = 5
success = 18(21824)
success(env) = 0
success = 0(1143)
edge collapsing done!
time = 3.42321s
#v = 66974
#t = 378398
max_energy = 974.658
avg_energy = 3.91224
edge swapping...
fixed 0 tangled element
success3 = 102
success4 = 1160
success5 = 42
success = 1304(223720)
edge swapping done!
time = 1.23691s
#v = 66974
#t = 378338
max_energy = 974.658
avg_energy = 3.91125
vertex smoothing...
success = 32093(63067)
vertex smoothing done!
time = 1.60738s
#v = 66974
#t = 378338
max_energy = 974.658
avg_energy = 3.9057
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 8 ////////////////
edge splitting...
fixed 0 tangled element
success = 517(517)
edge splitting done!
time = 0.214051s
#v = 67491
#t = 380682
max_energy = 974.658
avg_energy = 3.90774
edge collapsing...
fixed 0 tangled element
success(env) = 58
success = 607(734714)
success(env) = 1
success = 18(14040)
success(env) = 0
success = 0(1277)
edge collapsing done!
time = 3.48488s
#v = 66866
#t = 377726
max_energy = 974.658
avg_energy = 3.90714
edge swapping...
fixed 0 tangled element
success3 = 61
success4 = 657
success5 = 26
success = 744(222176)
edge swapping done!
time = 1.24192s
#v = 66866
#t = 377691
max_energy = 974.658
avg_energy = 3.9064
vertex smoothing...
success = 31044(62960)
vertex smoothing done!
time = 1.55286s
#v = 66866
#t = 377691
max_energy = 974.658
avg_energy = 3.90299
//////////////// pass 9 ////////////////
edge splitting...
fixed 0 tangled element
success = 537(537)
edge splitting done!
time = 0.211743s
#v = 67403
#t = 380119
max_energy = 974.658
avg_energy = 3.90543
edge collapsing...
fixed 0 tangled element
success(env) = 34
success = 562(733922)
success(env) = 2
success = 10(8484)
success(env) = 0
success = 0(589)
edge collapsing done!
time = 3.46146s
#v = 66831
#t = 377487
max_energy = 974.658
avg_energy = 3.9034
edge swapping...
fixed 0 tangled element
success3 = 31
success4 = 455
success5 = 20
success = 506(221613)
edge swapping done!
time = 1.25019s
#v = 66831
#t = 377476
max_energy = 974.658
avg_energy = 3.90287
vertex smoothing...
success = 30294(62926)
vertex smoothing done!
time = 1.51328s
#v = 66831
#t = 377476
max_energy = 974.658
avg_energy = 3.90085
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 10 ////////////////
edge splitting...
fixed 0 tangled element
success = 501(501)
edge splitting done!
time = 0.211472s
#v = 67332
#t = 379726
max_energy = 974.658
avg_energy = 3.90268
edge collapsing...
fixed 0 tangled element
success(env) = 8
success = 512(733676)
success(env) = 0
success = 8(5290)
success(env) = 0
success = 1(234)
success(env) = 0
success = 0(17)
edge collapsing done!
time = 3.5402s
#v = 66811
#t = 377374
max_energy = 974.658
avg_energy = 3.90131
edge swapping...
fixed 0 tangled element
success3 = 25
success4 = 285
success5 = 7
success = 317(221119)
edge swapping done!
time = 1.2212s
#v = 66811
#t = 377356
max_energy = 974.658
avg_energy = 3.90091
vertex smoothing...
success = 29735(62907)
vertex smoothing done!
time = 1.47351s
#v = 66811
#t = 377356
max_energy = 974.658
avg_energy = 3.89986
//////////////// pass 11 ////////////////
edge splitting...
fixed 0 tangled element
success = 500(500)
edge splitting done!
time = 0.21234s
#v = 67311
#t = 379589
max_energy = 974.658
avg_energy = 3.90149
edge collapsing...
fixed 0 tangled element
success(env) = 7
success = 505(733368)
success(env) = 3
success = 7(5352)
success(env) = 0
success = 0(263)
edge collapsing done!
time = 3.45176s
#v = 66799
#t = 377287
max_energy = 974.658
avg_energy = 3.9004
edge swapping...
fixed 0 tangled element
success3 = 10
success4 = 173
success5 = 5
success = 188(220841)
edge swapping done!
time = 1.2139s
#v = 66799
#t = 377282
max_energy = 974.658
avg_energy = 3.90016
vertex smoothing...
success = 29266(62890)
vertex smoothing done!
time = 1.45667s
#v = 66799
#t = 377282
max_energy = 974.658
avg_energy = 3.89944
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 12 ////////////////
edge splitting...
fixed 0 tangled element
success = 560(560)
edge splitting done!
time = 0.212336s
#v = 67359
#t = 379835
max_energy = 974.658
avg_energy = 3.90201
edge collapsing...
fixed 0 tangled element
success(env) = 3
success = 532(733847)
success(env) = 0
success = 15(3896)
success(env) = 0
success = 1(254)
success(env) = 0
success = 0(15)
edge collapsing done!
time = 3.56588s
#v = 66811
#t = 377352
max_energy = 974.658
avg_energy = 3.89977
edge swapping...
fixed 0 tangled element
success3 = 11
success4 = 159
success5 = 5
success = 175(220924)
edge swapping done!
time = 1.26678s
#v = 66811
#t = 377346
max_energy = 974.658
avg_energy = 3.89952
vertex smoothing...
success = 28955(62906)
vertex smoothing done!
time = 1.45251s
#v = 66811
#t = 377346
max_energy = 974.658
avg_energy = 3.899
//////////////// pass 13 ////////////////
edge splitting...
fixed 0 tangled element
success = 537(537)
edge splitting done!
time = 0.211352s
#v = 67348
#t = 379755
max_energy = 974.658
avg_energy = 3.90136
edge collapsing...
fixed 0 tangled element
success(env) = 1
success = 532(733707)
success(env) = 0
success = 15(3989)
success(env) = 0
success = 5(1046)
success(env) = 0
success = 0(511)
edge collapsing done!
time = 3.5762s
#v = 66796
#t = 377252
max_energy = 974.658
avg_energy = 3.90023
edge swapping...
fixed 0 tangled element
success3 = 6
success4 = 154
success5 = 2
success = 162(220815)
edge swapping done!
time = 1.20758s
#v = 66796
#t = 377248
max_energy = 974.658
avg_energy = 3.89987
vertex smoothing...
success = 28679(62889)
vertex smoothing done!
time = 1.44534s
#v = 66796
#t = 377248
max_energy = 974.658
avg_energy = 3.89903
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 14 ////////////////
edge splitting...
fixed 0 tangled element
success = 538(538)
edge splitting done!
time = 0.211869s
#v = 67334
#t = 379660
max_energy = 974.658
avg_energy = 3.90119
edge collapsing...
fixed 0 tangled element
success(env) = 2
success = 525(734030)
success(env) = 0
success = 6(3397)
success(env) = 0
success = 0(125)
edge collapsing done!
time = 3.37151s
#v = 66803
#t = 377286
max_energy = 974.658
avg_energy = 3.8994
edge swapping...
fixed 0 tangled element
success3 = 5
success4 = 117
success5 = 2
success = 124(220784)
edge swapping done!
time = 1.21721s
#v = 66803
#t = 377283
max_energy = 974.658
avg_energy = 3.89922
vertex smoothing...
success = 28437(62897)
vertex smoothing done!
time = 1.44227s
#v = 66803
#t = 377283
max_energy = 974.658
avg_energy = 3.89881
//////////////// pass 15 ////////////////
edge splitting...
fixed 0 tangled element
success = 532(532)
edge splitting done!
time = 0.209206s
#v = 67335
#t = 379654
max_energy = 974.658
avg_energy = 3.90086
edge collapsing...
fixed 0 tangled element
success(env) = 1
success = 528(734000)
success(env) = 0
success = 4(3449)
success(env) = 0
success = 0(34)
edge collapsing done!
time = 3.42997s
#v = 66803
#t = 377288
max_energy = 974.658
avg_energy = 3.89913
edge swapping...
fixed 0 tangled element
success3 = 9
success4 = 107
success5 = 4
success = 120(220782)
edge swapping done!
time = 1.20742s
#v = 66803
#t = 377283
max_energy = 974.658
avg_energy = 3.89895
vertex smoothing...
success = 28234(62901)
vertex smoothing done!
time = 1.43855s
#v = 66803
#t = 377283
max_energy = 974.658
avg_energy = 3.89858
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 16 ////////////////
edge splitting...
fixed 0 tangled element
success = 547(547)
edge splitting done!
time = 0.214451s
#v = 67350
#t = 379732
max_energy = 974.658
avg_energy = 3.90086
edge collapsing...
fixed 0 tangled element
success(env) = 1
success = 531(734134)
success(env) = 0
success = 7(3180)
success(env) = 0
success = 0(63)
edge collapsing done!
time = 3.40305s
#v = 66812
#t = 377329
max_energy = 974.658
avg_energy = 3.89895
edge swapping...
fixed 0 tangled element
success3 = 0
success4 = 90
success5 = 2
success = 92(220758)
edge swapping done!
time = 1.22439s
#v = 66812
#t = 377331
max_energy = 974.658
avg_energy = 3.89881
vertex smoothing...
success = 28060(62900)
vertex smoothing done!
time = 1.4332s
#v = 66812
#t = 377331
max_energy = 974.658
avg_energy = 3.89854
//////////////// pass 17 ////////////////
edge splitting...
fixed 0 tangled element
success = 537(537)
edge splitting done!
time = 0.211283s
#v = 67349
#t = 379716
max_energy = 974.658
avg_energy = 3.90088
edge collapsing...
fixed 0 tangled element
success(env) = 1
success = 537(734106)
success(env) = 0
success = 8(3427)
success(env) = 0
success = 0(158)
edge collapsing done!
time = 3.40603s
#v = 66804
#t = 377287
max_energy = 974.658
avg_energy = 3.89883
edge swapping...
fixed 0 tangled element
success3 = 1
success4 = 91
success5 = 1
success = 93(220719)
edge swapping done!
time = 1.21034s
#v = 66804
#t = 377287
max_energy = 974.658
avg_energy = 3.8987
vertex smoothing...
success = 27891(62896)
vertex smoothing done!
time = 1.43034s
#v = 66804
#t = 377287
max_energy = 974.658
avg_energy = 3.89841
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 18 ////////////////
edge splitting...
fixed 0 tangled element
success = 561(561)
edge splitting done!
time = 0.215108s
#v = 67365
#t = 379803
max_energy = 974.658
avg_energy = 3.90105
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 547(734225)
success(env) = 0
success = 3(3142)
success(env) = 0
success = 0(11)
edge collapsing done!
time = 3.46427s
#v = 66815
#t = 377343
max_energy = 974.658
avg_energy = 3.899
edge swapping...
fixed 0 tangled element
success3 = 0
success4 = 83
success5 = 3
success = 86(220761)
edge swapping done!
time = 1.24898s
#v = 66815
#t = 377346
max_energy = 974.658
avg_energy = 3.89885
vertex smoothing...
success = 27731(62908)
vertex smoothing done!
time = 1.42704s
#v = 66815
#t = 377346
max_energy = 974.658
avg_energy = 3.89854
//////////////// pass 19 ////////////////
edge splitting...
fixed 0 tangled element
success = 548(548)
edge splitting done!
time = 0.21084s
#v = 67363
#t = 379787
max_energy = 974.658
avg_energy = 3.9011
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 540(734217)
success(env) = 0
success = 13(3146)
success(env) = 0
success = 0(258)
edge collapsing done!
time = 3.45955s
#v = 66810
#t = 377318
max_energy = 974.658
avg_energy = 3.89883
edge swapping...
fixed 0 tangled element
success3 = 0
success4 = 82
success5 = 3
success = 85(220730)
edge swapping done!
time = 1.21697s
#v = 66810
#t = 377321
max_energy = 974.658
avg_energy = 3.89869
vertex smoothing...
success = 27537(62900)
vertex smoothing done!
time = 1.41512s
#v = 66810
#t = 377321
max_energy = 974.658
avg_energy = 3.89847
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 20 ////////////////
edge splitting...
fixed 0 tangled element
success = 546(546)
edge splitting done!
time = 0.210441s
#v = 67356
#t = 379751
max_energy = 974.658
avg_energy = 3.90088
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 535(734173)
success(env) = 0
success = 3(3291)
success(env) = 0
success = 0(11)
edge collapsing done!
time = 3.50897s
#v = 66818
#t = 377365
max_energy = 974.658
avg_energy = 3.89905
edge swapping...
fixed 0 tangled element
success3 = 5
success4 = 86
success5 = 2
success = 93(220794)
edge swapping done!
time = 1.25224s
#v = 66818
#t = 377362
max_energy = 974.658
avg_energy = 3.89888
vertex smoothing...
success = 27328(62908)
vertex smoothing done!
time = 1.41003s
#v = 66818
#t = 377362
max_energy = 974.658
avg_energy = 3.89859
//////////////// pass 21 ////////////////
edge splitting...
fixed 0 tangled element
success = 543(543)
edge splitting done!
time = 0.213542s
#v = 67361
#t = 379772
max_energy = 974.658
avg_energy = 3.90084
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 541(734200)
success(env) = 0
success = 6(3166)
success(env) = 0
success = 0(86)
edge collapsing done!
time = 3.39125s
#v = 66814
#t = 377342
max_energy = 974.658
avg_energy = 3.89889
edge swapping...
fixed 0 tangled element
success3 = 2
success4 = 82
success5 = 2
success = 86(220747)
edge swapping done!
time = 1.208s
#v = 66814
#t = 377342
max_energy = 974.658
avg_energy = 3.89871
vertex smoothing...
success = 27074(62910)
vertex smoothing done!
time = 1.42152s
#v = 66814
#t = 377342
max_energy = 974.658
avg_energy = 3.8985
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 22 ////////////////
edge splitting...
fixed 0 tangled element
success = 541(541)
edge splitting done!
time = 0.208438s
#v = 67355
#t = 379745
max_energy = 974.658
avg_energy = 3.9009
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 530(734194)
success(env) = 0
success = 10(3182)
success(env) = 0
success = 0(186)
edge collapsing done!
time = 3.39551s
#v = 66815
#t = 377346
max_energy = 974.658
avg_energy = 3.89914
edge swapping...
fixed 0 tangled element
success3 = 3
success4 = 90
success5 = 7
success = 100(220807)
edge swapping done!
time = 1.22283s
#v = 66815
#t = 377350
max_energy = 974.658
avg_energy = 3.89888
vertex smoothing...
success = 26817(62907)
vertex smoothing done!
time = 1.40049s
#v = 66815
#t = 377350
max_energy = 974.658
avg_energy = 3.89854
//////////////// pass 23 ////////////////
edge splitting...
fixed 0 tangled element
success = 557(557)
edge splitting done!
time = 0.210951s
#v = 67372
#t = 379827
max_energy = 974.658
avg_energy = 3.90137
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 553(734196)
success(env) = 0
success = 6(3219)
success(env) = 0
success = 0(119)
edge collapsing done!
time = 3.42688s
#v = 66813
#t = 377342
max_energy = 974.658
avg_energy = 3.89906
edge swapping...
fixed 0 tangled element
success3 = 5
success4 = 82
success5 = 1
success = 88(220765)
edge swapping done!
time = 1.2229s
#v = 66813
#t = 377338
max_energy = 974.658
avg_energy = 3.89885
vertex smoothing...
success = 26531(62908)
vertex smoothing done!
time = 1.39745s
#v = 66813
#t = 377338
max_energy = 974.658
avg_energy = 3.89847
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 24 ////////////////
edge splitting...
fixed 0 tangled element
success = 572(572)
edge splitting done!
time = 0.211653s
#v = 67385
#t = 379896
max_energy = 974.658
avg_energy = 3.90137
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 554(734302)
success(env) = 0
success = 8(3344)
success(env) = 0
success = 1(125)
success(env) = 0
success = 0(19)
edge collapsing done!
time = 3.4953s
#v = 66822
#t = 377387
max_energy = 974.658
avg_energy = 3.89914
edge swapping...
fixed 0 tangled element
success3 = 4
success4 = 91
success5 = 3
success = 98(220816)
edge swapping done!
time = 1.21202s
#v = 66822
#t = 377386
max_energy = 974.658
avg_energy = 3.89892
vertex smoothing...
success = 26120(62914)
vertex smoothing done!
time = 1.39364s
#v = 66822
#t = 377386
max_energy = 974.658
avg_energy = 3.89853
//////////////// pass 25 ////////////////
edge splitting...
fixed 0 tangled element
success = 565(565)
edge splitting done!
time = 0.210716s
#v = 67387
#t = 379904
max_energy = 974.658
avg_energy = 3.90124
edge collapsing...
fixed 0 tangled element
success(env) = 1
success = 555(734336)
success(env) = 0
success = 8(3525)
success(env) = 0
success = 2(125)
success(env) = 0
success = 0(18)
edge collapsing done!
time = 3.54823s
#v = 66822
#t = 377392
max_energy = 974.658
avg_energy = 3.89921
edge swapping...
fixed 0 tangled element
success3 = 6
success4 = 89
success5 = 4
success = 99(220825)
edge swapping done!
time = 1.2136s
#v = 66822
#t = 377390
max_energy = 974.658
avg_energy = 3.899
vertex smoothing...
success = 25748(62914)
vertex smoothing done!
time = 1.39154s
#v = 66822
#t = 377390
max_energy = 974.658
avg_energy = 3.89858
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 26 ////////////////
edge splitting...
fixed 0 tangled element
success = 561(561)
edge splitting done!
time = 0.211941s
#v = 67383
#t = 379888
max_energy = 974.658
avg_energy = 3.90121
edge collapsing...
fixed 0 tangled element
success(env) = 1
success = 550(734345)
success(env) = 0
success = 5(3402)
success(env) = 0
success = 0(33)
edge collapsing done!
time = 3.33875s
#v = 66828
#t = 377417
max_energy = 974.658
avg_energy = 3.89905
edge swapping...
fixed 0 tangled element
success3 = 3
success4 = 89
success5 = 6
success = 98(220861)
edge swapping done!
time = 1.21284s
#v = 66828
#t = 377420
max_energy = 974.658
avg_energy = 3.89887
vertex smoothing...
success = 25371(62923)
vertex smoothing done!
time = 1.38253s
#v = 66828
#t = 377420
max_energy = 974.658
avg_energy = 3.89852
//////////////// pass 27 ////////////////
edge splitting...
fixed 0 tangled element
success = 560(560)
edge splitting done!
time = 0.210534s
#v = 67388
#t = 379911
max_energy = 974.658
avg_energy = 3.90108
edge collapsing...
fixed 0 tangled element
success(env) = 1
success = 539(734321)
success(env) = 0
success = 7(3443)
success(env) = 0
success = 0(132)
edge collapsing done!
time = 3.36942s
#v = 66842
#t = 377492
max_energy = 974.658
avg_energy = 3.89893
edge swapping...
fixed 0 tangled element
success3 = 2
success4 = 79
success5 = 1
success = 82(220899)
edge swapping done!
time = 1.19437s
#v = 66842
#t = 377491
max_energy = 974.658
avg_energy = 3.89877
vertex smoothing...
success = 24856(62929)
vertex smoothing done!
time = 1.37562s
#v = 66842
#t = 377491
max_energy = 974.658
avg_energy = 3.89853
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// pass 28 ////////////////
edge splitting...
fixed 0 tangled element
success = 519(519)
edge splitting done!
time = 0.211897s
#v = 67361
#t = 379777
max_energy = 974.658
avg_energy = 3.90026
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 525(734235)
success(env) = 0
success = 5(3269)
success(env) = 0
success = 0(49)
edge collapsing done!
time = 3.379s
#v = 66831
#t = 377434
max_energy = 974.658
avg_energy = 3.89874
edge swapping...
fixed 0 tangled element
success3 = 0
success4 = 80
success5 = 4
success = 84(220851)
edge swapping done!
time = 1.19027s
#v = 66831
#t = 377438
max_energy = 974.658
avg_energy = 3.89858
vertex smoothing...
success = 24452(62922)
vertex smoothing done!
time = 1.37196s
#v = 66831
#t = 377438
max_energy = 974.658
avg_energy = 3.89845
//////////////// pass 29 ////////////////
edge splitting...
fixed 0 tangled element
success = 531(531)
edge splitting done!
time = 0.209948s
#v = 67362
#t = 379786
max_energy = 974.658
avg_energy = 3.90033
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 528(734241)
success(env) = 0
success = 3(3198)
success(env) = 0
success = 0(11)
edge collapsing done!
time = 3.41069s
#v = 66831
#t = 377439
max_energy = 974.658
avg_energy = 3.89871
edge swapping...
fixed 0 tangled element
success3 = 1
success4 = 80
success5 = 1
success = 82(220837)
edge swapping done!
time = 1.18639s
#v = 66831
#t = 377439
max_energy = 974.658
avg_energy = 3.89856
vertex smoothing...
success = 23995(62925)
vertex smoothing done!
time = 1.35901s
#v = 66831
#t = 377439
max_energy = 974.658
avg_energy = 3.89836
updating sclaing field ...
filter_energy = 9.74658
is_hit_min_edge_length = 1
//////////////// postprocessing ////////////////
edge collapsing...
fixed 0 tangled element
success(env) = 0
success = 11(733505)
success(env) = 0
success = 8(674)
success(env) = 0
success = 0(421)
edge collapsing done!
time = 3.52037s
#v = 66812
#t = 377339
max_energy = 974.658
avg_energy = 3.89837
edge collapsing...
fixed 0 tangled element
success(env) = 5459
success = 20627(527121)
success(env) = 369
success = 645(399412)
success(env) = 23
success = 36(68449)
success(env) = 3
success = 6(4196)
success(env) = 1
success = 2(823)
success(env) = 0
success = 0(235)
edge collapsing done!
time = 6.69823s
#v = 45496
#t = 255166
max_energy = 974.658
avg_energy = 4.79047
edge swapping...
fixed 0 tangled element
success3 = 4259
success4 = 17129
success5 = 1906
success = 23294(178366)
edge swapping done!
time = 1.15732s
#v = 45496
#t = 252813
max_energy = 974.658
avg_energy = 4.57965
edge collapsing...
fixed 0 tangled element
success(env) = 658
success = 1889(404061)
success(env) = 50
success = 84(163860)
success(env) = 3
success = 6(11285)
success(env) = 1
success = 2(866)
success(env) = 0
success = 0(125)
edge collapsing done!
time = 3.61327s
#v = 43515
#t = 240875
max_energy = 974.658
avg_energy = 4.72177
edge swapping...
fixed 0 tangled element
success3 = 946
success4 = 4023
success5 = 557
success = 5526(133307)
edge swapping done!
time = 0.806701s
#v = 43515
#t = 240486
max_energy = 974.658
avg_energy = 4.67336
edge collapsing...
fixed 0 tangled element
success(env) = 80
success = 250(385934)
success(env) = 5
success = 8(30278)
success(env) = 0
success = 0(1144)
edge collapsing done!
time = 2.51214s
#v = 43257
#t = 238923
max_energy = 974.658
avg_energy = 4.6951
edge swapping...
fixed 0 tangled element
success3 = 139
success4 = 662
success5 = 110
success = 911(122192)
edge swapping done!
time = 0.689117s
#v = 43257
#t = 238894
max_energy = 974.658
avg_energy = 4.68724
edge collapsing...
fixed 0 tangled element
success(env) = 8
success = 32(383485)
success(env) = 1
success = 2(4002)
success(env) = 0
success = 0(343)
edge collapsing done!
time = 2.46064s
#v = 43223
#t = 238684
max_energy = 974.658
avg_energy = 4.69098
edge swapping...
fixed 0 tangled element
success3 = 31
success4 = 122
success5 = 21
success = 174(120509)
edge swapping done!
time = 0.7153s
#v = 43223
#t = 238674
max_energy = 974.658
avg_energy = 4.68922
[2025-04-23 06:54:17.832] [float-tetwild] [info] mesh optimization 240.045s
[2025-04-23 06:54:17.832] [float-tetwild] [info] 
[2025-04-23 06:54:18.068] [float-tetwild] [info] correct_tracked_surface_orientation done

We want to use this mesh within FEniCSx, so we extract the cells, nodes and volume markers from the mesh instance in wildmeshing.

point_array, cell_array, marker = tetra.get_tet_mesh()
[2025-04-23 06:54:19.366] [float-tetwild] [info] after winding number
[2025-04-23 06:54:19.366] [float-tetwild] [info] #v = 43223
[2025-04-23 06:54:19.369] [float-tetwild] [info] #t = 187041
[2025-04-23 06:54:19.369] [float-tetwild] [info] winding number 1.74539e+09s
[2025-04-23 06:54:19.369] [float-tetwild] [info] 

We can pass these into the DOLFINx mesh constructor called create_mesh.

mesh = dolfinx.mesh.create_mesh(
    MPI.COMM_WORLD,
    cell_array.astype(np.int64),
    point_array,
    ufl.Mesh(basix.ufl.element("Lagrange", "tetrahedron", 1, shape=(3,))),
)

We add the cell markers as CellTags with the following lines of code

local_entities, local_values = dolfinx.io.gmshio.distribute_entity_data(
    mesh,
    mesh.topology.dim,
    cell_array.astype(np.int64),
    marker.flatten().astype(np.int32),
)
adj = dolfinx.graph.adjacencylist(local_entities)
ct = dolfinx.mesh.meshtags_from_entities(
    mesh,
    mesh.topology.dim,
    adj,
    local_values.astype(np.int32, copy=False),
)

We store the mesh and cell markers to file for usage in later tutorials

with dolfinx.io.XDMFFile(mesh.comm, folder / "brain.xdmf", "w") as xdmf:
    xdmf.write_mesh(mesh)
    xdmf.write_meshtags(ct, mesh.geometry)

We visualize the resulting mesh and cell tags with pyvista

pv_grid = pyvista.UnstructuredGrid(*dolfinx.plot.vtk_mesh(mesh))
pv_grid.cell_data["marker"] = ct.values
pyvista.start_xvfb(1.0)
plotter = pyvista.Plotter()
plotter.add_mesh(pv_grid)
plotter.show()

We also create som slices to see the various markers within the mesh

slices = pv_grid.slice_along_axis(n=7, axis="y")
slices.plot()