Behaviours generated by the generic interface of MFront are meant to be integrated in the material frame. This choice allows the calling solver to develop its own strategy and code for dealing with orthotropic behaviours.

However, dealing with orthotropic behaviours may be tedious and require an considerable amount of developments. Since version \(3.4\) of TFEL, the generic interface generates additional functions which performs those steps efficiently. The purpose of this document is to introduce some of the features provided by MGIS built around those functions.

Outline of this document

Section 1 provides an small introduction on the specific issues related to orthotropic behaviours from the behaviour integration point of view.

The functions generated by the MFrontgeneric interface are described in Section 2.

Section 3 is devoted to functionalities and wrappers around those functions and how to use them efficiently.

About orthotropic behaviours

Dealing with such behaviours implies, for the calling solver:

  1. Rotating the gradients from the global frame to the material frame before the behaviour integration. At each equilibrium iterations, only the gradients at the end of the time steps are generally rotated as the values of the gradients in the material frame at the beginning of the time step may be stored by the BehaviourData (for one integration point) or MaterialDataManager classes (for a set of integration points).
  2. Rotating the thermodynamic forces and tangent operator blocks from the material frame to the global frame after the behaviour integration for the assembly of the inner forces and stiffness matrices.

Functions generated by the generic interface

The first set of functions generated by the generic interface are named as follows:

Note that those functions are only generated for orthotropic behaviours.

All those functions take three arguments:

Note that, for efficiency reasons, those functions do no perform any consistency checks. In particular:

In place rotations is explicitly allowed, i.e. the first and second arguments can be a pointer to the same location.

The three previous functions works for an integration point. Three other functions are also generated:

Those functions takes an additional arguments which is the number of integration points to be treated.

Finite strain behaviours

Finite strain behaviours are a special case, because the returned stress measure and the returned tangent operator can be chosen at runtime time. A specific rotation function is generated for each supported stress measure and each supported tangent operator.

Here is the list of the generated functions:

Features of MGIS related to orthotropic behaviours

Additional members of the Behaviour class and some helper functions

The Behaviour class exposes the following data members:

which are functions pointers to the appropriate functions (i.e. the modelling hypothesis and, for finite strain behaviours, the choices of the stress measures and tangent operator kind, are automatically taken into account by the load functions.

Direct use of those functions is meant to be highly efficient but still don’t perform any consistency checks.

The following helper functions are thus provided:

Those functions takes memory blocks (wrapped in std::span objects) and checks that the size of those blocks are consistent. For each of those functions, two overloads exist, one for in place rotation, one for out of place rotations.

Those functions are perfectly fine for managing the rotations of the gradients of a set of integration points. However, their overhead is probably significant when dealing with only one integration point and their use shall probably be limited to debug builds.

Example of usage

The following example loads and orthotropic behaviour, defines a rotation matrix, and rotates a strain from the global frame to the material frame:

const auto b = load("src/libBehaviour.so", "OrthotropicElasticity",
                    Hypothesis::GENERALISEDPLANESTRAIN);
const std::array<real, 9> r = {0, 1, 0, //
                               1, 0, 0, //
                               0, 0, 1};
const std::array<real, 4> ge = {1e-3, 0, 0, 0};
std::array<real, 4> me;
rotateGradients(me, b, ge, r);

C bindings

The C bindings provides the following functions, declared in the MGIS/Behaviour/Behaviour.h header:

Fortran bindings

The mgis.behaviours bindings provides the following functions, declared in the mgis_behaviour header:

Example of usage

g = (/ 1, 0, 0, 0, 0, 0/)
r = (/ 0, 1, 0, &
       1, 0, 0, &
       0, 0, 1/)
call check_status(load_behaviour(b, &
     get_mfront_behaviour_test_library_path(), &
     'OrthotropicElasticity', 'Tridimensional'))
call check_status(rotate_gradients_in_place(g, b, r))

Python bindings

The Python module mgis.behaviour provide the following free functions:

Those functions takes NumPy arrays as arguments. Two overloads are available for in-place or out-of-place rotations.

For convenience, those functions are also available as members of the Behaviour class.

Example of usage

from mgis.behaviour import *
import numpy as np
r = np.array([0,1,0,1,0,0,0,0,1], dtype=np.double)
g1 = np.array([1,0,0,0,0,0], dtype=np.double)
g2 = np.array([1,0,0,0,0,0], dtype=np.double)
b = load('src/libBehaviour.so', 'OrthotropicElasticity', Hypothesis.TRIDIMENSIONAL)
b.rotateGradients(g1,r)
rotateGradients(g2, b, r)
print(g1)
print(g2)
g3 = np.array([1,0,0,0,0,0,1,0,0,0,0,0], dtype=np.double)
b.rotateGradients(g3, r)
print(g3)
g4 = np.array([1,0,0,0,0,0,1,0,0,0,0,0], dtype=np.double)
b.rotateGradients(g4[0:6], r)
print(g4)

Julia bindings

Not yet available.