OGS
SmallDeformationLocalAssemblerMatrix.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <optional>
7#include <vector>
8
16#include "SecondaryData.h"
18
19namespace ProcessLib
20{
21namespace LIE
22{
23namespace SmallDeformation
24{
25namespace MPL = MaterialPropertyLib;
26
27template <typename ShapeFunction, int DisplacementDim>
30{
31public:
39
44
49
51 MeshLib::Element const& e,
52 std::size_t const local_matrix_size,
53 NumLib::GenericIntegrationMethod const& integration_method,
54 bool const is_axially_symmetric,
56
57 void assemble(double const /*t*/, double const /*dt*/,
58 std::vector<double> const& /*local_x*/,
59 std::vector<double> const& /*local_x_prev*/,
60 std::vector<double>& /*local_M_data*/,
61 std::vector<double>& /*local_K_data*/,
62 std::vector<double>& /*local_b_data*/) override
63 {
65 "SmallDeformationLocalAssembler: assembly without jacobian is not "
66 "implemented.");
67 }
68
69 void assembleWithJacobian(double const t, double const dt,
70 std::vector<double> const& local_x,
71 std::vector<double> const& /*local_x_prev*/,
72 std::vector<double>& local_b_data,
73 std::vector<double>& local_Jac_data) override;
74
75 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
76 double const /*t*/,
77 double const /*delta_t*/) override
78 {
79 unsigned const n_integration_points =
80 _integration_method.getNumberOfPoints();
81
82 for (unsigned ip = 0; ip < n_integration_points; ip++)
83 {
84 _ip_data[ip].pushBackState();
85 }
86 }
87
89 double const t, Eigen::VectorXd const& local_x) override;
90
91 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
92 const unsigned integration_point) const override
93 {
94 auto const& N = _secondary_data.N[integration_point];
95
96 // assumes N is stored contiguously in memory
97 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
98 }
99
100 std::vector<double> const& getIntPtSigma(
101 const double t,
102 std::vector<GlobalVector*> const& x,
103 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
104 std::vector<double>& cache) const override;
105
106 std::vector<double> const& getIntPtEpsilon(
107 const double t,
108 std::vector<GlobalVector*> const& x,
109 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
110 std::vector<double>& cache) const override;
111
112 std::vector<double> const& getIntPtFractureStress(
113 const double /*t*/,
114 std::vector<GlobalVector*> const& /*x*/,
115 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
116 std::vector<double>& cache) const override
117 {
118 cache.resize(0);
119 return cache;
120 }
121
122 std::vector<double> const& getIntPtFractureAperture(
123 const double /*t*/,
124 std::vector<GlobalVector*> const& /*x*/,
125 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
126 std::vector<double>& cache) const override
127 {
128 cache.resize(0);
129 return cache;
130 }
131
132private:
134
137 DisplacementDim>;
138 std::vector<IntegrationPointDataType,
139 Eigen::aligned_allocator<IntegrationPointDataType>>
141
146
147 std::optional<BBarMatrixType> getDilatationalBBarMatrix() const
148 {
149 if (!(_process_data.use_b_bar))
150 {
151 return std::nullopt;
152 }
153
155 DisplacementDim, ShapeFunction::NPOINTS, ShapeFunction,
157 _ip_data, this->_element, this->_integration_method,
158 this->_is_axially_symmetric);
159 }
160};
161
162} // namespace SmallDeformation
163} // namespace LIE
164} // namespace ProcessLib
165
#define OGS_FATAL(...)
Definition Error.h:19
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
MatrixType< 3, ShapeFunction::NPOINTS > BBarMatrixType
MatrixType< _number_of_dof, _number_of_dof > StiffnessMatrixType
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
void computeSecondaryVariableConcreteWithVector(double const t, Eigen::VectorXd const &local_x) override
std::vector< double > const & getIntPtFractureStress(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFractureAperture(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
SmallDeformationLocalAssemblerMatrix(SmallDeformationLocalAssemblerMatrix &&)=delete
IntegrationPointDataMatrix< ShapeMatricesType, BMatricesType, DisplacementDim > IntegrationPointDataType
std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
SmallDeformationLocalAssemblerMatrix(SmallDeformationLocalAssemblerMatrix const &)=delete
BBarMatrixType computeDilatationalBbar(std::vector< IpData, Eigen::aligned_allocator< IpData > > const &ip_data, MeshLib::Element const &element, NumLib::GenericIntegrationMethod const &integration_method, const bool is_axially_symmetric)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType