Modules#
This section contains API documentation for the most important modules and interfaces that are necessary when using the library.
Thermomechanical properties through homogenization#
- fiberpy.mechanics.A2Eij(A)#
Calculate the orthotropic moduli from stiffness tensors written using the \((\phi_2, \phi)\) bases
- Parameters:
A (array_like of shape (6, 6)) – Stiffness tensor
- Returns:
\((E_1, E_2, E_3, \mu_{12}, \mu_{23}, \mu_{13}, \nu_{12}, \nu_{23}, \nu_{31})\)
- fiberpy.mechanics.AIsotropic(E, nu)#
Isotropic stiffness tensor given in the \((\phi, \phi)\) bases
- Parameters:
E (float) – Young’s modulus
nu (float) – Poisson coefficient
- Returns:
Stiffness tensor
- Return type:
array of shape (6, 6)
- class fiberpy.mechanics.FiberComposite(rve_data)#
Fiber-reinforced composite
- Parameters:
rve_data (dict) – Dictionary defining the microstructure
- ABar(a, model='TandonWeng', closure='orthotropic', recompute_UD=False)#
Homogenized Stiffness tensor in the principal frame
- Parameters:
a (array_like of shape (3,)) – Principal values of the 2nd fiber orientation tensor,
a[0] >= a[1] >= a[2]model (str) – Micromechanical model for the unidirectional RVE (“TandonWeng”, “MoriTanaka” or “Balanced”)
closure (str) – 4th-order fiber orientation closure model
A4_*, seefiberpy.closurerecompute_UD (bool) – Whether force recomputing elastic properties of the unidirectional RVE
- Returns:
Effective Stiffness tensor using the \((\phi_2, \phi)\) bases
- Return type:
array of shape (6, 6)
References
Advani, S. G. & Tucker III, C. L. The use of tensors to describe and predict fiber orientation in short fiber composites. Journal of Rheology, SOR, 1987, 31, 751-784
- Balanced()#
Unidirectional stiffness tensor using the balanced formulation
- Returns:
Stiffness tensor using the \((\phi, \phi)\) bases
- Return type:
array of shape (6, 6)
- Eshelby()#
Eshelby’s tensor
References
Tandon, G. P. & Weng, G. J. The effect of aspect ratio of inclusions on the elastic properties of unidirectionally aligned composites. Polymer Composites, Wiley Online Library, 1984, 5, 327-333
- MoriTanaka()#
Unidirectional stiffness tensor using the Mori-Tanaka formulation
- Returns:
Stiffness tensor using the \((\phi, \phi)\) bases
- Return type:
array of shape (6, 6)
- TandonWeng()#
Unidirectional stiffness tensor using Tandon-Weng’s equations
- Returns:
Stiffness tensor using the \((\phi, \phi)\) bases
- Return type:
array of shape (6, 6)
References
Tandon, G. P. & Weng, G. J. The effect of aspect ratio of inclusions on the elastic properties of unidirectionally aligned composites. Polymer Composites, Wiley Online Library, 1984, 5, 327-333
- alphaBar(ABar)#
Homogenized thermal expansion coefficients in the principal frame
- Parameters:
ABar (array_like of shape (6, 6)) – Stiffness tensor
- Returns:
Effective thermal dilatation coefficient matrix
- Return type:
array of shape (3, 3)
References
Rosen, B. W. & Hashin, Z. Effective thermal expansion coefficients and specific heats of composite materials. International Journal of Engineering Science, Elsevier BV, 1970, 8, 157-173
- get(variables)#
Retrieve the RVE variables
- read_rve_data()#
Parse a RVE data defining the microstructure
- vBar(x0, x1)#
Volume average
- fiberpy.mechanics.bulk_modulus(E, nu)#
Bulk modulus from \((E, \nu)\)
- fiberpy.mechanics.lmbda_mu(E, nu)#
Convert \((E, \nu)\) to \((\lambda, \mu)\)
FEA interface for integrative simulations#
Fiber orientation models#
- class fiberpy.orientation.Icosphere(radius=1, n_refinement=5)#
Discretization of a sphere using subdivision of a icosahedron
- Parameters:
radius (float) – Radius of the sphere
n_refinement (int) – Number of subdivision of the icosahedron
- centroid()#
Centroid coordinates of cells
- equal_earth_projection()#
Project node points from the sphere to 2-d plane using the Equal Earch projection
- integrate(fun, vectorized=False)#
Integrate a function defined on the icosphere
We assume that function is piecewisely constant in each cell
- Parameters:
f (callable) – Function depending on x
vectorized (bool) – Whether the given function is vectorized
- fiberpy.orientation.apply_Ax(v, A, x)#
Compute the contraction of a 4th-order tensor (with minor symmetry) given in the principal frame on a symmetric 2nd-order tensor given in the global frame
- Parameters:
v (ndarray of shape (3, 3)) – Principal directions along its columns
A (ndarray of shape (6, 6)) – 4th-order tensor (using
fiberpy.tensor.Mat4())x (ndarray of shape (3, 3)) – 2nd-order tensor
- Returns:
Result array
- Return type:
ndarray of shape (3, 3)
- fiberpy.orientation.distribution_function(a, n_refinement=5, return_mesh=False)#
Reconstruct the orientation distribution function (ODF) from the 2nd-order orientation tensor
The ODF is assumed to follow the Bingham distribution. Its probability density function is proportional to
\[\exp(\mathbf{x}^\mathsf{T}\mathbf{v}\mathbf{Z}\mathbf{v}^\mathsf{T}\mathbf{x})\]where \(\mathbf{v}\) are the principal directions and \(\mathbf{Z}\) is a diagonal matrix of trace 1
- Parameters:
a (array_like of shape (3, 3) or (3,)) – Fiber orientation tensor (or its principal values in decreasing order)
n_refinement (int) – Number of subdivision of the icosahedron
return_mesh (bool) – Also return the icosphere mesh containing the values on cells
- Returns:
Orientation distribution function
odf(x)defined for normal vectors on the unit sphere- Return type:
callable
- fiberpy.orientation.fiber_orientation(a0, t, L, ci, ar, kappa=1, D3=None, closure='orthotropic', method='RK45', debug=False, **kwargs)#
Compute fiber orientation tensor evolution using the Folgar-Tucker model as its variants
- Parameters:
a0 (ndarray of shape (3, 3)) – Initial fiber orientation tensor
t (ndarray of shape (1, )) – Time instants
L (ndarray of shape (3, 3)) – Velocity gradient
ci (float) – Interaction coefficient
ar (float) – Aspect ratio
kappa (float) – Reduction coefficient when using the RSC model,
0 < kappa <= 1D3 (ndarray of shape (3, )) – Coefficients \((D_1,D_2,D_3)\) when using the MRD model
closure (str) – 4th-order fiber orientation closure model
A4_*, seefiberpy.closuremethod (str) – Numerical method to integrate the IVP, see
scipy.integrate.solve_ivp(), orjuliadebug (bool) – Return instead the
solobject anddadt
- fiberpy.orientation.project_aij(a)#
Project fiber orientation tensors to the physical space
- Parameters:
a (array_like of shape (3,)) – Principal values of fiber orientation tensor in decreasing order
- fiberpy.orientation.shear_steady_state(ci, ar, closure='orthotropic', a0_isotropic='3d')#
Fiber orientation at steady state under simple shear
- Parameters:
ci (float) – Interaction coefficient
ar (float) – Aspect ratio
closure (str) – 4th-order fiber orientation closure model
A4_*, seefiberpy.closurea0_isotropic (str) – If
3d, using 3-d initial isotropic orientation; if2d, using planar initial isotropic orientation
Closure models for 4th-order fiber orientation tensor#
- fiberpy.closure.A4_exact(a)#
Compute the exact closure in the principal frame
- Parameters:
a (array_like of shape (3,)) – Fiber orientation principal values,
a[0] >= a[1] >= a[2]- Returns:
4th-order orientation tensor written using the \((\phi,\phi)\) bases
- Return type:
array of shape (6, 6)
- fiberpy.closure.A4_hybrid(a)#
Compute the hybrid closure
- Parameters:
a (array_like of shape (3, 3) or (3,)) – Fiber orientation tensor (or its principal values in decreasing order)
- Returns:
4th-order orientation tensor written using the \((\phi,\phi)\) bases
- Return type:
array of shape (6, 6)
References
Advani, S. G. & Tucker III, C. L. Closure approximations for three-dimensional structure tensors. Journal of Rheology, SOR, 1990, 34, 367-386
- fiberpy.closure.A4_invariants(a)#
Compute the IBOF closure
- Parameters:
a (array_like of shape (3, 3) or (3,)) – Fiber orientation tensor (or its principal values in decreasing order)
- Returns:
4th-order orientation tensor written using the \((\phi,\phi)\) bases
- Return type:
array of shape (6, 6)
- fiberpy.closure.A4_linear(a)#
Compute the linear closure
- Parameters:
a (array_like of shape (3, 3) or (3,)) – Fiber orientation tensor (or its principal values in decreasing order)
- Returns:
4th-order orientation tensor written using the \((\phi,\phi)\) bases
- Return type:
array of shape (6, 6)
References
Advani, S. G. & Tucker III, C. L. The use of tensors to describe and predict fiber orientation in short fiber composites. Journal of Rheology, SOR, 1987, 31, 751-784
- fiberpy.closure.A4_orthotropic(a)#
Compute the orthotropic closure in the principal frame
- Parameters:
a (array_like of shape (3,)) – Fiber orientation principal values,
a[0] >= a[1] >= a[2]- Returns:
4th-order orientation tensor using the \((\phi,\phi)\) bases
- Return type:
array of shape (6, 6)
References
VerWeyst, B. E. Numerical predictions of flow-induced fiber orientation in three-dimensional geometries. University of Illinois at Urbana-Champaign, 1998
- fiberpy.closure.A4_quadratic(a)#
Compute the quadratic closure
- Parameters:
a (array_like of shape (3, 3) or (3,)) – Fiber orientation tensor (or its principal values in decreasing order)
- Returns:
4th-order orientation tensor written using the \((\phi,\phi)\) bases
- Return type:
array of shape (6, 6)
References
Advani, S. G. & Tucker III, C. L. The use of tensors to describe and predict fiber orientation in short fiber composites. Journal of Rheology, SOR, 1987, 31, 751-784
Utility functions for manupulating 2nd/4th-order tensors#
- fiberpy.tensor.Mat2(sig)#
Bijection between a symmetric 2nd order tensor space and 6-dim vector space using the \(\phi\) basis
\[\begin{split}\begin{bmatrix} s_{11} & s_{12} & s_{13} \\ s_{12} & s_{22} & s_{23} \\ s_{13} & s_{23} & s_{33} \\ \end{bmatrix}\iff\begin{bmatrix} s_{11} \\ s_{22} \\ s_{33} \\ s_{12} \\ s_{23} \\ s_{13} \end{bmatrix}\end{split}\]- Parameters:
sig (ndarray of dim 1 or 2) – 2nd-order tensor
- fiberpy.tensor.Mat22(eps)#
Bijection between a symmetric 2nd order tensor space and 6-dim vector space using the \(\phi_2\) basis
\[\begin{split}\begin{bmatrix} e_{11} & e_{12} & e_{13} \\ e_{12} & e_{22} & e_{23} \\ e_{13} & e_{23} & e_{33} \\ \end{bmatrix}\iff\begin{bmatrix} e_{11} \\ e_{22} \\ e_{33} \\ 2e_{12} \\ 2e_{23} \\ 2e_{13} \end{bmatrix}\end{split}\]- Parameters:
eps (ndarray of dim 1 or 2) – 2nd-order tensor
- fiberpy.tensor.Mat2S(eps)#
Bijection between a symmetric 2nd order tensor space and 6-dim vector space using the \(\phi_\mathrm{S}\) basis
\[\begin{split}\begin{bmatrix} e_{11} & e_{12} & e_{13} \\ e_{12} & e_{22} & e_{23} \\ e_{13} & e_{23} & e_{33} \\ \end{bmatrix}\iff\begin{bmatrix} e_{11} \\ e_{22} \\ e_{33} \\ \sqrt{2}e_{12} \\ \sqrt{2}e_{23} \\ \sqrt{2}e_{13} \end{bmatrix}\end{split}\]- Parameters:
eps (ndarray of dim 1 or 2) – 2nd-order tensor
- fiberpy.tensor.Mat4(A)#
Matrix representation of a 4th order tensor with minor symmety using the \((\phi,\phi)\) bases
- Parameters:
A (ndarray of shape (3, 3, 3, 3)) – 4th-order tensor
- fiberpy.tensor.MatGP(v)#
Matrix that converts a 2nd-order tensor from the global frame (\(\phi\) basis) to the principal frame (\(\phi\) basis)
- Parameters:
v (ndarray of shape (3, 3)) – Principal directions along its columns
- fiberpy.tensor.MatGP2(v)#
Matrix that converts a 2nd-order tensor from the global frame (\(\phi_2\) basis) to the principal frame (\(\phi\) basis)
- Parameters:
v (ndarray of shape (3, 3)) – Principal directions along its columns
- fiberpy.tensor.MatPG(v)#
Matrix that converts a 2nd-order tensor from the principal frame (\(\phi\) basis) to the global frame (\(\phi\) basis)
- Parameters:
v (ndarray of shape (3, 3)) – Principal directions along its columns
- fiberpy.tensor.ij2M(ij)#
Convert (i, j) indices of a symmetric 2nd-order tensor to its vector index
- fiberpy.tensor.ijkl2MN(ijkl)#
Convert (i, j, k, l) indices of a symmetric 4nd-order tensor to its matrix index