This section describes the features available to analyse the failure of the integration of a behaviour at one integration point.
integrate_debug
functionsThe integrate_debug function are small wrappers around
the integrate functions which, in case of integration
failure, calls a user defined analyser or the default one.
The integrate_debug functions are direct substitutes for
the integrate functions and have the same prototypes, when
the default analyser is used. The following snippet shows how to call it
using a behaviour data view v and a behaviour
b:
const auto r = integrate_debug(v, b);See Section 1.2 for a description of the default analyser.
The user may pass his own analyser, the type of which must inherit
from the BehaviourIntegrationFailureAnalyser interface,
describe below, as an extra argument, as follows:
const auto a = CustomAnalyser{};
const auto r = integrate_debug(v, b, a);See Section 1.3 for details.
The default analyser generates one file per integration failure
containing the information required to re-run the simulation, using MTest
for instance.
By default, a markdown file is generated. Other formats may be supported in the future.
By default, the file name is composed by the name of the function implementing the behaviour (for the considered modelling hypothesis), a unique identifier (a number starting from zero and incremented for each integration failure report) and an extension. See Section 1.2.2 for details.
# Behaviour description
- library:/home/th202608/codes/mgis/master/src/build/tests/libBehaviourTest.so
- behaviour: NortonFailure
- function: NortonFailure_Tridimensional
# State at the beginning of the time step
## Gradients
- Strain (Stensor): {0, 0, 0, 0, 0, 0}
## Thermodynamic forces
- Stress (Stensor): {0, 0, 0, 0, 0, 0}
## Internal state variables
- ElasticStrain (Stensor): {0, 0, 0, 0, 0, 0}
- EquivalentViscoplasticStrain (Scalar): 0
## External state variables
- Temperature (Scalar): 293.15
# State at the end of the time step
## Gradients
- Strain (Stensor): {5e-05, 0, 0, 0, 0, 0}
## Thermodynamic forces
- Stress (Stensor): {0, 0, 0, 0, 0, 0}
## Internal state variables
- ElasticStrain (Stensor): {0, 0, 0, 0, 0, 0}
- EquivalentViscoplasticStrain (Scalar): 0
## External state variables
- Temperature (Scalar): 293.15The output file name is generated by the
getBehaviourIntegrationFailureAnalysisFileName function
which takes the behaviour’s function as first argument and the file
extension as second argument.
The behaviour of the
getBehaviourIntegrationFailureAnalysisFileName function can
be altered using the
setBehaviourIntegrationFailureAnalysisFileNameGenerator
which takes a std::function as first argument. This
function can be useful in MPI computations to distinguish the various
processes.
auto g = [](std::string_view b, std::string_view n, std::string_view e){
int r;
MPI_Comm_rank(MPI_COMM_WORLD, &r);
return std::string{b} + '-' + std::to_string(r) + '-' +
std::string(id) + "." + std::string{ext};
};
setBehaviourIntegrationFailureAnalysisFileNameGenerator(g);BehaviourIntegrationFailureAnalyser interfaceThe BehaviourIntegrationFailureAnalyser interface is
defined as follows:
struct MGIS_EXPORT BehaviourIntegrationFailureAnalyser {
virtual void analyse(const Behaviour&,
const BehaviourData&) const noexcept = 0;
virtual void analyse(const Behaviour&,
const BehaviourDataView&) const noexcept = 0;
virtual bool shallCopyBehaviourDataBeforeIntegration() const noexcept = 0;
virtual ~BehaviourIntegrationFailureAnalyser() noexcept;
};The two analyse functions do the actual work. The one
effectively called depends on the value returned by the
shallCopyBehaviourDataBeforeIntegration method.
If the shallCopyBehaviourDataBeforeIntegration return
true, the behaviour data view passed to the
integrate_debug function is copied before integrating the
behaviour using integrate:
integrate_debug shall be
almost as efficient than calling integrate.Note
The default analyser does a copy of the data, which indeed incurs an overhead.