OGS
SmallDeformationLocalAssemblerMatrix-impl.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 <Eigen/Core>
7#include <vector>
8
20#include "SecondaryData.h"
23
24namespace ProcessLib
25{
26namespace LIE
27{
28namespace SmallDeformation
29{
30template <typename ShapeFunction, int DisplacementDim>
33 MeshLib::Element const& e,
34 std::size_t const /*local_matrix_size*/,
35 NumLib::GenericIntegrationMethod const& integration_method,
36 bool const is_axially_symmetric,
38 : _process_data(process_data),
39 _integration_method(integration_method),
40 _element(e),
41 _is_axially_symmetric(is_axially_symmetric)
42{
43 unsigned const n_integration_points =
44 _integration_method.getNumberOfPoints();
45
46 _ip_data.reserve(n_integration_points);
47 _secondary_data.N.resize(n_integration_points);
48
49 auto const shape_matrices =
51 DisplacementDim>(e, is_axially_symmetric,
53
55 _process_data.solid_materials, _process_data.material_ids, e.getID());
56
57 for (unsigned ip = 0; ip < n_integration_points; ip++)
58 {
59 _ip_data.emplace_back(solid_material);
60 auto& ip_data = _ip_data[ip];
61 auto const& sm = shape_matrices[ip];
62 ip_data.N_u = sm.N;
63 ip_data.dNdx_u = sm.dNdx;
64 ip_data.integration_weight =
65 _integration_method.getWeightedPoint(ip).getWeight() *
66 sm.integralMeasure * sm.detJ;
67
68 // Initialize current time step values
69 static const int kelvin_vector_size =
71 ip_data._sigma.setZero(kelvin_vector_size);
72 ip_data._eps.setZero(kelvin_vector_size);
73
74 // Previous time step values are not initialized and are set later.
75 ip_data._sigma_prev.resize(kelvin_vector_size);
76 ip_data._eps_prev.resize(kelvin_vector_size);
77
78 ip_data._C.resize(kelvin_vector_size, kelvin_vector_size);
79
80 _secondary_data.N[ip] = sm.N;
81 }
82}
83
84template <typename ShapeFunction, int DisplacementDim>
86 assembleWithJacobian(double const t, double const dt,
87 std::vector<double> const& local_x,
88 std::vector<double> const& /*local_x_prev*/,
89 std::vector<double>& local_b_data,
90 std::vector<double>& local_Jac_data)
91{
92 assert(_element.getDimension() == DisplacementDim);
93
94 auto const local_matrix_size = local_x.size();
95
97 local_Jac_data, local_matrix_size, local_matrix_size);
98
100 local_b_data, local_matrix_size);
101
102 unsigned const n_integration_points =
103 _integration_method.getNumberOfPoints();
104
105 MPL::VariableArray variables;
106 MPL::VariableArray variables_prev;
107
108 auto const B_dil_bar = getDilatationalBBarMatrix();
109
110 for (unsigned ip = 0; ip < n_integration_points; ip++)
111 {
112 auto const& w = _ip_data[ip].integration_weight;
113
114 auto const& N = _ip_data[ip].N_u;
115 auto const& dNdx = _ip_data[ip].dNdx_u;
116
117 ParameterLib::SpatialPosition const x_position{
118 std::nullopt, this->_element.getID(),
121 this->_element, N))};
122
123 auto const x_coord =
124 x_position.getCoordinates().value()[0]; // r for axisymmetry
126 DisplacementDim, ShapeFunction::NPOINTS, BBarMatrixType,
127 typename BMatricesType::BMatrixType>(dNdx, N, B_dil_bar, x_coord,
129
130 auto const& eps_prev = _ip_data[ip]._eps_prev;
131 auto const& sigma_prev = _ip_data[ip]._sigma_prev;
132
133 auto& eps = _ip_data[ip]._eps;
134 auto& sigma = _ip_data[ip]._sigma;
135 auto& state = _ip_data[ip]._material_state_variables;
136
137 eps.noalias() =
138 B * Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
139 local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
140
141 variables.mechanical_strain
143 eps);
144 variables_prev.stress
146 sigma_prev);
147 variables_prev.mechanical_strain
149 eps_prev);
150 variables_prev.temperature = _process_data.reference_temperature;
151
152 auto&& solution = _ip_data[ip]._solid_material.integrateStress(
153 variables_prev, variables, t, x_position, dt, *state);
154
155 if (!solution)
156 {
157 OGS_FATAL("Computation of local constitutive relation failed.");
158 }
159
161 std::tie(sigma, state, C) = std::move(*solution);
162
163 local_b.noalias() -= B.transpose() * sigma * w;
164 local_Jac.noalias() += B.transpose() * C * B * w;
165 }
166}
167
168template <typename ShapeFunction, int DisplacementDim>
171 double const /*t*/, Eigen::VectorXd const& /*local_x*/)
172{
173 // Compute average value per element
175 KV sigma_avg = KV::Zero();
176 auto const e_id = _element.getID();
177
178 unsigned const n_integration_points =
179 _integration_method.getNumberOfPoints();
180 for (unsigned ip = 0; ip < n_integration_points; ip++)
181 {
182 sigma_avg += _ip_data[ip]._sigma;
183 }
184 sigma_avg /= n_integration_points;
185
186 Eigen::Map<KV>(
187 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
189}
190
191template <typename ShapeFunction, int DisplacementDim>
192std::vector<double> const&
195 const double /*t*/,
196 std::vector<GlobalVector*> const& /*x*/,
197 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
198 std::vector<double>& cache) const
199{
202}
203template <typename ShapeFunction, int DisplacementDim>
204std::vector<double> const&
207 const double /*t*/,
208 std::vector<GlobalVector*> const& /*x*/,
209 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
210 std::vector<double>& cache) const
211{
214}
215
216} // namespace SmallDeformation
217} // namespace LIE
218} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
std::optional< MathLib::Point3d > const getCoordinates() const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
void computeSecondaryVariableConcreteWithVector(double const t, Eigen::VectorXd const &local_x) override
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
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
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
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
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
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BMatrixType computeBMatrixPossiblyWithBbar(DNDX_Type const &dNdx, N_Type const &N, std::optional< BBarMatrixType > const &B_dil_bar, const double radius, const bool is_axially_symmetric)
Fills a B matrix, or a B bar matrix if required.
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)