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