OGS
MaterialLib::Solids::MFront Namespace Reference

Namespaces

namespace  detail

Classes

struct  DeformationGradient
 Meta data for deformation gradient. More...
struct  ForcesGradsCombinations
struct  GreenLagrangeStrain
 Meta data for Green-Lagrange strain. More...
struct  LiquidPressure
struct  MaterialStateVariablesMFront
class  MFront
struct  MFrontConfig
class  MFrontGeneric
class  OGSMFrontTangentOperatorBlocksView
struct  OGSMFrontTangentOperatorData
 Used for disambiguation with MFront's thermodynamic forces data. More...
struct  OGSMFrontThermodynamicForcesData
 Used for disambiguation with MFront's tangent operator blocks data. More...
class  OGSMFrontThermodynamicForcesView
struct  Saturation
struct  SecondPiolaKirchhoffStress
struct  Strain
 Meta data for strain. More...
struct  Stress
struct  Temperature
struct  Variable

Functions

template<int DisplacementDim>
std::unique_ptr< MechanicsBase< DisplacementDim > > createMFront (std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, BaseLib::ConfigTree const &config)
template std::unique_ptr< MechanicsBase< 2 > > createMFront< 2 > (std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, BaseLib::ConfigTree const &config)
template std::unique_ptr< MechanicsBase< 3 > > createMFront< 3 > (std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, BaseLib::ConfigTree const &config)
MFrontConfig createMFrontConfig (int const displacement_dim, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, BaseLib::ConfigTree const &config)
template<int DisplacementDim, typename Gradients, typename TDynForces, typename ExtStateVars>
std::unique_ptr< MFrontGeneric< DisplacementDim, Gradients, TDynForces, ExtStateVars > > createMFrontGeneric (std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, BaseLib::ConfigTree const &config)
const char * varTypeToString (int v)
int getEquivalentPlasticStrainOffset (mgis::behaviour::Behaviour const &b)
template<int DisplacementDim>
OGSMFrontTangentOperatorData tangentOperatorDataMFrontToOGS (std::vector< double > const &mfront_data, std::optional< MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > const &Q, mgis::behaviour::Behaviour const &behaviour)
template OGSMFrontTangentOperatorData tangentOperatorDataMFrontToOGS< 2 > (std::vector< double > const &mfront_data, std::optional< MathLib::KelvinVector::KelvinMatrixType< 2 > > const &Q, mgis::behaviour::Behaviour const &behaviour)
template OGSMFrontTangentOperatorData tangentOperatorDataMFrontToOGS< 3 > (std::vector< double > const &mfront_data, std::optional< MathLib::KelvinVector::KelvinMatrixType< 3 > > const &Q, mgis::behaviour::Behaviour const &behaviour)
template<typename Derived>
constexpr auto eigenSwap45View (Eigen::MatrixBase< Derived > const &matrix)
template<int DisplacementDim, typename Derived>
constexpr auto ogsTensorToMFrontTensor (Eigen::MatrixBase< Derived > const &matrix)

Variables

static constexpr Strain strain
 Instance that can be used for overload resolution/template type deduction.
static constexpr GreenLagrangeStrain green_lagrange_strain
 Instance that can be used for overload resolution/template type deduction.
static constexpr DeformationGradient deformation_gradient
 Instance that can be used for overload resolution/template type deduction.
static constexpr LiquidPressure liquid_pressure
static constexpr Stress stress
static constexpr SecondPiolaKirchhoffStress second_piola_kirchhoff_stress
static constexpr Saturation saturation
static constexpr Temperature temperature

Function Documentation

◆ createMFront()

template<int DisplacementDim>
std::unique_ptr< MechanicsBase< DisplacementDim > > MaterialLib::Solids::MFront::createMFront ( std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
BaseLib::ConfigTree const & config )

Definition at line 20 of file CreateMFront.cpp.

25{
26 auto conf = createMFrontConfig(DisplacementDim, parameters, config);
27
28 return std::make_unique<MFront<DisplacementDim>>(
29 std::move(conf.behaviour), std::move(conf.material_properties),
30 std::move(conf.state_variables_initial_properties),
31 local_coordinate_system);
32}
MFrontConfig createMFrontConfig(int const displacement_dim, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, BaseLib::ConfigTree const &config)

References createMFrontConfig().

Referenced by MaterialLib::Solids::createConstitutiveRelation().

◆ createMFront< 2 >()

template std::unique_ptr< MechanicsBase< 2 > > MaterialLib::Solids::MFront::createMFront< 2 > ( std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
BaseLib::ConfigTree const & config )

◆ createMFront< 3 >()

template std::unique_ptr< MechanicsBase< 3 > > MaterialLib::Solids::MFront::createMFront< 3 > ( std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
BaseLib::ConfigTree const & config )

◆ createMFrontConfig()

MFrontConfig MaterialLib::Solids::MFront::createMFrontConfig ( int const displacement_dim,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
BaseLib::ConfigTree const & config )
Input File Parameter
material__solid__constitutive_relation__type
Input File Parameter
material__solid__constitutive_relation__MFront__library
Input File Parameter
material__solid__constitutive_relation__MFront__library__path_is_relative_to_prj_file
Input File Parameter
material__solid__constitutive_relation__MFront__behaviour

Definition at line 318 of file CreateMFrontGeneric.cpp.

322{
323 INFO("### MFRONT ########################################################");
324
326 config.checkConfigParameter("type", "MFront");
327
329 auto const library = config.getConfigSubtreeOptional("library");
330
331 bool const library_path_is_relative_to_prj_file =
332 library /* If no library tag is specified in the prj file, the lib
333 shipped with OGS is used, whose path is not relative to the
334 project file. */
335 &&
337 library->getConfigAttribute("path_is_relative_to_prj_file", true);
338
339 std::string const library_name =
340 library ? library->getValue<std::string>() : "libOgsMFrontBehaviour";
341
342 auto const lib_path =
343 library_path_is_relative_to_prj_file
344 ? BaseLib::joinPaths(config.projectDirectory(), library_name)
345 : library_name;
346
347 mgis::behaviour::Hypothesis hypothesis;
348 if (displacement_dim == 2)
349 {
350 // TODO support the axial symmetry modelling hypothesis.
351 WARN(
352 "The model is defined in 2D. On the material level currently a "
353 "plane strain setting is used. In particular it is not checked if "
354 "axial symmetry or plane stress are assumed. Special material "
355 "behaviour for these settings is currently not supported.");
356 hypothesis = mgis::behaviour::Hypothesis::PLANESTRAIN;
357 }
358 else if (displacement_dim == 3)
359 {
360 hypothesis = mgis::behaviour::Hypothesis::TRIDIMENSIONAL;
361 }
362 else
363 {
364 OGS_FATAL("Displacement dim {} is not supported.", displacement_dim);
365 }
366
367 auto behaviour = loadBehaviour(
368 lib_path,
370 config.getConfigParameter<std::string>("behaviour"),
371 hypothesis);
372
373 INFO("Behaviour: `{:s}'.", behaviour.behaviour);
374 INFO("Hypothesis: `{:s}'.", mgis::behaviour::toString(hypothesis));
375 INFO("Source: `{:s}'.", behaviour.source);
376 INFO("TFEL version: `{:s}'.", behaviour.tfel_version);
377 INFO("Behaviour type: `{:s}'.", btypeToString(behaviour.btype));
378 INFO("Kinematic: `{:s}'.", toString(behaviour.kinematic));
379 INFO("Symmetry: `{:s}'.", toString(behaviour.symmetry));
380
381 varInfo("Mat. props.", behaviour.mps, hypothesis);
382 varInfo("Gradients", behaviour.gradients, hypothesis);
383 varInfo("Thdyn. forces", behaviour.thermodynamic_forces, hypothesis);
384 varInfo("Int. StVars.", behaviour.isvs, hypothesis);
385 varInfo("Ext. StVars.", behaviour.esvs, hypothesis);
386
387 // TODO read parameters from prj file, not yet (2018-11-05) supported by
388 // MGIS library.
389 varInfo("Real-valued parameters", behaviour.params);
390 varInfo("Integer parameters", behaviour.iparams);
391 varInfo("Unsigned parameters", behaviour.usparams);
392
393 INFO("#Tangent operator blocks: {}.", behaviour.to_blocks.size());
394 for (auto const& [var1, var2] : behaviour.to_blocks)
395 {
396 INFO(" --> ({}, {}).", var1.name, var2.name);
397 }
398
399 std::vector<ParameterLib::Parameter<double> const*> material_properties =
400 readMaterialProperties(behaviour, hypothesis, parameters, config);
401
402 std::map<std::string, ParameterLib::Parameter<double> const*>
403 state_variables_initial_properties =
404 readStateVariablesInitialValueProperties(behaviour, hypothesis,
405 parameters, config);
406
407 INFO("### MFRONT END ####################################################");
408
409 return {std::move(behaviour), std::move(material_properties),
410 std::move(state_variables_initial_properties)};
411}
#define OGS_FATAL(...)
Definition Error.h:19
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
std::string joinPaths(std::string const &pathA, std::string const &pathB)
std::vector< ParameterLib::Parameter< double > const * > readMaterialProperties(mgis::behaviour::Behaviour const &behaviour, mgis::behaviour::Hypothesis const &hypothesis, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, BaseLib::ConfigTree const &config)
std::map< std::string, ParameterLib::Parameter< double > const * > readStateVariablesInitialValueProperties(mgis::behaviour::Behaviour const &behaviour, mgis::behaviour::Hypothesis const &hypothesis, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, BaseLib::ConfigTree const &config)

References BaseLib::ConfigTree::checkConfigParameter(), BaseLib::ConfigTree::getConfigParameter(), BaseLib::ConfigTree::getConfigSubtreeOptional(), INFO(), BaseLib::joinPaths(), OGS_FATAL, BaseLib::ConfigTree::projectDirectory(), and WARN().

Referenced by createMFront(), and createMFrontGeneric().

◆ createMFrontGeneric()

template<int DisplacementDim, typename Gradients, typename TDynForces, typename ExtStateVars>
std::unique_ptr< MFrontGeneric< DisplacementDim, Gradients, TDynForces, ExtStateVars > > MaterialLib::Solids::MFront::createMFrontGeneric ( std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
BaseLib::ConfigTree const & config )

Definition at line 29 of file CreateMFrontGeneric.h.

34{
35 auto conf = createMFrontConfig(DisplacementDim, parameters, config);
36
37 return std::make_unique<
39 std::move(conf.behaviour), std::move(conf.material_properties),
40 std::move(conf.state_variables_initial_properties),
41 local_coordinate_system);
42}

References createMFrontConfig().

◆ eigenSwap45View()

template<typename Derived>
auto MaterialLib::Solids::MFront::eigenSwap45View ( Eigen::MatrixBase< Derived > const & matrix)
constexpr

Converts between OGS' and MFront's Kelvin vectors and matrices and vice versa. Numbering of the Kelvin vectors and matrices in the two cases is: MFront: 11 22 33 12 13 23 OGS: 11 22 33 12 23 13 The function swaps the 4th and 5th row and column of a matrix.

Definition at line 36 of file MFrontGeneric.h.

37{
38 using Matrix =
39 Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime,
40 Derived::ColsAtCompileTime>;
41
42 return Matrix::NullaryExpr(
43 matrix.rows(), matrix.cols(),
44 [&m = matrix.derived()](Eigen::Index const row, Eigen::Index const col)
45 {
46 constexpr std::ptrdiff_t result[6] = {0, 1, 2, 3, 5, 4};
47 return m(result[row], result[col]);
48 });
49}

Referenced by MaterialLib::Solids::MFront::detail::SetGradient< DisplacementDim >::operator()(), and tangentOperatorDataMFrontToOGS().

◆ getEquivalentPlasticStrainOffset()

int MaterialLib::Solids::MFront::getEquivalentPlasticStrainOffset ( mgis::behaviour::Behaviour const & b)

Definition at line 23 of file MFrontGeneric.cpp.

24{
25 return mgis::behaviour::contains(b.isvs, "EquivalentPlasticStrain")
26 ? mgis::behaviour::getVariableOffset(
27 b.isvs, "EquivalentPlasticStrain", b.hypothesis)
28 : -1;
29}

Referenced by MaterialLib::Solids::MFront::MFrontGeneric< DisplacementDim, boost::mp11::mp_list< Strain >, boost::mp11::mp_list< Stress >, boost::mp11::mp_list< Temperature > >::MFrontGeneric().

◆ ogsTensorToMFrontTensor()

template<int DisplacementDim, typename Derived>
auto MaterialLib::Solids::MFront::ogsTensorToMFrontTensor ( Eigen::MatrixBase< Derived > const & matrix)
constexpr

Converts between OGS' and MFront's tensors, which are represented as vectors. An OGS tensor 11 12 13 21 22 23 31 32 33 0 1 2 3 4 5 6 7 8 is converted to MFront tensor 11 22 33 12 21 13 31 23 32 0 4 8 1 3 2 6 5 7.

Definition at line 59 of file MFrontGeneric.h.

60{
61 using Matrix =
62 Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime,
63 Derived::ColsAtCompileTime>;
64
65 if constexpr (DisplacementDim == 2)
66 {
67 return Matrix::NullaryExpr(
68 matrix.rows(), matrix.cols(),
69 [&m = matrix.derived()](Eigen::Index const row,
70 Eigen::Index const col)
71 {
72 constexpr std::ptrdiff_t result[5] = {0, 3, 4, 1, 2};
73 return m(result[row], result[col]);
74 });
75 }
76 if constexpr (DisplacementDim == 3)
77 {
78 return Matrix::NullaryExpr(
79 matrix.rows(), matrix.cols(),
80 [&m = matrix.derived()](Eigen::Index const row,
81 Eigen::Index const col)
82 {
83 constexpr std::ptrdiff_t result[9] = {0, 4, 8, 1, 3,
84 2, 6, 5, 7};
85 return m(result[row], result[col]);
86 });
87 }
88}

References ogsTensorToMFrontTensor().

Referenced by ogsTensorToMFrontTensor(), and MaterialLib::Solids::MFront::detail::SetGradient< DisplacementDim >::operator()().

◆ tangentOperatorDataMFrontToOGS()

template<int DisplacementDim>
OGSMFrontTangentOperatorData MaterialLib::Solids::MFront::tangentOperatorDataMFrontToOGS ( std::vector< double > const & mfront_data,
std::optional< MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > const & Q,
mgis::behaviour::Behaviour const & behaviour )

Transforms MFront's to OGS's tangent operator data.

Essentially swaps some off-diagonal components of symmetric tensors and applies the given rotation of the frame of reference Q.

Definition at line 32 of file MFrontGeneric.cpp.

37{
38 using VT = mgis::behaviour::Variable::Type;
41
42 std::vector<double> ogs_data(mfront_data.size());
43
44 std::size_t offset = 0;
45 constexpr auto kv_size =
47
48 for (auto const& [var1, var2] : behaviour.to_blocks)
49 {
50 std::size_t size;
51
52 auto const* const d_in = mfront_data.data() + offset;
53 auto* const d_out = ogs_data.data() + offset;
54
55 if (var1.type == VT::SCALAR && var2.type == VT::SCALAR)
56 {
57 size = 1;
58 *d_out = *d_in;
59 }
60 else if (var1.type == VT::SCALAR && var2.type == VT::STENSOR)
61 {
62 assert(getVariableSize(var2, behaviour.hypothesis) == kv_size);
63 size = kv_size;
64
65 if (Q)
66 {
68 "Coordinate frame rotation not yet implemented for "
69 "dScalar/dSTensor.");
70 }
71
72 Eigen::Map<KV>{d_out} = eigenSwap45View(Eigen::Map<const KV>{d_in});
73 }
74 else if (var1.type == VT::STENSOR && var2.type == VT::SCALAR)
75 {
76 assert(getVariableSize(var1, behaviour.hypothesis) == kv_size);
77 size = kv_size;
78
79 if (Q)
80 {
82 "Coordinate frame rotation not yet implemented for "
83 "dSTensor/dScalar.");
84 }
85
86 Eigen::Map<KV>{d_out} = eigenSwap45View(Eigen::Map<const KV>{d_in});
87 }
88 else if (var1.type == VT::STENSOR && var2.type == VT::STENSOR)
89 {
90 assert(getVariableSize(var1, behaviour.hypothesis) == kv_size);
91 assert(getVariableSize(var2, behaviour.hypothesis) == kv_size);
92 size = kv_size * kv_size;
93
94 if (Q)
95 {
96 Eigen::Map<KM>{d_out} =
97 *Q * eigenSwap45View(Eigen::Map<const KM>{d_in}) *
98 Q->transpose();
99 }
100 else
101 {
102 Eigen::Map<KM>{d_out} =
103 eigenSwap45View(Eigen::Map<const KM>{d_in});
104 }
105 }
106 else
107 {
108 OGS_FATAL("unsupported variable type combination");
109 }
110
111 offset += size;
112 }
113
114 return {std::move(ogs_data)};
115}
constexpr auto eigenSwap45View(Eigen::MatrixBase< Derived > const &matrix)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
constexpr int size(int const displacement_dim)
Vectorized tensor size for given displacement dimension.

References eigenSwap45View(), MathLib::KelvinVector::kelvin_vector_dimensions(), and OGS_FATAL.

◆ tangentOperatorDataMFrontToOGS< 2 >()

template OGSMFrontTangentOperatorData MaterialLib::Solids::MFront::tangentOperatorDataMFrontToOGS< 2 > ( std::vector< double > const & mfront_data,
std::optional< MathLib::KelvinVector::KelvinMatrixType< 2 > > const & Q,
mgis::behaviour::Behaviour const & behaviour )

◆ tangentOperatorDataMFrontToOGS< 3 >()

template OGSMFrontTangentOperatorData MaterialLib::Solids::MFront::tangentOperatorDataMFrontToOGS< 3 > ( std::vector< double > const & mfront_data,
std::optional< MathLib::KelvinVector::KelvinMatrixType< 3 > > const & Q,
mgis::behaviour::Behaviour const & behaviour )

◆ varTypeToString()

const char * MaterialLib::Solids::MFront::varTypeToString ( int v)

Definition at line 8 of file MFrontGeneric.cpp.

9{
10 using V = mgis::behaviour::Variable;
11 if (v == V::SCALAR)
12 return "SCALAR";
13 if (v == V::VECTOR)
14 return "VECTOR";
15 if (v == V::STENSOR)
16 return "STENSOR";
17 if (v == V::TENSOR)
18 return "TENSOR";
19
20 OGS_FATAL("Unknown variable type {}.", v);
21}

References OGS_FATAL.

Referenced by MaterialLib::Solids::MFront::MFrontGeneric< DisplacementDim, boost::mp11::mp_list< Strain >, boost::mp11::mp_list< Stress >, boost::mp11::mp_list< Temperature > >::MFrontGeneric(), MaterialLib::Solids::MFront::detail::SetGradient< DisplacementDim >::operator()(), and anonymous_namespace{CreateMFrontGeneric.cpp}::varInfo().

Variable Documentation

◆ deformation_gradient

DeformationGradient MaterialLib::Solids::MFront::deformation_gradient
staticconstexpr

Instance that can be used for overload resolution/template type deduction.

Definition at line 151 of file Variable.h.

◆ green_lagrange_strain

GreenLagrangeStrain MaterialLib::Solids::MFront::green_lagrange_strain
staticconstexpr

Instance that can be used for overload resolution/template type deduction.

Definition at line 133 of file Variable.h.

◆ liquid_pressure

LiquidPressure MaterialLib::Solids::MFront::liquid_pressure
staticconstexpr

Definition at line 164 of file Variable.h.

◆ saturation

Saturation MaterialLib::Solids::MFront::saturation
staticconstexpr

Definition at line 201 of file Variable.h.

◆ second_piola_kirchhoff_stress

SecondPiolaKirchhoffStress MaterialLib::Solids::MFront::second_piola_kirchhoff_stress
staticconstexpr

Definition at line 188 of file Variable.h.

◆ strain

Strain MaterialLib::Solids::MFront::strain
staticconstexpr

Instance that can be used for overload resolution/template type deduction.

Definition at line 113 of file Variable.h.

Referenced by MaterialLib::Solids::MFront::MFront< DisplacementDim >::integrateStress().

◆ stress

Stress MaterialLib::Solids::MFront::stress
staticconstexpr

◆ temperature

Temperature MaterialLib::Solids::MFront::temperature
staticconstexpr

Definition at line 214 of file Variable.h.