cratepy.clustering.crve.CRVE

class CRVE(rve_dims, regular_grid, material_phases, strain_formulation, problem_type, global_data_matrix, clustering_type, phase_n_clusters, base_clustering_scheme, eff_elastic_properties=None, adaptive_clustering_scheme=None, adapt_criterion_data=None, adaptivity_type=None, adaptivity_control_feature=None)[source]

Bases: object

Cluster-Reduced Representative Volume Element.

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

_n_voxels_dims

Number of voxels in each dimension of the regular grid (spatial discretization of the RVE).

Type:

list[int]

_n_voxels

Total number of voxels of the regular grid (spatial discretization of the RVE).

Type:

int

_phase_voxel_flatidx

Flat (1D) voxels’ indexes (item, list[int]) associated with each material phase (key, str).

Type:

dict

_gop_X_dft_vox

Green operator material independent terms. Each term is stored in a dictionary, where each pair of strain/stress components (key, str) is associated with the Green operator material independent term evaluated in all spatial discrete points (item, numpy.ndarray (2d or 3d)).

Type:

list[dict]

_cluster_phases

Cluster-Reduced material phase instance (item, CRMP) associated with each material phase (key, str).

Type:

dict

_base_phase_n_clusters

Number of clusters (item, int) prescribed (base clustering) for each material phase (key, str).

Type:

dict

_adaptive_step

Counter of adaptive clustering steps, with 0 associated with the base clustering.

Type:

int

_voxels_clusters

Regular grid of voxels (spatial discretization of the RVE), where each entry contains the cluster label (int) assigned to the corresponding voxel.

Type:

numpy.ndarray (2d or 3d)

_phase_clusters

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

Type:

dict

_clusters_vf

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

Type:

dict

_cit_x_mf

Cluster interaction tensors associated with the Green operator material independent terms. Each term is stored in a dictionary (item, dict) for each pair of material phases (key, str), which in turn contains the corresponding matricial form (item, numpy.ndarray (2d)) associated with each pair of clusters (key, str).

Type:

list[dict]

_adapt_material_phases

RVE adaptive material phases labels (str).

Type:

list[str]

_adaptive_clustering_time

Total amount of time (s) spent in clustering adaptivity.

Type:

float

_adaptive_cit_time

Total amount of time (s) spent in clustering adaptivity cluster interaction tensors computation procedures.

Type:

float

perform_crve_base_clustering(self)[source]

Compute CRVE base clustering.

perform_crve_adaptivity(self, target_clusters, target_clusters_data)[source]

Perform CRVE clustering adaptivity.

get_rve_dims(self)[source]

Get RVE dimensions.

get_material_phases(self)[source]

Get RVE material phases.

get_n_voxels(self)[source]

Get number of voxels in each dimension and total number of voxels.

get_regular_grid(self)[source]

Get regular grid of voxels with material phase labels.

get_phase_n_clusters(self)[source]

Get number of clusters associated with each material phase.

get_phase_clusters(self)[source]

Get clusters associated with each material phase.

get_voxels_clusters(self)[source]

Get regular grid containing the cluster label of each voxel.

get_n_total_clusters(self)[source]

Get current total number of clusters.

get_cluster_phases(self)[source]

Get cluster-reduced material phase instance of each material phase.

get_clusters_vf(self)[source]

Get clusters volume fraction.

get_cit_x_mf(self)[source]

Get cluster interaction tensors (Green material independent terms).

get_eff_isotropic_elastic_constants(self)[source]

Get isotropic elastic constants from elastic tangent modulus.

get_adapt_material_phases(self)[source]

Get adaptive material phases labels.

get_adaptivity_control_feature(self)[source]

Get clustering adaptivity control feature of each material phase.

get_adapt_criterion_data(self)[source]

Get clustering adaptivity criterion data of each material phase.

get_voxels_array_variables(self)[source]

Get required variables to build a clusters state based voxels array.

get_adaptive_step(self)[source]

Get counter of adaptive clustering steps.

get_adaptive_clustering_time(self)[source]

Get total amount of time spent in clustering adaptivity.

get_adaptive_cit_time(self)[source]

Get total time spent in adaptivity cluster interaction tensors.

get_clustering_type(self)[source]

Get clustering type of each material phase.

get_adaptivity_output(self)[source]

Get required data for clustering adaptivity output file.

get_clustering_summary(self)[source]

Get summary of number of clusters of each material phase.

reset_adaptive_parameters(self)[source]

Reset CRVE adaptive progress parameters and set base clustering.

update_adaptive_parameters(self, adaptive_clustering_scheme, adapt_criterion_data, adaptivity_type, adaptivity_control_feature)[source]

Update CRVE clustering adaptivity attributes.

get_crmp_types()[source]

Get available cluster-reduced material phases types.

_get_features_indexes(clustering_scheme)[source]

Get unique clustering features indexes from clustering scheme.

_get_phase_idxs(regular_grid, material_phases)[source]

Get flat indexes of each material phase’s voxels.

_get_cluster_idxs(voxels_clusters, cluster_label)[source]

Get flat indexes of given cluster’s voxels.

_set_phase_clusters(self)[source]

Set CRVE cluster labels associated with each material phase.

_set_clusters_vf(self)[source]

Set CRVE clusters’ volume fractions.

_get_clusters_max_label(self)[source]

Get CRVE maximum cluster label.

_sort_cluster_labels(self)[source]

Reassign and sort CRVE cluster labels material phasewise.

compute_cit(self, mode='full', adaptive_clustering_map=None)[source]

Compute CRVE cluster interaction tensors.

_cluster_filter(self, cluster)[source]

Compute cluster discrete characteristic function.

_gop_convolution(self, cluster_filter_dft, gop_1_dft_vox, gop_2_dft_vox, gop_0_freq_dft_vox)[source]

Convolution of cluster characteristic function and Green operator.

_discrete_cit_integral(self, cluster_filter, gop_1_filt_vox, gop_2_filt_vox, gop_0_freq_filt_vox)[source]

Discrete integral over the spatial domain of material cluster.

_switch_pair(x, delimiter='_')[source]

Switch left and right sides of string with separating delimiter.

save_crve_file(crve, crve_file_path)[source]

Dump CRVE into file.

Constructor.

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

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

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

  • global_data_matrix (numpy.ndarray (2d)) – Data matrix (numpy.ndarray of shape (n_voxels, n_features_dims)) containing the required clustering features’ data to perform all the prescribed cluster analyses.

  • clustering_type (dict) – Clustering type (item, {‘static’, ‘adaptive’}) of each material phase (key, str).

  • phase_n_clusters (dict) – Number of clusters (item, int) associated with each material phase (key, str).

  • base_clustering_scheme (dict) – Prescribed base clustering scheme (item, numpy.ndarray of shape (n_clusterings, 3)) for each material phase (key, str). Each row is associated with a unique clustering characterized by a clustering algorithm (col 1, int), a list of features (col 2, list[int]) and a list of the features data matrix’ indexes (col 3, list[int]).

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

  • adaptive_clustering_scheme (dict) – Prescribed adaptive clustering scheme (item, numpy.ndarray of shape (n_clusterings, 3)) for each material phase (key, str). Each row is associated with a unique clustering characterized by a clustering algorithm (col 1, int), a list of features (col 2, list[int]) and a list of the features data matrix’ indexes (col 3, list[int]).

  • adapt_criterion_data (dict, default=None) – Clustering adaptivity criterion (item, dict) associated with each material phase (key, str). This dictionary contains the adaptivity criterion to be used and the required parameters.

  • adaptivity_type (dict, default=None) – Clustering adaptivity type (item, dict) associated with each material phase (key, str). This dictionary contains the adaptivity type to be used and the required parameters.

  • adaptivity_control_feature (dict, default=None) – Clustering adaptivity control feature (item, str) associated with each material phase (key, str).

List of Public Methods

compute_cit

Compute CRVE cluster interaction tensors.

get_adapt_criterion_data

Get clustering adaptivity criterion data of each material phase.

get_adapt_material_phases

Get adaptive material phases labels.

get_adaptive_cit_time

Get total time spent in adaptivity cluster interaction tensors.

get_adaptive_clustering_time

Get total amount of time spent in clustering adaptivity.

get_adaptive_step

Get counter of adaptive clustering steps.

get_adaptivity_control_feature

Get clustering adaptivity control feature of each material phase.

get_adaptivity_output

Get required data for clustering adaptivity output file.

get_cit_x_mf

Get cluster interaction tensors (Green material independent terms).

get_cluster_phases

Get cluster-reduced material phase instance of each material phase.

get_clustering_summary

Get summary of number of clusters of each material phase.

get_clustering_type

Get clustering type of each material phase.

get_clusters_vf

Get clusters volume fraction.

get_crmp_types

Get available cluster-reduced material phases types.

get_eff_isotropic_elastic_constants

Get isotropic elastic constants from elastic tangent modulus.

get_material_phases

Get RVE material phases.

get_n_total_clusters

Get current total number of clusters.

get_n_voxels

Get number of voxels in each dimension and total number of voxels.

get_phase_clusters

Get clusters associated with each material phase.

get_phase_n_clusters

Get number of clusters associated with each material phase.

get_regular_grid

Get regular grid of voxels with material phase labels.

get_rve_dims

Get RVE dimensions.

get_voxels_array_variables

Get required variables to build a clusters state based voxels array.

get_voxels_clusters

Get regular grid containing the cluster label of each voxel.

perform_crve_adaptivity

Perform CRVE clustering adaptivity.

perform_crve_base_clustering

Compute CRVE base clustering.

reset_adaptive_parameters

Reset CRVE adaptive progress parameters and set base clustering.

save_crve_file

Dump CRVE into file.

update_adaptive_parameters

Update CRVE clustering adaptivity attributes.

Methods

__init__(rve_dims, regular_grid, material_phases, strain_formulation, problem_type, global_data_matrix, clustering_type, phase_n_clusters, base_clustering_scheme, eff_elastic_properties=None, adaptive_clustering_scheme=None, adapt_criterion_data=None, adaptivity_type=None, adaptivity_control_feature=None)[source]

Constructor.

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

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

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

  • global_data_matrix (numpy.ndarray (2d)) – Data matrix (numpy.ndarray of shape (n_voxels, n_features_dims)) containing the required clustering features’ data to perform all the prescribed cluster analyses.

  • clustering_type (dict) – Clustering type (item, {‘static’, ‘adaptive’}) of each material phase (key, str).

  • phase_n_clusters (dict) – Number of clusters (item, int) associated with each material phase (key, str).

  • base_clustering_scheme (dict) – Prescribed base clustering scheme (item, numpy.ndarray of shape (n_clusterings, 3)) for each material phase (key, str). Each row is associated with a unique clustering characterized by a clustering algorithm (col 1, int), a list of features (col 2, list[int]) and a list of the features data matrix’ indexes (col 3, list[int]).

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

  • adaptive_clustering_scheme (dict) – Prescribed adaptive clustering scheme (item, numpy.ndarray of shape (n_clusterings, 3)) for each material phase (key, str). Each row is associated with a unique clustering characterized by a clustering algorithm (col 1, int), a list of features (col 2, list[int]) and a list of the features data matrix’ indexes (col 3, list[int]).

  • adapt_criterion_data (dict, default=None) – Clustering adaptivity criterion (item, dict) associated with each material phase (key, str). This dictionary contains the adaptivity criterion to be used and the required parameters.

  • adaptivity_type (dict, default=None) – Clustering adaptivity type (item, dict) associated with each material phase (key, str). This dictionary contains the adaptivity type to be used and the required parameters.

  • adaptivity_control_feature (dict, default=None) – Clustering adaptivity control feature (item, str) associated with each material phase (key, str).

_cluster_filter(cluster)[source]

Compute cluster discrete characteristic function.

\[\begin{split}\chi^{(I)} (\boldsymbol{Y}) = \begin{cases} 1 \quad \text{if} \; \; \; \boldsymbol{Y}\in \Omega^{(I)}_{\mu, \, 0} \\ 0 \quad \text{otherwise} \end{cases}\end{split}\]

where \(\chi^{(I)}\) is the characteristic function of the \(I\) th material cluster and \(\boldsymbol{Y}\) is a point of the microscale reference configuration (\(\Omega_{\mu,\,0}\)).

The detailed description of the cluster characteristic function can be found in Section 4.3.1 of Ferreira (2022) [1].


Parameters:

cluster (int) – Cluster label.

Returns:

  • cluster_filter (numpy.ndarray[bool] (2d or 3d)) – Cluster discrete characteristic function in spatial domain.

  • cluster_filter_dft (numpy.ndarray) – Cluster discrete characteristic function in frequency domain (discrete Fourier transform).

_discrete_cit_integral(cluster_filter, gop_1_filt_vox, gop_2_filt_vox, gop_0_freq_filt_vox)[source]

Discrete integral over the spatial domain of material cluster.

\[\int_{\Omega_{\mu, 0}} \chi^{(I)}(\boldsymbol{Y}) \left( \int_{\Omega_{\mu, 0}} \, \chi^{(J)}(\boldsymbol{Y}) \, \boldsymbol{\Phi}^{0}(\boldsymbol{Y}-\boldsymbol{Y}') \, \mathrm{d}v' \right) \mathrm{d}v \, ,\]

where \(\chi^{(I)}\) is the characteristic function of the \(I\) th material cluster, \(\boldsymbol{\mathsf{\Phi}}^{0}\) is the reference material Green operator (fourth-order tensor), and \(\boldsymbol{Y}\) and \(\boldsymbol{Y'}\) are points of the microscale reference configuration (\(\Omega_{\mu,\,0}\)).

Such a computation is required to compute the cluster interaction tensors. More information can be found in Ferreira (2022) [2] (see Equations (4.112) and surrounding text).


Parameters:
  • cluster_filter (numpy.ndarray) – Cluster discrete characteristic function in spatial domain.

  • gop_1_filt_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the convolution between the material cluster characteristic function and the first Green operator material independent term in the spatial domain (inverse discrete Fourier transform).

  • gop_2_filt_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the convolution between the material cluster characteristic function and the second Green operator material independent term in the spatial domain (inverse discrete Fourier transform).

  • gop_0_freq_filt_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the convolution between the material cluster characteristic function and the zero-frequency Green operator (material independent) term in the spatial domain (inverse discrete Fourier transform).

Returns:

  • cit_1_integral_mf (numpy.ndarray (2d)) – Discrete integral over the spatial domain of material cluster I of the discrete convolution between the material cluster J characteristic function and the first Green operator material independent term in the spatial domain (numpy.ndarray of shape (n_comps, n_comps)).

  • cit_2_integral_mf (numpy.ndarray (2d)) – Discrete integral over the spatial domain of material cluster I of the discrete convolution between the material cluster J characteristic function and the second Green operator material independent term in the spatial domain (numpy.ndarray of shape (n_comps, n_comps)).

  • cit_0_freq_integral_mf (numpy.ndarray (2d)) – Discrete integral over the spatial domain of material cluster I of the discrete convolution between the material cluster J characteristic function and the zero-frequency Green operator (material independent) term in the spatial domain.

static _get_cluster_idxs(voxels_clusters, cluster_label)[source]

Get flat indexes of given cluster’s voxels.

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

  • cluster_label (int) – Cluster label.

Returns:

cluster_voxel_flat_idx – Flat voxels’ indexes (int) associated with given material cluster.

Return type:

list[int]

_get_clusters_max_label()[source]

Get CRVE maximum cluster label.

Returns:

max_label – CRVE maximum cluster label.

Return type:

int

static _get_features_indexes(clustering_scheme)[source]

Get unique clustering features indexes from clustering scheme.

Parameters:

clustering_scheme (ndarry of shape (n_clusterings, 3)) – Clustering scheme (numpy.ndarray of shape (n_clusterings, 3)). Each row is associated with a unique clustering characterized by a clustering algorithm (col 1, int), a list of features (col 2, list[int]) and a list of the features data matrix’ indexes (col 3, list[int]).

Returns:

indexes – List of unique clustering features indexes.

Return type:

list[int]

static _get_phase_idxs(regular_grid, material_phases)[source]

Get flat indexes of each material phase’s voxels.

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

Returns:

phase_voxel_flat_idx – Flat voxels’ indexes (item, list[int]) associated with each material phase (key, str).

Return type:

dict

_gop_convolution(cluster_filter_dft, gop_1_dft_vox, gop_2_dft_vox, gop_0_freq_dft_vox)[source]

Convolution of cluster characteristic function and Green operator.

\[\int_{\Omega_{\mu, \, 0}} \chi^{(J)} (\boldsymbol{Y}') \, \boldsymbol{\mathsf{\Phi}}^{0} (\boldsymbol{Y} - \boldsymbol{Y}') \, \mathrm{d} v' = \mathscr{F}^{-1} \left( \breve{\chi}^{(J)}(\boldsymbol{\zeta}) \, \breve{\boldsymbol{\mathsf{\Phi}}}^{0} (\boldsymbol{\zeta}) \right) \, ,\]

where \(\chi^{(J)}\) is the characteristic function of the \(J\) th material cluster, \(\boldsymbol{\mathsf{\Phi}}^{0}\) is the reference material Green operator (fourth-order tensor), \(\boldsymbol{Y}\) and \(\boldsymbol{Y'}\) are points of the microscale reference configuration (\(\Omega_{\mu,\,0}\)), and \(\boldsymbol{\zeta}\) is the frequency wave vector. The operator \(\mathscr{F}^{-1}(\cdot)\) denotes the inverse Fourier transform and \(\breve{(\cdot)}\) denotes a field defined in the frequency domain.

Such a computation is required to compute the cluster interaction tensors. More information can be found in Ferreira (2022) [3] (see Equations (4.111) and surrounding text).


Parameters:
  • cluster_filter_dft (numpy.ndarray) – Cluster discrete characteristic function in frequency domain (discrete Fourier transform).

  • gop_1_dft_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the first Green operator material independent term in the frequency domain (discrete Fourier transform).

  • gop_2_dft_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the second Green operator material independent term in the frequency domain (discrete Fourier transform).

  • gop_0_freq_dft_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the Green operator zero-frequency (material independent) term in the frequency domain (discrete Fourier transform).

Returns:

  • gop_1_filt_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the convolution between the material cluster characteristic function and the first Green operator material independent term in the spatial domain (inverse discrete Fourier transform).

  • gop_2_filt_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the convolution between the material cluster characteristic function and the second Green operator materia independent term in the spatial domain (inverse discrete Fourier transform).

  • gop_0_freq_filt_vox (dict) – Regular grid shaped matrix (item, numpy.ndarray) containing each fourth-order matricial form component (key, str) of the convolution between the material cluster characteristic function and the zero-frequency Green operator (material independent) term in the spatial domain (inverse discrete Fourier transform).

_set_clusters_vf()[source]

Set CRVE clusters’ volume fractions.

_set_phase_clusters()[source]

Set CRVE cluster labels associated with each material phase.

_sort_cluster_labels()[source]

Reassign and sort CRVE cluster labels material phasewise.

Reassign CRVE cluster labels in the range (0, n_clusters) and sort them in ascending order of material phase’s labels.

static _switch_pair(x, delimiter='_')[source]

Switch left and right sides of string with separating delimiter.

Parameters:
  • x (str) – Target string.

  • delimiter (str, default='_') – Separating delimiter between target’s string left and right sides.

Returns:

y – Switched string.

Return type:

str

compute_cit(mode='full', adaptive_clustering_map=None)[source]

Compute CRVE cluster interaction tensors.

Cluster interaction tensors:

\[\boldsymbol{\mathsf{T}}^{(I)(J)} = \dfrac{1}{f^{(I)} v_{\mu}} \int_{\Omega_{\mu, \, 0}} \int_{\Omega_{\mu, \, 0}} \chi^{(I)}(\boldsymbol{Y}) \, \chi^{(J)} (\boldsymbol{Y}') \, \boldsymbol{\mathsf{\Phi}}^{0} (\boldsymbol{Y}-\boldsymbol{Y}') \, \mathrm{d} v' \mathrm{d} v \, ,\]
\[\quad I,J = 1,2,\dots, n_{\mathrm{c}}\]

where \(\boldsymbol{\mathsf{T}}^{(I)(J)}\) is the cluster interaction tensor (fourth-order tensor) between the \(I\) th and \(J\) th material clusters, \(f^{(I)}\) is the volume fraction of the \(I\) th material cluster, \(v_{\mu}\) is the volume of the CRVE, \(\boldsymbol{Y}\) and \(\boldsymbol{Y'}\) are points of the microscale reference configuration (\(\Omega_{\mu,\,0}\)), \(\chi^{(I)}\) is the characteristic function of the \(I\) th material cluster, \(\boldsymbol{\mathsf{\Phi}}^{0}\) is the reference material Green operator (fourth-order tensor), and \(n_{c}\) is the number of material clusters.

The detailed description of the custer interaction tensors can be found in Section 4.3 of Ferreira (2022) [4].


Adaptive update of cluster interaction matrix:

The cluster interaction tensors adaptive computation mode can only be performed after at least one full computation has been performed.

A detailed description of the adaptive update of the cluster interaction matrix can be found in Section 2.3.3 of Ferreira et. al (2022) [5].


Parameters:
  • mode ({'full', 'adaptive'}, default='full') – The default full mode performs the complete computation of all cluster interaction tensors. The ‘adaptive’ mode speeds up the computation of the new cluster interaction tensors resulting from an adaptive clustering characterized by adaptive_clustering_map.

  • 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). Required in adaptive mode, otherwise ignored.

get_adapt_criterion_data()[source]

Get clustering adaptivity criterion data of each material phase.

Returns:

adapt_criterion_data – Clustering adaptivity criterion (item, dict) associated with each material phase (key, str). This dictionary contains the adaptivity criterion to be used and the required parameters.

Return type:

dict, default=None

get_adapt_material_phases()[source]

Get adaptive material phases labels.

Returns:

adapt_material_phases – RVE adaptive material phases labels (str).

Return type:

list[str]

get_adaptive_cit_time()[source]

Get total time spent in adaptivity cluster interaction tensors.

Returns:

adaptive_cit_time – Total amount of time (s) spent in clustering adaptivity cluster interaction tensors computation procedures.

Return type:

float

get_adaptive_clustering_time()[source]

Get total amount of time spent in clustering adaptivity.

Returns:

adaptive_clustering_time – Total amount of time (s) spent in clustering adaptivity.

Return type:

float

get_adaptive_step()[source]

Get counter of adaptive clustering steps.

Returns:

adaptive_step – Counter of adaptive clustering steps, with 0 associated with the base clustering.

Return type:

int

get_adaptivity_control_feature()[source]

Get clustering adaptivity control feature of each material phase.

Returns:

adaptivity_control_feature – Clustering adaptivity control feature (item, str) associated with each material phase (key, str).

Return type:

dict, default=None

get_adaptivity_output()[source]

Get required data for clustering adaptivity output file.

Returns:

adaptivity_output – For each adaptive material phase (key, str), stores a list (item) containing the adaptivity metrics associated with the clustering adaptivity output file.

Return type:

dict

get_cit_x_mf()[source]

Get cluster interaction tensors (Green material independent terms).

Returns:

cit_x_mf – Cluster interaction tensors associated with the Green operator material independent terms. Each term is stored in a dictionary (item, dict) for each pair of material phases (key, str), which in turn contains the corresponding matricial form (item, numpy.ndarray (2d)) associated with each pair of clusters (key, str).

Return type:

list[dict]

get_cluster_phases()[source]

Get cluster-reduced material phase instance of each material phase.

Returns:

cluster_phases – Cluster-Reduced material phase instance (item, CRMP) associated with each material phase (key, str).

Return type:

dict

get_clustering_summary()[source]

Get summary of number of clusters of each material phase.

Returns:

clustering_summary – For each material phase (key, str), stores list (item) containing the associated type (‘static’ or ‘adaptive’), the base number of clusters (int) and the final number of clusters (int).

Return type:

dict

get_clustering_type()[source]

Get clustering type of each material phase.

Returns:

clustering_type – Clustering type (item, {‘static’, ‘adaptive’}) of each material phase (key, str).

Return type:

dict

get_clusters_vf()[source]

Get clusters volume fraction.

Returns:

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

Return type:

dict

static get_crmp_types()[source]

Get available cluster-reduced material phases types.

Returns:

available_crmp_types – Available cluster-reduced material phase classes (item, CRMP) and associated identifiers (key, str).

Return type:

dict

get_eff_isotropic_elastic_constants()[source]

Get isotropic elastic constants from elastic 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

get_material_phases()[source]

Get RVE material phases.

Returns:

material_phases – RVE material phases labels (str).

Return type:

list

get_n_total_clusters()[source]

Get current total number of clusters.

Returns:

n_total_clusters – Total number of clusters.

Return type:

int

get_n_voxels()[source]

Get number of voxels in each dimension and total number of voxels.

Returns:

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

  • n_voxels (int) – Total number of voxels of the regular grid (spatial discretization of the RVE).

get_phase_clusters()[source]

Get clusters associated with each material phase.

Returns:

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

Return type:

dict

get_phase_n_clusters()[source]

Get number of clusters associated with each material phase.

Returns:

phase_n_clusters – Number of clusters (item, int) associated with each material phase (key, str).

Return type:

dict

get_regular_grid()[source]

Get regular grid of voxels with material phase labels.

Returns:

regular_grid – Regular grid of voxels (spatial discretization of the RVE), where each entry contains the material phase label (int) assigned to the corresponding voxel.

Return type:

numpy.ndarray (2d or 3d)

get_rve_dims()[source]

Get RVE dimensions.

Returns:

rve_dims – RVE size in each dimension.

Return type:

list

get_voxels_array_variables()[source]

Get required variables to build a clusters state based voxels array.

Returns:

  • material_phases (list) – CRVE material phases labels (str).

  • phase_clusters (dict) – Clusters labels (item, list of int) associated with each material phase (key, str).

  • voxels_clusters (ndarray) – Regular grid of voxels (spatial discretization of the RVE), where each entry contains the cluster label (int) assigned to the corresponding voxel.

get_voxels_clusters()[source]

Get regular grid containing the cluster label of each voxel.

Returns:

voxels_clusters – Regular grid of voxels (spatial discretization of the RVE), where each entry contains the cluster label (int) assigned to the corresponding voxel.

Return type:

numpy.ndarray (2d or 3d)

perform_crve_adaptivity(target_clusters, target_clusters_data)[source]

Perform CRVE clustering adaptivity.

Parameters:
  • target_clusters (list[int]) – List with the labels (int) of clusters to be adapted.

  • target_clusters_data (dict) – For each target cluster (key, str), store dictionary (item, dict) containing cluster associated parameters relevant for the adaptive procedures.

Returns:

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

Return type:

dict

perform_crve_base_clustering()[source]

Compute CRVE base clustering.

reset_adaptive_parameters()[source]

Reset CRVE adaptive progress parameters and set base clustering.

static save_crve_file(crve, crve_file_path)[source]

Dump CRVE into file.

Parameters:
  • crve (CRVE) – Cluster-Reduced Representative Volume Element.

  • crve_file_path (str) – Path of file where the CRVE’s instance is dumped.

update_adaptive_parameters(adaptive_clustering_scheme, adapt_criterion_data, adaptivity_type, adaptivity_control_feature)[source]

Update CRVE clustering adaptivity attributes.

Parameters:
  • adaptive_clustering_scheme (dict) – Prescribed adaptive clustering scheme (item, numpy.ndarray of shape (n_clusterings, 3)) for each material phase (key, str). Each row is associated with a unique clustering characterized by a clustering algorithm (col 1, int), a list of features (col 2, list[int]) and a list of the features data matrix’ indexes (col 3, list[int]).

  • adapt_criterion_data (dict) – Clustering adaptivity criterion (item, dict) associated with each material phase (key, str). This dictionary contains the adaptivity criterion to be used and the required parameters.

  • adaptivity_type (dict) – Clustering adaptivity type (item, dict) associated with each material phase (key, str). This dictionary contains the adaptivity type to be used and the required parameters.

  • adaptivity_control_feature (dict) – Clustering adaptivity control feature (item, str) associated with each material phase (key, str).