Skip to content

File Conversion

fdtdx.core.conversion.load_stl.load_stl(stl, permittivity, target_shape, voxel_size, ambient_permittivity=1.0)

Loads an STL file and converts it to a voxelized permittivity array.

This function takes either an STL file path or a trimesh.Trimesh object, voxelizes it at the specified resolution, and converts it to a permittivity array where voxels inside the mesh have the specified permittivity value and voxels outside have the ambient permittivity value.

Parameters:

Name Type Description Default
stl str | Trimesh

Path to an STL file or a trimesh.Trimesh object containing the mesh

required
permittivity float

Relative permittivity value for the object

required
target_shape tuple[int, int, int]

Desired output shape as (nx, ny, nz)

required
voxel_size float

Size of each voxel in the initial voxelization

required
ambient_permittivity float

Relative permittivity value for the surrounding medium (default: 1.0)

1.0

Returns:

Type Description
Array

jax.Array: 4D array of shape (nx, ny, nz, 3) containing the permittivity values, where the last dimension represents the permittivity for each spatial direction

Source code in src/fdtdx/core/conversion/load_stl.py
def load_stl(
    stl: str | trimesh.Trimesh,
    permittivity: float,
    target_shape: tuple[int, int, int],
    voxel_size: float,
    ambient_permittivity: float = 1.0,
) -> jax.Array:
    """Loads an STL file and converts it to a voxelized permittivity array.

    This function takes either an STL file path or a trimesh.Trimesh object,
    voxelizes it at the specified resolution, and converts it to a permittivity
    array where voxels inside the mesh have the specified permittivity value
    and voxels outside have the ambient permittivity value.

    Args:
        stl: Path to an STL file or a trimesh.Trimesh object containing the mesh
        permittivity: Relative permittivity value for the object
        target_shape: Desired output shape as (nx, ny, nz)
        voxel_size: Size of each voxel in the initial voxelization
        ambient_permittivity: Relative permittivity value for the surrounding medium (default: 1.0)

    Returns:
        jax.Array: 4D array of shape (nx, ny, nz, 3) containing the permittivity values,
            where the last dimension represents the permittivity for each spatial direction
    """
    if isinstance(stl, str):
        your_mesh: trimesh.Trimesh = trimesh.load_mesh(stl)  # type: ignore
    else:
        your_mesh = stl

    mesh_voxel = your_mesh.voxelized(voxel_size).fill()

    occupancy_array = mesh_voxel.encoding.dense

    resized_occupancy_array = jax.image.resize(occupancy_array, target_shape, method="nearest")

    object_permittivity = jnp.where(resized_occupancy_array > 0, permittivity, ambient_permittivity)

    return jnp.tile(object_permittivity[..., jnp.newaxis], (1, 1, 1, 3))

fdtdx.core.conversion.export.export_stl(matrix, stl_filename, voxel_grid_size=(1, 1, 1))

Export a 3D boolean matrix to an STL file.

Converts a 3D boolean matrix into a mesh representation and saves it as an STL file. True values in the matrix are converted to solid voxels in the output mesh.

Parameters:

Name Type Description Default
matrix ndarray

3D boolean numpy array representing the voxel grid.

required
stl_filename Path | str

Output STL file path.

required
voxel_grid_size tuple[int, int, int]

Physical size of each voxel as (x, y, z) integers. Defaults to (1, 1, 1).

(1, 1, 1)

Raises:

Type Description
Exception

If input matrix is not 3-dimensional.

Source code in src/fdtdx/core/conversion/export.py
def export_stl(
    matrix: np.ndarray,
    stl_filename: Path | str,
    voxel_grid_size: tuple[int, int, int] = (1, 1, 1),
):
    """Export a 3D boolean matrix to an STL file.

    Converts a 3D boolean matrix into a mesh representation and saves it as an STL file.
    True values in the matrix are converted to solid voxels in the output mesh.

    Args:
        matrix: 3D boolean numpy array representing the voxel grid.
        stl_filename: Output STL file path.
        voxel_grid_size: Physical size of each voxel as (x, y, z) integers. Defaults to (1, 1, 1).

    Raises:
        Exception: If input matrix is not 3-dimensional.
    """
    if matrix.ndim != 3:
        raise Exception(f"Invalid matrix shape: {matrix.shape}")
    scaling_factor = np.asarray(voxel_grid_size, dtype=int)

    d0, d1, d2 = (
        matrix.shape[0] + 1,
        matrix.shape[1] + 1,
        matrix.shape[2] + 1,
    )
    vertex_shape = (d0, d1, d2)
    num_vertices = d0 * d1 * d2
    num_voxels = math.prod(matrix.shape)

    matrix_shape = (matrix.shape[0], matrix.shape[1], matrix.shape[2])
    x, y, z = idx_to_xyz(np.arange(num_voxels), matrix_shape)
    matrix_flat = matrix.reshape(-1)
    stacked_idx = np.stack([x, y, z], axis=-1)

    use_sides_arr = np.zeros((num_voxels, 6), dtype=bool)
    use_sides_arr[matrix_flat & (stacked_idx[..., 0] == 0), 0] = True  # left
    use_sides_arr[matrix_flat & (stacked_idx[..., 1] == 0), 1] = True  # front
    use_sides_arr[matrix_flat & (stacked_idx[..., 2] == 0), 2] = True  # bottom
    use_sides_arr[matrix_flat & (stacked_idx[..., 0] == matrix.shape[0] - 1), 3] = True  # right
    use_sides_arr[matrix_flat & (stacked_idx[..., 1] == matrix.shape[1] - 1), 4] = True  # back
    use_sides_arr[matrix_flat & (stacked_idx[..., 2] == matrix.shape[2] - 1), 5] = True  # top

    inv_matrix = ~matrix
    x_p1, y_p1, z_p1 = x + 1, y + 1, z + 1
    x_p1[x_p1 == matrix.shape[0]] = matrix.shape[0] - 1
    y_p1[y_p1 == matrix.shape[1]] = matrix.shape[1] - 1
    z_p1[z_p1 == matrix.shape[2]] = matrix.shape[2] - 1

    use_sides_arr[matrix_flat & inv_matrix[x - 1, y, z], 0] = True  # left
    use_sides_arr[matrix_flat & inv_matrix[x, y - 1, z], 1] = True  # front
    use_sides_arr[matrix_flat & inv_matrix[x, y, z - 1], 2] = True  # bottom
    use_sides_arr[matrix_flat & inv_matrix[x_p1, y, z], 3] = True  # right
    use_sides_arr[matrix_flat & inv_matrix[x, y_p1, z], 4] = True  # back
    use_sides_arr[matrix_flat & inv_matrix[x, y, z_p1], 5] = True  # top

    v_idx = np.asarray(
        [
            xyz_to_idx(x, y, z, shape=vertex_shape),
            xyz_to_idx(x, y, z + 1, shape=vertex_shape),
            xyz_to_idx(x, y + 1, z, shape=vertex_shape),
            xyz_to_idx(x, y + 1, z + 1, shape=vertex_shape),
            xyz_to_idx(x + 1, y, z, shape=vertex_shape),
            xyz_to_idx(x + 1, y, z + 1, shape=vertex_shape),
            xyz_to_idx(x + 1, y + 1, z, shape=vertex_shape),
            xyz_to_idx(x + 1, y + 1, z + 1, shape=vertex_shape),
        ]
    )

    faces_raw = np.asarray(
        [
            np.stack([[v_idx[0], v_idx[1], v_idx[2]], [v_idx[1], v_idx[3], v_idx[2]]], axis=-1),  # left
            np.stack([[v_idx[0], v_idx[4], v_idx[5]], [v_idx[5], v_idx[1], v_idx[0]]], axis=-1),  # front
            np.stack([[v_idx[0], v_idx[2], v_idx[6]], [v_idx[6], v_idx[4], v_idx[0]]], axis=-1),  # bottom
            np.stack([[v_idx[4], v_idx[6], v_idx[7]], [v_idx[7], v_idx[5], v_idx[4]]], axis=-1),  # right
            np.stack([[v_idx[2], v_idx[3], v_idx[7]], [v_idx[7], v_idx[6], v_idx[2]]], axis=-1),  # back
            np.stack([[v_idx[1], v_idx[5], v_idx[3]], [v_idx[3], v_idx[5], v_idx[7]]], axis=-1),  # top
        ]
    )
    faces_raw = faces_raw.transpose((2, 0, 3, 1))

    # select used faces only
    faces = faces_raw[use_sides_arr].reshape(-1, 3)
    vertex_arr = np.stack(
        idx_to_xyz(
            np.arange(num_vertices),
            shape=vertex_shape,
        ),
        axis=-1,
    )
    vertex_arr = vertex_arr * scaling_factor

    # export to trimesh
    mesh = trimesh.Trimesh(
        vertices=vertex_arr,
        faces=faces,
        validate=False,
    )
    mesh.export(stl_filename)

fdtdx.core.conversion.export.idx_to_xyz(idx, shape)

Convert flattened array indices to 3D coordinates.

Parameters:

Name Type Description Default
idx ndarray

Array of flattened indices.

required
shape tuple[int, int, int]

3D shape tuple (d0, d1, d2) of the original array.

required

Returns:

Name Type Description
tuple tuple[ndarray, ndarray, ndarray]

(x, y, z) coordinate arrays corresponding to the input indices.

Source code in src/fdtdx/core/conversion/export.py
def idx_to_xyz(idx: np.ndarray, shape: tuple[int, int, int]) -> tuple[np.ndarray, np.ndarray, np.ndarray]:
    """Convert flattened array indices to 3D coordinates.

    Args:
        idx: Array of flattened indices.
        shape: 3D shape tuple (d0, d1, d2) of the original array.

    Returns:
        tuple: (x, y, z) coordinate arrays corresponding to the input indices.
    """
    _, d1, d2 = shape
    x = idx // (d1 * d2)
    y = (idx // d2) % d1
    z = idx % d2
    return x, y, z

fdtdx.core.conversion.export.xyz_to_idx(x, y, z, shape)

Convert 3D coordinates to flattened array indices.

This is the inverse operation of idx_to_xyz(). Used for converting physical coordinates back to array indices.

Parameters:

Name Type Description Default
x ndarray

Array of x coordinates.

required
y ndarray

Array of y coordinates.

required
z ndarray

Array of z coordinates.

required
shape tuple[int, int, int]

3D shape tuple (d0, d1, d2) of the target array.

required

Returns:

Type Description
ndarray

np.ndarray: Flattened indices corresponding to the input coordinates.

Source code in src/fdtdx/core/conversion/export.py
def xyz_to_idx(
    x: np.ndarray,
    y: np.ndarray,
    z: np.ndarray,
    shape: tuple[int, int, int],
) -> np.ndarray:
    """Convert 3D coordinates to flattened array indices.

    This is the inverse operation of idx_to_xyz(). Used for converting physical coordinates
    back to array indices.

    Args:
        x: Array of x coordinates.
        y: Array of y coordinates.
        z: Array of z coordinates.
        shape: 3D shape tuple (d0, d1, d2) of the target array.

    Returns:
        np.ndarray: Flattened indices corresponding to the input coordinates.
    """
    """Convert 3D coordinates to flattened array indices.

    Args:
        x: Array of x coordinates.
        y: Array of y coordinates.
        z: Array of z coordinates.
        shape: 3D shape tuple (d0, d1, d2) of the target array.

    Returns:
        np.ndarray: Flattened indices corresponding to the input coordinates.
    """
    _, d1, d2 = shape
    return x * (d1 * d2) + y * (d2) + z

fdtdx.core.conversion.import_utils.gds_to_numpy(file_path, resolution, layer, datatype=None)

Converts GDSII geometry on a specific layer to a grid-based mask.

Parameters:

Name Type Description Default
file_path str | Path

Path to the GDSII file

required
resolution float

The size of each grid cell in the same units as the GDSII coordinates

required
layer int

The GDSII layer number to extract

required
datatype int | None

The GDSII datatype number to extract. Defaults to layer value.

None

Returns:

Type Description
ndarray[bool_]

np.ndarray[np.bool_]: A 2D NumPy array mask where True indicates geometry presence

Source code in src/fdtdx/core/conversion/import_utils.py
def gds_to_numpy(
    file_path: str | Path,
    resolution: float,
    layer: int,
    datatype: int | None = None,
) -> np.ndarray[np.bool_]:
    """
    Converts GDSII geometry on a specific layer to a grid-based mask.

    Args:
        file_path: Path to the GDSII file
        resolution: The size of each grid cell in the same units as the GDSII coordinates
        layer: The GDSII layer number to extract
        datatype: The GDSII datatype number to extract. Defaults to layer value.

    Returns:
        np.ndarray[np.bool_]: A 2D NumPy array mask where True indicates geometry presence
    """
    if datatype is None:
        datatype = layer
    gdsii_lib = gdspy.GdsLibrary()
    gdsii_lib.read_gds(file_path)

    cell = gdsii_lib.top_level()[0]  # Assuming one top-level cell

    # Determine bounding box for grid dimensions
    bbox = cell.get_bounding_box()
    x_min, y_min, x_max, y_max = bbox[0][0], bbox[0][1], bbox[1][0], bbox[1][1]

    # Create grid coordinates
    x_coords = np.arange(x_min, x_max + resolution, resolution)
    y_coords = np.arange(y_min, y_max + resolution, resolution)
    grid_x, grid_y = np.meshgrid(x_coords, y_coords)

    # Check each grid point for geometry presence
    polygons = cell.get_polygons(by_spec=(layer, datatype))
    points = np.stack([grid_x, grid_y], axis=-1).reshape(-1, 2)
    result_bool_list = gdspy.inside(
        points=points,
        polygons=polygons,
    )
    mask = np.asarray(result_bool_list).reshape(grid_x.shape).T
    return mask