The page declares the new functionalities of the 3.1 version of TFEL, MFront and MTest.

Highlights

Numerical stability

Tests with verrou

Enabling the -ffast-math with GCC

Single crystal behaviours

Introducing the TFELNUMODIS library

The SlipSystemsDescription class in TFEL/Material

TFEL

TFEL/Utilities

String algorithms

The starts_with string algorithm

The starts_with string algorithm is an helper function used to determine if a given string starts with another.

The ends_with string algorithm

The ends_with string algorithm is an helper function used to determine if a given string ends with another.

TFEL/System

The LibraryInformation class

This release introduces the LibraryInformation class that allow querying a library about exported symbols.

Note This class has been adapted from the boost/dll library version 1.63 and has been originally written by Antony Polukhin, Renato Tegon Forti and Antony Polukhin.

Improvements to the ExternalLibraryManager class

Completion of libraries names

If a library is not found, the ExternalLibraryManager class will try the following combinaisons:

The standard library suffix is:

Retrieving the path of a library

The getLibraryPath method returns the path to a shared library:

Better handling of behaviour parameters

The ExternalLibraryManager class has several new methods for better handling of behaviours' parameters:

Retrieving bounds values

The ExternalLibraryManager class has several new methods for better handling of a behaviour' variable bounds:

Retrieving physical bounds values

The ExternalLibraryManager class has several new methods for better handling of a behaviour' variable bounds:

Retrieving all mfront generated entry points and associated information

The getEntryPoints method returns a list containing all mfront generated entry points. Those can be functions or classes depending on the interface's needs.

The getMaterialKnowledgeType allows retrieving the material knowledge type associated with and entry point. The returned value has the following meaning:

The getInterface method allows retrieving the interface of used to generate an entry point.

The following code retrieves all the behaviours generated with the aster interface in the libAsterBehaviour.so library:

auto ab = std::vector<std::string>{};
const auto l = "AsterBehaviour";
auto& elm = ExternalLibraryManager::getExternalLibraryManager();
for(const auto& e : elm.getEntryPoints(l)){
  if((elm.getMaterialKnowledgeType(l,e)==1u)&&(elm.getInterface(l,e)=="Aster")){
    ab.push_back(e);
  }
}

Note that we did not mention the prefix and the suffix of the library. The library path is searched through the getLibraryPath method.

The equivalent python code is the following:

ab  = []
l   = 'AsterBehaviour';
elm = ExternalLibraryManager.getExternalLibraryManager();
for e in elm.getEntryPoints(l):
    if ((elm.getMaterialKnowledgeType(l,e)==1) and
        (elm.getInterface(l,e)=='Aster')):
        ab.append(e)

TFEL/Math

Symmetric tensor eigen values and eigen vectors

The computation of the eigen values and eigen vectors of a symmetric tensor has been improved in various ways:

New overloaded versions of the computeEigenValues, computeEigenVectors and computeEigenTensors methods

Return by value

Various overloaded versions of the computeEigenValues, computeEigenVectors and computeEigenTensors methods have been introduced for more readable usage and compatibility with structured bindings construct introduced in C++17: the results of the computations are returned by value.

For example:

tmatrix<3u,3u,real> m2;
tvector<3u,real>    vp2;
std::tie(vp,m)=s.computeEigenVectors();

Thanks to C++17 structured bindings construct, the previous code will be equivalent to this much shorter and more readable code:

auto [vp,m] = s.computeEigenVectors();

Even better, we could write:

const auto [vp,m] = s.computeEigenVectors();
Eigen values sorting

The computeEigenValues and computeEigenVectors methods now have an optional argument which specify if we want the eigen values to be sorted. Three options are available:

Here is how to use it:

tmatrix<3u,3u,real> m2;
tvector<3u,real>    vp2;
std::tie(vp,m)=s.computeEigenVectors(Stensor::ASCENDING);

New eigen solvers

The default eigen solver for symmetric tensors used in TFEL is based on analitical computations of the eigen values and eigen vectors. Such computations are more efficient but less accurate than the iterative Jacobi algorithm (see (Kopp 2008; Kopp 2017)).

With the courtesy of Joachim Kopp, we have created a C++11 compliant version of his routines that we gathered in header-only library called FSES (Fast Symmetric Eigen Solver). This library is included with TFEL and provides the following algorithms:

We have also introduced the Jacobi implementation of the Geometric Tools library (see (Eberly 2016; Eberly 2017)).

Those algorithms are available in 3D. For 2D symmetric tensors, we fall back to some default algorithm as described below.

Table 1: List of available eigen solvers.
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)

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::GTESYMMETRICQREIGENSOLVER>();
Table 2: Test on \(10^{6}\) random symmetric tensors in single precision (float).
Algorithm Failure ratio \(\Delta_{\infty}\) Times (ns) Time ratio
TFELEIGENSOLVER 0.000642 3.29e-05 250174564 1
GTESYMMETRICQREIGENSOLVER 0 1.10e-06 359854550 1.44
FSESJACOBIEIGENSOLVER 0 5.62e-07 473263841 1.89
FSESQLEIGENSOLVER 0.000397 1.69e-06 259080052 1.04
FSESCUPPENEIGENSOLVER 0.019469 3.49e-06 274547371 1.10
FSESHYBRIDEIGENSOLVER 0.076451 5.56e-03 126689850 0.51
FSESANALYTICALEIGENSOLVER 0.108877 1.58e-01 127236908 0.51
Table 3: Test on \(10^{6}\) random symmetric tensors in double precision (double).
Algorithm Failure ratio \(\Delta_{\infty}\) Times (ns) Time ratio
TFELEIGENSOLVER 0.000632 7.75e-14 252663338 1
GTESYMMETRICQREIGENSOLVER 0 2.06e-15 525845499 2.08
FSESJACOBIEIGENSOLVER 0 1.05e-15 489507133 1.94
FSESQLEIGENSOLVER 0.000422 3.30e-15 367599140 1.45
FSESCUPPENEIGENSOLVER 0.020174 5.79e-15 374190684 1.48
FSESHYBRIDEIGENSOLVER 0.090065 3.53e-10 154911762 0.61
FSESANALYTICALEIGENSOLVER 0.110399 1.09e-09 157613994 0.62
Table 4: Test on \(10^{6}\) random symmetric tensors in extended precision (long double).
Algorithm Failure ratio \(\Delta_{\infty}\) Times (ns) Time ratio
TFELEIGENSOLVER 0.000575 2.06e-17 428333721 1
GTESYMMETRICQREIGENSOLVER 0 1.00e-18 814990447 1.90
FSESJACOBIEIGENSOLVER 0 5.11e-19 748476926 1.75
FSESQLEIGENSOLVER 0.00045 1.83e-18 548604588 1.28
FSESCUPPENEIGENSOLVER 0.009134 3.23e-18 734707748 1.71
FSESHYBRIDEIGENSOLVER 0.99959 4.01e-10 272701749 0.64
FSESANALYTICALEIGENSOLVER 0.999669 1.36e-11 315243286 0.74

Some benchmarks

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:

Isotropic functions and Isotropic function derivatives of symmetric tensors

Given a scalar valuated function \(f\), one can define an associated isotropic function for symmetric tensors as: \[ f{\left({\underline{s}}\right)}=\sum_{i=1}^{3}f{\left(\lambda_{i}\right)}{\underline{n}}_{i} \]

where \(\left.\lambda_{i}\right|_{i\in[1,2,3]}\) are the eigen values of the symmetric tensor \({\underline{s}}\) and \(\left.{\underline{n}}_{i}\right|_{i\in[1,2,3]}\) the associated eigen tensors.

If \(f\) is \(\mathcal{C}^{1}\), then \(f\) is a differentiable function of \({\underline{s}}\).

\(f\) can be computed with the computeIsotropicFunction method of the stensor class. \({{\displaystyle\frac{\displaystyle \partial f}{\displaystyle \partial {\underline{s}}}}}\) can be computed with computeIsotropicFunctionDerivative. One can also compute \(f\) and \({{\displaystyle\frac{\displaystyle \partial f}{\displaystyle \partial {\underline{s}}}}}\) all at once by the computeIsotropicFunctionAndDerivative method. All those methods are templated by the name of the eigen solver (if no template parameter is given, the TFELEIGENSOLVER is used).

Various new overloaded versions of those methods have been introduced in TFEL-3.1. This overloaded methods are meant to:

To illustrate this new features, assuming that the structured bindings feature of the C++17 standard is available, one can now efficiently evaluate the exponential of a symmetric tensor and its derivative as follows:

const auto [vp,m] = s.computeEigenVectors();
const auto evp    = map([](const auto x){return exp(x)},vp);
const auto [f,df] = Stensor::computeIsotropicFunctionAndDerivative(evp,evp,vp,m,1.e-12);

Portable implementation of the fpclassify, isnan, isfinite functions

The C99 standard defines the fpclassify, isnan, isfinite functions to query some information about double precision floatting-point numbers (double):

The C++11 provides a set of overload for single precision (float) and extended precision (long double) floatting-point numbers.

Those functions are very handy to check the validity of a computation. However, those functions are not compatible with the use of the -ffast-math option of the GNU compiler which also implies the -ffinite-math-only option. This latter option allows optimizations for floating-point arithmetic that assume that arguments and results are finite numbers. As a consequence, when this option is enabled, the previous functions does not behave as expected. For example, isnan always returns false, whatever the value of its argument.

To overcome this issue, we have introduced in TFEL/Math the implementation of these functions provided by the musl library (see Musl development community (2017)). Those implementations are compatible with the -ffast-math option of the GNU compiler. Those implementations are defined in the TFEL/Math/General/IEEE754.hxx header file in the tfel::math::ieee754 namespace.

TFEL/Material

The SlipSystemsDescription class

New functionalities of the MFront code generator

The MFront gallery is meant to present well-written implementation of behaviours that will be updated to follow MFront latest evolutions.

Hyperelastic behaviours

The implementation of various hyperelastic behaviours can be found here:

Hyperviscoelastic behaviours

The following page describes how to implement standard hyperviscoelastic behaviours based on the development in Prony series:

http://tfel.sourceforge.net/hyperviscoelasticity.html

Plasticity

An example of a simple orthotropic behaviour is available here.

Viscoplasticity

This following article shows how to implement an isotropic viscoplastic behaviour combining isotropic hardening and multiple kinematic hardenings following an Armstrong-Frederic evolution of the back stress:

http://tfel.sourceforge.net/isotropicplasticityamstrongfrederickinematichardening.html

Behaviours interfaces

The Cast3M interface

The MieheApelLambrechtLogarithmic finite strain strategy

The pre- and post-computations performed by the MieheApelLambrechtLogarithmic finite strain strategy , which require the computation of the eigen values and eigen vectors of the right Cauchy strecth tensor, are now based the Jacobi algorithm from the FSES library for improved accuracy.

The Code_Aster interface

Support for the GROT_GDEP finite strain formulation

GROT_GDEP is the name in Code_Aster of a finite strain formulation based on the principle of virtual work in the reference configuration expressed in term of the Green-Lagrange strain and the second Piola-Kirchhoff stress. Such a formulation is also called Total Lagrangian in the litterature (see Belytschko (2000)) and in other finite element solvers.

Prior to this version, MFront behaviours were meant to be used with the SIMO_MIEHE finite strain formulation and could not be used with the GROT_GDEP finite strain formulation.

From the behaviour point of view, using SIMO_MIEHE or GROT_GDEP differs from the choice of the output stress and the definition of the consistent tangent operator.

The @AsterFiniteStrainFormulation keyword

The @AsterFiniteStrainFormulation keyword can now be used to choose one of these finite strain formulation.

This keyword must be followed by one of the following choice:

The choice SIMO_MIEHE remains the default for backward compatibility.

The Europlexus interface

The MieheApelLambrechtLogarithmic finite strain strategy

The pre- and post-computations performed by the MieheApelLambrechtLogarithmic finite strain strategy, which require the computation of the eigen values and eigen vectors of the right Cauchy strecth tensor, are now based the Jacobi algorithm from the FSES library for improved accuracy.

The Abaqus-Explicit interface

The MieheApelLambrechtLogarithmic finite strain strategy

The pre- and post-computations performed by the MieheApelLambrechtLogarithmic finite strain strategy, which require the computation of the eigen values and eigen vectors of the right Cauchy strecth tensor, are now based the Jacobi algorithm from the FSES library for improved accuracy.

New functionalities of MTest solver

Choice of the rounding mode from the command line

\(4\) rounding mode are defined in the IEEE754 standard. Changing the rounding mode is a gross way to check the numerical stability of the computations performed with MTest and MFront.

The rounding mode can be set using the --rounding-direction-mode option. Valid values for this option are:

Most unit-tests based on MTest are now executed five times, one for each available choice of the rounding mode.

Abritrary non linear constraints

Abritrary non linear constraints on driving variables and thermodynamic forces can now be added using the @NonLinearConstraint keyword.

Note

This keyword can also be used to define linear constraints, although the numerical treatment of such a constraint will be sub-optimal. A special treatment of such a constraint is planned.

Note

This development of this functionality highlighted the issue reported in Ticket #39. For more details, see: https://sourceforge.net/p/tfel/tickets/39/

Normalisation policy

This keyword must be followed by an option giving the normalisation policy. The normalisation policy can have one of the following values:

Examples

// ensure that the loading is isochoric in 1D
@NonLinearConstraint<Strain> 'FRR*FTT*FZZ-1';
// impose the first piola kirchoff stress
// in an uniaxial compression test
@Real 'Pi0' -40e6
@NonLinearConstraint<Stress> 'SXX*FYY*FZZ-Pi0';

The @Print and @Message keywords

The @Print keyword, or its alias named @Message, is used to display some informative message on the standard output.

This keyword is followed by floatting point values and/or strings.

Strings are first interpreted as formula. If the interpretation is successfull, the result is printed. Otherwise, the string is display witout interpretation.

All the following tokens are appended to the message up to a final semi-colon.

Example:

@Print "Complex computation result: " "12*5";

In this example, the first string can't be interpreted as a formula, so its contents is printed. The second part can be interpreted, so its result (\(60\)) is displayed. The message printed is thus:

Complex computation result: 60

The @Import keyword

Depending of the option used (given between '<' and '>'), the @Import keyword is meant to have various meanings.

In this version, the only option available is the castem option.

The castem option

The castem (or Castem or Cast3M) option let you import a function generated by MFront with the castem interface. This function can be used in every formula.

The keyword is followed by the library an function names.

Example

@Import<castem> 'CastemW' 'W_ThermalExpansion';
// height at 20°C
@Real 'h0' 16e-3;
// height at 1500°C
@Real 'h' 'h0*(1+W_ThermalExpansion(1723.15)*(1723.25-293.15))';

python bindings

The Behaviour class

The Behaviour class has been introduced in the mtest modules. This class can be used to determine at runtime time the material properties, internal state variables, parameters and external state variables required by a specific implementation.

Contrary the tfel.system.ExternalBehaviourDescription class, the information given by the Behaviour class takes into account the variables that are implicitly declared by the interface to match its specific (internal) requirements. For example:

The MTest class

In the python bindings, the setNonLinearConstraint method has been added to the MTest class.

This method takes two named arguments:

New functionalities of the mfront-query tool

New behaviours queries

Tickets fixed

Ticket #40: ImplicitDSL: Detect non finite values during resolution

During the resolution of the implicit system, invalid results may occur. In previous versions, no check were made leading to a propagation of those values and finally the failure of integration.

A test to check that the residual of the implicit system is finite have been added. If this test is not satisfied after the first iteration, the last increment of the unknowns is divided by two and the resolution is restarted with this guess. If this test is not satisfied at the first iteration, the behaviour integration can not be performed.

Ticket #41: MTest: check if the residual is finite and not NaN

In previous versions, if the behaviour integration returned a not-a-number value (NaN ), this value propagated throughout the computation.

This situation can be detected by checking that the convergence criteria are finite as defined by the IEEE754 standard.

For more details, see: https://sourceforge.net/p/tfel/tickets/41/

Ticket #42: Check for infinite and NaN values in material properties

In the previous versions of MFront, generated sources for material properties checked that the errno value to determine is something had gone wrong, but this check does not appear to portable nor reliable with the INTEL compiler or when the -ffast-math option of the GNU compiler is activated.

The current version now check that the return value is finite.

For more details, see: https://sourceforge.net/p/tfel/tickets/42/

Ticket #43: Add the list of parameters' names and types to generated library for the UMAT++ interface

In previous versions of MFront, the list of parameters' names and types were not exported to the generated library for the UMAT++ interface, i.e. the additional symbols defined in the generated shared libraries that can be read through the ExternalLibraryManager class.

For more details, see: https://sourceforge.net/p/tfel/tickets/43/

Ticket #45: Support for bounds on parameters

For more details, see: https://sourceforge.net/p/tfel/tickets/45/

Ticket #46: Improvements to the mfront python module

The following improvements to the mfront python module have been made:

For more details, see: https://sourceforge.net/p/tfel/tickets/46/

Ticket #47: Add python bindings for the mtest::Behaviour class

The mtest module now contains bindings for the mtest::Behaviour class. This class allow querying information about how to use a behaviour in a specific context (interface and modelling hypothesis): for example, if a behaviour has the requireStiffnessTensor attribute, the list of material properties is updated appropriately if required by the interface for the considered modelling hypothesis. The Behaviour class has the following useful methods:

For more details, see: https://sourceforge.net/p/tfel/tickets/47/

Here is an example of the usage of the Behaviour class in python:

import mtest
b= mtest.Behaviour('AsterBehaviour','asternorton','Tridimensional');
for p in b.getParametersNames():
    print('- '+p+': '+str(b.getRealParameterDefaultValue(p)))
for p in b.getIntegerParametersNames():
    print('- '+p+': '+str(b.getIntegerParameterDefaultValue(p)))
for p in b.getUnsignedShortParametersNames():
    print('- '+p+': '+str(b.getUnsignedShortParameterDefaultValue(p)))

Ticket #46: Improved python bindings for the mfront::BehaviourDescription class

The python bindings of the mfront::BehaviourDescription now gives access to the parameters default values, and information about a variable standard or physical bounds (type, range).

Here is an example of its usage:

from tfel.material import ModellingHypothesis
import mfront

def printBounds(n,b):
    print('Bounds of variable \''+n+'\':')
    if((b.boundsType==mfront.VariableBoundsTypes.LOWER) or
       (b.boundsType==mfront.VariableBoundsTypes.LOWERANDUPPER)):
        print('- lower bound: '+str(b.lowerBound))
    if((b.boundsType==mfront.VariableBoundsTypes.UPPER) or
       (b.boundsType==mfront.VariableBoundsTypes.LOWERANDUPPER)):
        print('- upper bound: '+str(b.upperBound))
        print('')

dsl = mfront.getDSL('Norton.mfront')
dsl.analyseFile('Norton.mfront',[])

# behaviour description
bd = dsl.getBehaviourDescription()

if(bd.getSymmetryType()==mfront.BehaviourSymmetryType.ISOTROPIC):
    print 'Isotropic behaviour\n'
else:
    print 'Orthropic behaviour\n'

if(bd.getElasticSymmetryType()==mfront.BehaviourSymmetryType.ISOTROPIC):
    print 'Isotropic elasticity\n'
else:
    print 'Orthropic elasticity\n'

# a deeper look at the 3D case
d = bd.getBehaviourData(ModellingHypothesis.TRIDIMENSIONAL)
for p in d.getParameters():
    if(p.arraySize==1):
        if(p.hasBounds()):
            printBounds(p.name,p.getBounds())
    else:
        for i in range(p.arraySize):
            if(p.hasBounds(i)):
                printBounds(p.name+'['+str(i)+']',p.getBounds(i))

Ticket #48: Add the ability to retrieve bounds for material properties and parameters from the mtest::Behaviour class

The following methods were added to the mtest.Behaviour class: - The hasBounds method returns true if the given variable has bounds. - The hasLowerBound method returns true if the given variable has a lower bound. - The hasUpperBound method hasUpperBound returns true if the given variable has an upper bound. - The getLowerBound method returns the lower bound of the given variable. - The getUpperBound method returns the uppert bound of the given variable. - The hasPhysicalBounds methodreturns true if the given variable has physical bounds. - The hasLowerPhysicalBound method returns true if the given variable has a physical lower bound. - The hasUpperPhysicalBound method returns true if the given variable has a physical upper bound. - The getLowerPhysicalBound method returns the lower bound of the given variable. - The getUpperPhysicalBound method returns the upper bound of the given variable.

Here is an example:

from mtest import Behaviour

b = Behaviour('AsterBehaviour','asternorton','Tridimensional')

for p in b.getParametersNames():
    if b.hasLowerBound(p):
        print(p+" lower bound: "+str(b.getLowerBound(p)))
    if b.hasUpperBound(p):
        print(p+" lower bound: "+str(b.getUpperBound(p)))

For more details, see: https://sourceforge.net/p/tfel/tickets/48/

Ticket #49: Add the ability to retrieve the symmetry of the behaviour and the symmetry of the elastic behaviour from mfront-query

The following queries are now available:

For more details, see: https://sourceforge.net/p/tfel/tickets/49/

Ticket #50: Add the ability to retrieve bounds values from mfront-query

The following queries are now available:

For more details, see: https://sourceforge.net/p/tfel/tickets/50/

Ticket #65: @ElasticMaterialProperties does not work for DSL describing isotropics behaviours

The @ElasticMaterialProperties is now available for domain specific languages (DSL) describing isotropics behaviours.

@DSL IsotropicStrainHardeningMisesCreep;
@Behaviour StrainHardeningCreep2;
@Author    Helfer Thomas;
@Date      23/11/06;

@ElasticMaterialProperties {"Inconel600_YoungModulus.mfront",0.3};

@MaterialProperty real A;
@MaterialProperty real Ns;
@MaterialProperty real Np;

@FlowRule{
  const real p0  = 1.e-6;
  const real tmp = A*pow(seq,Ns-1.)*pow(p+p0,-Np-1);
  f       = tmp*seq*(p+p0);
  df_dseq =  Ns*tmp*(p+p0);
  df_dp   = -Np*tmp*seq;
}

For more details, see: https://sourceforge.net/p/tfel/tickets/65/

Ticket 74 mfront-query: handle static variables

mfront-query now has the following options:

For more details, see: https://sourceforge.net/p/tfel/tickets/74/

Know regressions

Stricter rules on variable declarations

All variables, to the very exception of local variables, must be declared before the first user defined code block to allow appropriate analysis of those code blocks.

Some variables are automatically declared by keywords. For instance, the @Epsilon keyword defines implicitly a parameter named epsilon.

In previous versions of MFront, those rules were only partially enforced: it may happen that some keywords or variable declaration shall now be moved before the first user defined code block.

References

Belytschko, Ted. 2000. Nonlinear Finite Elements for Continua and Structures. Chichester ; New York: Wiley-Blackwell.

Eberly, David. 2016. “A Robust Eigensolver for 3 × 3 Symmetric Matrices.” https://www.geometrictools.com/Documentation/RobustEigenSymmetric3x3.pdf.

———. 2017. “Geometric Tools.” http://www.geometrictools.com/.

Kopp, Joachim. 2008. “Efficient Numerical Diagonalization of Hermitian 3x3 Matrices.” International Journal of Modern Physics C 19 (3): 523–48. doi:10.1142/S0129183108012303.

———. 2017. “Numerical Diagonalization of 3x3 Matrices.” https://www.mpi-hd.mpg.de/personalhomes/globes/3x3/.

Musl development community. 2017. “Musl Libc.” https://www.musl-libc.org/.