OGS
MaterialLib::Solids::MFront Namespace Reference

Classes

class  MFront
 

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)
 
const char * toString (mgis::behaviour::Behaviour::Kinematic)
 Converts MGIS kinematic to a string representation. More...
 
const char * toString (mgis::behaviour::Behaviour::Symmetry)
 Converts MGIS symmetry to a string representation. More...
 
const char * btypeToString (int)
 Converts MGIS behaviour type to a string representation. More...
 
const char * varTypeToString (int)
 Converts MGIS variable type to a string representation. More...
 
int getEquivalentPlasticStrainOffset (mgis::behaviour::Behaviour const &b)
 

Function Documentation

◆ btypeToString()

const char * MaterialLib::Solids::MFront::btypeToString ( int  btype)

Converts MGIS behaviour type to a string representation.

Definition at line 136 of file MFront.cpp.

137 {
138  using B = mgis::behaviour::Behaviour;
139  if (btype == B::GENERALBEHAVIOUR)
140  return "GENERALBEHAVIOUR";
141  if (btype == B::STANDARDSTRAINBASEDBEHAVIOUR)
142  return "STANDARDSTRAINBASEDBEHAVIOUR";
143  if (btype == B::STANDARDFINITESTRAINBEHAVIOUR)
144  return "STANDARDFINITESTRAINBEHAVIOUR";
145  if (btype == B::COHESIVEZONEMODEL)
146  return "COHESIVEZONEMODEL";
147 
148  OGS_FATAL("Unknown behaviour type {:d}.", btype);
149 }
#define OGS_FATAL(...)
Definition: Error.h:26

References OGS_FATAL.

Referenced by createMFront().

◆ 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 
)
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__behaviour
Input File Parameter:
material__solid__constitutive_relation__MFront__material_properties
Input File Parameter:
material__solid__constitutive_relation__MFront__material_properties__material_property
Input File Parameter:
material__solid__constitutive_relation__MFront__material_properties__material_property__name
Input File Parameter:
material__solid__constitutive_relation__MFront__material_properties__material_property__parameter

Definition at line 59 of file CreateMFront.cpp.

64 {
65  INFO("### MFRONT ########################################################");
66 
68  config.checkConfigParameter("type", "MFront");
69 
70  auto const library_name =
72  config.getConfigParameterOptional<std::string>("library");
73  auto const lib_path =
74  library_name
76  : "libOgsMFrontBehaviour";
77 
78  auto const behaviour_name =
80  config.getConfigParameter<std::string>("behaviour");
81 
82  static_assert(DisplacementDim == 2 || DisplacementDim == 3,
83  "Given DisplacementDim not supported.");
84 
85  mgis::behaviour::Hypothesis hypothesis;
86  if (DisplacementDim == 2)
87  {
88  // TODO support the axial symmetry modelling hypothesis.
89  WARN(
90  "The model is defined in 2D. On the material level currently a "
91  "plane strain setting is used. In particular it is not checked if "
92  "axial symmetry or plane stress are assumed. Special material "
93  "behaviour for these settings is currently not supported.");
94  hypothesis = mgis::behaviour::Hypothesis::PLANESTRAIN;
95  }
96  else if (DisplacementDim == 3)
97  {
98  hypothesis = mgis::behaviour::Hypothesis::TRIDIMENSIONAL;
99  }
100 
101  // Fix for https://gitlab.opengeosys.org/ogs/ogs/-/issues/3073
102  // Pre-load dependencies of mfront lib
103 #ifndef _WIN32
104  dlopen("libTFELNUMODIS.so", RTLD_NOW);
105  dlopen("libTFELUtilities.so", RTLD_NOW);
106  dlopen("libTFELException.so", RTLD_NOW);
107 #endif
108 
109  auto behaviour =
110  mgis::behaviour::load(lib_path, behaviour_name, hypothesis);
111 
112  INFO("Behaviour: `{:s}'.", behaviour.behaviour);
113  INFO("Hypothesis: `{:s}'.", mgis::behaviour::toString(hypothesis));
114  INFO("Source: `{:s}'.", behaviour.source);
115  INFO("TFEL version: `{:s}'.", behaviour.tfel_version);
116  INFO("Behaviour type: `{:s}'.", btypeToString(behaviour.btype));
117  INFO("Kinematic: `{:s}'.", toString(behaviour.kinematic));
118  INFO("Symmetry: `{:s}'.", toString(behaviour.symmetry));
119 
120  varInfo("Mat. props.", behaviour.mps, hypothesis);
121  varInfo("Gradients", behaviour.gradients, hypothesis);
122  varInfo("Thdyn. forces", behaviour.thermodynamic_forces, hypothesis);
123  varInfo("Int. StVars.", behaviour.isvs, hypothesis);
124  varInfo("Ext. StVars.", behaviour.esvs, hypothesis);
125 
126  // TODO read parameters from prj file, not yet (2018-11-05) supported by
127  // MGIS library.
128  varInfo("Real-valued parameters", behaviour.params);
129  varInfo("Integer parameters", behaviour.iparams);
130  varInfo("Unsigned parameters", behaviour.usparams);
131 
132  std::vector<ParameterLib::Parameter<double> const*> material_properties;
133 
134  if (!behaviour.mps.empty())
135  {
136  std::map<std::string, std::string> map_name_to_param;
137 
138  // gather material properties from the prj file
140  auto const mps_config = config.getConfigSubtree("material_properties");
141  for (
142  auto const mp_config :
144  mps_config.getConfigParameterList("material_property"))
145  {
147  auto name = mp_config.getConfigAttribute<std::string>("name");
148  auto const param_name =
150  mp_config.getConfigAttribute<std::string>("parameter");
151 
152  map_name_to_param.emplace(std::move(name), std::move(param_name));
153  }
154 
155  for (auto& mp : behaviour.mps)
156  {
157  auto const it = map_name_to_param.find(mp.name);
158  if (it == map_name_to_param.end())
159  OGS_FATAL(
160  "Material Property `{:s}' has not been configured in the "
161  "project file.",
162  mp.name);
163 
164  auto const param_name = it->second;
165  auto const num_comp =
166  mgis::behaviour::getVariableSize(mp, hypothesis);
167  auto const* param = &ParameterLib::findParameter<double>(
168  param_name, parameters, num_comp);
169 
170  INFO("Using OGS parameter `{:s}' for material property `{:s}'.",
171  param_name, mp.name);
172 
174  if (mp.type == V::STENSOR || mp.type == V::TENSOR)
175  {
176  WARN(
177  "Material property `{:s}' is a tensorial quantity. You, "
178  "the "
179  "user, have to make sure that the component order of "
180  "parameter `{:s}' matches the one required by MFront!",
181  mp.name, param_name);
182  }
183 
184  material_properties.push_back(param);
185  map_name_to_param.erase(it);
186  }
187 
188  if (!map_name_to_param.empty())
189  {
190  ERR("Some material parameters that were configured are not used by "
191  "the material model.");
192  ERR("These parameters are:");
193 
194  for (auto& e : map_name_to_param)
195  {
196  ERR(" name: `{:s}', parameter: `{:s}'.", e.first, e.second);
197  }
198 
199  OGS_FATAL(
200  "Configuration errors occurred. Please fix the project file.");
201  }
202  }
203 
204  INFO("### MFRONT END ####################################################");
205 
206  return std::make_unique<MFront<DisplacementDim>>(
207  std::move(behaviour), std::move(material_properties),
208  local_coordinate_system);
209 }
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
void WARN(char const *fmt, Args const &... args)
Definition: Logging.h:37
std::string const & getProjectDirectory()
Returns the directory where the prj file resides.
Definition: FileTools.cpp:217
std::string joinPaths(std::string const &pathA, std::string const &pathB)
Definition: FileTools.cpp:212
const char * btypeToString(int btype)
Converts MGIS behaviour type to a string representation.
Definition: MFront.cpp:136
const char * toString(mgis::behaviour::Behaviour::Kinematic kin)
Converts MGIS kinematic to a string representation.
Definition: MFront.cpp:103
void varInfo(std::string const &msg, std::vector< std::string > const &parameters)
Prints info about MFront parameters.

References btypeToString(), BaseLib::ConfigTree::checkConfigParameter(), ERR(), BaseLib::ConfigTree::getConfigParameter(), BaseLib::ConfigTree::getConfigParameterOptional(), BaseLib::ConfigTree::getConfigSubtree(), BaseLib::getProjectDirectory(), INFO(), BaseLib::joinPaths(), MaterialPropertyLib::name, OGS_FATAL, toString(), anonymous_namespace{CreateMFront.cpp}::varInfo(), and WARN().

◆ 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 
)

◆ getEquivalentPlasticStrainOffset()

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

Definition at line 165 of file MFront.cpp.

166 {
167  return mgis::behaviour::contains(b.isvs, "EquivalentPlasticStrain")
168  ? mgis::behaviour::getVariableOffset(
169  b.isvs, "EquivalentPlasticStrain", b.hypothesis)
170  : -1;
171 }
bool contains(Container const &container, typename Container::value_type const &element)
Definition: Algorithm.h:256

References BaseLib::contains().

◆ toString() [1/2]

const char * MaterialLib::Solids::MFront::toString ( mgis::behaviour::Behaviour::Kinematic  kin)

Converts MGIS kinematic to a string representation.

Definition at line 103 of file MFront.cpp.

104 {
105  using K = mgis::behaviour::Behaviour::Kinematic;
106  switch (kin)
107  {
108  case K::UNDEFINEDKINEMATIC:
109  return "UNDEFINEDKINEMATIC";
110  case K::SMALLSTRAINKINEMATIC:
111  return "SMALLSTRAINKINEMATIC";
112  case K::COHESIVEZONEKINEMATIC:
113  return "COHESIVEZONEKINEMATIC";
114  case K::FINITESTRAINKINEMATIC_F_CAUCHY:
115  return "FINITESTRAINKINEMATIC_F_CAUCHY";
116  case K::FINITESTRAINKINEMATIC_ETO_PK1:
117  return "FINITESTRAINKINEMATIC_ETO_PK1";
118  }
119 
120  OGS_FATAL("Unknown kinematic {:d}.", kin);
121 }

References OGS_FATAL.

Referenced by MainWindow::about(), ProcessLib::LIE::PostProcessTool::copyPropertyValues(), createMFront(), ProcessLib::LIE::PostProcessTool::createProperty(), MeshLib::Properties::getPropertyVector(), DiagramPrefsDialog::loadFile(), GeoTreeModel::removeGeoList(), StationTreeModel::vtkSource(), GeoTreeModel::vtkSource(), and ApplicationUtils::writeProperties().

◆ toString() [2/2]

const char * MaterialLib::Solids::MFront::toString ( mgis::behaviour::Behaviour::Symmetry  sym)

Converts MGIS symmetry to a string representation.

Definition at line 123 of file MFront.cpp.

124 {
125  using S = mgis::behaviour::Behaviour::Symmetry;
126  switch (sym)
127  {
128  case S::ISOTROPIC:
129  return "ISOTROPIC";
130  case S::ORTHOTROPIC:
131  return "ORTHOTROPIC";
132  }
133 
134  OGS_FATAL("Unknown symmetry {:d}.", sym);
135 }

References OGS_FATAL.

◆ varTypeToString()

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

Converts MGIS variable type to a string representation.

Definition at line 150 of file MFront.cpp.

151 {
153  if (v == V::SCALAR)
154  return "SCALAR";
155  if (v == V::VECTOR)
156  return "VECTOR";
157  if (v == V::STENSOR)
158  return "STENSOR";
159  if (v == V::TENSOR)
160  return "TENSOR";
161 
162  OGS_FATAL("Unknown variable type {:d}.", v);
163 }

References OGS_FATAL.

Referenced by anonymous_namespace{CreateMFront.cpp}::varInfo().