DAGMC H5M Files#
The primary output format. Creates triangular surface meshes for use with DAGMC-enabled transport codes (OpenMC, MCNP, FLUKA, etc.).
Basic Export#
import cadquery as cq
from cad_to_dagmc import CadToDagmc
assembly = cq.Assembly()
assembly.add(cq.Workplane("XY").sphere(10), name="sphere")
model = CadToDagmc()
model.add_cadquery_object(assembly, material_tags=["tungsten"])
model.export_dagmc_h5m_file(
filename="dagmc.h5m",
min_mesh_size=0.5,
max_mesh_size=10.0,
)
With Different Backends#
Choose between h5py (default) and pymoab backends:
import cadquery as cq
from cad_to_dagmc import CadToDagmc
sphere1 = cq.Workplane().sphere(5)
sphere2 = cq.Workplane().moveTo(10, 0).sphere(2)
assembly = cq.Assembly()
assembly.add(sphere1)
assembly.add(sphere2)
model = CadToDagmc()
model.add_cadquery_object(cadquery_object=assembly, material_tags=["mat1", "mat2"])
# Export using h5py backend (default, no MOAB needed)
model.export_dagmc_h5m_file(
filename="dagmc_h5py.h5m",
h5m_backend="h5py",
)
# Export using pymoab backend (requires MOAB installation)
model.export_dagmc_h5m_file(
filename="dagmc_pymoab.h5m",
h5m_backend="pymoab",
)
With Meshing Backends#
Choose between GMSH (default) and CadQuery meshing:
# GMSH backend (default) - full control over mesh parameters
model.export_dagmc_h5m_file(
filename="dagmc.h5m",
meshing_backend="gmsh",
min_mesh_size=0.5,
max_mesh_size=10.0,
mesh_algorithm=1,
)
# CadQuery backend - simpler, uses CadQuery's tessellation
model.export_dagmc_h5m_file(
filename="dagmc.h5m",
meshing_backend="cadquery",
tolerance=0.1,
angular_tolerance=0.1,
)
API Reference#
export_dagmc_h5m_file()#
Common Parameters:
Parameter |
Type |
Default |
Description |
|---|---|---|---|
|
str |
“dagmc.h5m” |
Output file path |
|
float |
1.0 |
Geometry scale factor |
|
bool |
True |
Imprint shared surfaces |
|
str |
None |
Void space material tag |
Backend Selection:
Parameter |
Type |
Default |
Description |
|---|---|---|---|
|
str |
auto |
|
|
str |
“h5py” |
|
GMSH Backend Parameters:
These parameters only apply when using meshing_backend="gmsh" (or when auto-selected):
Parameter |
Type |
Default |
Description |
|---|---|---|---|
|
float |
None |
Minimum mesh element size |
|
float |
None |
Maximum mesh element size |
|
int |
1 |
GMSH meshing algorithm (1-10) |
|
str |
“file” |
CAD transfer method: |
|
dict |
None |
Per-volume mesh sizes. Keys can be volume IDs (int) or material tag names (str). |
|
list |
None |
Volume IDs (int) or material tags (str) for conformal volume mesh |
|
str |
“umesh.vtk” |
Output filename for unstructured volume mesh |
CadQuery Backend Parameters:
These parameters only apply when using meshing_backend="cadquery":
Parameter |
Type |
Default |
Description |
|---|---|---|---|
|
float |
0.1 |
Linear tolerance for tessellation |
|
float |
0.1 |
Angular tolerance for tessellation |
Warning
Do not mix GMSH and CadQuery backend parameters in the same call. If you provide parameters from both backends without explicitly setting meshing_backend, an error will be raised.
Returns:
str- Path to the created h5m fileOr
tuple[str, str]- (h5m_path, vtk_path) ifunstructured_volumesis set
Using in OpenMC#
Load the DAGMC geometry in OpenMC:
import openmc
# Define materials (names must match material tags)
tungsten = openmc.Material(name="tungsten")
tungsten.add_element("W", 1.0)
tungsten.set_density("g/cm3", 19.3)
steel = openmc.Material(name="steel")
steel.add_element("Fe", 1.0)
steel.set_density("g/cm3", 7.8)
materials = openmc.Materials([tungsten, steel])
# Load the DAGMC geometry
dag_universe = openmc.DAGMCUniverse(filename="dagmc.h5m")
geometry = openmc.Geometry(root=dag_universe.bounded_universe())
# Set up settings
settings = openmc.Settings()
settings.batches = 10
settings.particles = 1000
settings.run_mode = "fixed source"
model = openmc.Model(geometry=geometry, materials=materials, settings=settings)
model.run()
See Also#
H5M Backends - h5py vs pymoab comparison
GMSH Backend - GMSH meshing options
CadQuery Backend - CadQuery meshing options
Conformal Meshes - Combined surface and volume output