MGIS
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.
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 MFront
’
generic
interface are described in Section 2.
Section 3 is devoted to functionalities and wrappers around those functions and how to use them efficiently.
Dealing with such behaviours implies, for the calling solver:
BehaviourData
(for one
integration point) or MaterialDataManager
classes (for a
set of integration points).generic
interfaceThe first set of functions generated by the generic
interface are named as follows:
<behaviour_function_name>_<hypothesis>_rotateGradients
<behaviour_function_name>_<hypothesis>_rotateThermodynamicForces
<behaviour_function_name>_<hypothesis>_rotateTangentOperatorBlocks
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:
<behaviour_function_name>_<hypothesis>_rotateArrayOfGradients
<behaviour_function_name>_<hypothesis>_rotateArrayOfThermodynamicForces
<behaviour_function_name>_<hypothesis>_rotateArrayOfTangentOperatorBlocks
Those functions takes an additional arguments which is the number of integration points to be treated.
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:
<behaviour_function_name>_<hypothesis>_rotateThermodynamicForces_CauchyStress
.
This function assumes that its first argument is the Cauchy stress in
the material frame.<behaviour_function_name>_<hypothesis>_rotateThermodynamicForces_PK1Stress
.
This function assumes that its first argument is the first
Piola-Kirchhoff stress in the material frame.<behaviour_function_name>_<hypothesis>_rotateThermodynamicForces_PK2Stress
.
This function assumes that its first argument is the second
Piola-Kirchhoff stress in the material frame.
<behaviour_function_name>_<hypothesis>_rotateTangentOperatorBlocks_dsig_dF
.
This function assumes that its first argument is the derivative of the
Cauchy stress with respect to the deformation gradient in the material
frame.<behaviour_function_name>_<hypothesis>_rotateTangentOperatorBlocks_dPK1_dF
.
This function assumes that its first argument is the derivative of the
first Piola-Kirchhoff stress with respect to the deformation gradient in
the material frame.<behaviour_function_name>_<hypothesis>_rotateTangentOperatorBlocks_PK2Stress
.
This function assumes that its first argument is the derivative of the
second Piola-Kirchhoff stress with respect to the Green-Lagrange strain
in the material frame.MGIS
related to orthotropic behavioursBehaviour
class and some helper
functionsThe Behaviour
class exposes the following data
members:
rotate_gradients_ptr
rotate_array_of_gradients_ptr
rotate_thermodynamic_forces_ptr
rotate_array_of_thermodynamic_forces_ptr
rotate_tangent_operator_blocks_ptr
rotate_array_of_tangent_operator_blocks_ptr
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:
rotateGradients
rotateThermodynamicForces
rotateTangentOperatorBlocks
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.
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",
::GENERALISEDPLANESTRAIN);
Hypothesisconst 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;
(me, b, ge, r); rotateGradients
C
bindingsThe C
bindings provides the following functions,
declared in the MGIS/Behaviour/Behaviour.h
header:
mgis_bv_behaviour_rotate_gradients_in_place
mgis_bv_behaviour_rotate_gradients_out_of_place
mgis_bv_behaviour_rotate_array_of_gradients_in_place
mgis_bv_behaviour_rotate_array_of_gradients_out_of_place
mgis_bv_behaviour_rotate_thermodynamic_forces_in_place
mgis_bv_behaviour_rotate_thermodynamic_forces_out_of_place
mgis_bv_behaviour_rotate_array_of_thermodynamic_forces_in_place
mgis_bv_behaviour_rotate_array_of_thermodynamic_forces_out_of_place
mgis_bv_behaviour_rotate_tangent_operator_blocks_in_place
mgis_bv_behaviour_rotate_tangent_operator_blocks_out_of_place
mgis_bv_behaviour_rotate_array_of_tangent_operator_blocks_in_place
mgis_bv_behaviour_rotate_array_of_tangent_operator_blocks_out_of_place
Fortran
bindingsThe mgis.behaviours
bindings provides the following
functions, declared in the mgis_behaviour
header:
rotate_gradients_in_place
rotate_gradients_out_of_place
rotate_array_of_gradients_in_place
rotate_array_of_gradients_out_of_place
rotate_thermodynamic_forces_in_place
rotate_thermodynamic_forces_out_of_place
rotate_array_of_thermodynamic_forces_in_place
rotate_array_of_thermodynamic_forces_out_of_place
rotate_tangent_operator_blocks_in_place
rotate_tangent_operator_blocks_out_of_place
rotate_array_of_tangent_operator_blocks_in_place
rotate_array_of_tangent_operator_blocks_out_of_place
= (/ 1, 0, 0, 0, 0, 0/)
g = (/ 0, 1, 0, &
r 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
bindingsThe Python
module mgis.behaviour
provide
the following free functions:
rotateGradients
rotateThermodynamicForces
rotateTangentOperatorBlocks
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.
from mgis.behaviour import *
import numpy as np
= np.array([0,1,0,1,0,0,0,0,1], dtype=np.double)
r = np.array([1,0,0,0,0,0], dtype=np.double)
g1 = np.array([1,0,0,0,0,0], dtype=np.double)
g2 = load('src/libBehaviour.so', 'OrthotropicElasticity', Hypothesis.TRIDIMENSIONAL)
b
b.rotateGradients(g1,r)
rotateGradients(g2, b, r)print(g1)
print(g2)
= np.array([1,0,0,0,0,0,1,0,0,0,0,0], dtype=np.double)
g3
b.rotateGradients(g3, r)print(g3)
= np.array([1,0,0,0,0,0,1,0,0,0,0,0], dtype=np.double)
g4 0:6], r)
b.rotateGradients(g4[print(g4)
Julia
bindingsNot yet available.