OGS
MaterialLib::Solids::MFront::MFront< DisplacementDim > Class Template Reference

Detailed Description

template<int DisplacementDim>
class MaterialLib::Solids::MFront::MFront< DisplacementDim >

Uses a material model provided by MFront (via MFront's generic interface and the MGIS library).

Definition at line 40 of file MFront.h.

#include <MFront.h>

Inheritance diagram for MaterialLib::Solids::MFront::MFront< DisplacementDim >:
[legend]
Collaboration diagram for MaterialLib::Solids::MFront::MFront< DisplacementDim >:
[legend]

Classes

struct  MaterialStateVariables
 

Public Types

using KelvinVector = MathLib::KelvinVector::KelvinVectorType< DisplacementDim >
 
using KelvinMatrix = MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >
 
- Public Types inherited from MaterialLib::Solids::MechanicsBase< DisplacementDim >
using KelvinVector = MathLib::KelvinVector::KelvinVectorType< DisplacementDim >
 
using KelvinMatrix = MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >
 

Public Member Functions

 MFront (mgis::behaviour::Behaviour &&behaviour, std::vector< ParameterLib::Parameter< double > const * > &&material_properties, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system)
 
std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariablescreateMaterialStateVariables () const override
 
std::optional< std::tuple< KelvinVector, std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables >, KelvinMatrix > > integrateStress (MaterialPropertyLib::VariableArray const &variable_array_prev, MaterialPropertyLib::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x, double const dt, typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &material_state_variables) const override
 
std::vector< typename MechanicsBase< DisplacementDim >::InternalVariablegetInternalVariables () const override
 
double getBulkModulus (double const, ParameterLib::SpatialPosition const &, KelvinMatrix const *const) const override
 
double computeFreeEnergyDensity (double const t, ParameterLib::SpatialPosition const &x, double const dt, KelvinVector const &eps, KelvinVector const &sigma, typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &material_state_variables) const override
 
- Public Member Functions inherited from MaterialLib::Solids::MechanicsBase< DisplacementDim >
virtual std::optional< std::tuple< KelvinVector, std::unique_ptr< MaterialStateVariables >, KelvinMatrix > > integrateStress (MaterialPropertyLib::VariableArray const &variable_array_prev, MaterialPropertyLib::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x, double const dt, MaterialStateVariables const &material_state_variables) const =0
 
virtual ConstitutiveModel getConstitutiveModel () const
 Gets the type of constitutive model. More...
 
virtual double getTemperatureRelatedCoefficient (double const, double const, ParameterLib::SpatialPosition const &, double const, double const) const
 
virtual double computeFreeEnergyDensity (double const t, ParameterLib::SpatialPosition const &x, double const dt, KelvinVector const &eps, KelvinVector const &sigma, MaterialStateVariables const &material_state_variables) const =0
 
virtual ~MechanicsBase ()=default
 

Private Attributes

mgis::behaviour::Behaviour _behaviour
 
int const equivalent_plastic_strain_offset_
 
std::vector< ParameterLib::Parameter< double > const * > _material_properties
 
ParameterLib::CoordinateSystem const *const _local_coordinate_system
 

Member Typedef Documentation

◆ KelvinMatrix

template<int DisplacementDim>
using MaterialLib::Solids::MFront::MFront< DisplacementDim >::KelvinMatrix = MathLib::KelvinVector::KelvinMatrixType<DisplacementDim>

Definition at line 71 of file MFront.h.

◆ KelvinVector

template<int DisplacementDim>
using MaterialLib::Solids::MFront::MFront< DisplacementDim >::KelvinVector = MathLib::KelvinVector::KelvinVectorType<DisplacementDim>

Definition at line 69 of file MFront.h.

Constructor & Destructor Documentation

◆ MFront()

template<int DisplacementDim>
MaterialLib::Solids::MFront::MFront< DisplacementDim >::MFront ( mgis::behaviour::Behaviour &&  behaviour,
std::vector< ParameterLib::Parameter< double > const * > &&  material_properties,
std::optional< ParameterLib::CoordinateSystem > const &  local_coordinate_system 
)

Definition at line 174 of file MFront.cpp.

179  : _behaviour(std::move(behaviour)),
182  _material_properties(std::move(material_properties)),
184  local_coordinate_system ? &local_coordinate_system.value() : nullptr)
185 {
186  auto const hypothesis = behaviour.hypothesis;
187 
188  if (_behaviour.gradients.size() != 1)
189  OGS_FATAL(
190  "The behaviour must have exactly a single gradient as input.");
191 
192  if (_behaviour.gradients[0].name != "Strain")
193  OGS_FATAL("The behaviour must be driven by strain.");
194 
195  if (_behaviour.gradients[0].type != mgis::behaviour::Variable::STENSOR)
196  OGS_FATAL("Strain must be a symmetric tensor.");
197 
198  if (mgis::behaviour::getVariableSize(_behaviour.gradients[0], hypothesis) !=
199  MFront<DisplacementDim>::KelvinVector::SizeAtCompileTime)
200  OGS_FATAL("Strain must have {:d} components instead of {:d}.",
201  MFront<DisplacementDim>::KelvinVector::SizeAtCompileTime,
202  mgis::behaviour::getVariableSize(_behaviour.gradients[0],
203  hypothesis));
204 
205  if (_behaviour.thermodynamic_forces.size() != 1)
206  OGS_FATAL(
207  "The behaviour must compute exactly one thermodynamic force.");
208 
209  if (_behaviour.thermodynamic_forces[0].name != "Stress")
210  OGS_FATAL("The behaviour must compute stress.");
211 
212  if (_behaviour.thermodynamic_forces[0].type !=
213  mgis::behaviour::Variable::STENSOR)
214  OGS_FATAL("Stress must be a symmetric tensor.");
215 
216  if (mgis::behaviour::getVariableSize(_behaviour.thermodynamic_forces[0],
217  hypothesis) !=
218  MFront<DisplacementDim>::KelvinVector::SizeAtCompileTime)
219  OGS_FATAL("Stress must have {:d} components instead of {:d}.",
220  MFront<DisplacementDim>::KelvinVector::SizeAtCompileTime,
221  mgis::behaviour::getVariableSize(
222  _behaviour.thermodynamic_forces[0], hypothesis));
223 
224  if (!_behaviour.esvs.empty())
225  {
226  if (_behaviour.esvs[0].name != "Temperature")
227  {
228  OGS_FATAL(
229  "Only temperature is supported as external state variable.");
230  }
231 
232  if (mgis::behaviour::getVariableSize(_behaviour.esvs[0], hypothesis) !=
233  1)
234  OGS_FATAL(
235  "Temperature must be a scalar instead of having {:d} "
236  "components.",
237  mgis::behaviour::getVariableSize(
238  _behaviour.thermodynamic_forces[0], hypothesis));
239  }
240 
241  if (_behaviour.mps.size() != _material_properties.size())
242  {
243  ERR("There are {:d} material properties in the loaded behaviour:",
244  _behaviour.mps.size());
245  for (auto const& mp : _behaviour.mps)
246  {
247  ERR("\t{:s}", mp.name);
248  }
249  OGS_FATAL("But the number of passed material properties is {:d}.",
250  _material_properties.size());
251  }
252 }
#define OGS_FATAL(...)
Definition: Error.h:26
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
mgis::behaviour::Behaviour _behaviour
Definition: MFront.h:114
std::vector< ParameterLib::Parameter< double > const * > _material_properties
Definition: MFront.h:116
int const equivalent_plastic_strain_offset_
Definition: MFront.h:115
ParameterLib::CoordinateSystem const *const _local_coordinate_system
Definition: MFront.h:117
int getEquivalentPlasticStrainOffset(mgis::behaviour::Behaviour const &b)
Definition: MFront.cpp:165

References MaterialLib::Solids::MFront::MFront< DisplacementDim >::_behaviour, MaterialLib::Solids::MFront::MFront< DisplacementDim >::_material_properties, ERR(), and OGS_FATAL.

Member Function Documentation

◆ computeFreeEnergyDensity()

template<int DisplacementDim>
double MaterialLib::Solids::MFront::MFront< DisplacementDim >::computeFreeEnergyDensity ( double const  t,
ParameterLib::SpatialPosition const &  x,
double const  dt,
KelvinVector const &  eps,
KelvinVector const &  sigma,
typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &  material_state_variables 
) const
override

Definition at line 493 of file MFront.cpp.

501 {
502  // TODO implement
503  return std::numeric_limits<double>::quiet_NaN();
504 }

◆ createMaterialStateVariables()

template<int DisplacementDim>
std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables > MaterialLib::Solids::MFront::MFront< DisplacementDim >::createMaterialStateVariables
overridevirtual

Polymorphic creator for MaterialStateVariables objects specific for a material model.

Reimplemented from MaterialLib::Solids::MechanicsBase< DisplacementDim >.

Definition at line 256 of file MFront.cpp.

257 {
258  return std::make_unique<MaterialStateVariables>(
260 }

◆ getBulkModulus()

template<int DisplacementDim>
double MaterialLib::Solids::MFront::MFront< DisplacementDim >::getBulkModulus ( double const  ,
ParameterLib::SpatialPosition const &  ,
KelvinMatrix const * const  C 
) const
overridevirtual

Reimplemented from MaterialLib::Solids::MechanicsBase< DisplacementDim >.

Definition at line 475 of file MFront.cpp.

479 {
480  if (C == nullptr)
481  {
482  OGS_FATAL(
483  "MFront::getBulkModulus() requires the tangent stiffness C input "
484  "argument to be valid.");
485  }
486  auto const& identity2 = MathLib::KelvinVector::Invariants<
488  DisplacementDim)>::identity2;
489  return 1. / 9. * identity2.transpose() * *C * identity2;
490 }
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23

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

◆ getInternalVariables()

template<int DisplacementDim>
std::vector< typename MechanicsBase< DisplacementDim >::InternalVariable > MaterialLib::Solids::MFront::MFront< DisplacementDim >::getInternalVariables
overridevirtual

Returns internal variables defined by the specific material model, if any.

Reimplemented from MaterialLib::Solids::MechanicsBase< DisplacementDim >.

Definition at line 420 of file MFront.cpp.

421 {
422  std::vector<typename MechanicsBase<DisplacementDim>::InternalVariable>
423  internal_variables;
424 
425  for (auto const& iv : _behaviour.isvs)
426  {
427  auto const name = iv.name;
428  auto const offset = mgis::behaviour::getVariableOffset(
429  _behaviour.isvs, name, _behaviour.hypothesis);
430  auto const size =
431  mgis::behaviour::getVariableSize(iv, _behaviour.hypothesis);
432 
433  // TODO (naumov): For orthotropic materials the internal variables
434  // should be rotated to the global coordinate system before output.
435  // MFront stores the variables in local coordinate system.
436  // The `size` variable could be used to find out the type of variable.
437  typename MechanicsBase<DisplacementDim>::InternalVariable new_variable{
438  name, static_cast<int>(size),
439  [offset, size](
440  typename MechanicsBase<
441  DisplacementDim>::MaterialStateVariables const& state,
442  std::vector<double>& cache) -> std::vector<double> const&
443  {
444  assert(dynamic_cast<MaterialStateVariables const*>(&state) !=
445  nullptr);
446  auto const& internal_state_variables =
447  static_cast<MaterialStateVariables const&>(state)
448  ._behaviour_data.s1.internal_state_variables;
449 
450  cache.resize(size);
451  std::copy_n(internal_state_variables.data() + offset,
452  size,
453  begin(cache));
454  return cache;
455  },
456  [offset, size](
457  typename MechanicsBase<DisplacementDim>::MaterialStateVariables&
459  {
460  assert(dynamic_cast<MaterialStateVariables const*>(&state) !=
461  nullptr);
462  auto& internal_state_variables =
463  static_cast<MaterialStateVariables&>(state)
464  ._behaviour_data.s1.internal_state_variables;
465 
466  return {internal_state_variables.data() + offset, size};
467  }};
468  internal_variables.push_back(new_variable);
469  }
470 
471  return internal_variables;
472 }

References MaterialPropertyLib::name.

◆ integrateStress()

template<int DisplacementDim>
std::optional< std::tuple< typename MFront< DisplacementDim >::KelvinVector, std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables >, typename MFront< DisplacementDim >::KelvinMatrix > > MaterialLib::Solids::MFront::MFront< DisplacementDim >::integrateStress ( MaterialPropertyLib::VariableArray const &  variable_array_prev,
MaterialPropertyLib::VariableArray const &  variable_array,
double const  t,
ParameterLib::SpatialPosition const &  x,
double const  dt,
typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &  material_state_variables 
) const
override

Definition at line 267 of file MFront.cpp.

275 {
276  using namespace MathLib::KelvinVector;
277 
278  assert(
279  dynamic_cast<MaterialStateVariables const*>(&material_state_variables));
280  // New state, copy of current one, packed in unique_ptr for return.
281  auto state = std::make_unique<MaterialStateVariables>(
282  static_cast<MaterialStateVariables const&>(material_state_variables));
283  auto& behaviour_data = state->_behaviour_data;
284 
285  // TODO add a test of material behaviour where the value of dt matters.
286  behaviour_data.dt = dt;
287  behaviour_data.rdt = 1.0;
288  behaviour_data.K[0] = 4.0; // if K[0] is greater than 3.5, the consistent
289  // tangent operator must be computed.
290 
291  // evaluate parameters at (t, x)
292  {
293  auto out = behaviour_data.s1.material_properties.begin();
294  for (auto* param : _material_properties)
295  {
296  auto const& vals = (*param)(t, x);
297  out = std::copy(vals.begin(), vals.end(), out);
298  }
299  assert(out == behaviour_data.s1.material_properties.end());
300  }
301 
302  if (!behaviour_data.s0.external_state_variables.empty())
303  {
304  // assuming that there is only temperature
305  if (!std::holds_alternative<double>(
306  variable_array_prev[static_cast<int>(
308  {
309  OGS_FATAL(
310  "MPL::Variable::temperature is not set for variable_array_prev "
311  "before calling MFront::integrateStress.");
312  }
313 
314  behaviour_data.s1.external_state_variables[0] = std::get<double>(
315  variable_array_prev[static_cast<int>(MPL::Variable::temperature)]);
316  }
317 
318  if (!behaviour_data.s1.external_state_variables.empty())
319  {
320  // assuming that there is only temperature
321  if (!std::holds_alternative<double>(
322  variable_array[static_cast<int>(MPL::Variable::temperature)]))
323  {
324  OGS_FATAL(
325  "MPL::Variable::temperature is not set for variable_array "
326  "before calling MFront::integrateStress");
327  }
328 
329  behaviour_data.s1.external_state_variables[0] = std::get<double>(
330  variable_array[static_cast<int>(MPL::Variable::temperature)]);
331  }
332 
333  // rotation tensor
334  auto const Q = [this, &x]() -> KelvinMatrixType<DisplacementDim>
335  {
337  {
339  }
341  _local_coordinate_system->transformation<DisplacementDim>(x));
342  }();
343 
344  auto const& eps_m_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
345  variable_array_prev[static_cast<int>(
346  MPL::Variable::mechanical_strain)]);
347  auto const eps_prev_MFront = OGSToMFront(
348  Q.transpose()
349  .template topLeftCorner<kelvin_vector_dimensions(DisplacementDim),
351  DisplacementDim)>() *
352  eps_m_prev);
353  std::copy_n(eps_prev_MFront.data(), KelvinVector::SizeAtCompileTime,
354  behaviour_data.s0.gradients.data());
355 
356  auto const& eps = std::get<MPL::SymmetricTensor<DisplacementDim>>(
357  variable_array[static_cast<int>(MPL::Variable::mechanical_strain)]);
358  auto const eps_MFront = OGSToMFront(
359  Q.transpose()
360  .template topLeftCorner<kelvin_vector_dimensions(DisplacementDim),
362  DisplacementDim)>() *
363  eps);
364  std::copy_n(eps_MFront.data(), KelvinVector::SizeAtCompileTime,
365  behaviour_data.s1.gradients.data());
366 
367  auto const& sigma_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
368  variable_array_prev[static_cast<int>(MPL::Variable::stress)]);
369  auto const sigma_prev_MFront = OGSToMFront(
370  Q.transpose()
371  .template topLeftCorner<kelvin_vector_dimensions(DisplacementDim),
373  DisplacementDim)>() *
374  sigma_prev);
375  std::copy_n(sigma_prev_MFront.data(), KelvinVector::SizeAtCompileTime,
376  behaviour_data.s0.thermodynamic_forces.data());
377  std::copy_n(sigma_prev_MFront.data(), KelvinVector::SizeAtCompileTime,
378  behaviour_data.s1.thermodynamic_forces.data());
379 
380  auto v = mgis::behaviour::make_view(behaviour_data);
381  auto const status = mgis::behaviour::integrate(v, _behaviour);
382  if (status != 1)
383  {
385  "MFront: integration failed with status" + std::to_string(status) +
386  ".");
387  }
388 
389  KelvinVector sigma;
390  std::copy_n(behaviour_data.s1.thermodynamic_forces.data(),
391  KelvinVector::SizeAtCompileTime, sigma.data());
392  sigma =
393  Q.template topLeftCorner<kelvin_vector_dimensions(DisplacementDim),
394  kelvin_vector_dimensions(DisplacementDim)>() *
395  MFrontToOGS(sigma);
396 
397  // TODO row- vs. column-major storage order. This should only matter for
398  // anisotropic materials.
399  if (behaviour_data.K.size() !=
400  KelvinMatrix::RowsAtCompileTime * KelvinMatrix::ColsAtCompileTime)
401  OGS_FATAL("Stiffness matrix has wrong size.");
402 
403  KelvinMatrix C =
404  Q * MFrontToOGS(Eigen::Map<KelvinMatrix>(behaviour_data.K.data())) *
405  Q.transpose()
406  .template topLeftCorner<kelvin_vector_dimensions(DisplacementDim),
408  DisplacementDim)>();
409 
410  return std::make_optional(
411  std::make_tuple<typename MFront<DisplacementDim>::KelvinVector,
412  std::unique_ptr<typename MechanicsBase<
413  DisplacementDim>::MaterialStateVariables>,
415  std::move(sigma), std::move(state), std::move(C)));
416 }
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > KelvinVector
Definition: MFront.h:70
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > KelvinMatrix
Definition: MFront.h:72
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Definition: KelvinVector.h:56
KelvinMatrixType< DisplacementDim > fourthOrderRotationMatrix(Eigen::Matrix< double, DisplacementDim, DisplacementDim, Eigen::ColMajor, DisplacementDim, DisplacementDim > const &transformation)
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
Derived::PlainObject OGSToMFront(Eigen::DenseBase< Derived > const &m)
Converts between OGSes and MFront's Kelvin vectors and matrices.
Definition: MFront.cpp:45
Derived::PlainObject MFrontToOGS(Eigen::DenseBase< Derived > const &m)
Converts between OGSes and MFront's Kelvin vectors and matrices.
Definition: MFront.cpp:71
Eigen::Matrix< double, Dimension, Dimension > transformation(SpatialPosition const &pos) const

References MathLib::LinAlg::copy(), MathLib::KelvinVector::fourthOrderRotationMatrix(), MathLib::KelvinVector::kelvin_vector_dimensions(), anonymous_namespace{MFront.cpp}::MFrontToOGS(), OGS_FATAL, anonymous_namespace{MFront.cpp}::OGSToMFront(), and MaterialPropertyLib::temperature.

Member Data Documentation

◆ _behaviour

template<int DisplacementDim>
mgis::behaviour::Behaviour MaterialLib::Solids::MFront::MFront< DisplacementDim >::_behaviour
private

◆ _local_coordinate_system

template<int DisplacementDim>
ParameterLib::CoordinateSystem const* const MaterialLib::Solids::MFront::MFront< DisplacementDim >::_local_coordinate_system
private

Definition at line 117 of file MFront.h.

◆ _material_properties

template<int DisplacementDim>
std::vector<ParameterLib::Parameter<double> const*> MaterialLib::Solids::MFront::MFront< DisplacementDim >::_material_properties
private

◆ equivalent_plastic_strain_offset_

template<int DisplacementDim>
int const MaterialLib::Solids::MFront::MFront< DisplacementDim >::equivalent_plastic_strain_offset_
private

Definition at line 115 of file MFront.h.


The documentation for this class was generated from the following files: