cratepy.clustering.rveelasticdatabase.FFTBasicScheme

class FFTBasicScheme(strain_formulation, problem_type, rve_dims, n_voxels_dims, regular_grid, material_phases, material_phases_properties)[source]

Bases: DNSHomogenizationMethod

FFT-based homogenization basic scheme.

FFT-based homogenization basic scheme proposed by Moulinec and Suquet (1998) [1] for the solution of micro-scale equilibrium problems of linear elastic heterogeneous materials. In particular, for a given RVE discretized in a regular grid of voxels, the method solves the microscale equilibrium problem when the RVE is subjected to a macroscale strain tensor and is constrained by periodic boundary conditions. Finite strain extension is also available, following the formulation of Kabel and coworkers (2022) [2] and adopting the Hencky hyperelastic constitutive model.

A detailed description of the computational implementation can be found in Appendix B (infinitesimal strains) and Section 4.6 (finite strains) of Ferreira (2022) [3].

_n_dim

Problem number of spatial dimensions.

Type:

int

_comp_order_sym

Strain/Stress components symmetric order.

Type:

list[str]

_comp_order_nsym

Strain/Stress components nonsymmetric order.

Type:

list[str]

_max_n_iterations

Maximum number of iterations to convergence.

Type:

int

_conv_criterion

Convergence criterion: ‘stress_div’ is the original convergence criterion based on the evaluation of the divergence of the stress tensor; ‘avg_stress_norm’ is based on the iterative change of the average stress norm.

Type:

{‘stress_div’, ‘avg_stress_norm’}

_conv_tol

Convergence tolerance.

Type:

float

_max_subinc_level

Maximum level of macroscale loading subincrementation.

Type:

int

_max_cinc_cuts

Maximum number of consecutive macroscale loading increment cuts.

Type:

int

_hom_stress_strain

Homogenized stress-strain material response. The homogenized strain and homogenized stress tensor components of the i-th loading increment are stored columnwise in the i-th row, sorted respectively. Infinitesimal strain tensor and Cauchy stress tensor (infinitesimal strains) or Deformation gradient and first Piola-Kirchhoff stress tensor (finite strains).

Type:

numpy.ndarray (2d)

compute_rve_local_response(self, mac_strain_id, mac_strain, verbose=False)[source]

Compute RVE local elastic strain response.

_elastic_constitutive_model(self, strain_vox, evar1, evar2, evar3, finite_strains_model='stvenant-kirchhoff', is_optimized=True)[source]

Elastic or hyperelastic material constitutive model.

stress_div_conv_criterion(self, freqs_dims, stress_DFT_vox)

Convergence criterion based on the divergence of the stress tensor.

compute_avg_state_vox(self, state_vox)

Compute average norm of strain or stress local field.

_compute_homogenized_field(self, state_vox)[source]

Perform homogenization over regular grid spatial discretization.

get_hom_stress_strain(self)[source]

Get homogenized strain-stress material response.

_display_greetings()[source]

Output greetings.

_display_increment_init(inc, subinc_level, total_lfact, inc_lfact)[source]

Output increment initial data.

_display_increment_end(strain_formulation, hom_strain, hom_stress, inc_time, total_time)[source]

Output increment end data.

_display_iteration(iter, iter_time, discrete_error)[source]

Output iteration data.

_display_increment_cut(cut_msg='')[source]

Output increment cut data.

Constructor.

Parameters:
  • strain_formulation ({'infinitesimal', 'finite'}) – Problem strain formulation.

  • problem_type (int) – Problem type: 2D plane strain (1), 2D plane stress (2), 2D axisymmetric (3) and 3D (4).

  • rve_dims (list[float]) – RVE size in each dimension.

  • n_voxels_dims (list[int]) – Number of voxels in each dimension of the regular grid (spatial discretization of the RVE).

  • regular_grid (numpy.ndarray (2d or 3d)) – Regular grid of voxels (spatial discretization of the RVE), where each entry contains the material phase label (int) assigned to the corresponding voxel.

  • material_phases (list[str]) – RVE material phases labels (str).

  • material_phases_properties (dict) – Constitutive model material properties (item, dict) associated to each material phase (key, str).

List of Public Methods

compute_rve_local_response

Compute RVE local elastic strain response.

get_hom_stress_strain

Get homogenized strain-stress material response.

Methods

__init__(strain_formulation, problem_type, rve_dims, n_voxels_dims, regular_grid, material_phases, material_phases_properties)[source]

Constructor.

Parameters:
  • strain_formulation ({'infinitesimal', 'finite'}) – Problem strain formulation.

  • problem_type (int) – Problem type: 2D plane strain (1), 2D plane stress (2), 2D axisymmetric (3) and 3D (4).

  • rve_dims (list[float]) – RVE size in each dimension.

  • n_voxels_dims (list[int]) – Number of voxels in each dimension of the regular grid (spatial discretization of the RVE).

  • regular_grid (numpy.ndarray (2d or 3d)) – Regular grid of voxels (spatial discretization of the RVE), where each entry contains the material phase label (int) assigned to the corresponding voxel.

  • material_phases (list[str]) – RVE material phases labels (str).

  • material_phases_properties (dict) – Constitutive model material properties (item, dict) associated to each material phase (key, str).

_compute_avg_state_vox(state_vox)[source]

Compute average norm of strain or stress local field.

Parameters:

state_vox (dict) – Local strain or stress response (item, numpy.ndarray of shape equal to RVE regular grid discretization) for each strain or stress component (key, str).

Returns:

avg_state_norm – Average norm of strain or stress local field.

Return type:

float

_compute_homogenized_field(state_vox)[source]

Perform homogenization over regular grid spatial discretization.

Parameters:

state_vox (dict) – Local strain or stress response (item, numpy.ndarray of shape equal to RVE regular grid discretization) for each strain or stress component (key, str).

Returns:

hom_state – Homogenized strain or stress tensor.

Return type:

numpy.ndarray (2d)

static _display_greetings()[source]

Output greetings.

static _display_increment_cut(cut_msg='')[source]

Output increment cut data.

Parameters:

cut_msg (str) – Increment cut output message.

static _display_increment_end(strain_formulation, hom_strain, hom_stress, inc_time, total_time)[source]

Output increment end data.

Parameters:
  • strain_formulation ({'infinitesimal', 'finite'}) – Problem strain formulation.

  • hom_strain (numpy.ndarray (2d)) – Homogenized strain tensor.

  • hom_stress (numpy.ndarray (2d)) – Homogenized stress tensor.

  • inc_time (float) – Increment running time.

  • total_time (float) – Total running time.

static _display_increment_init(inc, subinc_level, total_lfact, inc_lfact)[source]

Output increment initial data.

Parameters:
  • inc (int) – Increment number.

  • subinc_level (int) – Subincrementation level.

  • total_lfact (float) – Total load factor.

  • inc_lfact (float) – Incremental load factor.

static _display_iteration(iter, iter_time, discrete_error)[source]

Output iteration data.

Parameters:
  • iter (int) – Iteration number.

  • iter_time (float) – Iteration running time.

  • discrete_error (float) – Discrete error associated to convergence criterion.

_elastic_constitutive_model(strain_vox, evar1, evar2, evar3, finite_strains_model='stvenant-kirchhoff', is_optimized=True)[source]

Elastic or hyperelastic material constitutive model.

Available constitutive models:

Infinitesimal strains:

  • Isotropic linear elastic constitutive model:

    \[\boldsymbol{\sigma} = \boldsymbol{\mathsf{D}}^{e} : \boldsymbol{\varepsilon}\]

    where \(\boldsymbol{\sigma}\) is the Cauchyevar1 stress tensor, \(\boldsymbol{\mathsf{D}}^{e}\) is the elasticity tensor, and \(\boldsymbol{\varepsilon}\) is the infinitesimal strain tensor.


Finite strains:

  • Hencky hyperelastic (isotropic) constitutive model:

    \[\boldsymbol{\tau} = \boldsymbol{\mathsf{D}}^{e} : \boldsymbol{\varepsilon}\]

    where \(\boldsymbol{\tau}\) is the Kirchhoff stress tensor, \(\boldsymbol{\mathsf{D}}^{e}\) is the elasticity tensor, and \(\boldsymbol{\varepsilon}\) is the spatial logarithmic strain tensor.

  • St.Venant-Kirchhoff hyperelastic (isotropic) constitutive model:

    \[\boldsymbol{S} = \boldsymbol{\mathsf{D}}^{e} : \boldsymbol{E}^{(2)}\]

    where \(\boldsymbol{S}\) is the second Piola-Kirchhoff stress tensor, \(\boldsymbol{\mathsf{D}}^{e}\) is the elasticity tensor, and \(\boldsymbol{E}^{(2)}\) is the Green-Lagrange strain tensor.

A detailed description of the computational implementation based on Hadamard (element-wise) operations can be found in Section 4.6 of Ferreira (2022) [4].


Parameters:
  • strain_vox (dict) – Local strain response (item, numpy.ndarray of shape equal to RVE regular grid discretization) for each strain component (key, str). Infinitesimal strain tensor (infinitesimal strains) or deformation gradient (finite strains).

  • evar1 (numpy.ndarray (2d or 3d)) –

    Auxiliar elastic properties array (numpy.ndarray of shape equal to RVE regular grid discretization) containing an elastic properties-related quantity associated to each voxel,

    \[\dfrac{E \nu}{(1 + \nu)(1 - 2 \nu)} \, ,\]

    where \(E\) and \(\nu\) are the Young’s Modulus and Poisson’s ratio, respectively.

  • evar2 (numpy.ndarray (2d or 3d)) –

    Auxiliar elastic properties array (numpy.ndarray of shape equal to RVE regular grid discretization) containing an elastic properties-related quantity associated to each voxel,

    \[\dfrac{E}{1 + \nu} \, ,\]

    where \(E\) and \(\nu\) are the Young’s Modulus and Poisson’s ratio, respectively.

  • evar3 (numpy.ndarray (2d or 3d)) –

    Auxiliar elastic properties array (numpy.ndarray of shape equal to RVE regular grid discretization) containing an elastic properties-related quantity associated to each voxel,

    \[\dfrac{E \nu}{(1 + \nu)(1 - 2 \nu)} - \dfrac{E}{1 + \nu} \, ,\]

    where \(E\) and \(\nu\) are the Young’s Modulus and Poisson’s ratio, respectively.

  • finite_strains_model ({'hencky', 'stvenant-kirchhoff'}, default='hencky') – Finite strains hyperelastic isotropic constitutive model.

  • is_optimized (bool) – Optimization flag (minimizes loops over spatial discretization voxels).

Returns:

stress_vox – Local stress response (item, numpy.ndarray of shape equal to RVE regular grid discretization) for each stress component (key, str). Cauchy stress tensor (infinitesimal strains) or first Piola-Kirchhoff stress tensor (finite strains).

Return type:

dict

_stress_div_conv_criterion(freqs_dims, stress_DFT_vox)[source]

Convergence criterion based on the divergence of the stress tensor.

Convergence criterion proposed by Moulinec and Suquet (1998) [5].


Parameters:
  • freqs_dims (list[numpy.ndarray (1d)]) – List of discrete frequencies (numpy.ndarray (1d)) associated to each spatial dimension.

  • stress_DFT_vox (dict) – Discrete Fourier Transform of local stress response (item, numpy.ndarray of shape equal to RVE regular grid discretization) for each stress component (key, str). Cauchy stress tensor (infinitesimal strains) or First Piola-Kirchhoff stress tensor (finite strains).

Returns:

discrete_error – Discrete error associated to the convergence criterion.

Return type:

float

compute_rve_local_response(mac_strain_id, mac_strain, verbose=False)[source]

Compute RVE local elastic strain response.

Compute the RVE local strain response (solution of microscale equilibrium problem) when subjected to a given macroscale strain loading, namely a macroscale infinitesimal strain tensor (infinitesimal strains) or a macroscale deformation gradient (finite strains). It is assumed that the RVE is spatially discretized in a regular grid of voxels.


Parameters:
  • mac_strain_id (int) – Macroscale strain second-order tensor identifier.

  • mac_strain (numpy.ndarray (2d)) – Macroscale strain second-order tensor. Infinitesimal strain tensor (infinitesimal strains) or deformation gradient (finite strains).

  • verbose (bool, default=False) – Enable verbose output.

Returns:

strain_vox – RVE local strain response (item, numpy.ndarray of shape equal to RVE regular grid discretization) for each strain component (key, str). Infinitesimal strain tensor (infinitesimal strains) or material logarithmic strain tensor (finite strains).

Return type:

dict

get_hom_stress_strain()[source]

Get homogenized strain-stress material response.

Returns:

_hom_stress_strain – RVE homogenized stress-strain response (item, numpy.ndarray (2d)) for each macroscale strain loading identifier (key, int). The homogenized strain and homogenized stress tensor components of the i-th loading increment are stored columnwise in the i-th row, sorted respectively. Infinitesimal strain tensor and Cauchy stress tensor (infinitesimal strains) or Deformation gradient and first Piola-Kirchhoff stress tensor (finite strains).

Return type:

numpy.ndarray (2d)