Release notes of Version 1.1
This version inherits from all the features introduced in:
Version 1.0.1
Version 1.0.2
Version 1.0.3
Highlights
The name of the materials and boundaries are automatically retrieved from mfem’s mesh.
Many methods have been deprecated to have a consistent error handling scheme based on MGIS’s one. As such, many methods and functions now takes and MGIS’s
Contextas their first argument.The regularization proposed by Faltus et al. in the context of the third medium contact has been implemented for plane strain, plane stress and tridmensional hypotheses
New features
Prediction of the solution
By default, a nonlinear evolution problem uses the solution at the beginning of the time step, modified by applying Dirichlet boundary conditions, as the initial guess of the solution at the end of the time step, see below for details.
This can be changed by using the setPredictionPolicy
method, as follows:
// use the elastic operator by default
mechanics.setPredictionPolicy(
{.strategy = PredictionStrategy::BEGINNING_OF_TIME_STEP_PREDICTION});
Available strategies
Default prediction (PredictionStrategy::DEFAULT_PREDICTION)
By default, a nonlinear evolution problem uses the solution at the beginning of the time step, modified by applying Dirichlet boundary conditions, as the initial guess of the solution at the end of the time step.
Warning
In mechanics, this may lead to very high increments of the deformation gradients or the strain in the neighboring elements of boundaries where evolving displacements are imposed.
Prediction for the state at the beginning of the time step (PredictionStrategy::BEGINNING_OF_TIME_STEP_PREDICTION)
The prediction for the state at the beginning of the time step strategy determines the increment of the displacement \(\Delta\,\mathbb{u}\) by solving the following linear system:
where:
\(\mathbb{K}_{e}\) denotes one of the prediction operator (see below).
\(\bts{\mathbb{F}_{e}}\) denotes the external forces at the beginning of the time step.
\(\bts{\mathbb{F}_{i}}\) denotes the inner forces at the beginning of the time step.
\(\Delta\,\mathbb{u}\) is submitted to the the increment of the imposed Dirichlet boundary conditions.
Note
Although the wording explicitly refers to mechanics, this equation applies to all physics.
The following prediction operators can be chosen:
PredictionOperator::ELASTIC: the elastic operatorPredictionOperator::SECANT: the secant operator is typically defined by the elastic operatorPredictionOperator::TANGENT_PREDICTION: the tangent operator, defined by the time-continuous derivative of the thermodynamic force with respect to the gradients.PredictionOperator::LAST_ITERATE_OPERATOR: this operator reuses the one computed at the last iteration of the previous time step. At the first time step, the elastic operator is used.
Prediction based on a behaviour integration with constant gradients (PredictionStrategy::CONSTANT_GRADIENTS_INTEGRATION_PREDICTION)
The CONSTANT_GRADIENTS_INTEGRATION strategy determines the
increment of the unknown \(\Delta\,{u}\) by solving the following
linear system:
where:
\(\tilde{\mathbb{K}}\) denotes the operator computed at the end of the behaviour integration.
\(\bts{\mathbb{F}_{e}}\) denotes the external forces at the beginning of the time step.
\(\ets{\tilde{\mathbb{F}}_{i}}\) denotes an approximation inner forces at the end of the time step computed by assuming that the gradients are constant over the time (and thus egal to their values at the beginning of the time step).
\(\Delta\,\mathbb{u}\) is submitted to the the increment of the imposed Dirichlet boundary conditions.
The behaviour integration allows taking into account:
the evolution of stress-free strain over the time step (thermal expansion, swelling, etc..),
the viscoplastic relaxation of the stress. This relaxation can be discarded by integrating the behaviour with a null time step.
The following operators are available:
IntegrationOperator::ELASTIC: the elastic operator,IntegrationOperator::SECANT: the secant operator is typically defined by the elastic operator affected by damage,IntegrationOperator::TANGENT: the tangent operator, defined by the time-continuous derivative of the thermodynamic force with respect to the gradients,IntegrationOperator::CONSISTENT_TANGENT: the constistent tangent operator, defined by the derivative of the thermodynamic force with respect to the gradients at the end of the time step. See [ST85] for details.
Faltus 2026 regularization
The regularization proposed by Faltus et al. in the context of contact mechanics using a third medium [FAH]. This regularization only applies to finite strain behaviours. Currently, only this regularization is only available for isotropic behaviours.
This regularization adds a contribution to the standard variational operator in finite strain to can be derived from an energy \(W\) which penalizes the difference between the deformation gradient \(\underline{F}\) at a given quadrature point and its value \(\bar{\underline{F}}\) at the centroid of the element:
where \(\alpha\) is a penalization coefficient.
This regularization is enabled by passing an additional parameter to the
the Mechanics behaviour integrator, as follows:
const auto faltus_parameters = mfem_mgis::Parameters{
{"Regularization",
mfem_mgis::Parameters{
{"Faltus2026",
mfem_mgis::Parameters{{"PenalizationCoefficient", 1e11}}}}}};
mechanics.addBehaviourIntegrator(ctx, "Mechanics", "ThirdMedium", library,
behaviour2, faltus_parameters) | or_die;
The info function
The info function allows to display information about an object
in an output stream.
Retrieving information on a finite element discretization
Example of usage
const auto& fed = problem.getFiniteElementDiscretization();
const auto success = mfem_mgis::info(ctx, fed);
Retrieving information on a partial quadrature space
The PartialQuadratureSpaceInformation structure contains some
relevant information about a partial quadrature space:
the identifier and the name of the underlying material,
the total number of elements,
the total number of integration points,
the number of elements points per geometric type.
the number of quadrature points per geometric type.
This structure is created by:
getLocalInformation, which returns the information relative to the current process.getInformation, which returns the information gathered from all processes.
The PartialQuadratureSpaceInformation structure can be printed to
an output stream using the info function.
Example of usage
const auto& qspace =
problem.getBehaviourIntegrator(1).getPartialQuadratureSpace();
const auto success = mfem_mgis::info(ctx, std::cout, qspace);
Evaluators of quantities at integration points
Overview
The setMaterialProperty and setExternalStateVariable methods of behaviour integrators
Resolution of dependencies of a nonlinear evolution problemn using another
The Simulation class
Physical system, coupling schemes and models
Line-search-like handling of behaviour integration failures
Issues fixed
Issue 218: [performance] synchronize success of the setup methods at a higher level to minimize collective communications enhancement
Issue 213:  Add a simple way to resolve dependencies (material properties, external state variables) of a
NonLinearEvolutionProblemusing the gradients, thermodynamic forces and internal state variables of another one.Issue 211: Add coupling schemes and models
Issue 209: Introduce the
SimulationclassIssue 206: Line-search-like handling of behaviour integration failures
Issue 200: automatically assign materials and boundaries’s names from mfem’s attributes 
Issue 198: Add the ability to define multiple behaviour integrators on the same material
Issue 193: Add support for other types of search operators for computing a prediction of the solution at the end of the time step enhancement
Issue 192: Take external forces into account when computing the prediction of the solution
Issue 188: retrieve information about a quadrature space
Issue 149: work on the prediction of the solution