API reference

This section contains the API reference and usage information for ciclope.

ciclope core modules

Voxel Finite Elements

Ciclope module for voxel Finite Element model generation

ciclope.core.voxelFE.matpropdictionary(proplist)

Compose dictionary of material properties and property mapping GV ranges.

proplist

List of material property files followed by the corresponding Gray Value range for material mapping.

matpropdict

Dictionary of material properties for material property mapping: matprop = {

“file”: [“prop.inp”, “property_temp_bone.inp”, …], “range”: [[250, 255], [0, 250], …],

}

ciclope.core.voxelFE.mesh2voxelfe(mesh, templatefile, fileout='pippo.inp', matprop=None, keywords=['NSET', 'ELSET'], eltype='C3D8', matpropbits=8, refnode=None, verbose=False)

Generate ABAQUS voxel Finite Element (FE) input file from 3D Unstructured Grid mesh data. The file written is an input file (.INP) in ABAQUS syntax that can be solved using ABAQUS or CALCULIX. The user can define a material mapping strategy for the conversion of local GVs to local material properties in the FE model. Material mapping laws are defined in separate template file(s) (see “prop.inp” and “property_temp_bone.inp” for examples). Boundary conditions, analysis type and output requests are defined in a separate template file (see “tmp.inp” for an example). Info on analysis definition at: https://abaqus-docs.mit.edu/2017/English/SIMACAECAERefMap/simacae-m-Sim-sb.htm#simacae-m-Sim-sb

meshmeshio

Unstructured grid mesh.

templatefilestr

Analysis template file.

fileoutstr

Output .INP file.

matpropdict

Dictionary of material properties for material property mapping: matprop = {

“file”: [“prop.inp”, “property_temp_bone.inp”, …], “range”: [[250, 255], [0, 250], …],

}

keywordsstr

SUPPORTED ABAQUS KEYWORDS:

For a list of all Abaqus keywords and their description visit: https://abaqus-docs.mit.edu/2017/English/SIMACAECAERefMap/simacae-c-gen-kwbrowser.htm#simacae-c-gen-kwbrowser__simacae-gen-xsl-U

  • ‘NSET’:

    Create boundary node sets. (Default = ON) If ‘NSET’ is specified, the following node sets are created:

    • NODES_X0: Nodes on WEST (X-) surface of 3D model.

    • NODES_X1: Nodes on EAST (X+) surface of 3D model.

    • NODES_Y0: Nodes on SOUTH (Y-) surface of 3D model.

    • NODES_Y1: Nodes on NORTH (Y+) surface of 3D model.

    • NODES_Z0: Nodes on BOTTOM (Z-) surface of 3D model.

    • NODES_Z1: Nodes on TOP (Z+) surface of 3D model.

    • NODES_X0Y0Z0: 2 nodes on (0,0,0) model corner.

    • NODES_X0Y0Z1: 2 nodes on (0,0,1) model corner.

    These node sets are available for boundary conditions definition.

  • ‘ELSET’:

    Create boundary element sets. (Default = ON) If ‘ELSET’ is specified, the following element sets are created:

    • ELEMS_X0: Elements of WEST (X-) surface of 3D model.

    • ELEMS_X1: Elements of EAST (X+) surface of 3D model.

    • ELEMS_Y0: Elements of SOUTH (Y-) surface of 3D model.

    • ELEMS_Y1: Elements of NORTH (Y+) surface of 3D model.

    • ELEMS_Z0: Elements of BOTTOM (Z-) surface of 3D model.

    • ELEMS_Z1: Elements of TOP (Z+) surface of 3D model.

  • ‘PROPERTY’:

    Define an external material mapping law from template file. (Default = None) Use in combination with ‘matprop’ dictionary of material property files and corresponding GV ranges for the material mapping.

eltypestr

FE element type. The default is eight-node brick element (C3D8 and F3D8). See CalculiX node convention (Par. 6.2.1) at: http://www.dhondt.de/ccx_2.15.pdf

matpropbitsint

Bit depth for material mapping.

refnode

Reference node coordinates [REF_NODE_x, REF_NODE_y, REF_NODE_z] for kinematic coupling.

verbosebool

Activate verbose output.

ciclope.core.voxelFE.vol2ugrid(voldata, voxelsize=[1, 1, 1], GVmin=0, refnodes=False, verbose=False)

Generate unstructured grid mesh from 3D volume data.

voldatandarray

3D voxel data.

voxelsizefloat

3D model voxelsize.

matpropbitsint

Bit depth for material mapping.

GVmin

Minimum Grey Value considered for meshing. By default, GVmin=0: all zeroes are considered background.

refnodesbool

Return dictionary of reference nodes on the model boundaries. Even if this option is not activated, the returned mesh will contain the following nodes and elements sets:

  • NODES_X0: Nodes on WEST (X-) surface of 3D model.

  • NODES_X1: Nodes on EAST (X+) surface of 3D model.

  • NODES_Y0: Nodes on SOUTH (Y-) surface of 3D model.

  • NODES_Y1: Nodes on NORTH (Y+) surface of 3D model.

  • NODES_Z0: Nodes on BOTTOM (Z-) surface of 3D model.

  • NODES_Z1: Nodes on TOP (Z+) surface of 3D model.

  • NODES_X0Y0Z0: 2 nodes on (0,0,0) model corner.

  • NODES_X0Y0Z1: 2 nodes on (0,0,1) model corner.

  • ELEMS_X0: Elements of WEST (X-) surface of 3D model.

  • ELEMS_X1: Elements of EAST (X+) surface of 3D model.

  • ELEMS_Y0: Elements of SOUTH (Y-) surface of 3D model.

  • ELEMS_Y1: Elements of NORTH (Y+) surface of 3D model.

  • ELEMS_Z0: Elements of BOTTOM (Z-) surface of 3D model.

  • ELEMS_Z1: Elements of TOP (Z+) surface of 3D model.

verbosebool

Activate verbose output.

meshmeshio

Unstructured grid mesh.

refnodeslist

centroids on the model boundaries (X0, X1, Y0, Y1, Z0, Z1)

Tetrahedra Finite Elements

Ciclope module for tetrahedra Finite Element model generation

ciclope.core.tetraFE.cgal_mesh(bwimage, voxelsize, meshtype='both', max_facet_distance=0.0, max_cell_circumradius=0.0)

Generate mesh of from binary volume data using CGAL. The mesh is generated using the PyGalmesh module. For more info visit: https://github.com/nschloe/pygalmesh#volume-meshes-from-surface-meshes The pygalmesh.generate_from_array method returns a mesh containing both a cells set of tetrahedra (volume mesh) and a cells set of triangles (shell mesh). The parameter ‘meshtype’ is used to control which type of mesh is returned.

bwimage

Binary image.

voxelsizefloat

Image voxelsize.

meshtypestr

‘triangle’: Outer mesh (shell) of triangles. ‘tetra’: Volume mesh of tetrahedra. ‘both’: Both shell and volume cells sets.

max_facet_distancefloat

CGAL parameter.

max_cell_circumradiusfloat

CGAL parameter.

meshmeshio

Mesh data.

ciclope.core.tetraFE.check_cgal_params(max_facet_distance, max_cell_circumradius, voxelsize)

Check CGAL mesher parameters. # https://github.com/nschloe/pygalmesh#volume-meshes-from-surface-meshes

max_facet_distance max_cell_circumradius

voxelsizefloat

Image voxel size.

max_facet_distance : float max_cell_circumradius : float

ciclope.core.tetraFE.mesh2tetrafe(meshdata, templatefile, fileout, keywords=['NSET', 'ELSET'], float_fmt='.6e', verbose=False)

Generate ABAQUS tetrahedra Finite Element (FE) input file from 3D mesh. The output can be solved using ABAQUS or CalculiX.

Boundary conditions (BCs), simulation steps and associated output requests are defined in a separate template file. See the templates contained in the folder “input_templates” for examples.

meshdatameshio

Mesh data.

templatefilestr

Analysis template file.

fileoutstr

Output .INP file.

keywordsstr

SUPPORTED ABAQUS KEYWORDS:

For a list of all Abaqus keywords and their description visit: https://abaqus-docs.mit.edu/2017/English/SIMACAECAERefMap/simacae-c-gen-kwbrowser.htm#simacae-c-gen-kwbrowser__simacae-gen-xsl-U

  • ‘NSET’:

    Create boundary node sets. (Default = ON) If ‘NSET’ is specified, the following node sets are created:

    • NODES_X0: Nodes on WEST (X-) surface of 3D model.

    • NODES_X1: Nodes on EAST (X+) surface of 3D model.

    • NODES_Y0: Nodes on SOUTH (Y-) surface of 3D model.

    • NODES_Y1: Nodes on NORTH (Y+) surface of 3D model.

    • NODES_Z0: Nodes on BOTTOM (Z-) surface of 3D model.

    • NODES_Z1: Nodes on TOP (Z+) surface of 3D model.

    • NODES_X0Y0Z0: 2 nodes on (0,0,0) model corner.

    • NODES_X0Y0Z1: 2 nodes on (0,0,1) model corner.

    These node sets are available for boundary conditions definition.

  • ‘ELSET’:

    Create boundary element sets. (Default = ON) If ‘ELSET’ is specified, the following element sets are created:

    • ELEMS_X0: Elements of WEST (X-) surface of 3D model.

    • ELEMS_X1: Elements of EAST (X+) surface of 3D model.

    • ELEMS_Y0: Elements of SOUTH (Y-) surface of 3D model.

    • ELEMS_Y1: Elements of NORTH (Y+) surface of 3D model.

    • ELEMS_Z0: Elements of BOTTOM (Z-) surface of 3D model.

    • ELEMS_Z1: Elements of TOP (Z+) surface of 3D model.

float_fmtfloat

Precision for Abaqus input file writing.

verbosebool

Verbose output.

ciclope.core.tetraFE.shell_mesh(bwimage, method='pymcubes', voxelsize=[1.0, 1.0, 1.0], max_facet_distance=0.0, max_cell_circumradius=0.0)

Generate outer shell mesh of triangles from binary volume data. The mesh is generated using the PyMCubes module and the smooth function contained in it: https://github.com/pmneila/PyMCubes

Alternatively, the marching cube algorithm from the scikit-image python module can be used: https://scikit-image.org/docs/dev/api/skimage.measure.html?highlight=marching#skimage.measure.marching_cubes

bwimage

Binary image.

methodstr

‘pymcubes’: PyMCubes module. ‘marching_cubes’: scikit-image’s marching cube algorithm. ‘pygalmesh’: pygalmesh module (CGAL).

voxelsizefloat

Image voxelsize.

max_facet_distancefloat

CGAL parameter.

max_cell_circumradiusfloat

CGAL parameter.

vertices

Mesh vertices.

triangles

Mesh triangles.

shellmeshmeshio

Mesh data.

ciclope utilities

Data pre-processing

MicroCT reconstruction utilities

MicroCT image processing utilities.

ciclope.utils.recon_utils.bbox(bw, pad=0, dsize=None, verbose=None)

Bounding BOX limits of input binary image.

bwbool

Binary image.

padint

Add padding of given number of pixels to the BBOX limits.

dsizeint

perform image close with disk structuring element of radius ‘dsize’ before calculating the BBOX.

verbose

Activate verbose graphical output

bbox_origin: int

Origin [row col (slice)] of the BBOX inscribing True values in input image bw.

bbox_size: int

BBOX size [s_row s_col (s_slice)].

ciclope.utils.recon_utils.crop(data_3D, crop_origin, crop_size)

Crop 3D image given crop origin and size.

data_3D

Input data.

crop_origin[int, int, int]

Crop origin [Z,Y,X].

crop_size[int, int, int]

Crop size [Z,Y,X].

output

Cropped data.

ciclope.utils.recon_utils.plot_midplanes(data_3D, slice_x=- 1, slice_y=- 1, slice_z=- 1)

Plot orthogonal cross-sections through 3D dataset.

data_3D

Input 3D image data.

slice_xint

X-slice number.

slice_yint

Y-slice number.

slice_zint

Z-slice number.

ciclope.utils.recon_utils.plot_projections(data_3D, projection='max')

Plot orthogonal projections of 3D dataset.

data_3D

Input 3D image data.

projectionstr

Projection method. Available choices are ‘max’, ‘min’.

ciclope.utils.recon_utils.read_tiff_stack(filename, range=None, zfill=4)

Read stack of tiff files. Searches all files in parent folder and opens them as a stack of images.

filename

One of the stack images.

range[int, int]

Control load slices range.

zfillint

Number of leading zeros in file names.

  • check that folder contains only .TIFF files; skip the rest

ciclope.utils.recon_utils.to01(data_3D)

Normalize data to 0-1 range.

data_3D

Input data.

data_3Dfloat32

Normalized data.

ciclope.utils.recon_utils.touint8(data_3D, range=None, quantiles=None, numexpr=True)

Normalize and convert data to uint8.

data_3D

Input data.

range[float, float]

Control range for data normalization.

quantiles[float, float]

Define range for data normalization through input data quantiles. If range is given this input is ignored.

numexprbool

Use fast numerical expression evaluator for NumPy (memory expensive).

outputuint8

Normalized data.

ciclope.utils.recon_utils.writemidplanes(data_3D, fileout, slice_x=- 1, slice_y=- 1, slice_z=- 1)

Plot orthogonal mid-planes through 3D dataset and save them as images. Uses pypng for writing .PNG files.

data

Input 3D image data.

fileoutstr

Output .PNG image file name.

slice_xint

X-slice number.

slice_yint

Y-slice number.

slice_zint

Z-slice number.

Post-processing of FE results