Skip to content

Commit

Permalink
Merge pull request #116 from fusion-energy/adding_scaling_option_when…
Browse files Browse the repository at this point in the history
…_adding_cadquery_object

Adding scaling option when adding cadquery objects
  • Loading branch information
shimwell authored Jan 31, 2025
2 parents fd5c734 + f432de9 commit 943f698
Show file tree
Hide file tree
Showing 2 changed files with 133 additions and 9 deletions.
16 changes: 13 additions & 3 deletions src/cad_to_dagmc/core.py
Original file line number Diff line number Diff line change
Expand Up @@ -503,6 +503,7 @@ def add_cadquery_object(
cq.assembly.Assembly | cq.occ_impl.shapes.Compound | cq.occ_impl.shapes.Solid
),
material_tags: list[str] | None,
scale_factor: float = 1.0,
) -> int:
"""Loads the parts from CadQuery object into the model.
Expand All @@ -514,6 +515,10 @@ def add_cadquery_object(
same order as the volumes in the geometry added (STP file and
CadQuery objects) and match the material tags used in the
neutronics code (e.g. OpenMC).
scale_factor: a scaling factor to apply to the geometry that can be
used to increase the size or decrease the size of the geometry.
Useful when converting the geometry to cm for use in neutronics
simulations.
Returns:
int: number of volumes in the stp file.
Expand All @@ -527,12 +532,17 @@ def add_cadquery_object(
else:
iterable_solids = cadquery_object.val().Solids()

_check_material_tags(material_tags, iterable_solids)
if scale_factor == 1.0:
scaled_iterable_solids = iterable_solids
else:
scaled_iterable_solids = [part.scale(scale_factor) for part in iterable_solids]

_check_material_tags(material_tags, scaled_iterable_solids)
if material_tags:
self.material_tags = self.material_tags + material_tags
self.parts = self.parts + iterable_solids
self.parts = self.parts + scaled_iterable_solids

return len(iterable_solids)
return len(scaled_iterable_solids)

def export_unstructured_mesh_file(
self,
Expand Down
126 changes: 120 additions & 6 deletions tests/test_python_api.py
Original file line number Diff line number Diff line change
@@ -1,16 +1,16 @@
import os
from pathlib import Path
import pytest
import cadquery as cq
from cad_to_dagmc import CadToDagmc
import warnings
from cad_to_dagmc.core import _check_material_tags

from pathlib import Path

import cadquery as cq
import gmsh
import pymoab as mb
import pytest
from pymoab import core, types

from cad_to_dagmc import CadToDagmc
from cad_to_dagmc.core import _check_material_tags


def get_volumes_and_materials_from_h5m(filename: str) -> dict:
"""Reads in a DAGMC h5m file and uses PyMoab to find the volume ids with
Expand Down Expand Up @@ -219,3 +219,117 @@ def test_check_material_tags_too_long():
assert issubclass(w[-1].category, UserWarning)
assert "Material tag" in str(w[-1].message)
assert "a" * 29 in str(w[-1].message)


@pytest.mark.parametrize(
"scale_factor, expected_width",
[
(1, 10.0),
(2, 20.0),
(10, 100.0),
],
)
def test_scaling_factor_when_adding_stp(scale_factor, expected_width):

c2d = CadToDagmc()
c2d.add_stp_file("tests/single_cube.stp", scale_factor=scale_factor)
c2d.export_gmsh_mesh_file(f"st_test_scaling_factor_{scale_factor}.msh")

gmsh.initialize()
gmsh.open(f"st_test_scaling_factor_{scale_factor}.msh")
_, node_coords, _ = gmsh.model.mesh.getNodes()

# Reshape the node coordinates into a 2D array
node_coords = node_coords.reshape(-1, 3)

# Calculate the bounding box
min_coords = node_coords.min(axis=0)
max_coords = node_coords.max(axis=0)

width_x = max_coords[0] - min_coords[0]
width_y = max_coords[1] - min_coords[1]
width_z = max_coords[2] - min_coords[2]

gmsh.finalize()

assert width_x == expected_width
assert width_y == expected_width
assert width_z == expected_width


@pytest.mark.parametrize(
"scale_factor, expected_width",
[
(1, 10.0),
(2, 20.0),
(10, 100.0),
],
)
def test_scaling_factor_when_adding_cq_object(scale_factor, expected_width):

box = cq.Workplane("XY").box(10, 10, 10)
c2d = CadToDagmc()
c2d.add_cadquery_object(box, scale_factor=scale_factor, material_tags=["mat1"])
c2d.export_gmsh_mesh_file(f"cq_test_scaling_factor_{scale_factor}.msh")

gmsh.initialize()
gmsh.open(f"cq_test_scaling_factor_{scale_factor}.msh")
_, node_coords, _ = gmsh.model.mesh.getNodes()

# Reshape the node coordinates into a 2D array
node_coords = node_coords.reshape(-1, 3)

# Calculate the bounding box
min_coords = node_coords.min(axis=0)
max_coords = node_coords.max(axis=0)

width_x = max_coords[0] - min_coords[0]
width_y = max_coords[1] - min_coords[1]
width_z = max_coords[2] - min_coords[2]

gmsh.finalize()

assert width_x == expected_width
assert width_y == expected_width
assert width_z == expected_width


@pytest.mark.parametrize(
"scale_factor, expected_x_width, expected_y_width, expected_z_width",
[
(1, 20.0, 10.0, 10.0),
(2, 40.0, 20.0, 20.0),
(10, 200.0, 100.0, 100.0),
],
)
def test_two_box_scaling_factor_when_adding_cq_object(
scale_factor, expected_x_width, expected_y_width, expected_z_width
):

box = cq.Workplane("XY").box(10, 10, 10)
box2 = cq.Workplane("XY").moveTo(10, 0).box(10, 10, 10)
c2d = CadToDagmc()
c2d.add_cadquery_object(box, scale_factor=scale_factor, material_tags=["mat1"])
c2d.add_cadquery_object(box2, scale_factor=scale_factor, material_tags=["mat1"])
c2d.export_gmsh_mesh_file(f"cq_test_2_box_scaling_factor_{scale_factor}.msh")

gmsh.initialize()
gmsh.open(f"cq_test_2_box_scaling_factor_{scale_factor}.msh")
_, node_coords, _ = gmsh.model.mesh.getNodes()

# Reshape the node coordinates into a numpy 2D array
node_coords = node_coords.reshape(-1, 3)

# Calculate the bounding box
min_coords = node_coords.min(axis=0)
max_coords = node_coords.max(axis=0)

width_x = max_coords[0] - min_coords[0]
width_y = max_coords[1] - min_coords[1]
width_z = max_coords[2] - min_coords[2]

gmsh.finalize()

assert width_x == expected_x_width
assert width_y == expected_y_width
assert width_z == expected_z_width

0 comments on commit 943f698

Please sign in to comment.