pumapy.material_properties
pumapy.material_properties.anisotropic_radiation
- class pumapy.material_properties.anisotropic_radiation.AnisotropicRadiation(workspace, void_cutoff, sources_number, direction, bin_density, export_plot)[source]
Bases:
object
- pumapy.material_properties.anisotropic_radiation.compute_extinction_coefficients_anisotropic(ws, rays_distances, sources_number, 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)
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.material_properties.anisotropic_radiation.compute_radiation_anisotropic(workspace, void_cutoff, sources_number, bin_density=10000, export_pathname=None)[source]
Compute the anisotropic radiative thermal conductivity through ray tracing
- Parameters
workspace (pumapy.Workspace) – domain
void_cutoff ((int, int)) – specify the void phase
sources_number (int) – number of light sources spread randomly in the void space (i.e. 0)
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 and of ray distances
- Return type
(float, float, np.ndarray)
pumapy.material_properties.conductivity
Further explained in publication: Semeraro, F., Ferguson, J.C., Acin, M., Panerai, F. and Mansour, N.N., 2021. Anisotropic analysis of fibrous and woven materials part 2: Computation of effective conductivity. Computational Materials Science, 186, p.109956. https://www.sciencedirect.com/science/article/abs/pii/S092702562030447X
- pumapy.material_properties.conductivity.compute_electrical_conductivity(workspace, cond_map, direction, side_bc='p', prescribed_bc=None, tolerance=1e-05, maxiter=10000, solver_type='bicgstab', display_iter=True, method='fv', matrix_free=True)[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’,’z’, or ‘’ for prescribed_bc). If provided, prescribed_bc is ignored
side_bc (string) – side boundary conditions can be symmetric (‘s’) or periodic (‘p’). Only periodic available when method=’fe’
prescribed_bc (pumapy.IsotropicConductivityBC or pumapy.AnisotropicConductivityBC or None) – object holding dirichlet BC, only available for isotropic of MPFA implementations. Need to set direction=’’ in order to provide it. When prescribed_bc is provided, keff is not printed but it is computed as if direction==’x’ for testing purposes.
tolerance (float) – tolerance for iterative solver
maxiter (int) – maximum Iterations for solver
solver_type (string) – solver type, options: ‘bicgstab’ (default), ‘cg’, ‘gmres’, ‘minres’ (only for method=’fe’), ‘direct’
display_iter (bool) – display iterations and residual
method (string) – whether to use finite volume solver (‘fv’, either isotropic solver if IsotropicConductivityMap is passed, or mpfa if AnisotropicConductivityMap) or finite element Q1-Q1 EBE solver (‘fe’). For the latter method, it is recommended to use solver_type=’minres’ as lighter and faster
matrix_free (bool) – if True, then use matrix-free method if possible (only available for fv isotropic solver or for fe solver when the solver type is not ‘direct’)
- Returns
electrical conductivity, potential field, flux
- Return type
((float, float, float), numpy.ndarray, numpy.ndarray)
- pumapy.material_properties.conductivity.compute_thermal_conductivity(workspace, cond_map, direction, side_bc='p', prescribed_bc=None, tolerance=1e-05, maxiter=10000, solver_type='bicgstab', display_iter=True, method='fv', matrix_free=True)[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’,’z’, or ‘’ for prescribed_bc). If provided, prescribed_bc is ignored
side_bc (string) – side boundary conditions can be symmetric (‘s’) or periodic (‘p’). Only periodic available when method=’fe’
prescribed_bc (pumapy.IsotropicConductivityBC or pumapy.AnisotropicConductivityBC or None) – object holding dirichlet BC, only available for isotropic of MPFA implementations. Need to set direction=’’ in order to provide it. When prescribed_bc is provided, keff is not printed but it is computed as if direction==’x’ for testing purposes.
tolerance (float) – tolerance for iterative solver
maxiter (int) – maximum Iterations for solver
solver_type (string) – solver type, options: ‘bicgstab’ (default), ‘cg’, ‘gmres’, ‘minres’ (only for method=’fe’), ‘direct’
display_iter (bool) – display iterations and residual
method (string) – whether to use finite volume solver (‘fv’, either isotropic solver if IsotropicConductivityMap is passed, or mpfa if AnisotropicConductivityMap) or finite element Q1-Q1 EBE solver (‘fe’). For the latter method, it is recommended to use solver_type=’minres’ as lighter and faster
matrix_free (bool) – if True, then use matrix-free method if possible (only available for fv isotropic solver or for fe solver when the solver type is not ‘direct’)
- 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) Importing .../200_fiberform.tif ... Done >>> ws_fiberform.rescale(0.5, segmented=False) Rescaled workspace size: (100, 100, 100) >>> # Conductivity with Isotropic local phases >>> 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', matrix_free=True) Approximate memory requirement for simulation... >>> # Conductivity with Anisotropic local phases >>> puma.compute_orientation_st(ws_fiberform, cutoff=(90, 255)) First gradient computation... >>> cond_map = puma.AnisotropicConductivityMap() >>> # conductivity of the void phase to be 0.0257 (air at STP) >>> cond_map.add_isotropic_material((0, 89), 0.0257) >>> # axial fiber conductivity of 12, radial fiber conductivity of 0.7 >>> cond_map.add_material_to_orient((90, 255), 12., 0.7) >>> k_eff_z, T_z, q_z = puma.compute_thermal_conductivity(ws_fiberform, cond_map, 'z', 's') # method='fe' for finite element Approximate memory requirement for simulation... >>> # plot_conductivity_fields(ws, T_z, q_z) # to visualize fields >>> # export_conductivity_fields_vti("path/to/folder", ws, T_z, q_z) # to export fields
- pumapy.material_properties.conductivity.export_conductivity_fields_vti(filepath, workspace, T, q)[source]
Export conductivity fields, as output by the conductivity function
- Parameters
filepath (string) –
workspace (puma.Workspace or numpy.ndarray) – domain
T (numpy.ndarray) – temperature field
q (numpy.ndarray) – flux field
- pumapy.material_properties.conductivity.plot_conductivity_fields(workspace, T, q, show_cbar=True, show_edges=False, xy_view=False, rm_id=None, notebook=False)[source]
Plot the workspace colored by the temperature and flux fields, as output by the conductivity function
- Parameters
workspace (pumapy.Workspace) – domain
T (numpy.ndarray) – temperature field
q (numpy.ndarray) – flux field
show_cbar (bool) – show colorbar in each plot
show_edges (bool) – show edges in mesh
xy_view (bool) – show plot aligned with xy plane
rm_id (float or None) – remove a phase of the material from thresholded mesh
notebook (bool) – plotting interactively in a jupyter notebook (overwrites show_grid to False)
pumapy.material_properties.elasticity
- pumapy.material_properties.elasticity.compute_elasticity(workspace, elast_map, direction, side_bc='p', tolerance=1e-05, maxiter=100000, solver_type='bicgstab', display_iter=True, method='fv', matrix_free=True)[source]
Compute the effective elasticity coefficient
- 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’) or periodic (‘p’). Only periodic available when method=’fe’
tolerance – tolerance for iterative solver
maxiter (int) – maximum Iterations for solver
solver_type (string) – solver type, options: ‘bicgstab’ (default), ‘minres’ (only for method=’fe’), ‘gmres’, ‘direct’
display_iter (bool) – display iterations and residual
method (string) – whether to use finite volume solver (‘fv’, i.e. mpsa) or finite element Q1-Q1 EBE solver (‘fe’). For the latter method, it is recommended to use solver_type=’minres’ as lighter and faster
matrix_free (bool) – if True, then use matrix-free method if possible (only available for fe solver when the solver type is not ‘direct’)
- 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.experimental.ElasticityMap() >>> elast_map.add_isotropic_material((1, 1), 200, 0.3) >>> elast_map.add_isotropic_material((2, 2), 400, 0.1) >>> C = np.zeros((6, 6)) >>> C[:, 0], u, s, t = puma.experimental.compute_elasticity(ws, elast_map, direction='x', solver_type='direct') Approximate memory requirement for simulation: ... >>> C[:, 1], u, s, t = puma.experimental.compute_elasticity(ws, elast_map, direction='y', solver_type='direct') Approximate memory requirement for simulation: ... >>> C[:, 2], u, s, t = puma.experimental.compute_elasticity(ws, elast_map, direction='z', solver_type='direct') Approximate memory requirement for simulation: ... >>> C[:, 3], u, s, t = puma.experimental.compute_elasticity(ws, elast_map, direction='yz', solver_type='direct') Approximate memory requirement for simulation: ... >>> C[:, 4], u, s, t = puma.experimental.compute_elasticity(ws, elast_map, direction='xz', solver_type='direct') Approximate memory requirement for simulation: ... >>> C[:, 5], u, s, t = puma.experimental.compute_elasticity(ws, elast_map, direction='xy', solver_type='direct') Approximate memory requirement for simulation: ... >>> coeffs = puma.experimental.get_E_nu_from_elasticity(C) E1 ... >>> #puma.plot_elasticity_fields(ws, u, s, t) # to visualize fields >>> #puma.warp_elasticity_fields(ws, u, s, t, 5) # to visualize fields >>> #puma.export_elasticity_fields_vti("path/to/folder", ws, u, s, t) # to export fields >>> C1, u, s, t = puma.experimental.compute_elasticity(ws, elast_map, direction='x', solver_type='direct', method='fe') # finite element solver Approximate memory requirement for simulation: ...
- pumapy.material_properties.elasticity.compute_stress_analysis(workspace, elast_map, prescribed_bc, side_bc='p', tolerance=1e-05, maxiter=100000, solver_type='bicgstab', display_iter=True)[source]
Compute stress analysis based on specific user-input dirichlet boundary conditions
- Parameters
workspace (pumapy.Workspace) – domain
elast_map (pumapy.ElasticityMap) – local elasticity of the constituents
prescribed_bc (pumapy.ElasticityBC) – object holding the elasticity dirichlet boundary conditions
side_bc (string) – side boundary conditions can be symmetric (‘s’) or periodic (‘p’)
tolerance – tolerance for iterative solver
maxiter (int) – maximum Iterations for solver
solver_type (string) – solver type, options: ‘gmres’ (default), ‘bicgstab’, ‘direct’
display_iter (bool) – display iterations and residual
- Type
tolerance: float
- Returns
elasticity, 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, 22, 22), 1) >>> ws[ws.matrix.shape[0]//2:] = 2 >>> ws[:, [0, -1]] = 0 >>> ws[:, :, [0, -1]] = 0 >>> elast_map = puma.experimental.ElasticityMap() >>> elast_map.add_isotropic_material((0, 0), 1e-5, 0.3) >>> elast_map.add_isotropic_material((1, 1), 200, 0.3) >>> elast_map.add_isotropic_material((2, 2), 400, 0.1) >>> bc = puma.experimental.ElasticityBC(ws) >>> bc.xfaces[0, :, :, 0] = 0 >>> bc.xfaces[0, :, :, 1] = 0 >>> bc.xfaces[0, :, :, 2] = 0 >>> bc.xfaces[1, :, :, 0] = 0 >>> bc.xfaces[1, :, :, 1] = 1 >>> bc.xfaces[1, :, :, 2] = 0 >>> u, s, t = puma.experimental.compute_stress_analysis(ws, elast_map, bc) Approximate memory requirement for simulation: ... >>> # puma.warp_elasticity_fields(ws[:, 1:-1, 1:-1], u[:, 1:-1, 1:-1], s[:, 1:-1, 1:-1], t[:, 1:-1, 1:-1], 20, show_original=0., show_edges=True) # to visualize
- pumapy.material_properties.elasticity.export_elasticity_fields_vti(filepath, workspace, u, s, t)[source]
Export conductivity fields, as output by the conductivity function
- Parameters
filepath (string) –
workspace (puma.Workspace or numpy.ndarray) – domain
u (numpy.ndarray) – displacement field
s (numpy.ndarray) – direct stress field
t (numpy.ndarray) – shear stress field
- pumapy.material_properties.elasticity.get_E_nu_from_elasticity(C)[source]
Compute Young’s moduli E1, E2, E3, Shear moduli G23, G13, G12, and Poisson’s ratios nu23, nu13, nu12 for an orthotropic material from its symmetric elastic stiffness tensor
- Parameters
C (np.ndarray) – 6x6 elasticity tensor
- Returns
Young’s moduli E1, E2, E3, Shear moduli G23, G13, G12, and Poisson’s ratios nu23, nu13, nu12
- Return type
(float, float, float, float, float, float, float, float, float)
- pumapy.material_properties.elasticity.plot_elasticity_fields(workspace, u, s, t, show_cbar=True, show_edges=False, xy_view=False, rm_id=None, notebook=False)[source]
Plot the workspace according to the displacement and stress fields output by the elasticity functions
- Parameters
workspace (pumapy.Workspace) – domain
u (numpy.ndarray) – displacement field
s (numpy.ndarray) – direct stress field
t (numpy.ndarray) – shear stress field
show_cbar (bool) – show colorbar in each plot
show_edges (bool) – show edges in mesh
xy_view (bool) – show plot aligned with xy plane
rm_id (float or None) – remove a phase of the material from thresholded mesh
notebook (bool) – plotting interactively in a jupyter notebook (overwrites show_grid to False)
- pumapy.material_properties.elasticity.remove_rbms(workspace, void_id, direction)[source]
Remove Rigid Body Movements (RBMs), i.e. unconnected or “floating” voxels, in a segmented domain along a specified direction
- Parameters
workspace (pumapy.Workspace) – domain
void_id (int) – specify the void ID to discard from RBMs identification
direction (str) – Cartesian direction that has to be connected, options: ‘x’, ‘y’, ‘z’
- Returns
workspace without the possible RBMs determined by not being connected from side to side N.B. The output workspace is segmented into 0=void, 1=solid
- Return type
pumapy.Workspace
- Example
>>> import pumapy as puma >>> workspace = puma.import_3Dtiff(puma.path_to_example_file("100_fiberform.tif")) Importing ... >>> workspace.binarize(103) >>> new_ws = puma.experimental.remove_rbms(workspace, void_id=0, direction='y') >>> # puma.render_volume(workspace, (1, 1), solid_color=(255, 255, 255)) # to visualize before and after >>> # puma.render_volume(new_ws, (1, new_ws.max()), solid_color=(255, 255, 255))
- pumapy.material_properties.elasticity.warp_elasticity_fields(workspace, u, s, t, scale_factor=1, show_original=0.0, show_cbar=True, show_edges=False, xy_view=False, rm_id=None, show_axes=True, notebook=False)[source]
Warp the workspace according to the displacement field output by the elasticity functions, and colored by displacement and stress components
- Parameters
workspace (pumapy.Workspace) – domain
u (numpy.ndarray) – displacement field
s (numpy.ndarray) – direct stress field
t (numpy.ndarray) – shear stress field
scale_factor (float) – scaling factor for warp
show_original (float) – opacity of the original workspace before warp
show_cbar (bool) – show colorbar in each plot
show_edges (bool) – show edges in mesh
xy_view (bool) – show plot aligned with xy plane
rm_id (float or None) – remove a phase of the material from warped mesh (only works for 2D slice)
show_axes (float) – show the axes and side dimensions
notebook (bool) – plotting interactively in a jupyter notebook (overwrites show_grid to False)
pumapy.material_properties.mean_intercept_length
- pumapy.material_properties.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) Importing ... >>> mil = puma.compute_mean_intercept_length(ws, (0, 89))
pumapy.material_properties.orientation
Further explained in publication: Semeraro, F., Ferguson, J.C., Panerai, F., King, R.J. and Mansour, N.N., 2020. Anisotropic analysis of fibrous and woven materials part 1: Estimation of local orientation. Computational Materials Science, 178, p.109631. (https://www.sciencedirect.com/science/article/abs/pii/S0927025620301221)
- class pumapy.material_properties.orientation.OrientationST(ws, sigma, rho, cutoff, edt)[source]
Bases:
object
- pumapy.material_properties.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.material_properties.orientation.compute_orientation_st(ws, cutoff, sigma=0.7, rho=1.4, edt=False)[source]
Compute orientation of the material by the structure tensor algorithm
- Parameters
ws (pumapy.Workspace) – domain
cutoff ((int, int)) – which grayscales to consider
sigma (float) – kernel size parameter for Gaussian derivatives (should be lower than rho)
rho (float) – kernel size parameter for Gaussian filter of derivatives of grayscales (should be higher than sigma)
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 Importing ... >>> puma.compute_orientation_st(ws, cutoff=(90, 255), sigma=0.7, rho=1.4) # compute orientation First gradient computation ... >>> # p = pv.Plotter(shape=(1, 2)) # to visualize the orientation field >>> # p.subplot(0, 0) >>> # p.add_text("Microstructure") >>> # puma.render_contour(ws, (90, 255), add_to_plot=p, plot_directly=False) >>> # p.subplot(0, 1) >>> # p.add_text("Detected fiber orientation") >>> # puma.render_orientation(ws, add_to_plot=p, plot_directly=False) >>> # p.show() # to visualize it
pumapy.material_properties.permeability
The following FE numerical method and implementation are based on the following research paper:
Lopes, P.C., Vianna, R.S., Sapucaia, V.W., Semeraro, F., Leiderman, R. and Pereira, A.M., 2023. Simulation toolkit for digital material characterization of large image-based microstructures. Computational Materials Science, 219, p.112021. (https://www.sciencedirect.com/science/article/pii/S0927025623000150)
- pumapy.material_properties.permeability.compute_permeability(workspace, solid_cutoff, direction='xyz', tol=1e-08, maxiter=10000, solver_type='minres', display_iter=True, matrix_free=True, precondition=False, output_fields=True)[source]
Compute the permeability using first-order Q1-Q1 Finite Element solver and periodic BC on the sides.
Note: the iterative solvers have been observed to struggle converging for some cases. This is often connected to having the preconditioner option as True.
- Parameters
workspace (pumapy.Workspace) – domain
solid_cutoff ((int, int)) – specify the solid phase
direction (str) – direction for solve (‘xyz’,’x’,’y’,’z’). Note that if solver_type=”direct”, then the direction is considered as “xyz” automatically because there is no need to invert the A sparse matrix multiple times
tol (float) – tolerance for iterative solver
maxiter (int) – maximum Iterations for solver
solver_type (string) – solver type, options: ‘minres’ (default, only works if precondition=False), ‘cg’, ‘direct’
display_iter (bool) – display iteration in iterative solver
matrix_free (bool) – solve system using matrix-free method (True, recommended) or building the sparse A matrix (False)
precondition (bool) – solve system with Jacobi preconditioner (True) or without (False, default because it reduces memory, but more iterations. The default minres iterative solver does not support this kind of preconditioner)
output_fields (bool) – export velocity and pressure fields (True, default) or not (False, lower memory)
- Returns
effective permeability (3x3 matrix) and, if output_fields=True, the normalized velocity field for the corresponding direction (arranged as tuple of numpy.ndarrays, i.e. (ux, uy, uz). If output_fields=False, then (None, None, None) is output
- Return type
numpy.ndarray, (numpy.ndarray, numpy.ndarray, numpy.ndarray)
- Example
>>> import pumapy as puma >>> import pyvista as pv >>> ws = puma.generate_random_fibers_transverseisotropic(shape=(50, 50, 50), radius=2, porosity=0.7, direction='x', variation=15, length=200, allow_intersect=True, segmented=True) Fibers created... >>> keff, (ux, _, _) = puma.compute_permeability(ws, (1, ws.max()), direction='x', tol=1e-7) Approximate memory requirement for simulation: ... >>> # p = pv.Plotter() # to visualize it >>> # puma.render_orientation(ux, add_to_plot=p, scale_factor=2e12, plot_directly=False) >>> # ws.voxel_length = 1 # the voxel_length is converted to 1 for plotting the workspace together with the velocity >>> # puma.render_volume(ws, cutoff=(1, ws.max()), add_to_plot=p, plot_directly=False, cmap='jet') >>> # p.show()
pumapy.material_properties.radiation
- class pumapy.material_properties.radiation.Radiation(workspace, void_cutoff, sources_number, particles_number, boundary_behavior, bin_density, export_plot)[source]
Bases:
object
- pumapy.material_properties.radiation.compute_extinction_coefficients(ws, rays_distances, sources_number, particles_number, 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.material_properties.radiation.compute_radiation(workspace, void_cutoff, sources_number, particles_number, boundary_behavior=1, bin_density=10000, 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
void_cutoff ((int, int)) – specify the void phase
sources_number (int) – number of light sources spread randomly in the void space (i.e. 0)
particles_number (int) – number of particles emitted at each source point
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
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"), 0.65e-6) Importing ... >>> beta, beta_std, rays_distances = puma.experimental.compute_radiation(ws_fiberform, (0, 89), 100, 500) Number of particles in Ray Tracing simulation...
pumapy.material_properties.surface_area
- class pumapy.material_properties.surface_area.SurfaceArea(workspace, cutoff, flag_gaussian)[source]
Bases:
object
- pumapy.material_properties.surface_area.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 Importing ... >>> area, specific_area = puma.compute_surface_area(ws, cutoff=(90, 255)) # computing surface area
pumapy.material_properties.tortuosity
- pumapy.material_properties.tortuosity.compute_continuum_tortuosity(workspace, cutoff, direction, side_bc='p', prescribed_bc=None, tolerance=0.0001, maxiter=10000, solver_type='cg', display_iter=True, matrix_free=True)[source]
Compute the tortuosity modelling the local conductivity as isotropic
- Parameters
workspace (pumapy.Workspace) – domain
cutoff ((int, int)) – cutoff to binarize domain specifying the void phase
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.IsotropicConductivityBC or None) – 3D array holding dirichlet BC
tolerance (float) – tolerance for iterative solver
maxiter (int) – maximum Iterations for solver
solver_type (string) – solver type, options: ‘cg’ (default), ‘bicgstab’, ‘direct’
display_iter (bool) – display iterations and residual
matrix_free (bool) – if True, then use matrix-free method
- 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) Importing ... >>> n_eff_x, Deff_x, poro, C_x = puma.compute_continuum_tortuosity(ws_fiberform, (0, 89), 'x', side_bc='s', tolerance=1e-4) Approximate memory requirement for simulation...
pumapy.material_properties.volume_fraction
- class pumapy.material_properties.volume_fraction.VolumeFraction(workspace, cutoff, display)[source]
Bases:
object
- pumapy.material_properties.volume_fraction.compute_volume_fraction(workspace, cutoff, display=True)[source]
Compute the volume fraction
- Parameters
workspace (pumapy.Workspace or np.ndarray) – domain
cutoff ((int, int)) – to binarize domain
display (bool) – print result
- 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 Importing ... >>> vf = puma.compute_volume_fraction(ws, cutoff=(90, 255)) # compute volume fraction Volume Fraction for cutoff ...