1% New functionalities of the 3.2 version of `TFEL`, `MFront` and `MTest` 2% Thomas Helfer 3% 2017 4 5\newcommand{\absvalue}[1]{{\left|#1\right|}} 6\newcommand{\Frac}[2]{\displaystyle\frac{\displaystyle #1}{\displaystyle #2}} 7\newcommand{\paren}[1]{\left(#1\right)} 8\newcommand{\deriv}[2]{\Frac{\partial #1}{\partial #2}} 9\newcommand{\tenseur}[1]{\underline{#1}} 10\newcommand{\tenseurq}[1]{\underline{\underline{\mathbf{#1}}}} 11\newcommand{\sigmaeq}{\sigma_{\mathrm{eq}}} 12\newcommand{\tsigma}{\underline{\sigma}} 13\newcommand{\trace}[1]{{\mathrm{tr}\paren{#1}}} 14\newcommand{\sigmaH}{\sigma_{H}} 15\newcommand{\tepsilonto}{\underline{\epsilon}^{\mathrm{to}}} 16\newcommand{\tepsilonel}{\underline{\epsilon}^{\mathrm{el}}} 17\newcommand{\tepsilonp}{\underline{\epsilon}^{\mathrm{p}}} 18\newcommand{\tepsilonvp}{\underline{\epsilon}^{\mathrm{vp}}} 19\newcommand{\tepsilonth}{\underline{\epsilon}^{\mathrm{th}}} 20 21The page declares the new functionalities of the 3.2 version of 22the `TFEL` project. 23 24The `TFEL` project is a collaborative development of 25[CEA](http://www.cea.fr/english-portal "Commissariat à l'énergie 26atomique") and [EDF](http://www.edf.com/ "Électricité de France") 27dedicated to material knowledge manangement with special focus on 28mechanical behaviours. It provides a set of libraries (including 29`TFEL/Math` and `TFEL/Material`) and several executables, in 30particular `MFront` and `MTest`. 31 32`TFEL` is available on a wide variety of operating systems and 33compilers. 34 35# Documentation 36 37A new page dedicated to the `python` bindings of the `TFEL` libraries 38is available [here](tfel-python.html). 39 40# Updates in `TFEL` libraries 41 42The `TFEL` project provides several libraries. This paragraph is about 43updates made in those libraries. 44 45## TFEL/Math 46 47### Improvements to the `Evaluator` class 48 49The `Evaluator` class is used to interpret textual formula, as 50follows: 51 52~~~~{.cxx} 53Evalutator ev("sin(x)"); 54ev.setVariableValue("x",12); 55const auto s = ev.getValue(); 56~~~~ 57 58#### An overload for the `getValue` method 59 60In the previous example, each variable value had to be set using the 61`setVariableValue` method. The new overloaded version of the `getValue` 62method can take a map as argument as follows: 63 64~~~~{.cxx} 65Evalutator ev("sin(x)"); 66const auto s = ev.getValue({{"x",12}}); 67~~~~ 68 69#### `Operator()` 70 71Two overloaded versions of the `Evaluator::operator()` has been 72introduced as a synonyms for the `getValue` method: 73 74~~~~{.cxx} 75Evalutator ev("sin(x)"); 76const auto s = ev({{"x",12}}); 77~~~~ 78 79#### The `getCxxFormula` method 80 81The `getCxxFormula` method returns a string representing the evaluation 82of the formula in standard `C++`. This method takes a map as argument 83which describes how certain variables shall be represented. This method 84can be used, as follows: 85 86~~~~{.cxx} 87Evalutator ev("2*sin(x)"); 88std::cout << ev({"x":"this->x"}}) << '\n'; 89~~~~ 90 91The previous code displays: 92 93~~~~{.sh} 94sin((2)*(this->x)) 95~~~~ 96 97This function is the basis of a new functionality of the `MFront` code 98generator (inline material properties), see 99Section @sec:inline-mprops for details. 100 101#### New mathematical functions 102 103The following new mathematical functions have been introduced: 104 105- `exp2`: returns the base-2 exponential. 106- `expm1`: returns the base-e exponential minus one. 107- `cbrt`: returns the cubic root 108- `log2`: computes the base-2 logarithm of the argment. 109- `log1p`: computes the logarithm of the argment plus one. 110- `acosh`: computes the inverse of the hyperbolic cosine. 111- `asinh`: computes the inverse of the hyperbolic sine. 112- `atanh`: computes the inverse of the hyperbolic tangent. 113- `erf`: computes the error function. 114- `erfc`: computes the complementary error function. 115- `tgamma`: computes the Gamma function. 116- `lgamma`: compute the natural logarithm of the gamma function. 117- `hypot`: returns the hypotenuse of a right-angled triangle whose 118 legs are x and y, i.e. computes \(\sqrt{x^{2}+y^{2}}. 119- `atan2`: returns the principal value of the arc tangent of \(y/x\), 120 expressed in radians. 121 122### Efficient computations of the first and second derivatives of the invariants of the stress deviator tensor with respect to the stress {#sec:deviatoric:invariants} 123 124Let \(\tsigma\) be a stress tensor. Its deviatoric part 125\(\tenseur{s}\) is: 126 127\[ 128\tenseur{s}=\tsigma-\Frac{1}{3}\,\trace{\tsigma}\,\tenseur{I} 129=\paren{\tenseurq{I}-\Frac{1}{3}\,\tenseur{I}\,\otimes\,\tenseur{I}}\,\colon\,\tsigma 130\] 131 132The deviator of a tensor can be computed using the `deviator` 133function. 134 135As it is a second order tensor, the stress deviator tensor also has a 136set of invariants, which can be obtained using the same procedure used 137to calculate the invariants of the stress tensor. It can be shown that 138the principal directions of the stress deviator tensor \(s_{ij}\) are 139the same as the principal directions of the stress tensor 140\(\sigma_{ij}\). Thus, the characteristic equation is 141 142\[ 143\left| s_{ij}- \lambda\delta_{ij} \right| = -\lambda^3+J_1\lambda^2-J_2\lambda+J_3=0, 144\] 145 146where \(J_1\), \(J_2\) and \(J_3\) are the first, second, and third 147*deviatoric stress invariants*, respectively. Their values are the same 148(invariant) regardless of the orientation of the coordinate system 149chosen. These deviatoric stress invariants can be expressed as a 150function of the components of \(s_{ij}\) or its principal values \(s_1\), 151\(s_2\), and \(s_3\), or alternatively, as a function of \(\sigma_{ij}\) or 152its principal values \(\sigma_1\), \(\sigma_2\), and \(\sigma_3\). Thus, 153 154\[ 155\begin{aligned} 156J_1 &= s_{kk}=0,\, \\ 157J_2 &= \textstyle{\frac{1}{2}}s_{ij}s_{ji} = \Frac{1}{2}\trace{\tenseur{s}^2}\\ 158&= \Frac{1}{2}(s_1^2 + s_2^2 + s_3^2) \\ 159&= \Frac{1}{6}\left[(\sigma_{11} - \sigma_{22})^2 + (\sigma_{22} - \sigma_{33})^2 + (\sigma_{33} - \sigma_{11})^2 \right ] + \sigma_{12}^2 + \sigma_{23}^2 + \sigma_{31}^2 \\ 160&= \Frac{1}{6}\left[(\sigma_1 - \sigma_2)^2 + (\sigma_2 - \sigma_3)^2 + (\sigma_3 - \sigma_1)^2 \right ] \\ 161&= \Frac{1}{3}I_1^2-I_2 = \frac{1}{2}\left[\trace{\tenseur{\sigma}^2} - \frac{1}{3}\trace{\tenseur{\sigma}}^2\right],\,\\ 162J_3 &= \det\paren{\tenseur{s}} \\ 163&= \Frac{1}{3}s_{ij}s_{jk}s_{ki} = \Frac{1}{3} \trace{\tenseur{s}^3}\\ 164&= \Frac{1}{3}(s_1^3 + s_2^3 + s_3^3) \\ 165&= s_1s_2s_3 \\ 166&= \Frac{2}{27}I_1^3 - \Frac{1}{3}I_1 I_2 + I_3 = \Frac{1}{3}\left[\trace{\tenseur{\sigma}^3} - \trace{\tenseur{\sigma}^2}\trace{\tenseur{\sigma}} +\Frac{2}{9}\trace{\tenseur{\sigma}}^3\right]. 167\end{aligned} 168\] 169 170where \(I_{1}\), \(I_{2}\) and \(I_{3}\) are the invariants of 171\(\tsigma\). 172 173\(J_{2}\) and \(J_{3}\) are building blocks for many isotropic yield 174critera. Classically, \(J_{2}\) is directly related to the von Mises 175stress \(\sigmaeq\): 176 177\[ 178\sigmaeq=\sqrt{\Frac{3}{2}\,\tenseur{s}\,\colon\,\tenseur{s}}=\sqrt{3\,J_{2}} 179\] 180 181The first and second derivatives of \(J_{2}\) with respect to 182\(\sigma\) can be trivially implemented, as follows: 183 184~~~~{.cpp} 185constexpr const auto id = stensor<N,real>::Id(); 186constexpr const auto id4 = st2tost2<N,real>::Id(); 187// first derivative of J2 188const auto dJ2 = deviator(sig); 189// second derivative of J2 190const auto d2J2 = eval(id4-(id^id)/3); 191~~~~ 192 193In comparison, the computation of the first and second derivatives of 194\(J_{3}\) with respect to \(\sigma\) are more cumbersome. In previous 195versions `TFEL`, one had to write: 196 197~~~~{.cpp} 198constexpr const auto id = stensor<N,real>::Id(); 199constexpr const auto id4 = st2tost2<N,real>::Id(); 200const auto I1 = trace(sig); 201const auto I2 = (I1*I1-trace(square(sig)))/2; 202const auto dI2 = I1*id-sig; 203const auto dI3 = computeDeterminantDerivative(sig); 204const auto d2I2 = (id^id)-id4; 205const auto d2I3 = computeDeterminantSecondDerivative(sig); 206// first derivative of J3 207const auto dJ3 = eval((2*I1*I1/9)*id-(I2*id+I1*dI2)/3+dI3); 208// second derivative of J3 209const auto d2J3 = eval((4*I1/9)*(id^id)-((id^dI2)+(dI2^id)+i1*d2I2)/3+d2I3); 210~~~~ 211 212More efficient implementations are now available using the 213`computeDeviatorDeterminantDerivative` and 214`computeDeviatorDeterminantSecondDerivative` functions: 215 216~~~~{.cpp} 217// first derivative of J3 218const auto dJ3 = computeDeviatorDeterminantDerivative(sig); 219// second derivative of J3 220const auto d2J3 = computeDeviatorDeterminantSecondDerivative(sig); 221~~~~ 222 223## TFEL/Material 224 225### Isotropic Plasticity 226 227By definition, \(J_{2}\) and \(J_{3}\) are the second and third 228invariants of the deviatoric part \(\tenseur{s}\) of the stress tensor 229\(\tsigma\) (see also Section @sec:deviatoric:invariants): 230 231\[ 232\left\{ 233\begin{aligned} 234J_2 &= \Frac{1}{2}\trace{\tenseur{s}^2}\\ 235J_3 &= \det(\tenseur{s}) \\ 236\end{aligned} 237\right. 238\] 239 240The first and second derivatives of \(J_{2}\) with respect to the 241stress tensor \(\tsigma\) are trivially computed and implemented (see 242Section @sec:deviatoric:invariants). 243 244The first and second derivatives of \(J_{2}\) with respect to the 245stress tensor \(\tsigma\) can be computed respectively by: 246 247- The `computesJ3Derivative` function, which is a synonym for the 248 `computeDeviatorDeterminantDerivative` function defined in the 249 `tfel::math` namespace (see 250 Section @sec:deviatoric:invariants for details). 251- The `computeJ3SecondDerivative` function, which is a synonym for the 252 `computeDeviatorDeterminantSecondDerivative` function defined in the 253 `tfel::math` namespace (see Section @sec:deviatoric:invariants 254 for details). 255 256### Orthotropic plasticity 257 258Within the framework of the theory of representation, generalizations 259to anisotropic conditions of the invariants of the deviatoric stress 260have been proposed by Cazacu and Barlat (see 261@cazacu_generalization_2001): 262 263- The generalization of \(J_{2}\) is denoted \(J_{2}^{O}\). It is 264 defined by: 265 \[ 266 J_{2}^{O}= a_6\,s_{yz}^2+a_5\,s_{xz}^2+a_4\,s_{xy}^2+\frac{a_2}{6}\,(s_{yy}-s_{zz})^2+\frac{a_3}{6}\,(s_{xx}-s_{zz})^2+\frac{a_1}{6}\,(s_{xx}-s_{yy})^2 267 \] 268 where the \(\left.a_{i}\right|_{i\in[1:6]}\) are six coefficients 269 describing the orthotropy of the material. 270- The generalization of \(J_{3}\) is denoted \(J_{3}^{O}\). It is 271 defined by: 272 \[ 273 \begin{aligned} 274 J_{3}^{O}= 275 &\frac{1}{27}\,(b_1+b_2)\,s_{xx}^3+\frac{1}{27}\,(b_3+b_4)\,s_{yy}^3+\frac{1}{27}\,(2\,(b_1+b_4)-b_2-b_3)\,s_{zz}^3\\ 276 &-\frac{1}{9}\,(b_1\,s_{yy}+b_2s_{zz})\,s_{xx}^2\\ 277 &-\frac{1}{9}\,(b_3\,s_{zz}+b_4\,s_{xx})\,s_{yy}^2\\ 278 &-\frac{1}{9}\,((b_1-b_2+b_4)\,s_{xx}+(b_1-b3+b_4)\,s_{yy})\,s_{zz}^3\\ 279 &+\frac{2}{9}\,(b_1+b_4)\,s_{xx}\,s_{yy}\,s_{zz}\\ 280 &-\frac{s_{xz}^2}{3}\,(2\,b_9\,s_{yy}-b_8\,s_{zz}-(2\,b_9-b_8)\,s_{xx})\\ 281 &-\frac{s_{xy}^2}{3}\,(2\,b_{10}\,s_{zz}-b_5\,s_{yy}-(2\,b_{10}-b_5)\,s_{xx})\\ 282 &-\frac{s_{yz}^2}{3}\,((b_6+b_7)\,s_{xx}-b_6\,s_{yy}-b_7\,s_{zz})\\ 283 &+2\,b_{11}\,s_{xy}\,s_{xz}\,s_{yz} 284 \end{aligned} 285 \] 286 where the \(\left.b_{i}\right|_{i\in[1:11]}\) are eleven coefficients 287 describing the orthotropy of the material. 288 289Those invariants may be used to generalize isotropic yield criteria 290based on \(J_{2}\) and \(J_{3}\) invariants to orthotropy. 291 292The following functions 293 294\(J_{2}^{0}\), \(J_{3}^{0}\) and their first and second derivatives 295with respect to the stress tensor \(\tsigma\) can be computed 296by the following functions: 297 298- `computesJ2O`, `computesJ2ODerivative` and 299 `computesJ2OSecondDerivative`. 300- `computesJ3O`, `computesJ3ODerivative` and 301 `computesJ3OSecondDerivative`. 302 303Those functions take the stress tensor as first argument and each 304orthotropic coefficients. Each of those functions has an overload 305taking the stress tensor as its firs arguments and a tiny vector 306(`tfel::math::tvector`) containing the orthotropic coefficients. 307 308### \(\pi\)-plane 309 310![Comparison of the isosurfaces of various equivalent stresses (Tresca, von Mises, Hosford \(a=8\)) in the \(\pi\)-plane](img/IsotropicEquivalentStressesInPiPlane.svg 311 "Comparison of the isosurfaces of various equivalent stresses 312 (Tresca, von Mises, Hosford \(a=8\)) in the 313 \(\pi\)-plane"){width=70%} 314 315The \(\pi\)-plane is defined in the space defined by the three 316eigenvalues \(S_{0}\), \(S_{1}\) and \(S_{2}\) of the stress by the 317following equations: 318\[ 319S_{0}+S_{1}+S_{2}=0 320\] 321 322This plane contains deviatoric stress states and is perpendicular to 323the hydrostatic axis. A basis of this plane is given by the following 324vectors: 325\[ 326\vec{n}_{0}= 327\frac{1}{\sqrt{2}}\, 328\begin{pmatrix} 3291 \\ 330-1 \\ 3310 332\end{pmatrix} 333\quad\text{and}\quad 334\vec{n}_{1}= 335\frac{1}{\sqrt{6}}\, 336\begin{pmatrix} 337-1 \\ 338-1 \\ 339 2 340\end{pmatrix} 341\] 342 343This plane is used to characterize the iso-values of equivalent 344stresses which are not sensitive to the hydrostatic pression. 345 346Various functions are available: 347 348- `projectOnPiPlane`: this function projects a stress state on the 349 \(\pi\)-plane. 350- `buildFromPiPlane`: this function builds a stress state, defined by 351 its three eigenvalues, from its coordinate in the \(\pi\)-plane. 352 353## `python` bindings 354 355### `tfel.math` module 356 357#### Updated bindings for the `stensor` class 358 359The following operations are supported: 360 361- addition of two symmetric tensors. 362- substraction of two symmetric tensors. 363- multiplication by scalar. 364- in-place addition by a symmetric tensor. 365- in-place substraction by a symmetric tensor. 366- in-place multiplication by scalar. 367- in-place division by scalar. 368 369The following functions have been introduced: 370 371- `makeStensor1D`: builds a \(1D\) symmetric tensor from a tuple of 372 three values. 373- `makeStensor2D`: builds a \(2D\) symmetric tensor from a tuple of 374 three values. 375- `makeStensor3D`: builds a \(3D\) symmetric tensor from a tuple of 376 three values. 377 378### `tfel.material` module 379 380#### Bindings related to the \(\pi\)-plane 381 382The following functions are available: 383 384- `buildFromPiPlane`: returns a tuple containing the three eigenvalues 385 of the stress corresponding to the given point in the \(\pi\)-plane. 386- `projectOnPiPlane`: projects a stress state, defined its three 387 eigenvalues or by a symmetric tensor, on the \(\pi\)-plane. 388 389The following script shows how to build an isosurface of the von Mises 390equivalent stress in the \(\pi\)-plane: 391 392~~~~{.python} 393from math import pi,cos,sin 394import tfel.math as tmath 395import tfel.material as tmaterial 396 397nmax = 100 398for a in [pi*(-1.+(2.*i)/(nmax-1)) for i in range(0,nmax)]: 399 s = tmath.makeStensor1D(tmaterial.buildFromPiPlane(cos(a),sin(a))) 400 seq = tmath.sigmaeq(s); 401 s *= 1/seq; 402 s1,s2 = tmaterial.projectOnPiPlane(s); 403 print(s1,s2); 404~~~~ 405 406#### Bindings related to the Hosford equivalent stress 407 408The `computeHosfordStress` function, which compute the Hosford 409equivalent stress, is available. 410 411The following script shows how to print an iso-surface of the Hosford 412equivalent stress in the \(\pi\)-plane: 413 414~~~~{.python} 415from math import pi,cos,sin 416import tfel.math as tmath 417import tfel.material as tmaterial 418 419nmax = 100 420for a in [pi*(-1.+(2.*i)/(nmax-1)) for i in range(0,nmax)]: 421 s = tmath.makeStensor1D(tmaterial.buildFromPiPlane(cos(a),sin(a))) 422 seq = tmaterial.computeHosfordStress(s,8,1.e-12); 423 s /= seq; 424 s1,s2 = tmaterial.projectOnPiPlane(s); 425 print(s1,s2); 426~~~~ 427 428#### Bindings related to the Barlat equivalent stress 429 430The following functions are available: 431 432- `makeBarlatLinearTransformation1D`: builds a \(1D\) linear 433 transformation of the stress tensor. 434- `makeBarlatLinearTransformation2D`: builds a \(2D\) linear 435 transformation of the stress tensor. 436- `makeBarlatLinearTransformation3D`: builds a \(3D\) linear 437 transformation of the stress tensor. 438- `computeBarlatStress`: computes the Barlat equivalent Barlat stress. 439 440The following script shows how to print an iso-surface of the Barlat 441equivalent stress for the 2090-T3 aluminum alloy in the \(\pi\)-plane 442(see @barlat_linear_2005): 443 444~~~~{.python} 445from math import pi,cos,sin 446import tfel.math as tmath 447import tfel.material as tmaterial 448 449nmax = 100 450l1 = tmaterial.makeBarlatLinearTransformation1D(-0.069888,0.079143,0.936408, 451 0.524741,1.00306,1.36318, 452 0.954322,1.06906,1.02377); 453l2 = tmaterial.makeBarlatLinearTransformation1D(0.981171,0.575316,0.476741, 454 1.14501,0.866827,-0.079294, 455 1.40462,1.1471,1.05166); 456for a in [pi*(-1.+(2.*i)/(nmax-1)) for i in range(0,nmax)]: 457 s = tmath.makeStensor1D(tmaterial.buildFromPiPlane(cos(a),sin(a))) 458 seq = tmaterial.computeBarlatStress(s,l1,l2,8,1.e-12); 459 s *= 1/seq; 460 s1,s2 = tmaterial.projectOnPiPlane(s); 461 print(s1,s2); 462~~~~ 463 464# New functionalities in `MFront` 465 466## The `StandardElastoViscoPlasticity` brick 467 468This brick is used to describe a specific class of strain based 469behaviours based on an additive split of the total strain 470\(\tepsilonto\) into an elastic part \(\tepsilonel\) and an one or 471several inelastic strains describing plastic (time-independent) flows 472and/or viscoplastic (time-dependent) flows: 473\[ 474 \tepsilonto=\tepsilonel 475+\sum_{i_{\mathrm{p}}=0}^{n_{\mathrm{p}}}\tepsilonp_{i_{\mathrm{p}}} 476+\sum_{i_{\mathrm{vp}}=0}^{n_{\mathrm{vp}}}\tepsilonvp_{i_{\mathrm{vp}}} 477\] 478 479This equation defines the equation associated with the elastic strain 480\(\tepsilonel\). 481 482The brick decomposes the behaviour into two components: 483 484- the stress potential which defines the relation between the elastic 485 strain \(\tepsilonel\) and possibly some damage variables and the 486 stress measure \(\tsigma\). As the definition of the elastic 487 properties can be part of the definition of the stress potential, the 488 thermal expansion coefficients can also be defined in the block 489 corresponding to the stress potential. 490- a list of inelastic flows. Inelastic flows are defined by: 491 - a criterion 492 - a flow criterion 493 - a set of isotropic hardening rules 494 - a set of kinematic hardening rules 495 496Here is a complete usage: 497 498~~~~{.cpp} 499@Brick "StandardElastoViscoPlasticity" { 500 // Here the stress potential is given by the Hooke law. We define: 501 // - the elastic properties (Young modulus and Poisson ratio). 502 // Here the Young modulus is a function of the temperature. 503 // The Poisson ratio is constant. 504 // - the thermal expansion coefficient 505 // - the reference temperature for the thermal expansion 506 stress_potential : "Hooke" { 507 young_modulus : "2.e5 - (1.e5*((T - 100.)/960.)**2)", 508 poisson_ratio : 0.3, 509 thermal_expansion : "1.e-5 + (1.e-5 * ((T - 100.)/960.) ** 4)", 510 thermal_expansion_reference_temperature : 0 511 }, 512 // Here we define only one viscplastic flow defined by the Norton law, 513 // which is based: 514 // - the von Mises stress criterion 515 // - one isotorpic hardening rule based on Voce formalism 516 // - one kinematic hardening rule following the Armstrong-Frederick law 517 inelastic_flow : "Norton" { 518 criterion : "Mises", 519 isotropic_hardening : "Voce" {R0 : 200, Rinf : 100, b : 20}, 520 kinematic_hardening : "Armstrong-Frederick" { 521 C : "1.e6 - 98500 * (T - 100) / 96", 522 D : "5000 - 5* (T - 100)" 523 }, 524 K : "(4200. * (T + 20.) - 3. * (T + 20.0)**2)/4900.", 525 n : "7. - (T - 100.) / 160.", 526 Ksf : 3 527 } 528}; 529~~~~ 530 531## Logarithmic strain framework support in the `cyrano` interface 532 533`Cyrano` is the state of the art fuel performance code developed by EDF 534(See [@thouvenin_edf_2010;@petry_advanced_2015]). 535 536### Extension of `Cyrano` to finite strain analysis 537 538As Cyrano, provides a mono-dimensional description of the fuel rod, its 539extension to finite strain follows the same treatment used for: 540 541- the extension of `MTest` to pipes. 542- the extension of `Alcyone1D` to finite strain analysis. 543 544### Logarithmic strain framework 545 546Miehe et al. have introduced a finite strain framework based on the 547Hencky strain measure (logarithmic strain) and its dual stress (see 548@miehe_anisotropic_2002). This framework has been introduced in 549`Code_Aster` in 2011 (see @edf_modeles_2013). 550 551A behaviour is declared to follow this framework by using the 552`@StrainMeasure` keyword: 553 554~~~~{.cpp} 555@StrainMeasure Hencky; 556~~~~ 557 558The `cyrano` interface now provides support for behaviours based on this 559strain measure. 560 561The interface handles (See @helfer_extension_2015 for details): 562 563- the pre-processing stage (computation of the logarithmic strain from 564 the linearised strain). 565- the post-processing stage (conversion of the dual stress to the first 566 Piola-Kirchoff stress, conversion of the consistent tangent operator), 567 568## The generic behaviour interface 569 570The generic behaviour interface has been introduced to provide an 571interface suitable for most needs. 572 573Whereas other interfaces target a specific solver and thus are 574restricted by choices made by this specific solver, this interface has 575been created for developers of homebrew solvers who are able to modify 576their code to take full advantage of `MFront` behaviours. 577 578This interface is tightly linked with the 579`MFrontGenericInterfaceSupport` project which is available on `github`: 580<https://github.com/thelfer/MFrontGenericInterfaceSupport>. This project 581has a more liberal licence which allows it to be included in both 582commercial and open-source solvers/library. This licensing choice 583explains why this project is not part of the `TFEL`. 584 585## Various improvements 586 587### Inline material properties in mechanical behaviours {#sec:inline-mprops} 588 589Various keywords (such as `@ElasticMaterialProperties`, 590`@ComputeThermalExpansion`, `@HillTensor`, etc.) expects one or more 591material properties. In previous versions, those material properties 592were constants or defined by an external `MFront` file. 593 594This new version allows those material properties to be defined by 595formulae, as follows: 596 597~~~~{.cxx} 598@Parameter E0 =2.1421e11,E1 = -3.8654e7,E2 = -3.1636e4; 599@ElasticMaterialProperties {"E0+(T-273.15)*(E1+E2*(T-273.15))",0.3} 600~~~~ 601 602As for material properties defined in external `MFront` files, the 603material properties evaluated by formulae will be computed for updated 604values of their parameters. For example, if the previous lines were used 605in the `Implicit` DSL, two variables `young` and `young_tdt` will be 606automatically made available: 607 608- the `young` variable will be computed for the temperature 609 \(T+\theta\,\Delta\,T\). 610- the `young_tdt` variable will be computed for the temperature 611 \(T+\Delta\,T\). 612 613## New command line arguments 614 615### The `--list-behaviour-bricks` argument 616 617The `--list-behaviour-bricks` argument returns the list of all available 618behaviour bricks. 619 620~~~~{.sh} 621$ mfront --list-behaviour-bricks 622DDIF2 (undocumented) 623FiniteStrainSingleCrystal (undocumented) 624StandardElasticity (documented) 625StandardElastoViscoPlasticity (undocumented) 626~~~~ 627 628### The `--help-behaviour-brick` argument 629 630The `--help-behaviour-brick` argument returns the help associated with a 631given behaviour brick. 632 633 634~~~~{.sh} 635$ mfront --help-behaviour-brick=StandardElasticity 636The `StandardElasticity` brick describes the linear elastic part of 637the behaviour of an isotropic or orthotropic material. 638~~~~ 639 640### Arguments related to the `StandardElastoViscoplasticity` brick 641 642The `StandardElastoViscoplasticity` brick is based on the concepts of: 643 644- stress-potential 645- criterion 646- isotropic hardening rule 647- kinematic hardening rule 648 649For each concept, one can query the list of concrete implementations of 650this concept and the associated help, as in the following example: 651 652~~~~{.sh} 653$ mfront --list-stress-potentials 654- DDIF2 (documented) 655- Hooke (documented) 656- IsotropicDamage (undocumented) 657$ mfront --help-stress-potential=Hooke 658The `Hooke` brick describes the linear elastic part of 659the behaviour of an isotropic or orthotropic material. 660..... 661~~~~ 662 663The following arguments are thus avaiable: 664 665- `--list-isotropic-hardening-rules` 666- `--list-kinematic-hardening-rules` 667- `--list-stress-criteria` 668- `--list-stress-potentials` 669- `--help-isotropic-hardening-rule` 670- `--help-kinematic-hardening-rule` 671- `--help-stress-criterion` 672- `--help-stress-potential` 673 674 675# New functionalities in `MTest` 676 677## Events, activation and desactivation of constraints 678 679Defining an event is a way to active/desactivate a constraint, (see 680`@ImposedStrain`, `@ImposedCohesiveForce`, 681`@ImposedOpeningDisplacement`, `@ImposedStrain`, 682`@ImposedDeformationGradient`, `@NonLinearConstraint`). 683 684### Defining a new event 685 686The `@Event` keyword is followed by the name of the event and a time or 687a list of times (given as an array of values). 688 689#### Example 690 691~~~~ {.cpp} 692@Event 'Stop' 1; 693~~~~~~~~ 694 695### Activation and desactivation of constraints 696 697At the end of the definition of a constraint, one may now optionnally 698set options on the active state of this contraint at the beginning of 699the computation, the activating events and desactivating events using a 700JSON-like structure. This structure starts with an opening curly brace 701(`{`) and ends with a closing curly brace (`}`). An option is given by 702its name, an double-dot character (`:`) and the value of the option. 703Consecutive options are separated by a comma `,` (see below for an 704example). The following options are available: 705 706- `active`: states if the constraint is active at the beginning of the 707 computation. This is a boolean, so the expected value are either 708 `true` or `false`. 709- `activating_event`: gives the name of the event which activate the 710 constraint. An array of string is expected. 711- `activating_events`: gives the list of events which activate the 712 constraint. An array of string is expected. 713- `desactivating_event`: gives the name of the event which desactivate 714 the constraint. An array of string is expected. 715- `desactivating_events`: gives the list of events which desactivate the 716 constraint. An array of string is expected. 717 718#### Example 719 720~~~~ {.cpp} 721@ImposedStrain<evolution> 'EXX' {0.:0.,1:1e-3}{ 722 desactivating_event : 'Stop' 723}; 724~~~~~~~~ 725 726### User defined post-processings 727 728The `@UserDefinedPostProcessing` lets the user define is own 729post-processings. 730 731This keywords is followed by: 732 733- the name of the output file 734- the list of post-processings given by a string or and an array of 735 strings. Those strings defines formulae which are evaluated at the end 736 of the time step. Those formulae may depend on: 737 - the behaviour' driving variable 738 - the behaviour' thermodynamic forces 739 - the behaviour' internal state variables 740 - the behaviour' external state variables 741 - any evolution defined in the input file 742 743#### Example 744 745~~~~{.cpp} 746@UserDefinedPostProcessing 'myoutput.txt' {'SXX','EquivalentPlasticStrain'}; 747~~~~ 748 749# Tickets solved during the development of this version 750 751## Ticket #137: Compiling error fix `Mfront/Python` under `Windows` 752 753When compiling `MFront` with python bindings support using `Visual 754Studio 2015 Community`, the following error message is triggered: 755 756~~~~ 757C:\Codes\tfel\sources\tfel\include\TFEL/Material/PiPlane.ixx(60): error C2131: l'expression n'a pas été évaluée en co 758nstante [C:\Codes\tfel\sources\build-vs-2015\bindings\python\tfel\py_tfel_material.vcxproj] 759 7600 Avertissement(s) 7611 Erreur(s) 762~~~~ 763 764The problem is that `constexpr` variables are poorly handled by `Visual 765Studio 2015`. The `TFEL_CONSTEXPR` macro was introduced to handle this 766compiler (it boils down to `constexpr` for others compilers). The issue 767can be solved by changing the faulty line as as follows: 768 769~~~~ 770TFEL_CONSTEXPR const auto isqrt6 = isqrt2 * isqrt3; 771~~~~ 772 773For more details, see: <https://sourceforge.net/p/tfel/tickets/137/> 774 775 776## Ticket #112: MTEST outputs 777 778The `@UserDefinedPostProcessing` lets the user define is own 779post-processings. 780 781This keywords is followed by: 782 783- the name of the output file 784- the list of post-processings given by a string or and an array of 785 strings. Those strings defines formulae which are evaluated at the end 786 of the time step. Those formulae may depend on: 787 - the behaviour' driving variable 788 - the behaviour' thermodynamic forces 789 - the behaviour' internal state variables 790 - the behaviour' external state variables 791 - any evolution defined in the input file 792 793For more details, see: <https://sourceforge.net/p/tfel/tickets/112/> 794 795# Known incompatibilities 796 797## Header files 798 799The following header files have been renamed: 800 801- `Hosford.hxx` has been moved in `Hosford1972YieldCriterion.hxx`. 802- `Barlat.hxx` has been moved in `Barlat2004YieldCriterion.hxx`. 803 804# References 805 806<!-- Local IspellDict: english --> 807