cratepy.clustering.clusteringdata.RVEElasticDatabase

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

Bases: object

RVE local elastic response database class.

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

Type:

dict

_elastic_eff_modulus_matrix

RVE’s elastic effective tangent modulus in matrix form, where each entry contains a given elastic modulus (similar to Voigt notation). The elastic effective tangent modulus is computed from the stress conjugate to the infinitesimal strain tensor (infinitesimal strains) or to the material logarithmic strain tensor (finite strains).

Type:

numpy.ndarray (2d)

_eff_elastic_properties

Elastic properties (key, str) and their values (item, float) estimated from the RVE’s elastic effective tangent modulus.

Type:

dict

rve_global_response

RVE local elastic strain response for a given set of macroscale loadings, where each macroscale loading is associated with a set of independent strain components (numpy.ndarray of shape (n_voxels, n_mac_strains*n_strain_comps)). Each column is associated associated with a independent strain component of the infinitesimal strain tensor (infinitesimal strains) or material logarithmic strain tensor (finite strains).

Type:

numpy.ndarray (2d)

compute_rve_response_database(self, dns_method, dns_method_data, mac_strains, is_strain_sym)[source]

Compute RVE’s local elastic strain response database.

compute_rve_elastic_tangent_modulus(self, strain_magnitude_factor=1.0)[source]

Compute RVE’s elastic effective tangent modulus.

set_eff_isotropic_elastic_constants(self):

Set isotropic elastic constants from effective tangent modulus.

get_eff_isotropic_elastic_constants(self):

Get isotropic elastic constants from effective tangent modulus.

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_elastic_tangent_modulus

Compute RVE's elastic effective tangent modulus.

compute_rve_response_database

Compute RVE's local elastic strain response database.

get_eff_isotropic_elastic_constants

Get isotropic elastic constants from effective tangent modulus.

set_eff_isotropic_elastic_constants

Set isotropic elastic constants from effective tangent modulus.

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_rve_elastic_tangent_modulus(strain_magnitude_factor=1.0)[source]

Compute RVE’s elastic effective tangent modulus.

Parameters:

strain_magnitude_factor (float, default=1.0) – Macroscale strain magnitude factor.

compute_rve_response_database(dns_method, dns_method_data, mac_strains, is_strain_sym)[source]

Compute RVE’s local elastic strain response database.

Build a RVE’s local elastic strain response database by solving one or more microscale equilibrium problems (each associated with a given macroscale strain loading) through a given homogenization-based multi-scale method.


Parameters:
  • dns_method ({'fft-basic'}) – DNS homogenization-based multi-scale method.

  • dns_method_data (dict) – Parameters of DNS homogenization-based multi-scale method.

  • mac_strains (list[numpy.ndarray (2d)]) – List of macroscale strain loadings (numpy.ndarray (2d)). Infinitesimal strain tensor (infinitesimal strains) or deformation gradient (finite strains).

  • is_strain_sym (bool) – True if the macroscale strain second-order tensor associated with the RVE’s local elastic strain response database is symmetric by definition, False otherwise.

get_eff_isotropic_elastic_constants()[source]

Get isotropic elastic constants from effective tangent modulus.

Returns:

eff_elastic_properties – Elastic properties (key, str) and their values (item, float) estimated from the RVE’s elastic effective tangent modulus.

Return type:

dict

set_eff_isotropic_elastic_constants()[source]

Set isotropic elastic constants from effective tangent modulus.