pumapy.materialproperties

pumapy.materialproperties.conductivity

pumapy.materialproperties.conductivity.compute_electrical_conductivity(workspace, cond_map, direction, side_bc='p', prescribed_bc=None, tolerance=0.0001, maxiter=10000, solver_type='bicgstab', display_iter=True, print_matrices=(0, 0, 0, 0, 0))[source]

Compute the electrical conductivity

Parameters
  • workspace (pumapy.Workspace) – domain

  • cond_map (IsotropicConductivityMap or AnisotropicConductivityMap) – local constituents electrical conductivities

  • direction (string) – direction for solve (‘x’,’y’, or ‘z’)

  • side_bc (string) – side boundary conditions can be symmetric (‘s’), periodic (‘p’) or dirichlet (‘d’)

  • prescribed_bc (pumapy.ConductivityBC) – 3D array holding dirichlet BC

  • tolerance (float) – tolerance for iterative solver

  • maxiter (int) – maximum Iterations for solver

  • solver_type (string) – solver type, options: ‘bicgstab’, ‘cg’, ‘gmres’, ‘direct’

  • display_iter (bool) – display iterations and residual

  • print_matrices ((bool, bool, bool, bool, bool)) – corresponding to E, A, b, T, q decimal places. If 0, they are not printed

Returns

electrical conductivity, potential field, flux

Return type

((float, float, float), numpy.ndarray, numpy.ndarray)

Example

>>> import pumapy as puma
>>> ws_fiberform = puma.import_3Dtiff(puma.path_to_example_file("200_fiberform.tif"), 1.3e-6)
>>> ws_fiberform.matrix = ws_fiberform.matrix[50:150, 50:150, 50:150]
>>> cond_map = puma.IsotropicConductivityMap()
>>> cond_map.add_material((0, 89), 0.0257)
>>> cond_map.add_material((90, 255), 12)
>>> k_eff_x, P_x, q_x = puma.compute_electrical_conductivity(ws_fiberform, cond_map, 'x', 's', tolerance=1e-2, solver_type='cg')
>>> print("Effective electrical conductivity tensor:")
>>> print(k_eff_x)
pumapy.materialproperties.conductivity.compute_thermal_conductivity(workspace, cond_map, direction, side_bc='s', prescribed_bc=None, tolerance=0.0001, maxiter=10000, solver_type='bicgstab', display_iter=True, print_matrices=(0, 0, 0, 0, 0))[source]

Compute the thermal conductivity

Parameters
  • workspace (pumapy.Workspace) – domain

  • cond_map (IsotropicConductivityMap or AnisotropicConductivityMap) – local constituents themal conductivities

  • direction (string) – direction for solve (‘x’,’y’, or ‘z’)

  • side_bc (string) – side boundary conditions can be symmetric (‘s’), periodic (‘p’) or dirichlet (‘d’)

  • prescribed_bc (pumapy.ConductivityBC) – 3D array holding dirichlet BC.

  • tolerance (float) – tolerance for iterative solver

  • maxiter (int) – maximum Iterations for solver

  • solver_type (string) – solver type, options: ‘bicgstab’, ‘cg’, ‘gmres’, ‘direct’

  • display_iter (bool) – display iterations and residual

  • print_matrices ((bool, bool, bool, bool, bool)) – corresponding to b, E, A, T, q decimal places. If 0, they are not printed

Returns

thermal conductivity, temperature field, flux

Return type

((float, float, float), numpy.ndarray, numpy.ndarray)

Example

>>> import pumapy as puma
>>> ws_fiberform = puma.import_3Dtiff(puma.path_to_example_file("200_fiberform.tif"), 1.3e-6)
>>> ws_fiberform.matrix = ws_fiberform.matrix[50:150, 50:150, 50:150]
>>> cond_map = puma.IsotropicConductivityMap()
>>> cond_map.add_material((0, 89), 0.0257)
>>> cond_map.add_material((90, 255), 12)
>>> k_eff_x, T_x, q_x = puma.compute_thermal_conductivity(ws_fiberform, cond_map, 'x', 's', tolerance=1e-2, solver_type='cg')
>>> print("Effective thermal conductivity tensor:")
>>> print(k_eff_x)

pumapy.materialproperties.elasticity

pumapy.materialproperties.elasticity.compute_elasticity(workspace, elast_map, direction, side_bc='p', prescribed_bc=None, tolerance=0.0001, maxiter=10000, solver_type='bicgstab', display_iter=True, print_matrices=(0, 0, 0, 0, 0))[source]

Compute the thermal conductivity (N.B. 0 material ID in workspace refers to air unless otherwise specified)

Parameters
  • workspace (pumapy.Workspace) – domain

  • elast_map (pumapy.ElasticityMap) – local elasticity of the constituents

  • direction (string) – direction for solve (‘x’,’y’, ‘z’, ‘yz’, ‘xz’, ‘xy’)

  • side_bc (string) – side boundary conditions can be symmetric (‘s’), periodic (‘p’), dirichlet (‘d’) or free (‘f’)

  • prescribed_bc (pumapy.ElasticityBC) – 3D array holding dirichlet BC

  • tolerance – tolerance for iterative solver

  • maxiter (int) – maximum Iterations for solver

  • solver_type (string) – solver type, options: ‘bicgstab’, ‘cg’, ‘gmres’, ‘direct’

  • display_iter (bool) – display iterations and residual

  • print_matrices ((int, int, int, int, int)) – corresponding to b, E, A, u, s decimal places. If 0, they are not printed

Type

tolerance: float

Returns

elasticity, displacement field, direct stresses (sigma xx, yy, zz), shear stresses (tau yz, xz, xy)

Return type

((float, float, float, float, float, float), numpy.ndarray, numpy.ndarray, numpy.ndarray)

Example

>>> import pumapy as puma
>>> ws = puma.Workspace.from_shape_value((20, 20, 20), 1)
>>> ws[int(ws.matrix.shape[0] / 2):] = 2
>>> elast_map = puma.ElasticityMap()
>>> elast_map.add_isotropic_material((1, 1), 200, 0.3)
>>> elast_map.add_isotropic_material((2, 2), 400, 0.1)
>>> C, u, s, t = puma.compute_elasticity(ws, elast_map, direction='x', side_bc='f', solver_type="direct")
>>> print(C)
pumapy.materialproperties.elasticity.compute_stress_analysis(workspace, elast_map, prescribed_bc, side_bc='p', tolerance=0.0001, maxiter=10000, solver_type='bicgstab', display_iter=True, print_matrices=(0, 0, 0, 0, 0))[source]

Compute the thermal conductivity (N.B. 0 material ID in workspace refers to air unless otherwise specified)

Parameters
  • workspace (pumapy.Workspace) – domain

  • elast_map (pumapy.ElasticityMap) – local elasticity of the constituents

  • prescribed_bc (pumapy.ElasticityBC) – 3D array holding dirichlet BC

  • side_bc (string) – side boundary conditions can be symmetric (‘s’), periodic (‘p’), dirichlet (‘d’) or free (‘f’)

  • tolerance (float) – tolerance for iterative solver

  • maxiter (int) – maximum Iterations for solver

  • solver_type (string) – solver type, options: ‘bicgstab’, ‘cg’, ‘gmres’, ‘direct’

  • display_iter (bool) – display iterations and residual

  • print_matrices ((int, int, int, int, int)) – corresponding to b, E, A, u, s decimal places. If 0, they are not printed

Returns

displacement field, direct stresses (sigma xx, yy, zz), shear stresses (tau yz, xz, xy)

Return type

(numpy.ndarray, numpy.ndarray, numpy.ndarray)

Example

>>> import pumapy as puma
>>> ws = puma.Workspace.from_shape_value((20, 20, 20), 1)
>>> ws[ws.matrix.shape[0]//2:] = 2
>>> elast_map = puma.ElasticityMap()
>>> elast_map.add_isotropic_material((1, 1), 200, 0.3)
>>> elast_map.add_isotropic_material((2, 2), 400, 0.1)
>>> bc = puma.ElasticityBC.from_workspace(ws)
>>> bc[0] = 0  # hold x -ve face
>>> bc[-1, :, :, 0] = 1   # displace x +ve face by 1 in x direction
>>> bc[-1, :, :, 1:] = 0  # hold x +ve face in y and z directions
>>> u, s, t = puma.compute_stress_analysis(ws, elast_map, bc, side_bc='f', solver_type="direct")
>>> puma.render_volume(u[:, :, :, 1], cmap='jet')  # displacement magnitude in y direction

pumapy.materialproperties.mean_intercept_length

pumapy.materialproperties.mean_intercept_length.compute_mean_intercept_length(workspace, void_cutoff)[source]

Computation of the mean intercept length

Parameters
  • workspace (pumapy.Workspace) – domain

  • void_cutoff ((int, int)) – specify the void or gaseous phase of the domain

Returns

mean intercept length in x,y,z

Return type

(float, float, float)

Example

>>> import pumapy as puma
>>> ws = puma.import_3Dtiff(puma.path_to_example_file("200_fiberform.tif"), 1.3e-6)
>>> mil = puma.compute_mean_intercept_length(ws.matrix, (0, 89))
>>> print(mil)

pumapy.materialproperties.orientation

class pumapy.materialproperties.orientation.OrientationST(ws, sigma, rho, cutoff, edt)[source]

Bases: object

compute()[source]
error_check()[source]
pumapy.materialproperties.orientation.compute_angular_differences(matrix, orientation1, orientation2, cutoff)[source]

Compute angular difference between two orientation ndarrays

Parameters
  • matrix (np.ndarray) – domain matrix

  • orientation1 (np.ndarray) – orientation as (x, y, z, 3)

  • orientation2 (np.ndarray) – orientation as (x, y, z, 3)

  • cutoff ((int, int)) – to binarize domain

Returns

angle_errors in degrees, mean, std

Return type

(np.ndarray, float, float)

pumapy.materialproperties.orientation.compute_orientation_st(ws, sigma, rho, cutoff, edt=False)[source]

Compute orientation of the material by the structure tensor algorithm

Parameters
  • ws (pumapy.Workspace) – domain

  • sigma (float) – kernel size parameter for Gaussian derivatives (should be smaller than rho)

  • rho (float) – kernel size parameter for Gaussian filter (should be bigger than sigma)

  • cutoff ((int, int)) – which grayscales to consider

  • edt (bool) – indicating if we need to apply Euclidean Distance Transform before computing ST

Returns

True if successful, False otherwise.

Return type

bool

Example

>>> import pumapy as puma
>>> import pyvista as pv
>>> ws = puma.import_3Dtiff(puma.path_to_example_file("100_fiberform.tif"), 1.3e-6) # import example file
>>> puma.compute_orientation_st(ws, sigma=1.4, rho=0.7, cutoff=(90, 255)) # compute orientation
>>> p = pv.Plotter(shape=(1, 2))
>>> p.subplot(0, 0)
>>> p.add_text("Microstructure")
>>> puma.render_contour(ws, (90, 255), add_to_plot=p, plot_directly=False) # visualize the workspace
>>> p.subplot(0, 1)
>>> p.add_text("Detected fiber orientation")
>>> puma.render_orientation(ws, add_to_plot=p, plot_directly=False)
>>> p.show()

pumapy.materialproperties.permeability

pumapy.materialproperties.permeability.compute_permeability(workspace, solid_cutoff, tol=1e-08, maxiter=10000, solver_type='minres', display_iter=True)[source]

Compute the permeability using first order Q1-Q1 Finite Element solver and periodic BC on the sides

Parameters
  • workspace (pumapy.Workspace) – domain

  • solid_cutoff ((int, int)) – specify the solid phase

  • tol (float) – tolerance for iterative solver

  • maxiter (int) – maximum Iterations for solver

  • solver_type (string) – solver type, options: ‘minres’ (default), ‘direct’, ‘cg’, ‘bicgstab’

Returns

permeability, velocity field

Return type

(numpy.ndarray, numpy.ndarray)

Example

>>> import pumapy as puma
>>> ws = puma.generate_random_fibers(shape=(100, 100, 100), radius=3, porosity=0.7, phi=90, theta=90, length=200)
>>> keff, u_x, p_x, u_y, p_y, u_z, p_z = puma.compute_permeability(ws, (1, ws.max()))
>>> puma.render_orientation(u_x, scale_factor=5e11, solid_color=None)

pumapy.materialproperties.radiation

class pumapy.materialproperties.radiation.Radiation(workspace, cutoff, sources_number, degree_accuracy, void_phase, boundary_behavior, bin_density, rayexport_filepathname, export_plot)[source]

Bases: object

compute()[source]
error_check()[source]
generate_sources()[source]
log_input()[source]
log_output()[source]
pumapy.materialproperties.radiation.compute_extinction_coefficients(ws, rays_distances, sources_number, degree_accuracy, bin_density=10000, export_pathname=None)[source]

Compute the extinction coefficient based on the ray casting radiation simulation (this is normally a step inside the compute_radiation function)

Parameters
  • ws (pumapy.Workspace) – domain

  • rays_distances (np.ndarray) – rays distances, as output by compute_radiation function

  • sources_number (int) – number of light sources spread randomly in the void space (i.e. 0)

  • degree_accuracy (int) – angle difference between rays emitted in degrees (has to be an exact divider of 180°)

  • bin_density (int) – number of bins used to create histogram of ray distances

  • export_pathname (str) – path to save curve plot of ray distance distribution

Returns

extinction coefficient (beta), its standard deviation

Return type

(float, float)

pumapy.materialproperties.radiation.compute_radiation(workspace, solid_cutoff, sources_number, degree_accuracy, void_phase=0, boundary_behavior=1, bin_density=10000, exportparticles_filepathname='', export_pathname=None)[source]

Compute the radiative thermal conductivity through ray tracing (N.B. 0 material ID in workspace refers to gas phases unless otherwise specified)

Parameters
  • workspace (pumapy.Workspace) – domain

  • solid_cutoff ((int, int)) – specify the solid phase

  • sources_number (int) – number of light sources spread randomly in the void space (i.e. 0)

  • degree_accuracy (int) – angle difference between rays emitted in degrees (has to be an exact divider of 180°)

  • void_phase (int) – ID of the void phase, defaulted as 0

  • boundary_behavior (int) – how to treat particles exiting the domain: 0=kill, 1=periodic (default)

  • bin_density (int) – number of bins used to create histogram of ray distances

  • exportparticles_filepathname (string) – path and name of the particle files exported as VTK points

  • export_pathname (str) – path to save curve plot of ray distance distribution

Returns

extinction coefficient (beta), its standard deviation and of ray distances

Return type

(float, float, np.ndarray)

Example

>>> import pumapy as puma
>>> ws_fiberform = puma.import_3Dtiff(puma.path_to_example_file("200_fiberform.tif"), 1.3e-6)
>>> beta, beta_std, rays_distances = puma.compute_radiation(ws_fiberform, (104, 255), 100, 15)

pumapy.materialproperties.surfacearea

class pumapy.materialproperties.surfacearea.SurfaceArea(workspace, cutoff, flag_gaussian)[source]

Bases: object

compute()[source]
log_input()[source]
log_output()[source]
pumapy.materialproperties.surfacearea.compute_surface_area(workspace, cutoff, flag_gaussian=False)[source]

Computation of the surface area based on isosurface

Parameters
  • workspace (pumapy.Workspace) – domain

  • cutoff ((int, int)) – specify the solid phase

  • flag_gaussian (bool) – apply Gaussian filter before generating surface

Returns

area, specific_area

Return type

(float, float)

Example

>>> import pumapy as puma
>>> ws = puma.import_3Dtiff(puma.path_to_example_file("200_fiberform.tif"), 1.3e-6) # import workspace
>>> area, specific_area = puma.compute_surface_area(ws, cutoff=(90, 255)) # computing surface area
>>> print("Surface area: ", area)
>>> print("Specific surface area: ", specific_area)

pumapy.materialproperties.tortuosity

pumapy.materialproperties.tortuosity.compute_continuum_tortuosity(workspace, cutoff, direction, side_bc='p', prescribed_bc=None, tolerance=0.0001, maxiter=10000, solver_type='cg', display_iter=True)[source]

Compute the tortuosity modelling the local conductivity as isotropic

Parameters
  • workspace (pumapy.Workspace) – domain

  • cutoff ((int, int)) – to binarize domain

  • direction (string) – direction for solve (‘x’,’y’, or ‘z’)

  • side_bc (string) – side boundary conditions (string) can be symmetric (‘s’), periodic (‘p’) or dirichlet (‘d’)

  • prescribed_bc (pumapy.ConductivityBC) – 3D array holding dirichlet BC

  • tolerance (float) – tolerance for iterative solver

  • maxiter (int) – maximum Iterations for solver

  • solver_type (string) – solver type, options: ‘cg’, ‘bicgstab’, ‘direct’

  • display_iter (bool) – display iterations and residual

Returns

tortuosity, diffusivity, porosity, concentration field

Return type

((float, float, float), float, float, numpy.ndarray)

Example

>>> import pumapy as puma
>>> ws_fiberform = puma.import_3Dtiff(puma.path_to_example_file("200_fiberform.tif"), 1.3e-6) # import image
>>> n_eff_x, Deff_x, poro, C_x = puma.compute_continuum_tortuosity(ws_fiberform, (0,89), 'x', side_bc='s', tolerance=1e-2, solver_type='cg')
>>> print('Effective tortuosity factors:', n_eff_x)
>>> print("Porosity of the material:", poro)

pumapy.materialproperties.volumefraction

class pumapy.materialproperties.volumefraction.VolumeFraction(workspace, cutoff)[source]

Bases: object

compute()[source]
error_check()[source]
log_input()[source]
log_output()[source]
pumapy.materialproperties.volumefraction.compute_volume_fraction(workspace, cutoff)[source]

Compute the volume fraction

Parameters
  • workspace (pumapy.Workspace) – domain

  • cutoff ((int, int)) – to binarize domain

Returns

volume fraction

Return type

float

Example

>>> import pumapy as puma
>>> ws = puma.import_3Dtiff(puma.path_to_example_file("100_fiberform.tif"), 1.3e-6)) # import example file
>>> vf = puma.compute_volume_fraction(ws, cutoff=(90, 255)) # compute volume fraction
>>> print ("Volume Fraction of workspace:", vf)