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.

The module implements various fatigue damage accumulation rules used to estimate the total damage under variable amplitude loading. They include both classical linear models, such as Palmgren–Miner’s rule, and more advanced non-linear approaches that account for load sequence effects, mean stress influence, and material-specific behavior for improved fatigue life prediction.

decompositions

Decompositions module.

Methods for breaking down complex load signals to cycles. Contains both the uniaxial and the multiaxial procedures.

multiaxial

Multiaxial decompositions methods.

These routines provide tools for decomposing multiaxial stress and strain states into their individual components, such as normal and shear parts, principal values, and others. They are essential for analyzing complex loading conditions, enabling accurate evaluation of material behavior under multiaxial fatigue.

load_path_analysis_5d

Load Path Analysis 5D methods of multiaxial decompositions.

These methods compute the amplitude (J_2,a) and mean (J_2,m) values of the second invariant of the deviatoric stress tensor in five-dimensional deviatoric space. This evaluation is particularly useful in advanced multiaxial fatigue models, where accurate representation of cyclic loading paths and stress states in the deviatoric space is critical for predicting material response.

shear_path_analysis_2d

Shear Path Analysis 2D methods of multiaxial decompositions.

These methods are designed to evaluate the mean shear stress and shear stress amplitude acting on a specified material plane, typically used in multiaxial fatigue analysis. By projecting the stress tensor onto the plane of interest, they provide key parameters for assessing crack initiation risks and fatigue life under complex loading paths.

uniaxial

Uniaxial decompositions methods.

Implements algorithms designed to decompose one-dimensional signals into individual load cycles, enabling detailed analysis of variable amplitude loading. By identifying turning points and extracting cycles based on established criteria (e.g., rainflow counting), they facilitate fatigue analysis, damage assessment, and other time-domain signal evaluations.

energy_life

Energy-based fatigue life analysis methods.

Fatigue analysis methods based on the relationship between the strain energy (usually strain energy density) and number of cycles to failure.

base_methods

Base methods for the energy-life.

correction_methods

Correction methods for the energy-life approach.

decompositions

Decomposition methods for the energy-life.

demage_parameters

Damage parameters calculation methods for the energy-life.

plane_based_methods

Plane-based methods module.

Methods for processing stress tensor path on a material plane. Provides basic infrastructure for prediction methods based on critical-plane and integral approaches.

Plane Search methods.

Algorithms to identify a subset of potentially critical planes for further detailed analysis and processing by critical-plane and integral methods.

strain_life

Strain-based fatigue life analysis methods.

These methods perform fatigue analysis using the relationship between strain amplitude and the number of cycles to failure, as expressed by strain-life (ε-N) approaches such as the Coffin-Manson and Basquin laws. They are particularly suited for low-cycle and transitional fatigue regimes, where plastic strain plays a significant role, enabling life prediction under variable amplitude or multiaxial loading when combined with appropriate mean stress corrections and cycle counting techniques.

base_methods

Base methods for the strain-life.

correction_methods

Correction methods for the strain-life approach.

Methods to estimate notch stress/strain with local yielding and apply mean stress correction.

elastic_plastic_conversion

Elastic-Plastic Conversion methods for Strain-Life analysis.

These methods convert purely elastic stress or strain states into their elastic-plastic equivalents to account for local yielding effects in fatigue or fracture assessments. Using approaches such as Neuber's rule or Glinka's strain energy density method, they provide corrected stress-strain values that better represent actual material behavior under high local stresses, especially at notches or other stress concentrators.

mean_stress_effect

Mean Stress Effect correction methods for Strain-Life Analysis.

The methods for incorporating the effect of mean stress in strain-life (ε-N) analysis.

damage_params

Damage parameters calculation methods for the strain-life approach.

This module includes both uniaxial and multiaxial parameters for low-cycle fatigue assessment, covering a range of strain-life-based damage models suitable for high plastic strain regimes. It also incorporates methods based on the FKM Guideline, enabling standardized evaluation of fatigue strength under complex loading conditions.

fkm_nonlin

fkm non-linear damage parameter calculations for strain-life analysis.

The module contains damage parameters based on FKM-Guideline Non-Linear.

multiaxial

Multiaxial fatigue criteria methods for the strain-life approach.

These methods perform low-cycle fatigue assessment using multiaxial damage parameters that combine stress, strain, and mean stress effects into a single fatigue-relevant quantity. By incorporating critical plane, invariant-based, and energy-based criteria, they enable life prediction for components subjected to complex non-proportional loading paths where significant plastic deformation occurs.

uniaxial

Uniaxial fatigue criteria methods for the strain-life approach.

These methods compute uniaxial damage parameters used in strain-life (ε-N) analysis to correlate cyclic strain amplitudes with fatigue life. They include parameters such as Coffin-Manson, Morrow, and Smith-Watson-Topper (SWT), which combine elastic and plastic strain components—and, when applicable, mean stress effects—into a single scalar value for life prediction under uniaxial loading conditions.

stress_life

Stress-based fatigue life analysis methods.

Stress-based fatigue assessment methods. Contains approach to estimate fatigue damage and their correction on various fatigue effects (stress concentration, size, surface quality, etc.)

base_methods

Base methods for the stress-life.

correction_methods

Correction methods for the stress-life.

These methods project how various factors—such as surface finish, size effect, mean stress, residual stresses, temperature, manufacturing defects, and notch geometry—modify the fatigue strength of a real component by quantifying shifts in the S-N curve or fatigue limit.

eq_stress_at_crit_dist

Equivalent stress at critical distance correction methods for stress-life approach.

These modules implement various Theory of Critical Distances (TCD)-based approaches for stress-gradient correction, which improve fatigue strength predictions near stress concentrators such as notches or cracks. By incorporating characteristic length parameters, the methods account for the spatial distribution of stress and mitigate the over-conservatism of local stress-based criteria, offering more accurate fatigue assessments under both uniaxial and multiaxial loading conditions.

fatigue_limit

Fatigue limit correction methods for the stress-life approach.

These routines calculate correction factors that modify the fatigue limit in an S-N curve to account for real-world influences such as surface finish, size effect, mean stress, residual stresses, temperature, and notches. By combining these empirical correction models, adjusted endurance limits may be produced for component-level fatigue prediction.

power_law

Power law correction methods for the stress-life approach.

These methods are used to project the influence of various factors—such as surface finish, size effect, mean stress, or residual stresses—on the slope of the S-N curve for a real component. The fatigue behavior is modeled using a power law formulation.

damage_params

Damage parameters calculation methods for the stress-life approach.

These methods provide parameters for evaluating stress-based fatigue damage under uniaxial or multiaxial loading conditions. They include critical plane and invariant-based approaches, such as maximum principal stress, Findley's and Dang Van's criteria, which relate the applied stress state to fatigue life by identifying the most damaging stress components or orientations.

multiaxial

Multiaxial fatigue criteria methods for the stress-life approach.

These modules implement a range of criteria for assessing fatigue strength under multiaxial stress states, addressing both proportional and non-proportional loading paths. The implemented models include critical plane approaches, energy-based methods, and stress invariants, enabling prediction of fatigue limits/strengths for complex loading scenarios across various materials and geometries.

uniaxial

Uniaxial fatigue criteria methods for the stress-life approach.

Contains criteria that address uniaxial high-cycle fatigue by incorporating the mean stress effect through an equivalent stress amplitude approach. By adjusting the stress amplitude to account for mean stress influences—using models such as Goodman, Gerber, or Soderberg—they enable more accurate fatigue life predictions where mean stresses significantly affect material endurance.

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 classes and functions that define material laws and enable their calibration for the material under consideration.

cyclic_stress_strain_curve

Cyclic stress–strain curve approximation models.

Implements various formulations, such as the Ramberg–Osgood model, describing the relationship between elastic–plastic stress and strain.

en_curve

Strain-life curve approximations.

Implements strain-life (ε-N) curve models such as the Manson-Coffin and Basquin relations, and provides conversions between strain amplitude and fatigue life in both directions.

hookes_law

Hooke's law.

Converts between elastic stress and strain with varying levels of complexity depending on dimensionality.

regression_analysis

Regression analysis methods of material laws.

Implements linear and non-linear regression methods for calibrating material laws.

sn_curve

Stress-life curve methods of material laws.

Provides implementations of Wöhler (S-N) curve models along with methods for converting between stress amplitude and fatigue life in both directions.

struct_mech

Structural mechanics module.

Contains structural mechanics related classes and methods.

strain

Calculate fundamental strain metrics and invariants.

These functions provide principal strains, hydrostatic strain, von Mises equivalent strain, and invariants of the strain 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 strains 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_strains_and_directions(strain_vector_voigt)

Calculate principal strains and principal directions for each state.

Math Equations

Principal strains and directions are found by solving the eigenvalue problem for the strain tensor:

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

Parameters:

Name Type Description Default
strain_vector_voigt NDArray[float64]

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

required

Returns:

Name Type Description
Tuple (eigvals, eigvecs)
  • eigvals: Array of shape (..., 3). Principal strains (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:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/strain.py
def calc_principal_strains_and_directions(
    strain_vector_voigt: NDArray[np.float64],
) -> tuple[NDArray[np.float64], NDArray[np.float64]]:
    r"""Calculate principal strains and principal directions for each state.

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

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

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

    Returns:
        Tuple (eigvals, eigvecs):
            - eigvals: Array of shape (..., 3). Principal strains
            (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(strain_vector_voigt)

    tensor = voigt.voigt_to_tensor(strain_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_strains(strain_vector_voigt)

Calculate principal strains for each strain state.

Math Equations

Principal strains are found by solving the eigenvalue problem for the strain tensor:

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

Parameters:

Name Type Description Default
strain_vector_voigt NDArray[float64]

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

required

Returns:

Type Description
NDArray[float64]

Array of shape (..., 3). Principal strains (descending: ε1 ≥ ε2 ≥ ε3).

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/strain.py
def calc_principal_strains(
    strain_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate principal strains for each strain state.

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

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

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

    Returns:
        Array of shape (..., 3). Principal strains (descending: ε1 ≥ ε2 ≥ ε3).

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

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

    return np.sort(eigvals, axis=-1)[..., ::-1]

calc_strain_invariants(strain_vector_voigt)

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

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

Parameters:

Name Type Description Default
strain_vector_voigt NDArray[float64]

Array of shape (..., 6). The last dimension contains the Voigt strain 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/strain.py
def calc_strain_invariants(
    strain_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate the first, second, and third invariants for each strain state.

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

    Args:
        strain_vector_voigt: Array of shape (..., 6). The last dimension contains the
            Voigt strain 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(strain_vector_voigt)

    tensor = voigt.voigt_to_tensor(strain_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_volumetric_strain(strain_vector_voigt)

Calculate the volumetric (mean normal)strain for each strain state.

Math Equations
\[ \varepsilon_{vol} = \frac{1}{3} \, tr(\varepsilon) = \frac{1}{3}(\varepsilon_{11} + \varepsilon_{22} + \varepsilon_{33}) \]

Parameters:

Name Type Description Default
strain_vector_voigt NDArray[float64]

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

required

Returns:

Type Description
NDArray[float64]

Array of shape (...). Volumetric (mean normal) strain 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/strain.py
def calc_volumetric_strain(
    strain_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate the volumetric (mean normal)strain for each strain state.

    ??? abstract "Math Equations"
        $$ \varepsilon_{vol} = \frac{1}{3} \, tr(\varepsilon) =
           \frac{1}{3}(\varepsilon_{11} + \varepsilon_{22} + \varepsilon_{33})
        $$

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

    Returns:
        Array of shape (...). Volumetric (mean normal) strain for each input state.
            Tensor rank is reduced by one.

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

    return (
        strain_vector_voigt[..., 0]
        + strain_vector_voigt[..., 1]
        + strain_vector_voigt[..., 2]
    ) / 3.0

calc_deviatoric_strain(strain_vector_voigt)

Calculate the deviatoric strain for each strain state.

Math Equations

The strain tensor decomposes as:

\[ \varepsilon = \varepsilon_{dev} + \varepsilon_{vol} \mathbf{I} \]

where the deviatoric part is traceless and obtained by subtracting the volumetric part from the normal components.

\[ \varepsilon_{dev} = \varepsilon - \frac{1}{3} tr(\varepsilon) \]

Parameters:

Name Type Description Default
strain_vector_voigt NDArray[float64]

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

required

Returns:

Type Description
NDArray[float64]

Array of shape (..., 6). Deviatoric strain for each input state.

Raises:

Type Description
ValueError

If the last dimension is not of size 6.

Source code in src/fatpy/struct_mech/strain.py
def calc_deviatoric_strain(
    strain_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Calculate the deviatoric strain for each strain state.

    ??? abstract "Math Equations"
        The strain tensor decomposes as:

        $$ \varepsilon = \varepsilon_{dev} + \varepsilon_{vol} \mathbf{I} $$

        where the deviatoric part is traceless and obtained by subtracting the
        volumetric part from the normal components.

        $$ \varepsilon_{dev} = \varepsilon - \frac{1}{3} tr(\varepsilon) $$

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

    Returns:
        Array of shape (..., 6). Deviatoric strain for each input state.

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

    volumetric = calc_volumetric_strain(strain_vector_voigt)
    deviatoric = strain_vector_voigt.copy()
    deviatoric[..., :3] = deviatoric[..., :3] - volumetric[..., None]

    return deviatoric

calc_von_mises_strain(strain_vector_voigt)

Von Mises equivalent strain computed directly from Voigt components.

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

Parameters:

Name Type Description Default
strain_vector_voigt NDArray[float64]

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

required

Returns:

Type Description
NDArray[float64]

Array of shape (...). Von Mises equivalent strain 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/strain.py
def calc_von_mises_strain(
    strain_vector_voigt: NDArray[np.float64],
) -> NDArray[np.float64]:
    r"""Von Mises equivalent strain computed directly from Voigt components.

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

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

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

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

    e11 = strain_vector_voigt[..., 0]
    e22 = strain_vector_voigt[..., 1]
    e33 = strain_vector_voigt[..., 2]
    e23 = strain_vector_voigt[..., 3]  # epsilon_23
    e13 = strain_vector_voigt[..., 4]  # epsilon_13
    e12 = strain_vector_voigt[..., 5]  # epsilon_12
    return np.sqrt(
        (2.0 / 9.0)
        * (
            (e11 - e22) ** 2
            + (e22 - e33) ** 2
            + (e33 - e11) ** 2
            + 6.0 * (e12**2 + e23**2 + e13**2)
        )
    )

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

Calculate signed von Mises equivalent strain for each strain state.

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

Sign Convention

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
\[ \varepsilon_{SvM} = \begin{cases} +\varepsilon_{vM} & \text{if } \frac{\varepsilon_1 + \varepsilon_3}{2} \geq 0 \\ -\varepsilon_{vM} & \text{if } \frac{\varepsilon_1 + \varepsilon_3}{2} < 0 \end{cases} \]

Parameters:

Name Type Description Default
strain_vector_voigt NDArray[float64]

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

required
rtol float

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

1e-05
atol float

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

1e-08

Returns:

Type Description
NDArray[float64]

Array of shape (...). Signed von Mises equivalent strain 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/strain.py
def calc_signed_von_mises_by_max_abs_principal(
    strain_vector_voigt: NDArray[np.float64],
    rtol: float = 1e-5,
    atol: float = 1e-8,
) -> NDArray[np.float64]:
    r"""Calculate signed von Mises equivalent strain for each strain state.

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

    ??? note "Sign Convention"
        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"
        $$
        \varepsilon_{SvM} = \begin{cases}
        +\varepsilon_{vM} & \text{if } \frac{\varepsilon_1 + \varepsilon_3}{2} \geq 0 \\
        -\varepsilon_{vM} & \text{if } \frac{\varepsilon_1 + \varepsilon_3}{2} < 0
        \end{cases}
        $$

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

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

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

    von_mises = calc_von_mises_strain(strain_vector_voigt)
    principals = calc_principal_strains(strain_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 * von_mises

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)
  • 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:

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

The stress tensor decomposes as:

\[ \sigma = \mathbf{s} + \sigma_H \mathbf{I} \]

where the deviatoric part \(\mathbf{s}\) is traceless and obtained by subtracting the hydrostatic part from the normal components.

\[ \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"
        The stress tensor decomposes as:

        $$ \sigma = \mathbf{s} + \sigma_H \mathbf{I} $$

        where the deviatoric part $\mathbf{s}$ is traceless and obtained by subtracting
        the hydrostatic part from the normal components.

        $$ \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.

Sign Convention

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.

    ??? note "Sign Convention"
        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_vector_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.

Sign Convention

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_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 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_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 average of the maximum and minimum principal stresses.

    ??? note "Sign Convention"
        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_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 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)
    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 * 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.

Sign Convention

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.

    ??? note "Sign Convention"
        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.

Sign Convention

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.

    ??? note "Sign Convention"
        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.

Sign Convention

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.

    ??? note "Sign Convention"
        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

Transformations methods of structural mechanics.

These routines perform mathematical operations to transform stress and strain tensors between different coordinate systems or reference frames. They include tensor rotation, transformation to principal axes, and projection onto specified planes, enabling detailed analysis of material response under complex loading conditions and supporting criteria based on critical plane, invariant, or principal stress/strain approaches. Includes 5D deviatoric stress reduction.

utils

Utilities module.

A collection of general utility functions.

mesh

Mesh processing utilities module.

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

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