cratepy.ioput.readinputdata.MaterialState

class MaterialState(strain_formulation, problem_type, material_phases, material_phases_properties, material_phases_vf)[source]

Bases: object

CRVE material constitutive state.

_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]

_material_phases_models

Material constitutive model (item, ConstitutiveModel) associated to each material phase (key, str).

Type:

dict

_phase_clusters

Clusters labels (item, list[int]) associated to each material phase (key, str).

Type:

dict

_clusters_vf

Volume fraction (item, float) associated to each material cluster (key, str).

Type:

dict

_clusters_def_gradient_mf

Deformation gradient (item, numpy.darray (1d)) associated to each material cluster (key, str), stored in matricial form.

Type:

dict

_clusters_def_gradient_old_mf

Last converged deformation gradient (item, numpy.darray (1d)) associated to each material cluster (key, str), stored in matricial form.

Type:

dict

_clusters_state

Material constitutive model state variables (item, dict) associated to each material cluster (key, str).

Type:

dict

_clusters_state_old

Last converged material constitutive model state variables (item, dict) associated to each material cluster (key, str).

Type:

dict

_clusters_tangent_mf

Material consistent tangent modulus (item, numpy.ndarray) associated to each material cluster (key, str), stored in matricial form.

Type:

dict

_hom_strain_mf

Homogenized strain tensor stored in matricial form: infinitesimal strain tensor (infinitesimal strains) or deformation gradient (finite strains).

Type:

numpy.ndarray (1d)

_hom_strain_33

Homogenized strain tensor out-of-plane component: infinitesimal strain tensor (infinitesimal strains) or deformation gradient (finite strains).

Type:

float

_hom_stress_mf

Homogenized stress tensor stored in matricial form: Cauchy stress tensor (infinitesimal strains) or first Piola-Kirchhoff stress tensor (finite strains).

Type:

numpy.ndarray (1d)

_hom_stress_33

Homogenized stress tensor out-of-plane component: Cauchy stress tensor (infinitesimal strains) or first Piola-Kirchhoff stress tensor (finite strains).

Type:

float

_hom_strain_old_mf

Last converged homogenized strain tensor stored in matricial form: infinitesimal strain tensor (infinitesimal strains) or deformation gradient (finite strains).

Type:

numpy.ndarray (1d)

_hom_stress_old_mf

Last converged homogenized stress tensor stored in matricial form: Cauchy stress tensor (infinitesimal strains) or first Piola-Kirchhoff stress tensor (finite strains).

Type:

numpy.ndarray (1d)

init_constitutive_model(self, mat_phase, model_keyword, model_source='crate')[source]

Initialize material phase constitutive model.

init_clusters_state(self)[source]

Initialize clusters state variables.

set_phase_clusters(self, phase_clusters, clusters_vf)[source]

Set CRVE cluster labels and volume fractions of each material phase.

get_clusters_inc_strain_mf(self, global_strain_mf)[source]

Get clusters incremental strain in matricial form.

update_clusters_state(self, clusters_inc_strain_mf)[source]

Update clusters state variables and consistent tangent modulus.

update_state_homogenization(self)[source]

Update homogenized strain and stress tensors.

update_converged_state(self)[source]

Update last converged material state variables.

set_rewind_state_updated_clustering(self, phase_clusters, clusters_vf, clusters_state, clusters_def_gradient_mf)[source]

Set rewind state variables according to updated clustering.

get_hom_strain_mf(self)[source]

Get homogenized strain tensor (matricial form).

get_hom_stress_mf(self)[source]

Get homogenized stress tensor (matricial form).

get_oop_hom_comp(self)[source]

Get homogenized strain or stress tensor out-of-plane component.

get_inc_hom_strain_mf(self)[source]

Get incremental homogenized strain tensor (matricial form).

get_inc_hom_stress_mf(self)[source]

Get incremental homogenized stress tensor (matricial form).

get_hom_strain_old_mf(self)[source]

Get last converged homogenized strain tensor (matricial form).

get_hom_stress_old_mf(self)[source]

Get last converged homogenized stress tensor (matricial form).

get_material_phases(self)[source]

Get RVE material phases.

get_material_phases_properties(self)[source]

Get RVE material phases constitutive properties.

get_material_phases_models(self)[source]

Get RVE material phases constitutive models.

get_material_phases_vf(self)[source]

Get RVE material phases volume fraction.

get_clusters_def_gradient_mf(self)[source]

Get deformation gradient of each material cluster.

get_clusters_def_gradient_old_mf(self)[source]

Get last converged deformation gradient of each material cluster.

get_clusters_state(self)[source]

Get material state variables of each material cluster.

get_clusters_state_old(self)[source]

Get last converged material state variables of each material cluster.

get_clusters_tangent_mf(self)[source]

Get material consistent tangent modulus of each material cluster.

_material_su_interface(strain_formulation, problem_type, constitutive_model, def_gradient_old, inc_strain, state_variables_old)[source]

Material constitutive state update interface.

compute_inc_log_strain(e_log_strain_old, inc_def_gradient)[source]

Compute incremental spatial logarithmic strain.

compute_spatial_tangent_modulus(e_log_strain_old, def_gradient_old, inc_def_gradient, cauchy_stress, inf_consistent_tangent)[source]

Compute finite strain spatial consistent tangent modulus.

clustering_adaptivity_update(self, phase_clusters, clusters_vf, adaptive_clustering_map)[source]

Update cluster variables according to clustering adaptivity step.

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).

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

  • material_phases_vf (dict) – Volume fraction (item, float) associated to each material phase (key, str).

List of Public Methods

clustering_adaptivity_update

Update cluster variables according to clustering adaptivity step.

compute_inc_log_strain

Compute incremental spatial logarithmic strain.

compute_spatial_tangent_modulus

Compute finite strain spatial consistent tangent modulus.

get_clusters_def_gradient_mf

Get deformation gradient of each material cluster.

get_clusters_def_gradient_old_mf

Get last converged deformation gradient of each material cluster.

get_clusters_inc_strain_mf

Get clusters incremental strain in matricial form.

get_clusters_state

Get material state variables of each material cluster.

get_clusters_state_old

Get last converged material state variables of each mat.

get_clusters_tangent_mf

Get material consistent tangent modulus of each material cluster.

get_hom_strain_mf

Get homogenized strain tensor (matricial form).

get_hom_strain_old_mf

Get last converged homogenized strain tensor (matricial form).

get_hom_stress_mf

Get homogenized stress tensor (matricial form).

get_hom_stress_old_mf

Get last converged homogenized stress tensor (matricial form).

get_inc_hom_strain_mf

Get incremental homogenized strain tensor (matricial form).

get_inc_hom_stress_mf

Get incremental homogenized stress tensor (matricial form).

get_material_phases

Get RVE material phases.

get_material_phases_models

Get RVE material phases constitutive models.

get_material_phases_properties

Get RVE material phases constitutive properties.

get_material_phases_vf

Get RVE material phases volume fraction.

get_oop_hom_comp

Get homogenized strain or stress tensor out-of-plane component.

init_clusters_state

Initialize clusters state variables.

init_constitutive_model

Initialize material phase constitutive model.

set_phase_clusters

Set CRVE cluster labels and volume fractions of each material phase.

set_rewind_state_updated_clustering

Set rewind state variables according to updated clustering.

update_clusters_state

Update clusters state variables and consistent tangent modulus.

update_converged_state

Update last converged material state variables.

update_state_homogenization

Update homogenized strain and stress tensors.

Methods

__init__(strain_formulation, problem_type, material_phases, material_phases_properties, material_phases_vf)[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).

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

  • material_phases_vf (dict) – Volume fraction (item, float) associated to each material phase (key, str).

static _material_su_interface(strain_formulation, problem_type, constitutive_model, def_gradient_old, inc_strain, state_variables_old)[source]

Material constitutive state update interface.

This material constitutive state update interface contemplates three different families of constitutive models: (1) infinitesimal strains constitutive models, (2) finite strains constitutive models, and (3) isotropic hyperelastic-based finite strain constitutive models whose finite strain extension (from infinitesimal counterpart) is purely kinematical.

This interface is schematically illustrated in Figure 5.3 of Ferreira (2022) [1], and the last family of constitutive models is described on Appendix F.


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).

  • constitutive_model (ConstitutiveModel) – Material constitutive model.

  • def_gradient_old (numpy.ndarray (2d)) – Last converged deformation gradient.

  • inc_strain (numpy.ndarray (2d)) – Incremental strain second-order tensor.

  • state_variables_old (dict) – Last converged material constitutive model state variables.

Returns:

  • state_variables (dict) – Material constitutive model state variables.

  • consistent_tangent_mf (ndarray) – Material constitutive model material consistent tangent modulus in matricial form.

clustering_adaptivity_update(phase_clusters, clusters_vf, adaptive_clustering_map)[source]

Update cluster variables according to clustering adaptivity step.

Parameters:
  • phase_clusters (dict) – Clusters labels (item, list[int]) associated to each material phase (key, str).

  • clusters_vf (dict) – Volume fraction (item, float) associated to each material cluster (key, str).

  • adaptive_clustering_map (dict) – Adaptive clustering map (item, dict with list of new cluster labels (item, list[int]) resulting from the refinement of each target cluster (key, str)) for each material phase (key, str).

static compute_inc_log_strain(e_log_strain_old, inc_def_gradient)[source]

Compute incremental spatial logarithmic strain.

Incremental spatial logarithmic strain:

\[\Delta \boldsymbol{\varepsilon}_{n + 1} = \boldsymbol{\varepsilon}_{n + 1}^{e, \, \text{trial}} - \boldsymbol{\varepsilon}_{n}^{e}\]

where \(\Delta \boldsymbol{\varepsilon}\) is the incremental spatial logarithmic strain tensor, \(\boldsymbol{\varepsilon}^{e, \, \text{trial}}\) is the elastic trial spatial logarithmic strain tensor, \(\boldsymbol{\varepsilon}^{e}\) is the elastic spatial logarithmic strain tensor, \(n+1\) denotes the current increment, and \(n\) denotes the last converged increment.


Elastic trial left Cauchy-Green strain tensor:

\[\boldsymbol{\varepsilon}_{n + 1}^{e, \, \text{trial}} = \dfrac{1}{2} \ln ( \boldsymbol{B}^{e, \, \text{trial}}_{n+1} ) = \dfrac{1}{2} \ln \Big( (\boldsymbol{F}_{\Delta})_{n+1} \boldsymbol{B}^{e}_{n} (\boldsymbol{F}_{\Delta})_{n+1}^{T} \Big)\]

where \(\boldsymbol{\varepsilon}^{e, \, \text{trial}}\) is the elastic trial spatial logarithmic strain tensor, \(\boldsymbol{B}^{e, \, \text{trial}}\) is the elastic trial left Cauchy-Green strain tensor, \(\boldsymbol{F}_{\Delta}\) is the incremental deformation gradient, \(\boldsymbol{B}^{e}\) is the elastic left Cauchy-Green strain tensor, \(n+1\) denotes the current increment, and \(n\) denotes the last converged increment.

The definition of the elastic trial spatial logarithmic strain tensor can be found in Appendix F.4 of Ferreira (2022) [2] (see Equations (F.14) and (F.15)).


Parameters:
  • e_log_strain_old (numpy.ndarray (2d)) – Last converged elastic spatial logarithmic strain tensor.

  • inc_def_gradient (numpy.ndarray (2d)) – Incremental deformation gradient.

Returns:

inc_log_strain – Incremental spatial logarithmic strain.

Return type:

numpy.ndarray (2d)

static compute_spatial_tangent_modulus(e_log_strain_old, def_gradient_old, inc_def_gradient, cauchy_stress, inf_consistent_tangent)[source]

Compute finite strain spatial consistent tangent modulus.

\[\mathsf{a}_{ijkl} = \dfrac{1}{2 \det (\boldsymbol{F})} \, \left[ \mathsf{D} : \mathsf{L} : \mathsf{B} \right]_{ijkl} - \sigma_{il} \delta_{jk}\]

where \(\mathbf{\mathsf{a}}\) is the spatial consistent tangent modulus, \(\mathbf{\mathsf{D}}\) is the derivative of the Kirchhoff stress tensor with respect to the spatial logarithmic strain tensor, \(\mathbf{\mathsf{L}}\) is the derivative of the tensor logarithm function evaluated at the elastic trial left Cauchy-Green strain tensor, \(\mathbf{\mathsf{B}}\) is computed from the elastic trial left Cauchy-Green strain tensor components, \(\boldsymbol{\sigma}\) is the Cauchy stress tensor, and \(\delta_{ij}\) is the Kronecker delta.

The detailed definition of the finite strain spatial consistent tangent modulus of isotropic hyperelastic-based finite strain elastoplastic constitutive models whose finite strain formalism is purely kinematical can be found in Appendix F.4 of Ferreira (2022) [3] (see Equation (F.18)).


Parameters:
  • e_log_strain_old (numpy.ndarray (2d)) – Last converged elastic spatial logarithmic strain tensor.

  • def_gradient_old (numpy.ndarray (2d)) – Last converged deformation gradient.

  • inc_def_gradient (numpy.ndarray (2d)) – Incremental deformation gradient.

  • cauchy_stress (numpy.ndarray (2d)) – Cauchy stress tensor.

  • inf_consistent_tangent (numpy.ndarray (4d)) – Infinitesimal consistent tangent modulus.

Returns:

spatial_consistent_tangent – Spatial consistent tangent modulus.

Return type:

numpy.ndarray (4d)

get_clusters_def_gradient_mf()[source]

Get deformation gradient of each material cluster.

Returns:

clusters_def_gradient_mf – Deformation gradient (item, numpy.ndarray (1d)) associated to each material cluster (key, str), stored in matricial form.

Return type:

dict

get_clusters_def_gradient_old_mf()[source]

Get last converged deformation gradient of each material cluster.

Returns:

clusters_def_gradient_old_mf – Last converged deformation gradient (item, numpy.ndarray (1d)) associated to each material cluster (key, str), stored in matricial form.

Return type:

dict

get_clusters_inc_strain_mf(global_strain_mf)[source]

Get clusters incremental strain in matricial form.

Infinitesimal strains:

\[\Delta \boldsymbol{\varepsilon}_{\mu, n + 1}^{(I)} = \boldsymbol{\varepsilon}_{\mu, n + 1}^{(I)} - \boldsymbol{\varepsilon}_{\mu, n}^{(I)} \, , \quad I=1,\dots, n_{c}\]

where \(\Delta \boldsymbol{\varepsilon}_{\mu}^{(I)}\) is the \(I\) th material cluster incremental infinitesimal strain tensor, \(\boldsymbol{\varepsilon}_{\mu}^{(I)}\) is the \(I\) th material cluster infinitesimal strain tensor, \(n_{c}\) is the number of material clusters, \(n+1\) denotes the current increment, and \(n\) denotes the last converged increment.


Finite strains:

\[(\boldsymbol{F}_{\Delta})_{\mu, n + 1}^{(I)} = \boldsymbol{F}_{\mu, n + 1}^{(I)} ( \boldsymbol{F}_{\mu, n}^{(I)})^{-1} \, , \quad I=1,\dots, n_{c}\]

where \(\Delta \boldsymbol{F}_{\mu}^{(I)}\) is the \(I\) th material cluster incremental deformation gradient \(\boldsymbol{F}_{\mu}^{(I)}\) is the \(I\) th material cluster deformation gradient, \(n_{c}\) is the number of material clusters, \(n+1\) denotes the current increment, and \(n\) denotes the last converged increment.


Parameters:

global_strain_mf (numpy.ndarray (1d)) – Global vector of clusters strains stored in matricial form.

Returns:

clusters_inc_strain_mf – Incremental strain (item, dict) associated to each material cluster (key, str), stored in matricial form.

Return type:

dict

get_clusters_state()[source]

Get material state variables of each material cluster.

Returns:

clusters_state – Material constitutive model state variables (item, dict) associated to each material cluster (key, str).

Return type:

dict

get_clusters_state_old()[source]

Get last converged material state variables of each mat. cluster.

Returns:

clusters_state_old – Last converged material constitutive model state variables (item, dict) associated to each material cluster (key, str).

Return type:

dict

get_clusters_tangent_mf()[source]

Get material consistent tangent modulus of each material cluster.

Returns:

clusters_tangent_mf – Material consistent tangent modulus (item, numpy.ndarray) associated to each material cluster (key, str), stored in matricial form.

Return type:

dict

get_hom_strain_mf()[source]

Get homogenized strain tensor (matricial form).

Returns:

hom_strain_mf – Homogenized strain tensor stored in matricial form.

Return type:

numpy.ndarray (1d)

get_hom_strain_old_mf()[source]

Get last converged homogenized strain tensor (matricial form).

Returns:

hom_strain_old_mf – Last converged homogenized strain tensor stored in matricial form: infinitesimal strain tensor (infinitesimal strains) or deformation gradient (finite strains).

Return type:

numpy.ndarray (1d)

get_hom_stress_mf()[source]

Get homogenized stress tensor (matricial form).

Returns:

hom_stress_mf – Homogenized stress tensor stored in matricial form.

Return type:

numpy.ndarray (1d)

get_hom_stress_old_mf()[source]

Get last converged homogenized stress tensor (matricial form).

Returns:

hom_stress_old_mf – Last converged homogenized stress tensor stored in matricial form: Cauchy stress tensor (infinitesimal strains) or first Piola-Kirchhoff stress tensor (finite strains).

Return type:

numpy.ndarray (1d)

get_inc_hom_strain_mf()[source]

Get incremental homogenized strain tensor (matricial form).

Infinitesimal strains:

\[\Delta \boldsymbol{\varepsilon}_{n + 1} = \boldsymbol{\varepsilon}_{n + 1} - \boldsymbol{\varepsilon}_{n}\]

where \(\Delta \boldsymbol{\varepsilon}\) is the incremental homogenized infinitesimal strain tensor, \(\boldsymbol{\varepsilon}\) is the homogenized infinitesimal strain tensor, \(n+1\) denotes the current increment, and \(n\) denotes the last converged increment.


Finite strains:

\[(\boldsymbol{F}_{\Delta})_{n + 1} = \boldsymbol{F}_{n + 1} (\boldsymbol{F}_{n})^{-1}\]

where \(\boldsymbol{F}_{\Delta}\) is the homogenized incremental deformation gradient, \(\boldsymbol{F}\) is the homogenized deformation gradient, \(n+1\) denotes the current increment, and \(n\) denotes the last converged increment.


Returns:

inc_hom_strain_mf – Incremental homogenized strain tensor stored in matricial form: infinitesimal strain tensor (infinitesimal strains) or deformation gradient (finite strains).

Return type:

numpy.ndarray (1d)

get_inc_hom_stress_mf()[source]

Get incremental homogenized stress tensor (matricial form).

Infinitesimal strains:

\[\Delta \boldsymbol{\sigma}_{n + 1} = \boldsymbol{\sigma}_{n + 1} - \boldsymbol{\sigma}_{n}\]

where \(\Delta \boldsymbol{\sigma}\) is the incremental homogenized Cauchy stress tensor, \(\boldsymbol{\sigma}\) is the homogenized Cauchy stress tensor, \(n+1\) denotes the current increment, and \(n\) denotes the last converged increment.


Finite strains:

\[\Delta \boldsymbol{P}_{n + 1} = \boldsymbol{P}_{n + 1} - \boldsymbol{P}_{n}\]

where \(\Delta \boldsymbol{P}\) is the incremental homogenized first Piola-Kirchhoff stress tensor, \(\boldsymbol{P}\) is the homogenized first Piola-Kirchhoff stress tensor, \(n+1\) denotes the current increment, and \(n\) denotes the last converged increment.


Returns:

inc_hom_stress_mf – Incremental homogenized stress tensor stored in matricial form: Cauchy stress tensor (infinitesimal strains) or first Piola-Kirchhoff stress tensor (finite strains).

Return type:

numpy.ndarray (1d)

get_material_phases()[source]

Get RVE material phases.

Returns:

material_phases – RVE material phases labels (str).

Return type:

list[str]

get_material_phases_models()[source]

Get RVE material phases constitutive models.

Returns:

_material_phases_models – Material constitutive model (item, ConstitutiveModel) associated to each material phase (key, str).

Return type:

dict

get_material_phases_properties()[source]

Get RVE material phases constitutive properties.

Returns:

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

Return type:

dict

get_material_phases_vf()[source]

Get RVE material phases volume fraction.

Returns:

material_phases_vf – Volume fraction (item, float) associated to each material phase (key, str).

Return type:

dict

get_oop_hom_comp()[source]

Get homogenized strain or stress tensor out-of-plane component.

Returns:

oop_hom_comp – Homogenized strain or stress tensor out-of-plane component.

Return type:

float

init_clusters_state()[source]

Initialize clusters state variables.

init_constitutive_model(mat_phase, model_keyword, model_source='crate')[source]

Initialize material phase constitutive model.

Parameters:
  • mat_phase (str) – Material phase label.

  • model_keyword (str) – Material constitutive model input data file keyword.

  • model_source ({'crate',}, default='crate') – Material constitutive model source.

set_phase_clusters(phase_clusters, clusters_vf)[source]

Set CRVE cluster labels and volume fractions of each material phase.

Parameters:
  • phase_clusters (dict) – Clusters labels (item, list[int]) associated to each material phase (key, str).

  • clusters_vf (dict) – Volume fraction (item, float) associated to each material cluster (key, str).

set_rewind_state_updated_clustering(phase_clusters, clusters_vf, clusters_state, clusters_def_gradient_mf)[source]

Set rewind state variables according to updated clustering.

Parameters:
  • phase_clusters (dict) – Clusters labels (item, list[int]) associated to each material phase (key, str).

  • clusters_vf (dict) – Volume fraction (item, float) associated to each material cluster (key, str).

  • clusters_state (dict) – Material constitutive model state variables (item, dict) associated to each material cluster (key, str).

  • clusters_def_gradient_mf (dict) – Deformation gradient (item, numpy.darray (1d)) associated to each material cluster (key, str), stored in matricial form.

update_clusters_state(clusters_inc_strain_mf)[source]

Update clusters state variables and consistent tangent modulus.

Parameters:

clusters_inc_strain_mf (dict) – Incremental strain (item, numpy.ndarray) associated to each material cluster (key, str), stored in matricial form.

Returns:

su_fail_state – State update failure state.

Return type:

dict

update_converged_state()[source]

Update last converged material state variables.

update_state_homogenization()[source]

Update homogenized strain and stress tensors.

Infinitesimal strains:

\[\boldsymbol{\varepsilon}_{n + 1} = \sum_{I=1}^{n_{c}} f^{(I)} \boldsymbol{\varepsilon}_{\mu, n + 1}^{(I)}\]

where \(\boldsymbol{\varepsilon}\) is the homogenized infinitesimal strain tensor, \(f^{(I)}\) is the \(I\) th material cluster volume fraction, \(\boldsymbol{\varepsilon}_{\mu}^{(I)}\) is the \(I\) th material cluster infinitesimal strain tensor, \(n_{c}\) is the number of material clusters, and \(n+1\) denotes the current increment.

\[\boldsymbol{\sigma}_{n + 1} = \sum_{I=1}^{n_{c}} f^{(I)} \boldsymbol{\sigma}_{\mu, n + 1}^{(I)}\]

where \(\boldsymbol{\sigma}\) is the homogenized Cauchy stress tensor, \(f^{(I)}\) is the \(I\) th material cluster volume fraction, \(\boldsymbol{\sigma}_{\mu}^{(I)}\) is the \(I\) th material cluster Cauchy stress tensor, \(n_{c}\) is the number of material clusters, and \(n+1\) denotes the current increment.


Finite strains:

\[\boldsymbol{F}_{n + 1} = \sum_{I=1}^{n_{c}} f^{(I)} \boldsymbol{F}_{\mu, n + 1}^{(I)}\]

where \(\boldsymbol{F}\) is the homogenized deformation gradient, \(f^{(I)}\) is the \(I\) th material cluster volume fraction, \(\boldsymbol{F}_{\mu}^{(I)}\) is the \(I\) th material cluster deformation gradient, \(n_{c}\) is the number of material clusters, and \(n+1\) denotes the current increment.

\[\boldsymbol{P}_{n + 1} = \sum_{I=1}^{n_{c}} f^{(I)} \boldsymbol{P}_{\mu, n + 1}^{(I)}\]

where \(\boldsymbol{P}\) is the homogenized first Piola-Kirchhoff stress tensor, \(f^{(I)}\) is the \(I\) th material cluster volume fraction, \(\boldsymbol{P}_{\mu}^{(I)}\) is the \(I\) th material cluster first Piola-Kirchhoff stress tensor, \(n_{c}\) is the number of material clusters, and \(n+1\) denotes the current increment.