ThermoNTFA Documentation
The python package thermontfa with our reference implementation of the thermo-plastic NTFA contains following classes:
- class thermontfa.TabularInterpolation(temps: ndarray = None, data: ndarray = None, const_extrapolate: bool = False)
Bases:
object
Tabular interpolation for the thermo-mechanical NTFA
Performs a linear interpolation of the NTFA system matrices for a given temperature given tabular data at sufficiently many temperature points. It can be initialized with given tabular data or based on a HDF5 file (
*.h5
).- __init__(temps: ndarray = None, data: ndarray = None, const_extrapolate: bool = False) None
Initialize the tabular interpolator for given
data
at prescribed temperaturestemps
.- Parameters:
temps (np.ndarray) – temperature points on which tabular data is available. The shape of the
numpy
array is expected to be(N_t,)
.data (np.ndarray) – tabular data, e.g., a batch of NTFA system matrices with shape
(N_t, ...)
.const_extrapolate – If true, a constant extrapolation instead of a linear extrapolation is performed. The default value is false.
- temps: ndarray = array([], dtype=float64)
- temp_min: float = 0.0
- temp_max: float = 0.0
- dim: Tuple[int, ...] = ()
- const_extrapolate: bool = False
- data: ndarray = array([], dtype=float64)
- classmethod from_h5(file_name: str, dset_temps: str, dset_data: str, transpose_dims: Tuple[int, ...] | None = None, const_extrapolate: bool = False) Self
Initialize the tabular interpolator based on tabular data stored in a HDF5 file (
*.h5
).This is a factory method and returns a new instance of the
TabularInterpolation
class. It is expected that the HDF5 file contains a data set with pathdset_temps
that contains a list of the temperature points on which tabular data is available. The shape of this dataset is expected to be(N_t,)
. Additionaly, the HDF5 file must contain a data set with pathdset_data
that contains the tabular data, e.g., a batch of NTFA system matrices with shape(N_t, ...)
. The order of axes/dimensions of the data set with path dset_data can be changed by transposing to the axis order given intranspose_dims
.- Parameters:
file_name (str) – path of the HDF5 file
dset_temps (str) – path to the desired dataset in the HDF5 file
dset_data
const_extrapolate (bool) – “linear” or “constant”
transpose_dims (Tuple[int, ...], optional) – axis order for transposition
- Returns:
new instance of the
TabularInterpolation
class- Return type:
- interpolate(temp: float) ndarray
Perform a linear interpolation based on the available tabular data at a given temperature
temp
- Parameters:
temp (float) – temperature point for interpolation
- Returns:
interpolated quantity
- Return type:
np.ndarray
- class thermontfa.ThermoMechNTFA(file_name: str, group_name: str, sig_y: Callable[[float, float, bool], float], N_max: int | None = None, tol: float = 0.0001, verbose: bool = False)
Bases:
object
Thermo-mechanical NTFA
Represents a material routine that describes the effective behavior of a thermo-elasto-plastic composite material with temperature-dependent material parameters in both phases.
- __init__(file_name: str, group_name: str, sig_y: Callable[[float, float, bool], float], N_max: int | None = None, tol: float = 0.0001, verbose: bool = False) None
Initialize the thermo-mechanical NTFA from an HDF5 file (
*.h5
)Seek the data in HDF5 file named
file_name
within the groupgroup_name
. The following data sets containing tabular data for the NTFA are expected in the group: - temperatures: list of temperature points of the tabular data, shape:(N_temp,)
-A_bar
: shape:(N_temp, N_modes, 6)
-C_bar
: shape:(N_temp, 6, 6)
-D_xi
: shape:(N_temp, N_modes, N_modes)
-tau_theta
: shape:(N_temp, 6)
-tau_xi
: shape:(N_temp, N_modes)
In addition, the group in the HDF5 file must contain the following data sets: -
v_frac
: volume fraction of the different phases -SIG_phases
: stress data different phases- Parameters:
file_name (str) – path to the HDF5 file
group_name (str) – group in the HDF5 file that contains the NTFA tabular data
sig_y (Callable[[float, float, bool], float]) – function/callable that returns the yield stress
sig_y(theta, q_n, derivative)
given the temperaturetheta
and the current isotropic hardening variableq_n
. Ifderivative = True
, the derivative should also be returned.N_max (int) – maximum number of NTFA modes that should be used. If None, all available modes are used.
verbose (bool) – If debug information should be printed.
- file_name: str
- group_name: str
- sig_y: Callable[[float, float, bool], float]
- tol: float
- verbose: bool = False
- A_bar: TabularInterpolation
- C_bar: TabularInterpolation
- D_xi: TabularInterpolation
- tau_theta: TabularInterpolation
- tau_xi: TabularInterpolation
- sig_phases: List
- interpolate(theta: float) None
Interpolate NTFA matrices to current temperature
theta
if the given tolerance is exceeded- Parameters:
theta – Temperature
- stress(eps: ndarray, theta: float, xi: ndarray, i_phase: int | None = None)
Compute the stress given strain
eps
, plastic mode activationsxi
. Ifi_phase
is given, the stress is computed only for the phase with indexi_phase
.- Parameters:
eps (int, optional) – Strain
theta – Temperature
xi – Plastic mode activations
i_phase – Phase index for the stress computation. If
None
, overall stress is computed.
- Returns:
Computed stress
- UMAT_mixed(eps_idx: ndarray, eps_n: ndarray, deps: ndarray, sig_bc: ndarray, theta: float, q_n: float, xi_n: ndarray) Tuple[ndarray, ndarray, ndarray, float, ndarray]
Run the UMAT using partial eps-BC, e.g., uniaxial tension test.
- Parameters:
eps_idx (nd-array, dtype=int) – Indices of eps that are prescribed. If
None
, then all components ofsig_bc
are prescribed.eps_n (nd-array, dtype=float, shape=[6]) – Strain at the beginning of the increment.
deps (nd-array, dtype=float, shape=[6]) – Strain increment. Only deps[eps_idx] is used.
sig_bc (nd-array, dtype=float, shape=[6]) – Stress at the end of the increment. Only non
eps_idx
components are used.theta (float) – Temperature at the end of the time increment.
q_n (float) – Hardening variable at the beginning of the time increment.
xi_n (nd-array, dtype=float) – Reduced coefficients at the beginning of the time increment.
- Returns:
E - Full strain tensor at the end of the increment, i.e. the entries not within eps_idx are set
- Return type:
np.ndarray, dtype=float, shape=[6]
- Returns:
S - Stress at the end of the increment. Only S[eps_idx] is non-zero.
- Return type:
np.ndarray, dtype=float
- Returns:
C - Stiffness tensor at the end of the increment
- Return type:
np.ndarray, dtype=float, shape=[6, 6]
- Returns:
q - Hardening variable at the end of the time increment
- Return type:
float
- Returns:
xi - Reduced coefficients at the end of the time increment
- Return type:
nd.ndarray, dtype=float
- solve(eps: ndarray, deps: ndarray, theta: float, q_n: float, xi_n: ndarray) Tuple[ndarray, float, ndarray, ndarray]
Solve for stress
S
, hardening variableq
, reduced coefficientsxi
, and stiffnessC
given the straineps
, strain incrementdeps
, temperaturetheta
, hardening variableq_n
, and reduced coefficientsxi_n
.- Parameters:
eps (np.ndarray) – Strain
deps (np.ndarray) – Strain increment
theta (float) – Temperature
q_n (float) – Hardening variable at the beginning of the time increment
xi_n (np.ndarray) – Reduced coefficients at the beginning of the time increment
- Returns:
S - Stress
- Return type:
np.ndarray
- Returns:
q - Hardening variable at the end of the time increment
- Return type:
float
- Returns:
xi - Reduced coefficients at the end of the time increment
- Return type:
np.ndarray
- Returns:
C - Stiffness
- Return type:
np.ndarray