OGS
SmallDeformationLocalAssemblerMatrixNearFracture-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 <range/v3/range/conversion.hpp>
8#include <range/v3/view/transform.hpp>
9#include <vector>
10
15#include "MathLib/Point3d.h"
16#include "MeshLib/Node.h"
25#include "SecondaryData.h"
28
29namespace ProcessLib
30{
31namespace LIE
32{
33namespace SmallDeformation
34{
35template <typename ShapeFunction,
36
37 int DisplacementDim>
39
40 DisplacementDim>::
41 SmallDeformationLocalAssemblerMatrixNearFracture(
42 MeshLib::Element const& e,
43 std::size_t const n_variables,
44 std::size_t const /*local_matrix_size*/,
45 std::vector<unsigned> const& dofIndex_to_localIndex,
46 NumLib::GenericIntegrationMethod const& integration_method,
47 bool const is_axially_symmetric,
50 n_variables * ShapeFunction::NPOINTS * DisplacementDim,
51 dofIndex_to_localIndex),
52 _process_data(process_data),
53 _integration_method(integration_method),
54 _element(e),
55 _is_axially_symmetric(is_axially_symmetric)
56{
57 std::vector<ShapeMatrices, Eigen::aligned_allocator<
59 shape_matrices =
61 DisplacementDim>(e, is_axially_symmetric,
63
64 unsigned const n_integration_points =
65 _integration_method.getNumberOfPoints();
66
67 _ip_data.reserve(n_integration_points);
68 _secondary_data.N.resize(n_integration_points);
69
71 _process_data.solid_materials, _process_data.material_ids, e.getID());
72
73 for (unsigned ip = 0; ip < n_integration_points; ip++)
74 {
75 _ip_data.emplace_back(solid_material);
76 auto& ip_data = _ip_data[ip];
77 auto const& sm = shape_matrices[ip];
78 ip_data.N_u = sm.N;
79 ip_data.dNdx_u = sm.dNdx;
80 ip_data.integration_weight =
81 _integration_method.getWeightedPoint(ip).getWeight() *
82 sm.integralMeasure * sm.detJ;
83
84 // Initialize current time step values
85 static const int kelvin_vector_size =
87 ip_data._sigma.setZero(kelvin_vector_size);
88 ip_data._eps.setZero(kelvin_vector_size);
89
90 // Previous time step values are not initialized and are set later.
91 ip_data._sigma_prev.resize(kelvin_vector_size);
92 ip_data._eps_prev.resize(kelvin_vector_size);
93
94 ip_data._C.resize(kelvin_vector_size, kelvin_vector_size);
95
96 _secondary_data.N[ip] = sm.N;
97 }
98
99 for (auto fid : process_data.vec_ele_connected_fractureIDs[e.getID()])
100 {
101 _fracID_to_local.insert({fid, _fracture_props.size()});
102 _fracture_props.push_back(&_process_data.fracture_properties[fid]);
103 }
104
106 ranges::views::transform(
107 [&](auto const jid)
108 { return &_process_data.junction_properties[jid]; }) |
109 ranges::to<std::vector>;
110}
111
112template <typename ShapeFunction, int DisplacementDim>
115 DisplacementDim>::assembleWithJacobian(double const t, double const dt,
116 Eigen::VectorXd const& local_u,
117 Eigen::VectorXd& local_b,
118 Eigen::MatrixXd& local_J)
119{
120 assert(_element.getDimension() == DisplacementDim);
121
122 auto const N_DOF_PER_VAR = ShapeFunction::NPOINTS * DisplacementDim;
123 auto const n_fractures = _fracture_props.size();
124 auto const n_junctions = _junction_props.size();
125 auto const n_enrich_var = n_fractures + n_junctions;
126
127 using BlockVectorType =
128 typename Eigen::VectorXd::FixedSegmentReturnType<N_DOF_PER_VAR>::Type;
129 using BlockMatrixType =
130 Eigen::Block<Eigen::MatrixXd, N_DOF_PER_VAR, N_DOF_PER_VAR>;
131
132 //--------------------------------------------------------------------------------------
133 // prepare sub vectors, matrices for regular displacement (u) and
134 // displacement jumps (g)
135 //
136 // example with two fractures with one intersection:
137 // |b(u)|
138 // b = |b(g1)|
139 // |b(g2)|
140 // |b(j1)|
141 //
142 // |J(u,u) J(u,g1) J(u,g2) J(u,j1) |
143 // J = |J(g1,u) J(g1,g1) J(g1,g2) J(g1,j1)|
144 // |J(g2,u) J(g2,g1) J(g2,g2) J(g2,j1)|
145 // |J(j1,u) J(j1,g1) J(j1,g2) J(j1,j1)|
146 //--------------------------------------------------------------------------------------
147 auto local_b_u = local_b.segment<N_DOF_PER_VAR>(0);
148 std::vector<BlockVectorType> vec_local_b_g;
149 for (unsigned i = 0; i < n_enrich_var; i++)
150 {
151 vec_local_b_g.push_back(
152 local_b.segment<N_DOF_PER_VAR>(N_DOF_PER_VAR * (i + 1)));
153 }
154
155 auto local_J_uu = local_J.block<N_DOF_PER_VAR, N_DOF_PER_VAR>(0, 0);
156 std::vector<BlockMatrixType> vec_local_J_ug;
157 std::vector<BlockMatrixType> vec_local_J_gu;
158 std::vector<std::vector<BlockMatrixType>> vec_local_J_gg(n_enrich_var);
159 for (unsigned i = 0; i < n_enrich_var; i++)
160 {
161 auto sub_ug = local_J.block<N_DOF_PER_VAR, N_DOF_PER_VAR>(
162 0, N_DOF_PER_VAR * (i + 1));
163 vec_local_J_ug.push_back(sub_ug);
164
165 auto sub_gu = local_J.block<N_DOF_PER_VAR, N_DOF_PER_VAR>(
166 N_DOF_PER_VAR * (i + 1), 0);
167 vec_local_J_gu.push_back(sub_gu);
168
169 for (unsigned j = 0; j < n_enrich_var; j++)
170 {
171 auto sub_gg = local_J.block<N_DOF_PER_VAR, N_DOF_PER_VAR>(
172 N_DOF_PER_VAR * (i + 1), N_DOF_PER_VAR * (j + 1));
173 vec_local_J_gg[i].push_back(sub_gg);
174 }
175 }
176
177 auto const nodal_u = local_u.segment<N_DOF_PER_VAR>(0);
178 std::vector<BlockVectorType> vec_nodal_g;
179 for (unsigned i = 0; i < n_enrich_var; i++)
180 {
181 auto sub = const_cast<Eigen::VectorXd&>(local_u).segment<N_DOF_PER_VAR>(
182 N_DOF_PER_VAR * (i + 1));
183 vec_nodal_g.push_back(sub);
184 }
185
186 //------------------------------------------------
187 // integration
188 //------------------------------------------------
189 unsigned const n_integration_points =
190 _integration_method.getNumberOfPoints();
191
192 MPL::VariableArray variables;
193 MPL::VariableArray variables_prev;
194
195 auto const B_dil_bar = getDilatationalBBarMatrix();
196
197 for (unsigned ip = 0; ip < n_integration_points; ip++)
198 {
199 auto& ip_data = _ip_data[ip];
200 auto const& w = _ip_data[ip].integration_weight;
201
202 auto const& N = ip_data.N_u;
203 auto const& dNdx = ip_data.dNdx_u;
204
205 // levelset functions
206 Eigen::Vector3d const ip_physical_coords(
208 std::vector<double> const levelsets(
210 _fracID_to_local, ip_physical_coords));
211
212 // u = u^hat + sum_i(enrich^br_i(x) * [u]_i) + sum_i(enrich^junc_i(x) *
213 // [u]_i)
214 NodalDisplacementVectorType nodal_total_u = nodal_u;
215 for (unsigned i = 0; i < n_enrich_var; i++)
216 {
217 nodal_total_u += levelsets[i] * vec_nodal_g[i];
218 }
219
220 ParameterLib::SpatialPosition const x_position{
221 std::nullopt, this->_element.getID(),
224 this->_element, N))};
225
226 auto const x_coord =
227 x_position.getCoordinates().value()[0]; // r for axisymmetry
229 DisplacementDim, ShapeFunction::NPOINTS, BBarMatrixType,
230 typename BMatricesType::BMatrixType>(dNdx, N, B_dil_bar, x_coord,
232
233 // strain, stress
234 auto const& eps_prev = ip_data._eps_prev;
235 auto const& sigma_prev = ip_data._sigma_prev;
236
237 auto& eps = ip_data._eps;
238 auto& sigma = ip_data._sigma;
239 auto& state = ip_data._material_state_variables;
240
241 eps.noalias() = B * nodal_total_u;
242
243 variables.mechanical_strain
245 eps);
246
247 variables_prev.stress
249 sigma_prev);
250 variables_prev.mechanical_strain
252 eps_prev);
253 variables_prev.temperature = _process_data.reference_temperature;
254
255 auto&& solution = _ip_data[ip]._solid_material.integrateStress(
256 variables_prev, variables, t, x_position, dt, *state);
257
258 if (!solution)
259 {
260 OGS_FATAL("Computation of local constitutive relation failed.");
261 }
262
264 std::tie(sigma, state, C) = std::move(*solution);
265
266 // r_u = B^T * Sigma = B^T * C * B * (u+phi*[u])
267 // r_[u] = (phi*B)^T * Sigma = (phi*B)^T * C * B * (u+phi*[u])
268 local_b_u.noalias() -= B.transpose() * sigma * w;
269 for (unsigned i = 0; i < n_enrich_var; i++)
270 {
271 vec_local_b_g[i].noalias() -=
272 levelsets[i] * B.transpose() * sigma * w;
273 }
274
275 // J_uu += B^T * C * B
276 local_J_uu.noalias() += B.transpose() * C * B * w;
277
278 for (unsigned i = 0; i < n_enrich_var; i++)
279 {
280 // J_u[u] += B^T * C * (levelset * B)
281 vec_local_J_ug[i].noalias() +=
282 B.transpose() * C * (levelsets[i] * B) * w;
283
284 // J_[u]u += (levelset * B)^T * C * B
285 vec_local_J_gu[i].noalias() +=
286 (levelsets[i] * B.transpose()) * C * B * w;
287
288 for (unsigned j = 0; j < n_enrich_var; j++)
289 {
290 // J_[u][u] += (levelset * B)^T * C * (levelset * B)
291 vec_local_J_gg[i][j].noalias() +=
292 (levelsets[i] * B.transpose()) * C * (levelsets[j] * B) * w;
293 }
294 }
295 }
296}
297
298template <typename ShapeFunction,
299
300 int DisplacementDim>
302
303 DisplacementDim>::
304 computeSecondaryVariableConcreteWithVector(
305 double const /*t*/, Eigen::VectorXd const& /*local_x*/)
306{
307 // Compute average value per element
309 KV sigma_avg = KV::Zero();
310 auto const e_id = _element.getID();
311
312 unsigned const n_integration_points =
313 _integration_method.getNumberOfPoints();
314 for (unsigned ip = 0; ip < n_integration_points; ip++)
315 {
316 sigma_avg += _ip_data[ip]._sigma;
317 }
318 sigma_avg /= n_integration_points;
319
320 Eigen::Map<KV>(
321 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
323}
324
325template <typename ShapeFunction, int DisplacementDim>
326std::vector<double> const& SmallDeformationLocalAssemblerMatrixNearFracture<
327 ShapeFunction, DisplacementDim>::
328 getIntPtSigma(
329 const double /*t*/,
330 std::vector<GlobalVector*> const& /*x*/,
331 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
332 std::vector<double>& cache) const
333{
336}
337template <typename ShapeFunction, int DisplacementDim>
338std::vector<double> const& SmallDeformationLocalAssemblerMatrixNearFracture<
339 ShapeFunction, DisplacementDim>::
340 getIntPtEpsilon(
341 const double /*t*/,
342 std::vector<GlobalVector*> const& /*x*/,
343 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
344 std::vector<double>& cache) const
345{
348}
349
350} // namespace SmallDeformation
351} // namespace LIE
352} // 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 assembleWithJacobian(double const t, double const dt, Eigen::VectorXd const &local_u, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
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
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)
Eigen::Vector3d computePhysicalCoordinates(MeshLib::Element const &e, Eigen::MatrixBase< Derived > const &shape)
std::vector< double > uGlobalEnrichments(std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)
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)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices