OGS
ProcessLib::LargeDeformation Namespace Reference

Namespaces

namespace  ConstitutiveRelations
namespace  detail

Classes

struct  DeformationGradientData
struct  GravityModel
struct  IntegrationPointData
class  LargeDeformationLocalAssembler
struct  LargeDeformationLocalAssemblerInterface
class  LargeDeformationProcess
struct  LargeDeformationProcessData
class  MaterialStateData
struct  MediaData
struct  SecondaryData
struct  SolidDensityModel

Typedefs

template<int DisplacementDim>
using KelvinVector = KV::KelvinVectorType<DisplacementDim>
template<int DisplacementDim>
using KelvinMatrix = KV::KelvinMatrixType<DisplacementDim>
template<int DisplacementDim>
using GlobalDimVector = Eigen::Vector<double, DisplacementDim>
template<int DisplacementDim>
using GlobalDimMatrix
using Temperature = BaseLib::StrongType<double, struct TemperatureTag>
template<int DisplacementDim>
using VolumetricBodyForce
using SolidDensity = BaseLib::StrongType<double, struct SolidDensityTag>

Functions

template<int D>
constexpr GlobalDimVector< D > DVnan ()
 Used to set a D dimensional vector to all not-a-number.
template<int D>
constexpr GlobalDimMatrix< D > DMnan ()
 Used to set a D x D matrix to all not-a-number.
template<typename T>
constexpr bool areEvalArgumentTypesUnique ()
template<typename Model>
constexpr void assertEvalArgsUnique (Model const &)
template<int DisplacementDim>
std::unique_ptr< ProcesscreateLargeDeformationProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< ProcesscreateLargeDeformationProcess< 2 > (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< ProcesscreateLargeDeformationProcess< 3 > (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template<int GlobalDim, template< typename, int > class LocalAssemblerImplementation, typename LocalAssemblerInterface, IntegrationMethodProviderOrIntegrationOrder ProviderOrOrder, typename... ExtraCtorArgs>
void createLocalAssemblers (std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
template<int DisplacementDim, typename ShapeMatricesType>
Eigen::Matrix< double, MPL::tensorSize(DisplacementDim), MPL::tensorSize(DisplacementDim)> computeSigmaGeom (Eigen::Matrix3d const &sigma_tensor)

Typedef Documentation

◆ GlobalDimMatrix

template<int DisplacementDim>
using ProcessLib::LargeDeformation::GlobalDimMatrix
Initial value:
Eigen::Matrix<double, DisplacementDim, DisplacementDim, Eigen::RowMajor>

Definition at line 36 of file LargeDeformation/ConstitutiveRelations/Base.h.

◆ GlobalDimVector

template<int DisplacementDim>
using ProcessLib::LargeDeformation::GlobalDimVector = Eigen::Vector<double, DisplacementDim>

◆ KelvinMatrix

template<int DisplacementDim>
using ProcessLib::LargeDeformation::KelvinMatrix = KV::KelvinMatrixType<DisplacementDim>

◆ KelvinVector

template<int DisplacementDim>
using ProcessLib::LargeDeformation::KelvinVector = KV::KelvinVectorType<DisplacementDim>

◆ SolidDensity

◆ Temperature

using ProcessLib::LargeDeformation::Temperature = BaseLib::StrongType<double, struct TemperatureTag>

◆ VolumetricBodyForce

Function Documentation

◆ areEvalArgumentTypesUnique()

template<typename T>
bool ProcessLib::LargeDeformation::areEvalArgumentTypesUnique ( )
constexpr

Checks whether the argument types of the eval() method of the given type T are unique.

Argument types differing only in constness, reference or volatility are considered equal.

Definition at line 41 of file LargeDeformation/ConstitutiveRelations/Invoke.h.

42{
44}
constexpr bool areEvalArgumentTypesUnique(Result(Class::*)(Args...))

References ProcessLib::LargeDeformation::detail::areEvalArgumentTypesUnique().

Referenced by assertEvalArgsUnique().

◆ assertEvalArgsUnique()

template<typename Model>
void ProcessLib::LargeDeformation::assertEvalArgsUnique ( Model const & )
constexpr

Statically asserts that the argument types of the passed Model's eval() method are unique.

Definition at line 49 of file LargeDeformation/ConstitutiveRelations/Invoke.h.

References areEvalArgumentTypesUnique().

Referenced by ProcessLib::LargeDeformation::ConstitutiveRelations::ConstitutiveSetting< DisplacementDim >::eval().

◆ computeSigmaGeom()

template<int DisplacementDim, typename ShapeMatricesType>
Eigen::Matrix< double, MPL::tensorSize(DisplacementDim), MPL::tensorSize(DisplacementDim)> ProcessLib::LargeDeformation::computeSigmaGeom ( Eigen::Matrix3d const & sigma_tensor)

Definition at line 56 of file LargeDeformationFEM.h.

57{
58 static constexpr auto& sigma_geom_op = MathLib::eigenBlockMatrixView<
59 DisplacementDim,
60 Eigen::Matrix<double, DisplacementDim, DisplacementDim>>;
61
62 using SigmaGeom = Eigen::Matrix<double, MPL::tensorSize(DisplacementDim),
63 MPL::tensorSize(DisplacementDim)>;
64 if constexpr (DisplacementDim == 2)
65 {
66 SigmaGeom sigma_geom = SigmaGeom::Zero(5, 5);
67 sigma_geom.template block<4, 4>(0, 0) =
68 sigma_geom_op(sigma_tensor.template block<2, 2>(0, 0).eval());
69 sigma_geom(4, 4) = sigma_tensor(2, 2);
70
71 return sigma_geom;
72 }
73
74 if constexpr (DisplacementDim == 3)
75 {
76 return sigma_geom_op(sigma_tensor);
77 }
78}
constexpr int tensorSize(int dim)
See Tensor type for details.
Definition Tensor.h:19
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)

References MathLib::eigenBlockMatrixView(), and MaterialPropertyLib::tensorSize().

Referenced by ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian().

◆ createLargeDeformationProcess()

template<int DisplacementDim>
std::unique_ptr< Process > ProcessLib::LargeDeformation::createLargeDeformationProcess ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< ProcessVariable > const & variables,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
unsigned const integration_order,
BaseLib::ConfigTree const & config,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )
Input File Parameter
prj__processes__process__type

Process Variables

Input File Parameter
prj__processes__process__LARGE_DEFORMATION__process_variables

Primary process variables as they appear in the global component vector:

Input File Parameter
prj__processes__process__LARGE_DEFORMATION__process_variables__process_variable

Process Parameters

Input File Parameter
prj__processes__process__LARGE_DEFORMATION__specific_body_force
Input File Parameter
prj__processes__process__LARGE_DEFORMATION__reference_temperature
Input File Parameter
prj__processes__process__LARGE_DEFORMATION__initial_stress
Input File Parameter
prj__processes__process__LARGE_DEFORMATION__f_bar

Definition at line 28 of file CreateLargeDeformationProcess.cpp.

39{
41 config.checkConfigParameter("type", "LARGE_DEFORMATION");
42 DBUG("Create LargeDeformationProcess.");
43
45
47 auto const pv_config = config.getConfigSubtree("process_variables");
48
50 auto per_process_variables = findProcessVariables(
51 variables, pv_config,
52 {
53 "process_variable"});
54
55 DBUG("Associate displacement with process variable '{:s}'.",
56 per_process_variables.back().get().getName());
57
58 if (per_process_variables.back().get().getNumberOfGlobalComponents() !=
59 DisplacementDim)
60 {
62 "Number of components of the process variable '{:s}' is different "
63 "from the displacement dimension: got {:d}, expected {:d}",
64 per_process_variables.back().get().getName(),
65 per_process_variables.back().get().getNumberOfGlobalComponents(),
66 DisplacementDim);
67 }
68 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
69 process_variables;
70 process_variables.push_back(std::move(per_process_variables));
71
73 auto solid_constitutive_relations =
76 parameters, local_coordinate_system, materialIDs(mesh), config);
77
78 // Specific body force
79 Eigen::Matrix<double, DisplacementDim, 1> specific_body_force;
80 {
81 std::vector<double> const b =
83 config.getConfigParameter<std::vector<double>>(
84 "specific_body_force");
85 if (b.size() != DisplacementDim)
86 {
88 "The size of the specific body force vector does not match the "
89 "displacement dimension. Vector size is {:d}, displacement "
90 "dimension is {:d}",
91 b.size(), DisplacementDim);
92 }
93
94 std::copy_n(b.data(), b.size(), specific_body_force.data());
95 }
96
97 auto media_map =
99
100 // Reference temperature
102 double>(
104 config, "reference_temperature", parameters, 1, &mesh);
105 if (reference_temperature)
106 {
107 DBUG("Use '{:s}' as reference temperature parameter.",
108 (*reference_temperature).name);
109 }
110
111 // Initial stress conditions
112 auto const initial_stress = ParameterLib::findOptionalTagParameter<double>(
114 config, "initial_stress", parameters,
115 // Symmetric tensor size, 4 or 6, not a Kelvin vector.
117 &mesh);
118
119 std::string const f_bar_info =
121 config.getConfigParameter<std::string>("f_bar", "");
122
123 if (mesh.isAxiallySymmetric() && f_bar_info != "")
124 {
125 OGS_FATAL(
126 "The F-bar method is not available for axisymmetric problems.");
127 }
128
130 materialIDs(mesh),
131 std::move(media_map),
132 std::move(solid_constitutive_relations),
133 initial_stress,
134 specific_body_force,
137
138 SecondaryVariableCollection secondary_variables;
139
140 ProcessLib::createSecondaryVariables(config, secondary_variables);
141
142 return std::make_unique<LargeDeformationProcess<DisplacementDim>>(
143 std::move(name), mesh, std::move(jacobian_assembler), parameters,
144 integration_order, std::move(process_variables),
145 std::move(process_data), std::move(secondary_variables));
146}
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
bool isAxiallySymmetric() const
Definition Mesh.h:139
MaterialSpatialDistributionMap createMaterialSpatialDistributionMap(std::map< int, std::shared_ptr< Medium > > const &media, MeshLib::Mesh const &mesh)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition Mesh.cpp:269
Parameter< ParameterDataType > * findOptionalTagParameter(BaseLib::ConfigTree const &process_config, std::string const &tag, std::vector< std::unique_ptr< ParameterBase > > const &parameters, int const num_components, MeshLib::Mesh const *const mesh=nullptr)
BarDetFType convertStringToDetFBarType(std::string_view const bar_det_f_type_name)
std::vector< std::reference_wrapper< ProcessVariable > > findProcessVariables(std::vector< ProcessVariable > const &variables, BaseLib::ConfigTree const &pv_config, std::initializer_list< std::string > tags)
void createSecondaryVariables(BaseLib::ConfigTree const &config, SecondaryVariableCollection &secondary_variables)
static std::map< int, std::shared_ptr< SolidConstitutiveRelation< DisplacementDim > > > createSolidConstitutiveRelations(std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, MeshLib::PropertyVector< int > const *const material_ids, BaseLib::ConfigTree const &config)

References BaseLib::ConfigTree::checkConfigParameter(), ProcessLib::NonLinearFbar::convertStringToDetFBarType(), MaterialPropertyLib::createMaterialSpatialDistributionMap(), ProcessLib::createSecondaryVariables(), ProcessLib::LargeDeformation::ConstitutiveRelations::CreateConstitutiveSetting< DisplacementDim >::createSolidConstitutiveRelations(), DBUG(), ParameterLib::findOptionalTagParameter(), ProcessLib::findProcessVariables(), BaseLib::ConfigTree::getConfigParameter(), BaseLib::ConfigTree::getConfigSubtree(), MeshLib::Mesh::isAxiallySymmetric(), MathLib::KelvinVector::kelvin_vector_dimensions(), and OGS_FATAL.

◆ createLargeDeformationProcess< 2 >()

template std::unique_ptr< Process > ProcessLib::LargeDeformation::createLargeDeformationProcess< 2 > ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< ProcessVariable > const & variables,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
unsigned const integration_order,
BaseLib::ConfigTree const & config,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )

◆ createLargeDeformationProcess< 3 >()

template std::unique_ptr< Process > ProcessLib::LargeDeformation::createLargeDeformationProcess< 3 > ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< ProcessVariable > const & variables,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
std::optional< ParameterLib::CoordinateSystem > const & local_coordinate_system,
unsigned const integration_order,
BaseLib::ConfigTree const & config,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )

◆ createLocalAssemblers()

template<int GlobalDim, template< typename, int > class LocalAssemblerImplementation, typename LocalAssemblerInterface, IntegrationMethodProviderOrIntegrationOrder ProviderOrOrder, typename... ExtraCtorArgs>
void ProcessLib::LargeDeformation::createLocalAssemblers ( std::vector< MeshLib::Element * > const & mesh_elements,
NumLib::LocalToGlobalIndexMap const & dof_table,
std::vector< std::unique_ptr< LocalAssemblerInterface > > & local_assemblers,
ProviderOrOrder const & provider_or_order,
ExtraCtorArgs &&... extra_ctor_args )

Creates local assemblers for each element of the given mesh.

Template Parameters
LocalAssemblerImplementationthe individual local assembler type
LocalAssemblerInterfacethe general local assembler interface
ExtraCtorArgstypes of additional constructor arguments. Those arguments will be passed to the constructor of LocalAssemblerImplementation.

The first two template parameters cannot be deduced from the arguments. Therefore they always have to be provided manually.

Definition at line 76 of file LargeDeformation/CreateLocalAssemblers.h.

82{
83 DBUG("Create local assemblers.");
84
86 dof_table, mesh_elements, local_assemblers, provider_or_order,
87 std::forward<ExtraCtorArgs>(extra_ctor_args)...);
88}
void createLocalAssemblers(NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< MeshLib::Element * > const &mesh_elements, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)

References createLocalAssemblers(), ProcessLib::LargeDeformation::detail::createLocalAssemblers(), and DBUG().

Referenced by createLocalAssemblers(), and ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >::initializeConcreteProcess().

◆ DMnan()

template<int D>
GlobalDimMatrix< D > ProcessLib::LargeDeformation::DMnan ( )
constexpr

Used to set a D x D matrix to all not-a-number.

Definition at line 48 of file LargeDeformation/ConstitutiveRelations/Base.h.

49{
51}
static constexpr double nan
Convenience alias for not a number.
Eigen::Matrix< double, DisplacementDim, DisplacementDim, Eigen::RowMajor > GlobalDimMatrix

References ProcessLib::ConstitutiveRelations::nan.

◆ DVnan()

template<int D>
GlobalDimVector< D > ProcessLib::LargeDeformation::DVnan ( )
constexpr

Used to set a D dimensional vector to all not-a-number.

Definition at line 41 of file LargeDeformation/ConstitutiveRelations/Base.h.

42{
44}
Eigen::Vector< double, DisplacementDim > GlobalDimVector

References ProcessLib::ConstitutiveRelations::nan.