Skip to content

fatpy

FatPy - A Python package for fatigue life evaluation of materials.

core

Core fatigue analysis methods and functionality.

damage_cumulation

Damage cumulation module.

This module contains the classes and functions for damage cumulation methods.

decompositions

Decompositions module.

This module contains the classes and functions for decomposition methods.

energy_life

Energy-based fatigue life analysis methods.

This package provides implementations of energy-based approaches for fatigue life prediction and analysis.

base_methods

Base methods for the energy-life.

correction_methods

Correction methods for the energy-life.

decompositions

Decomposition methods for the energy-life.

demage_parameters

Damage parameters calculation methods for the energy-life.

plane_based_methods

Plane-based methods module.

This module contains the classes and functions for plane-based methods.

strain_life

Strain-based fatigue life analysis methods.

This package provides implementations of strain-based approaches for fatigue life prediction and analysis.

base_methods

Base methods for the strain-life.

correction_methods

Correction methods for the strain-life.

demage_parameters

Damage parameters calculation methods for the strain-life.

stress_life

Stress-based fatigue life analysis methods.

This package provides implementations of stress-based approaches for fatigue life prediction and analysis.

base_methods

Base methods for the stress-life.

correction_methods

Correction methods for the stress-life.

demage_parameters

Damage parameters calculation methods for the stress-life.

data_parsing

Data parsing utilities for the fatigue analysis package.

fe_model

Finite element model parsing module.

loads

Load data parsing module.

material

Material properties parsing module.

user_input

User input data parsing module.

examples

This module contains an example function with a detailed docstring.

docstring_example_tmp

This module contains an example function with a detailed docstring.

ExampleClass

An example class to demonstrate docstring formatting.

Source code in src/fatpy/examples/docstring_example_tmp.py
class ExampleClass:
    """An example class to demonstrate docstring formatting."""

    def __init__(self, value: int):
        """Initialize the ExampleClass with a value.

        Args:
            value (int): The initial value for the class instance.
        """
        self.value = value

    def increment(self, amount: int) -> int:
        """Increment the stored value by a specified amount.

        Args:
            amount (int): The amount to increment the value by.

        Returns:
            int: The new value after incrementing.
        """
        self.value += amount
        return self.value

    def example_method_with_docstring(self, a: int, b: int) -> int:
        # ruff: noqa: E501
        """Docstring with a cross-reference to the example function.

        This method demonstrates how to reference another function's docstring.
        It calls `example_function_with_docstring` with sample arguments.

        Cross-reference:
            1. [fatpy.examples.docstring_example_tmp.example_function_with_docstring][]
            2. [`Reference with title`][fatpy.examples.docstring_example_tmp.example_function_with_docstring]

        Args:
            a (int): The first integer value.
            b (int): The second integer value.

        Returns:
            int: The result of the example function.
        """
        return example_function_with_docstring(a, b)
__init__(value)

Initialize the ExampleClass with a value.

Parameters:

Name Type Description Default
value int

The initial value for the class instance.

required
Source code in src/fatpy/examples/docstring_example_tmp.py
def __init__(self, value: int):
    """Initialize the ExampleClass with a value.

    Args:
        value (int): The initial value for the class instance.
    """
    self.value = value
increment(amount)

Increment the stored value by a specified amount.

Parameters:

Name Type Description Default
amount int

The amount to increment the value by.

required

Returns:

Name Type Description
int int

The new value after incrementing.

Source code in src/fatpy/examples/docstring_example_tmp.py
def increment(self, amount: int) -> int:
    """Increment the stored value by a specified amount.

    Args:
        amount (int): The amount to increment the value by.

    Returns:
        int: The new value after incrementing.
    """
    self.value += amount
    return self.value
example_method_with_docstring(a, b)

Docstring with a cross-reference to the example function.

This method demonstrates how to reference another function's docstring. It calls example_function_with_docstring with sample arguments.

Cross-reference
  1. fatpy.examples.docstring_example_tmp.example_function_with_docstring
  2. Reference with title

Parameters:

Name Type Description Default
a int

The first integer value.

required
b int

The second integer value.

required

Returns:

Name Type Description
int int

The result of the example function.

Source code in src/fatpy/examples/docstring_example_tmp.py
def example_method_with_docstring(self, a: int, b: int) -> int:
    # ruff: noqa: E501
    """Docstring with a cross-reference to the example function.

    This method demonstrates how to reference another function's docstring.
    It calls `example_function_with_docstring` with sample arguments.

    Cross-reference:
        1. [fatpy.examples.docstring_example_tmp.example_function_with_docstring][]
        2. [`Reference with title`][fatpy.examples.docstring_example_tmp.example_function_with_docstring]

    Args:
        a (int): The first integer value.
        b (int): The second integer value.

    Returns:
        int: The result of the example function.
    """
    return example_function_with_docstring(a, b)

example_function_with_docstring(a, b)

This docstring highlights Mermaid diagrams, math expressions, and code examples.

This function takes two integers and returns their sum. It demonstrates various mkdocs-material features:

Information

This is an informational admonition block.

Important Note

Make sure both parameters are integers.

Mermaid Diagram
graph TD
    A[Start] --> B{Is it a number?}
    B -- Yes --> C[Process the number]
    B -- No --> D[Throw an error]
    C --> E[End]
    D --> E
Mathematical Expression
\[ \varphi = \frac{1 + \sqrt{5}}{2} \]

Inline math: \(E = mc^2\)

Code Examples
def add(a: int, b: int) -> int:
    return a + b
function add(a, b) {
    return a + b;
}
Expandable Example

This is a collapsible content section that shows more detailed usage:

result = example_function_with_docstring(5, 3)
print(f"Result: {result}")  # Output: Result: 8
Table Example
Input A Input B Result
1 2 3
5 7 12
0 0 0
Args:
  • a (int): The first integer.
  • b (int): The second integer.
Returns:
  • int: The sum of the two integers.
Raises:
  • TypeError: If inputs are not integers.
Examples:
>>> example_function_with_docstring(2, 3)
5
>>> example_function_with_docstring(0, 0)
0
Source code in src/fatpy/examples/docstring_example_tmp.py
def example_function_with_docstring(a: int, b: int) -> int:
    r"""This docstring highlights Mermaid diagrams, math expressions, and code examples.

    This function takes two integers and returns their sum.
    It demonstrates various mkdocs-material features:

    !!! info "Information"
        This is an informational admonition block.

    !!! warning "Important Note"
        Make sure both parameters are integers.

    ### Mermaid Diagram
    ```mermaid
    graph TD
        A[Start] --> B{Is it a number?}
        B -- Yes --> C[Process the number]
        B -- No --> D[Throw an error]
        C --> E[End]
        D --> E
    ```

    ### Mathematical Expression
    $$ \varphi = \frac{1 + \sqrt{5}}{2} $$

    Inline math: $E = mc^2$

    ### Code Examples
    === "Python"
        ```python
        def add(a: int, b: int) -> int:
            return a + b
        ```

    === "JavaScript"
        ```javascript
        function add(a, b) {
            return a + b;
        }
        ```

    ??? example "Expandable Example"
        This is a collapsible content section that shows more detailed usage:

        ```python
        result = example_function_with_docstring(5, 3)
        print(f"Result: {result}")  # Output: Result: 8
        ```

    ### Table Example
    | Input A | Input B | Result |
    |---------|---------|--------|
    | 1       | 2       | 3      |
    | 5       | 7       | 12     |
    | 0       | 0       | 0      |

    #### Args:
    - `a` (*int*): The first integer.
    - `b` (*int*): The second integer.

    #### Returns:
    - *int*: The sum of the two integers.

    #### Raises:
    - `TypeError`: If inputs are not integers.

    #### Examples:
    ```python
    >>> example_function_with_docstring(2, 3)
    5
    >>> example_function_with_docstring(0, 0)
    0
    ```
    """
    return a + b

material_laws

Material laws module.

Contains the classes and functions for defining and working with material laws.

struct_mech

Structural mechanics module.

This module contains the classes and functions for structural mechanics analysis.

stress

Calculate fundamental stress metrics and invariants.

These functions provide principal stresses, maximum shear (Tresca) stress, hydrostatic stress, von Mises equivalent stress, and invariants of the stress tensor. They are essential for strength, fatigue, and fracture analyses under both uniaxial and multiaxial loading conditions.

Conventions: - Vectors use Voigt notation with shape (..., 6), where the last dimension contains the six Voigt components and leading dimensions are preserved:

    (σ_11, σ_22, σ_33, σ_23, σ_13, σ_12)
    (σ_xx, σ_yy, σ_zz, σ_yz, σ_xz, σ_xy)
  • Principal stresses are ordered in descending order throughout the module (σ_1 ≥ σ_2 ≥ σ_3).

  • Principal directions (eigenvectors) are aligned to this ordering (columns correspond to σ_1, σ_2, σ_3).

calc_principal_stresses_and_directions(stress_vector_voigt)

Calculate principal stresses and principal directions for each state.

Math Equations

Principal stresses and directions are found by solving the eigenvalue problem for the stress tensor:

\[ \sigma \mathbf{v} = \lambda \mathbf{v} \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required

Returns:

Name Type Description
Tuple (eigvals, eigvecs)
NDArray[float64]
  • eigvals: Array of shape (..., 3). Principal stresses
tuple[NDArray[float64], NDArray[float64]]

(descending: σ_1 ≥ σ_2 ≥ σ_3) with leading dimensions preserved.

tuple[NDArray[float64], NDArray[float64]]
  • eigvecs: Array of shape (..., 3, 3). Principal directions (columns are
tuple[NDArray[float64], NDArray[float64]]

eigenvectors) aligned with eigvals in the same order. The last two

tuple[NDArray[float64], NDArray[float64]]

dimensions are the 3x3 eigenvector matrix for each input.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_principal_stresses_and_directions(
    stress_vector_voigt: NDArray[np.float64],
) -> tuple[NDArray[np.float64], NDArray[np.float64]]:
    r"""Calculate principal stresses and principal directions for each state.

    ??? abstract "Math Equations"
        Principal stresses and directions are found by solving the eigenvalue problem
        for the stress tensor:

        $$ \sigma \mathbf{v} = \lambda \mathbf{v} $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.

    Returns:
        Tuple (eigvals, eigvecs):
        - eigvals: Array of shape (..., 3). Principal stresses
        (descending: σ_1 ≥ σ_2 ≥ σ_3) with leading dimensions preserved.
        - eigvecs: Array of shape (..., 3, 3). Principal directions (columns are
        eigenvectors) aligned with eigvals in the same order. The last two
        dimensions are the 3x3 eigenvector matrix for each input.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    voigt.check_shape(stress_vector_voigt)

    tensor = voigt.voigt_to_tensor(stress_vector_voigt)
    eigvals, eigvecs = np.linalg.eigh(tensor)
    sorted_indices = np.argsort(eigvals, axis=-1)[..., ::-1]
    eigvals_sorted = np.take_along_axis(eigvals, sorted_indices, axis=-1)
    eigvecs_sorted = np.take_along_axis(
        eigvecs, np.expand_dims(sorted_indices, axis=-2), axis=-1
    )

    return eigvals_sorted, eigvecs_sorted

calc_principal_stresses(stress_vector_voigt)

Calculate principal stresses for each stress state.

Math Equations

Principal stresses are found by solving the eigenvalue problem for the stress tensor:

\[ \sigma \mathbf{v} = \lambda \mathbf{v} \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required

Returns:

Type Description
NDArray[float64]

Array of shape (..., 3). Principal stresses (descending: σ1 ≥ σ2 ≥ σ3).

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_principal_stresses(
    stress_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate principal stresses for each stress state.

    ??? abstract "Math Equations"
        Principal stresses are found by solving the eigenvalue problem
        for the stress tensor:

        $$ \sigma \mathbf{v} = \lambda \mathbf{v} $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.

    Returns:
        Array of shape (..., 3). Principal stresses (descending: σ1 ≥ σ2 ≥ σ3).

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    voigt.check_shape(stress_vector_voigt)

    tensor = voigt.voigt_to_tensor(stress_vector_voigt)
    eigvals = np.linalg.eigvalsh(tensor)

    # eigvals shape: (n, ..., 3). Return sorted descending along last axis.
    return np.sort(eigvals, axis=-1)[..., ::-1]

calc_stress_invariants(stress_vector_voigt)

Calculate the first, second, and third invariants for each stress state.

Math Equations
\[ \begin{align*} I_1 &= tr(\sigma), \\ I_2 &= \frac{1}{2}(I_1^{2} - tr(\sigma^{2})), \\ I_3 &= \det(\sigma) \end{align*} \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required

Returns:

Type Description
NDArray[float64]

Array of shape (..., 3). The last dimension contains (I1, I2, I3) for each entry.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_stress_invariants(
    stress_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate the first, second, and third invariants for each stress state.

    ??? abstract "Math Equations"
        $$
        \begin{align*}
        I_1 &=  tr(\sigma), \\
        I_2 &= \frac{1}{2}(I_1^{2} - tr(\sigma^{2})), \\
        I_3 &= \det(\sigma)
        \end{align*}
        $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.

    Returns:
        Array of shape (..., 3). The last dimension contains (I1, I2, I3) for
            each entry.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    voigt.check_shape(stress_vector_voigt)

    tensor = voigt.voigt_to_tensor(stress_vector_voigt)
    invariant_1 = np.trace(tensor, axis1=-2, axis2=-1)
    invariant_2 = 0.5 * (
        invariant_1**2 - np.trace(np.matmul(tensor, tensor), axis1=-2, axis2=-1)
    )
    invariant_3 = np.linalg.det(tensor)

    return np.stack((invariant_1, invariant_2, invariant_3), axis=-1)

calc_hydrostatic_stress(stress_vector_voigt)

Calculate the hydrostatic (mean normal) stress for each stress state.

Math Equations
\[ \sigma_H = \frac{1}{3} tr(\sigma) = \frac{1}{3} (\sigma_{11} + \sigma_{22} + \sigma_{33}) \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required

Returns:

Type Description
NDArray[float64]

Array of shape (...). Hydrostatic stress for each input state. Tensor rank is reduced by one.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_hydrostatic_stress(
    stress_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate the hydrostatic (mean normal) stress for each stress state.

    ??? abstract "Math Equations"
        $$
        \sigma_H = \frac{1}{3} tr(\sigma) =
        \frac{1}{3} (\sigma_{11} + \sigma_{22} + \sigma_{33})
        $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.

    Returns:
        Array of shape (...). Hydrostatic stress for each input state.
            Tensor rank is reduced by one.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    voigt.check_shape(stress_vector_voigt)

    # Voigt normal components are at last axis indices 0,1,2
    return (
        stress_vector_voigt[..., 0]
        + stress_vector_voigt[..., 1]
        + stress_vector_voigt[..., 2]
    ) / 3.0

calc_stress_deviator(stress_vector_voigt)

Calculate stress deviator for each stress state.

Math Equations
\[ \mathbf{s} = \sigma - \frac{1}{3} tr(\sigma) \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required

Returns:

Type Description
NDArray[float64]

Array of shape (..., 6). Stress deviator for each entry.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_stress_deviator(
    stress_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate stress deviator for each stress state.

    ??? abstract "Math Equations"
        $$ \mathbf{s} = \sigma - \frac{1}{3} tr(\sigma) $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.

    Returns:
        Array of shape (..., 6). Stress deviator for each entry.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    voigt.check_shape(stress_vector_voigt)
    hydrostatic = calc_hydrostatic_stress(stress_vector_voigt)

    deviator = stress_vector_voigt.copy()
    # Subtract hydrostatic from the first three Voigt components (last axis)
    deviator[..., :3] = deviator[..., :3] - hydrostatic[..., None]

    return deviator

calc_von_mises_stress(stress_vector_voigt)

Calculate von Mises equivalent stress for each stress state.

Math Equations
\[ \sigma_{vM} = \tfrac{\sqrt{2}}{2}\sqrt{ (\sigma_{11}-\sigma_{22})^2 +(\sigma_{22}-\sigma_{33})^2 +(\sigma_{33}-\sigma_{11})^2 + 3(\sigma_{12}^2+\sigma_{23}^2+\sigma_{13}^2)} \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required

Returns:

Type Description
NDArray[float64]

Array of shape (...). Von Mises equivalent stress for each entry. Tensor rank is reduced by one.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_von_mises_stress(
    stress_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate von Mises equivalent stress for each stress state.

    ??? abstract "Math Equations"
        $$
        \sigma_{vM} = \tfrac{\sqrt{2}}{2}\sqrt{
        (\sigma_{11}-\sigma_{22})^2
        +(\sigma_{22}-\sigma_{33})^2
        +(\sigma_{33}-\sigma_{11})^2
        + 3(\sigma_{12}^2+\sigma_{23}^2+\sigma_{13}^2)}
        $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.

    Returns:
        Array of shape (...). Von Mises equivalent stress for each entry.
            Tensor rank is reduced by one.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    voigt.check_shape(stress_vector_voigt)

    sx = stress_vector_voigt[..., 0]
    sy = stress_vector_voigt[..., 1]
    sz = stress_vector_voigt[..., 2]
    syz = stress_vector_voigt[..., 3]
    sxz = stress_vector_voigt[..., 4]
    sxy = stress_vector_voigt[..., 5]

    # Von Mises formula expanded to simplify computation
    return np.sqrt(
        0.5
        * (
            (sx - sy) ** 2
            + (sy - sz) ** 2
            + (sz - sx) ** 2
            + 6 * (sxy**2 + syz**2 + sxz**2)
        )
    )

calc_signed_von_mises_by_hydrostatic(stress_vector_voigt, rtol=1e-05, atol=1e-08)

Calculate signed von Mises stress for each stress state.

Sign is determined by the hydrostatic stress.

The sign assignment follows these rules: - Positive (+): When hydrostatic stress > 0 (tensile dominant state) - Negative (-): When hydrostatic stress < 0 (compressive dominant state) - Positive (+): When hydrostatic stress ≈ 0 (within tolerance, default fallback)

Tolerance parameters ensure numerical stability in edge cases where the determining value is very close to zero, preventing erratic sign changes that could occur due to floating-point precision limitations.

Math Equations
\[ \sigma_{SvM} = \begin{cases} +\sigma_{vM} & \text{if } \sigma_H \geq 0 \\ -\sigma_{vM} & \text{if } \sigma_H < 0 \end{cases} \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required
rtol float

Relative tolerance for comparing hydrostatic stress to zero. Default is 1e-5.

1e-05
atol float

Absolute tolerance for comparing hydrostatic stress to zero. Default is 1e-8.

1e-08

Returns:

Type Description
NDArray[float64]

Array of shape (...). Signed von Mises stress for each entry. Tensor rank is reduced by one.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_signed_von_mises_by_hydrostatic(
    stress_vector_voigt: NDArray[np.float64],
    rtol: float = 1e-5,
    atol: float = 1e-8,
) -> NDArray[np.float64]:
    r"""Calculate signed von Mises stress for each stress state.

    Sign is determined by the hydrostatic stress.

    The sign assignment follows these rules:
    - **Positive (+)**: When hydrostatic stress > 0 (tensile dominant state)
    - **Negative (-)**: When hydrostatic stress < 0 (compressive dominant state)
    - **Positive (+)**: When hydrostatic stress ≈ 0 (within tolerance, default fallback)

    Tolerance parameters ensure numerical stability in edge cases where the
    determining value is very close to zero, preventing erratic sign changes
    that could occur due to floating-point precision limitations.

    ??? abstract "Math Equations"
        $$
        \sigma_{SvM} = \begin{cases}
        +\sigma_{vM} & \text{if } \sigma_H \geq 0 \\
        -\sigma_{vM} & \text{if } \sigma_H < 0
        \end{cases}
        $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.
        rtol: Relative tolerance for comparing hydrostatic stress to zero.
            Default is 1e-5.
        atol: Absolute tolerance for comparing hydrostatic stress to zero.
            Default is 1e-8.

    Returns:
        Array of shape (...). Signed von Mises stress for each entry.
            Tensor rank is reduced by one.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    von_mises = calc_von_mises_stress(stress_vector_voigt)
    hydrostatic_stress = calc_hydrostatic_stress(stress_vector_voigt)

    sign = np.sign(hydrostatic_stress).astype(np.float64, copy=False)
    sign = np.where(np.isclose(hydrostatic_stress, 0, rtol=rtol, atol=atol), 1.0, sign)

    return sign * von_mises

calc_signed_von_mises_by_max_abs_principal(stress_voigt, rtol=1e-05, atol=1e-08)

Calculate signed von Mises stress for each stress state.

Sign is determined by average of the maximum and minimum principal stresses.

The sign assignment follows these rules: - Positive (+): When (σ₁ + σ₃)/2 > 0 (tension dominant) - Negative (-): When (σ₁ + σ₃)/2 < 0 (compression dominant) - Positive (+): When (σ₁ + σ₃)/2 ≈ 0 (within tolerance, default fallback)

Tolerance parameters ensure numerical stability in edge cases where the determining value is very close to zero, preventing erratic sign changes that could occur due to floating-point precision limitations.

Math Equations
\[ \sigma_{SvM} = \begin{cases} +\sigma_{vM} & \text{if } \frac{\sigma_1 + \sigma_3}{2} \geq 0 \\ -\sigma_{vM} & \text{if } \frac{\sigma_1 + \sigma_3}{2} < 0 \end{cases} \]

Parameters:

Name Type Description Default
stress_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required
rtol float

Relative tolerance for comparing the average of maximum and minimum principal stresses to zero. Default is 1e-5.

1e-05
atol float

Absolute tolerance for comparing the average of maximum and minimum principal stresses to zero. Default is 1e-8.

1e-08

Returns:

Type Description
NDArray[float64]

Array of shape (...). Signed von Mises stress for each entry. Tensor rank is reduced by one.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_signed_von_mises_by_max_abs_principal(
    stress_voigt: NDArray[np.float64],
    rtol: float = 1e-5,
    atol: float = 1e-8,
) -> NDArray[np.float64]:
    r"""Calculate signed von Mises stress for each stress state.

    Sign is determined by average of the maximum and minimum principal stresses.

    The sign assignment follows these rules:
    - **Positive (+)**: When (σ₁ + σ₃)/2 > 0 (tension dominant)
    - **Negative (-)**: When (σ₁ + σ₃)/2 < 0 (compression dominant)
    - **Positive (+)**: When (σ₁ + σ₃)/2 ≈ 0 (within tolerance, default fallback)

    Tolerance parameters ensure numerical stability in edge cases where the
    determining value is very close to zero, preventing erratic sign changes
    that could occur due to floating-point precision limitations.

    ??? abstract "Math Equations"
        $$
        \sigma_{SvM} = \begin{cases}
        +\sigma_{vM} & \text{if } \frac{\sigma_1 + \sigma_3}{2} \geq 0 \\
        -\sigma_{vM} & \text{if } \frac{\sigma_1 + \sigma_3}{2} < 0
        \end{cases}
        $$

    Args:
        stress_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.
        rtol: Relative tolerance for comparing the average of maximum and minimum
            principal stresses to zero. Default is 1e-5.
        atol: Absolute tolerance for comparing the average of maximum and minimum
            principal stresses to zero. Default is 1e-8.

    Returns:
        Array of shape (...). Signed von Mises stress for each entry.
            Tensor rank is reduced by one.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    von_mises = calc_von_mises_stress(stress_voigt)
    principals = calc_principal_stresses(stress_voigt)

    avg_13 = 0.5 * (principals[..., 0] + principals[..., 2])
    sign = np.sign(avg_13).astype(np.float64, copy=False)
    sign = np.where(np.isclose(avg_13, 0, rtol=rtol, atol=atol), 1.0, sign)

    return sign * von_mises

calc_signed_von_mises_by_first_invariant(stress_vector_voigt, rtol=1e-05, atol=1e-08)

Calculate signed von Mises stress for each stress state.

Sign is determined by the first invariant of the stress tensor.

The sign assignment follows these rules: - Positive (+): When tr(σ) > 0 (tensile volumetric stress) - Negative (-): When tr(σ) < 0 (compressive volumetric stress) - Positive (+): When tr(σ) ≈ 0 (within tolerance, default fallback)

Tolerance parameters ensure numerical stability in edge cases where the determining value is very close to zero, preventing erratic sign changes that could occur due to floating-point precision limitations.

Math Equations
\[ \sigma_{SvM} = \begin{cases} +\sigma_{vM} & \text{if } tr(\sigma) \geq 0 \\ -\sigma_{vM} & \text{if } tr(\sigma) < 0 \end{cases} \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required
rtol float

Relative tolerance for comparing the first invariant to zero. Default is 1e-5.

1e-05
atol float

Absolute tolerance for comparing the first invariant to zero. Default is 1e-8.

1e-08

Returns:

Type Description
NDArray[float64]

Array of shape (...). Signed von Mises stress for each entry. Tensor rank is reduced by one.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_signed_von_mises_by_first_invariant(
    stress_vector_voigt: NDArray[np.float64],
    rtol: float = 1e-5,
    atol: float = 1e-8,
) -> NDArray[np.float64]:
    r"""Calculate signed von Mises stress for each stress state.

    Sign is determined by the first invariant of the stress tensor.

    The sign assignment follows these rules:
    - **Positive (+)**: When tr(σ) > 0 (tensile volumetric stress)
    - **Negative (-)**: When tr(σ) < 0 (compressive volumetric stress)
    - **Positive (+)**: When tr(σ) ≈ 0 (within tolerance, default fallback)

    Tolerance parameters ensure numerical stability in edge cases where the
    determining value is very close to zero, preventing erratic sign changes
    that could occur due to floating-point precision limitations.

    ??? abstract "Math Equations"
        $$
        \sigma_{SvM} = \begin{cases}
        +\sigma_{vM} & \text{if } tr(\sigma) \geq 0 \\
        -\sigma_{vM} & \text{if } tr(\sigma) < 0
        \end{cases}
        $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.
        rtol: Relative tolerance for comparing the first invariant to zero.
            Default is 1e-5.
        atol: Absolute tolerance for comparing the first invariant to zero.
            Default is 1e-8.

    Returns:
        Array of shape (...). Signed von Mises stress for each entry.
            Tensor rank is reduced by one.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    von_mises = calc_von_mises_stress(stress_vector_voigt)
    invariant_1 = (
        stress_vector_voigt[..., 0]
        + stress_vector_voigt[..., 1]
        + stress_vector_voigt[..., 2]
    )

    sign = np.sign(invariant_1).astype(np.float64, copy=False)
    sign = np.where(np.isclose(invariant_1, 0, rtol=rtol, atol=atol), 1.0, sign)

    return sign * von_mises

calc_tresca_stress(stress_vector_voigt)

Calculate Tresca (maximum shear) stress for each stress state.

Math Equations
\[\sigma_{\tau_{max}} = \frac{\sigma_{1} - \sigma_{3}}{2}\]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required

Returns:

Type Description
NDArray[float64]

Array of shape (...). Tresca stress for each entry. For principal stresses σ1 ≥ σ2 ≥ σ3, Tresca = (σ1 − σ3)/2. Tensor rank is reduced by one.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_tresca_stress(stress_vector_voigt: NDArray[np.float64]) -> NDArray[np.float64]:
    r"""Calculate Tresca (maximum shear) stress for each stress state.

    ??? abstract "Math Equations"
        $$\sigma_{\tau_{max}} = \frac{\sigma_{1} - \sigma_{3}}{2}$$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.

    Returns:
        Array of shape (...). Tresca stress for each entry. For principal
            stresses σ1 ≥ σ2 ≥ σ3, Tresca = (σ1 − σ3)/2.
            Tensor rank is reduced by one.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    principals = calc_principal_stresses(stress_vector_voigt)
    return 0.5 * (principals[..., 0] - principals[..., 2])

calc_signed_tresca_by_hydrostatic(stress_vector_voigt, rtol=1e-05, atol=1e-08)

Calculate signed Tresca stress for each stress state.

Sign is determined by the hydrostatic stress.

The sign assignment follows these rules: - Positive (+): When hydrostatic stress > 0 (tensile dominant state) - Negative (-): When hydrostatic stress < 0 (compressive dominant state) - Positive (+): When hydrostatic stress ≈ 0 (within tolerance, default fallback)

Tolerance parameters ensure numerical stability in edge cases where the determining value is very close to zero, preventing erratic sign changes that could occur due to floating-point precision limitations.

Math Equations
\[ \sigma_{S\tau_{max}} = \begin{cases} +\sigma_{\tau_{max}} & \text{if } \sigma_H \geq 0 \\ -\sigma_{\tau_{max}} & \text{if } \sigma_H < 0 \end{cases} \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required
rtol float

Relative tolerance for comparing hydrostatic stress to zero. Default is 1e-5.

1e-05
atol float

Absolute tolerance for comparing hydrostatic stress to zero. Default is 1e-8.

1e-08

Returns:

Type Description
NDArray[float64]

Array of shape (...). Signed Tresca stress for each entry. Tensor rank is reduced by one.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_signed_tresca_by_hydrostatic(
    stress_vector_voigt: NDArray[np.float64],
    rtol: float = 1e-5,
    atol: float = 1e-8,
) -> NDArray[np.float64]:
    r"""Calculate signed Tresca stress for each stress state.

    Sign is determined by the hydrostatic stress.

    The sign assignment follows these rules:
    - **Positive (+)**: When hydrostatic stress > 0 (tensile dominant state)
    - **Negative (-)**: When hydrostatic stress < 0 (compressive dominant state)
    - **Positive (+)**: When hydrostatic stress ≈ 0 (within tolerance, default fallback)

    Tolerance parameters ensure numerical stability in edge cases where the
    determining value is very close to zero, preventing erratic sign changes
    that could occur due to floating-point precision limitations.

    ??? abstract "Math Equations"
        $$
        \sigma_{S\tau_{max}} = \begin{cases}
        +\sigma_{\tau_{max}} & \text{if } \sigma_H \geq 0 \\
        -\sigma_{\tau_{max}} & \text{if } \sigma_H < 0
        \end{cases}
        $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.
        rtol: Relative tolerance for comparing hydrostatic stress to zero.
            Default is 1e-5.
        atol: Absolute tolerance for comparing hydrostatic stress to zero.
            Default is 1e-8.

    Returns:
        Array of shape (...). Signed Tresca stress for each entry.
            Tensor rank is reduced by one.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    tresca = calc_tresca_stress(stress_vector_voigt)
    hydrostatic_stress = calc_hydrostatic_stress(stress_vector_voigt)

    sign = np.sign(hydrostatic_stress)
    sign = np.where(np.isclose(hydrostatic_stress, 0, rtol=rtol, atol=atol), 1.0, sign)

    return sign * tresca

calc_signed_tresca_by_max_abs_principal(stress_vector_voigt, rtol=1e-05, atol=1e-08)

Calculate signed Tresca stress for each stress state.

Sign is determined by the maximum absolute principal stress value.

The sign assignment follows these rules: - Positive (+): When (σ₁ + σ₃)/2 > 0 (tension dominant) - Negative (-): When (σ₁ + σ₃)/2 < 0 (compression dominant) - Positive (+): When (σ₁ + σ₃)/2 ≈ 0 (within tolerance, default fallback)

Tolerance parameters ensure numerical stability in edge cases where the determining value is very close to zero, preventing erratic sign changes that could occur due to floating-point precision limitations.

Math Equations
\[ \sigma_{S\tau_{max}} = \begin{cases} +\sigma_{\tau_{max}} & \text{if } \frac{\sigma_1 + \sigma_3}{2} \geq 0 \\ -\sigma_{\tau_{max}} & \text{if } \frac{\sigma_1 + \sigma_3}{2} < 0 \end{cases} \]

Parameters:

Name Type Description Default
stress_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt stress components. Leading dimensions are preserved.

required
rtol float

Relative tolerance for comparing the average of maximum and minimum principal stresses to zero. Default is 1e-5.

1e-05
atol float

Absolute tolerance for comparing the average of maximum and minimum principal stresses to zero. Default is 1e-8.

1e-08

Returns:

Type Description
NDArray[float64]

Array of shape (...). Signed Tresca stress for each entry. Tensor rank is reduced by one.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/stress.py
def calc_signed_tresca_by_max_abs_principal(
    stress_vector_voigt: NDArray[np.float64],
    rtol: float = 1e-5,
    atol: float = 1e-8,
) -> NDArray[np.float64]:
    r"""Calculate signed Tresca stress for each stress state.

    Sign is determined by the maximum absolute principal stress value.

    The sign assignment follows these rules:
    - **Positive (+)**: When (σ₁ + σ₃)/2 > 0 (tension dominant)
    - **Negative (-)**: When (σ₁ + σ₃)/2 < 0 (compression dominant)
    - **Positive (+)**: When (σ₁ + σ₃)/2 ≈ 0 (within tolerance, default fallback)

    Tolerance parameters ensure numerical stability in edge cases where the
    determining value is very close to zero, preventing erratic sign changes
    that could occur due to floating-point precision limitations.

    ??? abstract "Math Equations"
        $$
        \sigma_{S\tau_{max}} = \begin{cases}
        +\sigma_{\tau_{max}} & \text{if } \frac{\sigma_1 + \sigma_3}{2} \geq 0 \\
        -\sigma_{\tau_{max}} & \text{if } \frac{\sigma_1 + \sigma_3}{2} < 0
        \end{cases}
        $$

    Args:
        stress_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt stress components. Leading dimensions are preserved.
        rtol: Relative tolerance for comparing the average of maximum and minimum
            principal stresses to zero. Default is 1e-5.
        atol: Absolute tolerance for comparing the average of maximum and minimum
            principal stresses to zero. Default is 1e-8.

    Returns:
        Array of shape (...). Signed Tresca stress for each entry.
            Tensor rank is reduced by one.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    tresca = calc_tresca_stress(stress_vector_voigt)
    principals = calc_principal_stresses(stress_vector_voigt)

    avg_13 = 0.5 * (principals[..., 0] + principals[..., 2])
    sign = np.sign(avg_13).astype(np.float64, copy=False)
    sign = np.where(np.isclose(avg_13, 0, rtol=rtol, atol=atol), 1.0, sign)

    return sign * tresca

transformations

Contains the functions for tensor and other transformations.

utils

Utility functions for the fatigue analysis package.

mesh

Mesh processing utilities module.

This module provides various functions and classes for mesh processing tasks.

post_processing

Utilities for post-processing fatigue analysis results.

pre_processing

Utilities for pre-processing fatigue analysis results.

signal

Signal processing utilities module.

This module provides various functions and classes for signal processing tasks.

voigt

Tools for working with Voigt notation in continuum mechanics.

Overview

This module provides general utilities for handling vectors and tensors represented in Voigt notation. Voigt notation is widely used to express symmetric 3x3 tensors, such as stress and strain, in a compact vector form.

Input Shape Convention
  • Arrays are expected to have shape (..., 6), where the last dimension contains the six Voigt components: (σ_11, σ_22, σ_33, σ_23, σ_13, σ_12) (σ_xx, σ_yy, σ_zz, σ_yz, σ_xz, σ_xy) Leading dimensions (the '...' part) are preserved and can be arbitrary, allowing batch processing of multiple vectors or tensors.

check_shape(vector)

Validate the Voigt vector shape.

Parameters:

Name Type Description Default
vector NDArray[float64]

Array with shape (..., 6) where the last dimension has length 6.

required

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/utils/voigt.py
def check_shape(vector: NDArray[np.float64]) -> None:
    """Validate the Voigt vector shape.

    Args:
        vector: Array with shape (..., 6) where the last dimension has length 6.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    if vector.shape[-1] != VOIGT_COMPONENTS_COUNT:
        raise ValueError(
            "Last dimension must correspond to 6 Voigt "
            "stress/strain components (..., 6)."
        )

voigt_to_tensor(vector)

Convert Voigt vectors to symmetric 3x3 tensors.

Parameters:

Name Type Description Default
vector NDArray[float64]

Array of shape (..., 6) where the last dimension contains the stress/strain components in order according to Voigt notation: (σ_11, σ_22, σ_33, σ_23, σ_13, σ_12) (σ_xx, σ_yy, σ_zz, σ_yz, σ_xz, σ_xy)

required

Returns:

Type Description
NDArray[float64]

Array with shape (..., 3, 3). The last two axes are the symmetric

NDArray[float64]

3x3 tensor for each input index. Leading dimensions are preserved.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/utils/voigt.py
def voigt_to_tensor(vector: NDArray[np.float64]) -> NDArray[np.float64]:
    """Convert Voigt vectors to symmetric 3x3 tensors.

    Args:
        vector: Array of shape (..., 6) where the last dimension contains the
            stress/strain components in order according to Voigt notation:
                (σ_11, σ_22, σ_33, σ_23, σ_13, σ_12)
                (σ_xx, σ_yy, σ_zz, σ_yz, σ_xz, σ_xy)

    Returns:
        Array with shape (..., 3, 3). The last two axes are the symmetric
        3x3 tensor for each input index. Leading dimensions are preserved.

    Raises:
        ValueError: If the last dimension is not of size 6.
    """
    check_shape(vector)

    # (1e6, 50, 8, 6) -> (1e6, 50, 8) + (3, 3) -> (1e6, 50, 8, 3, 3)
    shape = vector.shape[:-1] + (3, 3)
    tensor = np.zeros(shape, dtype=vector.dtype)

    # Normal components
    tensor[(..., 0, 0)] = vector[..., 0]  # xx
    tensor[(..., 1, 1)] = vector[..., 1]  # yy
    tensor[(..., 2, 2)] = vector[..., 2]  # zz

    # Shear components
    tensor[(..., [1, 2], [2, 1])] = vector[..., [3]]  # yz
    tensor[(..., [0, 2], [2, 0])] = vector[..., [4]]  # xz
    tensor[(..., [0, 1], [1, 0])] = vector[..., [5]]  # xy

    return tensor

tensor_to_voigt(tensor)

Convert symmetric 3x3 tensors to Voigt vectors.

Parameters:

Name Type Description Default
tensor NDArray[float64]

Array of shape (..., 3, 3) where the last two dimensions represent symmetric 3x3 tensors.

required

Returns:

Type Description
NDArray[float64]

Array with shape (..., 6). The last dimension contains the stress/strain

NDArray[float64]

components in order according to Voigt notation: (σ_11, σ_22, σ_33, σ_23, σ_13, σ_12) (σ_xx, σ_yy, σ_zz, σ_yz, σ_xz, σ_xy)

NDArray[float64]

Leading dimensions are preserved.

Raises:

Type Description
ValueError

If the last two dimensions are not of size (3, 3).

Source code in src/fatpy/utils/voigt.py
def tensor_to_voigt(tensor: NDArray[np.float64]) -> NDArray[np.float64]:
    """Convert symmetric 3x3 tensors to Voigt vectors.

    Args:
        tensor: Array of shape (..., 3, 3) where the last two dimensions
            represent symmetric 3x3 tensors.

    Returns:
        Array with shape (..., 6). The last dimension contains the stress/strain
        components in order according to Voigt notation:
            (σ_11, σ_22, σ_33, σ_23, σ_13, σ_12)
            (σ_xx, σ_yy, σ_zz, σ_yz, σ_xz, σ_xy)
        Leading dimensions are preserved.

    Raises:
        ValueError: If the last two dimensions are not of size (3, 3).
    """
    if tensor.shape[-2:] != (3, 3):
        raise ValueError("Last two dimensions must correspond to (3, 3) tensors.")

    shape = tensor.shape[:-2] + (VOIGT_COMPONENTS_COUNT,)
    array = np.zeros(shape, dtype=tensor.dtype)

    # Normal components
    array[..., 0] = tensor[..., 0, 0]  # xx
    array[..., 1] = tensor[..., 1, 1]  # yy
    array[..., 2] = tensor[..., 2, 2]  # zz

    # Shear components
    array[..., 3] = tensor[..., 1, 2]  # yz
    array[..., 4] = tensor[..., 0, 2]  # xz
    array[..., 5] = tensor[..., 0, 1]  # xy

    return array