TFEL,
MFront and MTestTFEL/Utilities
features
TFEL/Math features
TFEL/Material
features
MTestCurrentState.copy() produces a shallow copy in python
bindings[mfront]
allow to specify dsl options in MFrontBase::getDSL[cmake]
export compile flags in dedicated variables[cmake]
better handling of dependencies in exported cmake
files@UseQt
keyword is not mentioned in the MaterialLaw’s keywords help
pageC++ digit separatorThe page describes the new functionalities of Version 5.0.0 of the
TFEL project.
Version 5.0.0 has been released on December, 17th 2024. This version inherits from the bug fixes of:
There are currently \(20\,576\) unit tests.
| Operating system | Compiler | Configuration | Successful Tests |
|---|---|---|---|
| Ubuntu 22.04 LTS | gcc 11.4.0 | Release | 100% |
| gcc 11.4.0 | Release + fast-math | 100% | |
| gcc 11.4.0 | Debug | 100% | |
| clang 18.1.8 | Release | 100% | |
| clang 19.1.0 | Release | 100% | |
| clang 19.1.0 | Debug | 100% | |
| intel 2024.2.1 | Release | 99% | |
| Fedora 41 | gcc 14.2.1 | Release | 100% |
| clang 17.0.6 | Release | 100% | |
| Debian 12 | gcc 12.2.0 | Release | 100% |
| Ubuntu 24.04 LTS | gcc 13.2.0 | Release | 100% |
| Mac Os | gcc | Release | 100% |
| Windows 10 | Visual Studio 2022 | Release | 95% |
Versions of clang prior to Version \(17\) are known not to work.
TFEL/Utilities
featuresC++
digit separator in CxxTokenizerThe CxxTokenizer class now supports C++
digit separators in numbers, (integer of floatting point numbers), as
illustrated in the following example:
@DSL IsotropicMisesCreep;
@Behaviour DigitSeparatorTest;
@UseQt true;
@Epsilon 0.000'000'000'1;
@IterMax 1'000;
@ElasticMaterialProperties{150'000'000'000, 0.3};
@Parameter stress s0 = 10'000'000;
@Parameter strainrate de0 = 8.230'512e-67;
@Parameter strainrate E = 8.2;
@FlowRule{
f = de0 * pow(seq / s0, E);
}TFEL/Math
featuresThe product of two tiny matrices has been implemented:
const auto m1 = tmatrix<2u, 2u, int>{0, 1, //
2, 3};
const auto m2 = tmatrix<2u, 2u, int>{4, -1, //
5, 2};
const auto m3 = m1 * m2;A new eigen solver based on Harari’s analytical solution have been
introduced for symmetric tensors. The computation of eigen values is
done with Harari’s algorithm [1]
and the computation of eigen vectors is done with the default eigen
solver for symmetric tensors of TFEL. Such computations are
more efficient and more accurate than the default TFEL
algorithm.
Those algorithms are available in 3D. For 2D symmetric tensors, we fall back to some default algorithm as described below.
| Name | Algorithm in 3D | Algorithm in 2D |
|---|---|---|
TFELEIGENSOLVER |
Analytical (TFEL) | Analytical (TFEL) |
FSESJACOBIEIGENSOLVER |
Jacobi | Analytical (FSES) |
FSESQLEIGENSOLVER |
QL with implicit shifts | Analytical (FSES) |
FSESCUPPENEIGENSOLVER |
Cuppen’s Divide & Conquer | Analytical (FSES) |
FSESANALYTICALEIGENSOLVER |
Analytical | Analytical (FSES) |
FSESHYBRIDEIGENSOLVER |
Hybrid | Analytical (FSES) |
GTESYMMETRICQREIGENSOLVER |
Symmetric QR | Analytical (TFEL) |
HARARIEIGENSOLVER |
Analytical (Harari) | Analytical (TFEL) |
The various eigen solvers available are enumerated in Table 1.
The eigen solver is passed as a template argument of the
computeEigenValues or the computeEigenVectors
methods as illustrated in the code below:
tmatrix<3u,3u,real> m2;
tvector<3u,real> vp2;
std::tie(vp,m)=s.computeEigenVectors<stensor::HARARIEIGENSOLVER>();| Algorithm | Failure ratio | \(\Delta_{\infty}\) | Times (ns) | Time ratio |
|---|---|---|---|---|
TFELEIGENSOLVER |
0.000557 | 2.37e-05 | 129452335 | 1 |
GTESYMMETRICQREIGENSOLVER |
0 | 9.57e-07 | 236544828 | 1.83 |
FSESJACOBIEIGENSOLVER |
0 | 4.61e-07 | 241131631 | 1.86 |
FSESQLEIGENSOLVER |
0.000173 | 1.67e-06 | 155435496 | 1.20 |
FSESCUPPENEIGENSOLVER |
0.018223 | 2.87e-06 | 151727321 | 1.17 |
FSESHYBRIDEIGENSOLVER |
0.068411 | 3.90e-03 | 80039266 | 0.62 |
FSESANALYTICALEIGENSOLVER |
0.102626 | 6.21e-02 | 76865961 | 0.59 |
HARARIEIGENSOLVER |
0.000018 | 2.46e-06 | 110028802 | 0.85 |
| Algorithm | Failure ratio | \(\Delta_{\infty}\) | Times (ns) | Time ratio |
|---|---|---|---|---|
TFELEIGENSOLVER |
0.00058 | 6.94e-14 | 137752068 | 1 |
GTESYMMETRICQREIGENSOLVER |
1e-06 | 2.30e-15 | 315593552 | 2.29 |
FSESJACOBIEIGENSOLVER |
0 | 9.08e-16 | 256285090 | 1.86 |
FSESQLEIGENSOLVER |
0.000202 | 3.04e-15 | 214537012 | 1.56 |
FSESCUPPENEIGENSOLVER |
0.019251 | 5.58e-15 | 219113965 | 1.59 |
FSESHYBRIDEIGENSOLVER |
0.081586 | 1.29e-10 | 81861668 | 0.59 |
FSESANALYTICALEIGENSOLVER |
0.103935 | 4.11e-10 | 79701256 | 0.58 |
HARARIEIGENSOLVER |
0.000037 | 2.27e-14 | 116977683 | 0.85 |
| Algorithm | Failure ratio | \(\Delta_{\infty}\) | Times (ns) | Time ratio |
|---|---|---|---|---|
TFELEIGENSOLVER |
0.000578 | 1.76e-17 | 304165023 | 1 |
GTESYMMETRICQREIGENSOLVER |
0 | 1.01e-18 | 558605772 | 1.84 |
FSESJACOBIEIGENSOLVER |
0 | 5.11e-19 | 408584703 | 1.34 |
FSESQLEIGENSOLVER |
0.00045 | 1.83e-18 | 311028180 | 1.02 |
FSESCUPPENEIGENSOLVER |
0.009134 | 3.23e-18 | 485590150 | 1.60 |
FSESHYBRIDEIGENSOLVER |
0.99959 | 4.01e-10 | 187308886 | 0.62 |
FSESANALYTICALEIGENSOLVER |
0.999669 | 1.36e-11 | 211710377 | 0.70 |
HARARIEIGENSOLVER |
0.000046 | 4.62e-18 | 338589049 | 1.11 |
We have compared the available algorithm on \(10^{6}\) random symmetric tensors whose components are in \([-1:1]\).
For a given symmetric tensor, we consider that the computation of the eigenvalues and eigenvectors failed if: \[ \Delta_{\infty}=\max_{i\in[1,2,3]}\left\|\underline{s}\,\cdot\,\vec{v}_{i}-\lambda_{i}\,\vec{v}_{i}\right\|>10\,\varepsilon \] where \(\varepsilon\) is the accuracy of the floatting point considered.
The results of those tests are reported on Tables 2, 3 and 4. The
Harari eigen solver offers a better compromise between accuracy and
numerical efficiency than the default TFEL solver.
TFEL/Material
featuresThe homogenization functions are part of the namespace
tfel::material::homogenization. A specialization for
elasticity is defined:
tfel::material::homogenization::elasticity, in prevision of
further homogenization developments in other physics.
The function computeEshelbyTensor computes the Eshelby
tensor of an ellipsoid whose semi-axis lengths are a,
b, c, embedded in an isotropic matrix. There
is also computeSphereEshelbyTensor,
computeAxisymmetricalEshelbyTensor, and also
computeCircularCylinderEshelbyTensor and
computeEllipticCylinderEshelbyTensor for plane strain
elasticity.
Three functions also compute the strain localisation tensor of an
ellipsoid embedded in an isotropic matrix and submitted to an external
uniform strain field : computeEllipsoidLocalisationTensor,
computeAxisymmetricalEllipsoidLocalisationTensor and
computeSphereLocalisationTensor.
Different schemes are implemented and return the homogenized
stiffness of the material. The available functions are
computeMoriTanakaScheme, computeDiluteScheme,
computeSphereDiluteScheme, computeSphereMoriTanakaScheme. If a
distribution of ellipsoids is considered, three types of distributions
are considered. The corresponding functions are
computeIsotropicDiluteScheme,
computeTransverseIsotropicDiluteScheme,
computeOrientedDiluteScheme,
computeIsotropicMoriTanakaScheme,
computeTransverseIsotropicMoriTanakaScheme and
computeOrientedMoriTanakaScheme.
MaterialProperty DSL@Data
keywordThe @Data keyword allow the definition of a material
properties using the interpolation of a set of values.
The @Data keyword is followed by a set of options
defining:
The precise syntax depends on the number of inputs of the material property.
The only option accepted is value.
The values option is required. It must be a map
associating values of the input and values of the output.
The interpolation option is optional. It must be a
string. The values linear and cubic_spline are
accepted.
The extrapolation option is optional. It must be a
boolean or a string:
true means that a linear
extrapolation will be performed.false or the string value
bound_to_last_value or the string value
constant means that a constant value will be used for
extrapolation.@DSL MaterialProperty;
@Law LinearDataInterpolation;
@UseQt true;
@UnitSystem SI;
@Output stress E;
E.setGlossaryName("YoungModulus");
@StateVariable temperature T;
T.setGlossaryName("Temperature");
@Data {
values: { 293.15 : 240e9, 693.15 : 180e9, 893.15 : 170e9 },
interpolation : "linear"
};StandardElastoViscoPlasticity brickThe StrainRateSensitive isotropic hardening rule is
defined by: \[
R{\left(p, \dot{p}\right)} =
R_{0}{\left(p\right)}\,R_{rs}{\left(\dot{p}\right)},
\]
where:
This isotropic hardening rule can be parametrised using the following options:
rate_independent_isotropic_hardening: this option
introduces a contribution to \(R_{0}{\left(p\right)}\). This option can be
repeated multiple times.rate_sensitivity_factor: this option introduces the
rate sensivity factor \(R_{rs}{\left(\dot{p}\right)}\). The list of
available rate sensivity factors is given in section 5.2.1.2.@Brick StandardElastoViscoPlasticity{
stress_potential : "Hooke" {young_modulus : 210e9, poisson_ratio : 0.3},
inelastic_flow : "Plastic" {
criterion : "Mises",
isotropic_hardening : "StrainRateSensitive" {
rate_independent_isotropic_hardening :
"Swift" {R0 : "alpha * Ks * (p0s ** ns)", p0 : "p0s", n : "ns"},
rate_independent_isotropic_hardening : "Voce" {
R0 : "(1 - alpha) * Q1 * (1 - Q2)",
Rinf : "(1 - alpha) * Q1",
b : "Q3"
},
rate_sensitivity_factor :
"CowperSymonds" {dp0 : "dp0cs", n : "ncs", Rs_eps : 1e-8}
}
}
};Cowper-Symonds’s rate sensitivity factor: \[ R_{rs}{\left(\dot{p}\right)}=1+A\,\left(\frac{\dot{p}}{\dot{p}_{0}}\right)^{n} \] When the exponent \(n\) is lower than one, the following regularised version, based on the user defined value \(R_{\epsilon}\) is available: \[ R_{rs}{\left(\dot{p}\right)}=1+ \left\{ \begin{aligned} R_{\epsilon}\,\left(\frac{\dot{p}}{\dot{p}_{\epsilon}}\right) &&\text{if}&&\dot{p}<\dot{p}_{\epsilon}\\ A\,\left(\frac{\dot{p}}{\dot{p}_{0}}\right)^{n} &&\text{if}&&\dot{p}\geq\dot{p}_{\epsilon} \end{aligned} \right. \] where \(\dot{p}_{\epsilon}\) is such that \(R_{\epsilon} = A\,\left(\frac{\dot{p}_{\epsilon}}{\dot{p}_{0}}\right)^{n}\)
Johnson-Cook’s rate sensitivity factor: \[ R_{rs}{\left(\dot{p}\right)}=1+ \left\{ \begin{aligned} 0 &&\text{if}&&\dot{p}<\dot{p}_{0}\\ A\,\log\left(\frac{\dot{p}}{\dot{p}_{0}}\right)&&\text{if}&&\dot{p}\geq\dot{p}_{0} \end{aligned} \right. \]
disable_runtime_checks optionIf this option is set to true, interfaces may disable as
much runtime checks as possible. Those runtime checks include checking
standard bounds and physical bounds for instance.
generic interface
improvements@SelectedModellingHypothesis and
@SelectedModellingHypotheses keywordsThe @SelectedModellingHypothesis and
@SelectedModellingHypotheses keywords allows to select
which modelling hypotheses will be generated.
Their syntaxes are similar to the keywords
@ModellingHypothesis and @ModellingHypotheses,
respectively.
The selected modelling hypotheses must be a sub-set of the modelling hypotheses supported by the behaviour.
$ mfront --obuild --interface=generic --@SelectedModellingHypothesis=PlaneStrain Plasticity.mfront
Treating target : all
The following library has been built :
- libBehaviour.so : Plasticity_PlaneStrainThe page Libaries usage in C++
describe how to use the TFEL libraries in C++
projects, using either the tfel-config utility or
cmake packages.
Allows mfront-doc to be used with the MaterialLaw DSL.
Example of use:
$ mfront-doc T91MartensiticSteel_gamma1_ROUX2007.mfrontFor more details, see https://github.com/thelfer/tfel/issues/535
Adds the correspondence between variable names and the 4-letter names
used in castem to the example .dgibi file, generated by
MFront with castem interfaces. The mapping is written in
the file as a comment, for example in the form :
** List of material properties:
**
** - YOUN: YoungModulus
** - NU: PoissonRatio
** - RHO: MassDensity
** - H: HardeningSlope
** - SO: Yield strength
For more details, see https://github.com/thelfer/tfel/issues/582
$ mfront-query --state-variables Inconel600_YoungModulus.mfront
- Temperature (TK): the temperature
$ mfront-query --parameters Inconel600_YoungModulus.mfront$ mfront-query --state-variables --dsl-option='overriding_parameters:{"TK": 400}}' Inconel600_YoungModulus.mfront
$ mfront-query --parameters --dsl-option='overriding_parameters:{"TK": 400}}' Inconel600_YoungModulus.mfront
- Temperature (TK): the temperatureFor more details, see https://github.com/thelfer/tfel/issues/585
MTestCurrentState.copy() produces a shallow copy in python
bindingsThis behaviour is consistent with the copy constructor of the
StudyCurrentState class.
To make a deep copy of this object, the makeDeepCopy
method has been introduced.
For more details, see https://github.com/thelfer/tfel/issues/567
[mfront] allow to specify dsl options in
MFrontBase::getDSLFor more details, see https://github.com/thelfer/tfel/issues/557
[cmake]
export compile flags in dedicated variablesFor more details, see https://github.com/thelfer/tfel/issues/556
[cmake]
better handling of dependencies in exported cmake
filesFor more details, see https://github.com/thelfer/tfel/issues/556
@UseQt keyword is not mentioned in the
MaterialLaw’s keywords help pageFor more details, see https://github.com/thelfer/tfel/issues/526
For more details, see https://github.com/thelfer/tfel/issues/476
C++ digit separatorFor more details, see https://github.com/thelfer/tfel/issues/370