OGS
HydroMechanicsLocalAssemblerMatrix-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 <limits>
7
21
22namespace ProcessLib
23{
24namespace LIE
25{
26namespace HydroMechanics
27{
28namespace MPL = MaterialPropertyLib;
29
30template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
31 int DisplacementDim>
32HydroMechanicsLocalAssemblerMatrix<ShapeFunctionDisplacement,
33 ShapeFunctionPressure, DisplacementDim>::
34 HydroMechanicsLocalAssemblerMatrix(
35 MeshLib::Element const& e,
36 std::size_t const n_variables,
37 std::size_t const /*local_matrix_size*/,
38 std::vector<unsigned> const& dofIndex_to_localIndex,
39 NumLib::GenericIntegrationMethod const& integration_method,
40 bool const is_axially_symmetric,
43 e, is_axially_symmetric, integration_method,
44 (n_variables - 1) * ShapeFunctionDisplacement::NPOINTS *
45 DisplacementDim +
46 ShapeFunctionPressure::NPOINTS,
47 dofIndex_to_localIndex),
48 _process_data(process_data)
49{
50 unsigned const n_integration_points =
51 integration_method.getNumberOfPoints();
52
53 _ip_data.reserve(n_integration_points);
54 _secondary_data.N.resize(n_integration_points);
55
56 auto const shape_matrices_u =
57 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
59 DisplacementDim>(e, is_axially_symmetric,
60 integration_method);
61
62 auto const shape_matrices_p =
63 NumLib::initShapeMatrices<ShapeFunctionPressure,
64 ShapeMatricesTypePressure, DisplacementDim>(
65 e, is_axially_symmetric, integration_method);
66
68 _process_data.solid_materials, _process_data.material_ids, e.getID());
69
71 x_position.setElementID(e.getID());
72 for (unsigned ip = 0; ip < n_integration_points; ip++)
73 {
74 _ip_data.emplace_back(solid_material);
75 auto& ip_data = _ip_data[ip];
76 auto const& sm_u = shape_matrices_u[ip];
77 auto const& sm_p = shape_matrices_p[ip];
78
79 ip_data.integration_weight =
80 sm_u.detJ * sm_u.integralMeasure *
81 integration_method.getWeightedPoint(ip).getWeight();
82 ip_data.darcy_velocity.setZero();
83
84 ip_data.N_u = sm_u.N;
85 ip_data.dNdx_u = sm_u.dNdx;
86 ip_data.H_u.setZero(DisplacementDim, displacement_size);
87 for (int i = 0; i < DisplacementDim; ++i)
88 {
89 ip_data.H_u
90 .template block<1, displacement_size / DisplacementDim>(
91 i, i * displacement_size / DisplacementDim)
92 .noalias() = ip_data.N_u;
93 }
94
95 ip_data.N_p = sm_p.N;
96 ip_data.dNdx_p = sm_p.dNdx;
97
98 _secondary_data.N[ip] = sm_u.N;
99
100 ip_data.sigma_eff.setZero(kelvin_vector_size);
101 ip_data.eps.setZero(kelvin_vector_size);
102
103 ip_data.sigma_eff_prev.resize(kelvin_vector_size);
104 ip_data.eps_prev.resize(kelvin_vector_size);
105 ip_data.C.resize(kelvin_vector_size, kelvin_vector_size);
106
107 auto const initial_effective_stress =
108 _process_data.initial_effective_stress(0, x_position);
109 for (unsigned i = 0; i < kelvin_vector_size; i++)
110 {
111 ip_data.sigma_eff[i] = initial_effective_stress[i];
112 ip_data.sigma_eff_prev[i] = initial_effective_stress[i];
113 }
114 }
115}
116
117template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
118 int DisplacementDim>
120 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
121 assembleWithJacobianConcrete(double const t, double const dt,
122 Eigen::VectorXd const& local_x,
123 Eigen::VectorXd const& local_x_prev,
124 Eigen::VectorXd& local_rhs,
125 Eigen::MatrixXd& local_Jac)
126{
127 auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
129 auto const p_prev = local_x_prev.segment(pressure_index, pressure_size);
130
131 if (_process_data.deactivate_matrix_in_flow)
132 {
134 }
135
136 auto u = local_x.segment(displacement_index, displacement_size);
137 auto u_prev = local_x_prev.segment(displacement_index, displacement_size);
138
139 auto rhs_p = local_rhs.template segment<pressure_size>(pressure_index);
140 auto rhs_u =
141 local_rhs.template segment<displacement_size>(displacement_index);
142
143 auto J_pp = local_Jac.template block<pressure_size, pressure_size>(
145 auto J_pu = local_Jac.template block<pressure_size, displacement_size>(
147 auto J_uu = local_Jac.template block<displacement_size, displacement_size>(
149 auto J_up = local_Jac.template block<displacement_size, pressure_size>(
151
152 assembleBlockMatricesWithJacobian(t, dt, p, p_prev, u, u_prev, rhs_p, rhs_u,
153 J_pp, J_pu, J_uu, J_up);
154}
155
156template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
157 int DisplacementDim>
159 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
160 assembleBlockMatricesWithJacobian(
161 double const t, double const dt,
162 Eigen::Ref<const Eigen::VectorXd> const& p,
163 Eigen::Ref<const Eigen::VectorXd> const& p_prev,
164 Eigen::Ref<const Eigen::VectorXd> const& u,
165 Eigen::Ref<const Eigen::VectorXd> const& u_prev,
166 Eigen::Ref<Eigen::VectorXd> rhs_p, Eigen::Ref<Eigen::VectorXd> rhs_u,
167 Eigen::Ref<Eigen::MatrixXd> J_pp, Eigen::Ref<Eigen::MatrixXd> J_pu,
168 Eigen::Ref<Eigen::MatrixXd> J_uu, Eigen::Ref<Eigen::MatrixXd> J_up)
169{
170 assert(this->_element.getDimension() == DisplacementDim);
171
173 ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
175
177 ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
179
180 typename ShapeMatricesTypeDisplacement::template MatrixType<
182 Kup = ShapeMatricesTypeDisplacement::template MatrixType<
185
186 auto const& gravity_vec = _process_data.specific_body_force;
187
188 MPL::VariableArray variables;
189 MPL::VariableArray variables_prev;
191 x_position.setElementID(_element.getID());
192
193 auto const& medium = _process_data.media_map.getMedium(_element.getID());
194 auto const& liquid_phase =
196 auto const& solid_phase =
198
199 auto const T_ref =
201 .template value<double>(variables, x_position, t, dt);
202 variables.temperature = T_ref;
203 variables_prev.temperature = T_ref;
204
205 bool const has_storage_property =
206 medium->hasProperty(MPL::PropertyType::storage);
207
208 auto const B_dil_bar = getDilatationalBBarMatrix();
209
210 unsigned const n_integration_points = _ip_data.size();
211 for (unsigned ip = 0; ip < n_integration_points; ip++)
212 {
213 auto& ip_data = _ip_data[ip];
214 auto const& ip_w = ip_data.integration_weight;
215 auto const& N_u = ip_data.N_u;
216 auto const& dNdx_u = ip_data.dNdx_u;
217 auto const& N_p = ip_data.N_p;
218 auto const& dNdx_p = ip_data.dNdx_p;
219 auto const& H_u = ip_data.H_u;
220
221 variables.liquid_phase_pressure = N_p.dot(p);
222
223 x_position = {
224 std::nullopt, _element.getID(),
226 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
228 _element, N_u))};
229
230 auto const x_coord = x_position.getCoordinates().value()[0];
231
232 auto const B =
234 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
236 dNdx_u, N_u, B_dil_bar, x_coord, this->_is_axially_symmetric)
237 .eval();
238
239 auto const& eps_prev = ip_data.eps_prev;
240 auto const& sigma_eff_prev = ip_data.sigma_eff_prev;
241 auto& sigma_eff = ip_data.sigma_eff;
242
243 auto& eps = ip_data.eps;
244 auto& state = ip_data.material_state_variables;
245
246 auto const rho_sr =
247 solid_phase.property(MPL::PropertyType::density)
248 .template value<double>(variables, x_position, t, dt);
249 auto const rho_fr =
250 liquid_phase.property(MPL::PropertyType::density)
251 .template value<double>(variables, x_position, t, dt);
252
253 auto const alpha =
255 .template value<double>(variables, x_position, t, dt);
256 auto const porosity =
257 medium->property(MPL::PropertyType::porosity)
258 .template value<double>(variables, x_position, t, dt);
259
260 double const rho = rho_sr * (1 - porosity) + porosity * rho_fr;
261 auto const& identity2 =
263
264 eps.noalias() = B * u;
265
266 variables.mechanical_strain
268 eps);
269
270 variables_prev.stress
272 sigma_eff_prev);
273 variables_prev.mechanical_strain
275 eps_prev);
276 variables_prev.temperature = T_ref;
277
278 auto&& solution = _ip_data[ip].solid_material.integrateStress(
279 variables_prev, variables, t, x_position, dt, *state);
280
281 if (!solution)
282 {
283 OGS_FATAL("Computation of local constitutive relation failed.");
284 }
285
287 std::tie(sigma_eff, state, C) = std::move(*solution);
288
289 J_uu.noalias() += (B.transpose() * C * B) * ip_w;
290
291 rhs_u.noalias() -= B.transpose() * sigma_eff * ip_w;
292 rhs_u.noalias() -= -H_u.transpose() * rho * gravity_vec * ip_w;
293
294 //
295 // pressure equation, pressure part and displacement equation, pressure
296 // part
297 //
298 if (!_process_data.deactivate_matrix_in_flow) // Only for hydraulically
299 // active matrix
300 {
301 Kup.noalias() += B.transpose() * alpha * identity2 * N_p * ip_w;
302
303 variables.density = rho_fr;
304 auto const mu =
305 liquid_phase.property(MPL::PropertyType::viscosity)
306 .template value<double>(variables, x_position, t, dt);
307
308 auto const k_over_mu =
310 medium->property(MPL::PropertyType::permeability)
311 .value(variables, x_position, t, dt)) /
312 mu)
313 .eval();
314
315 double const S =
316 has_storage_property
317 ? medium->property(MPL::PropertyType::storage)
318 .template value<double>(variables, x_position, t, dt)
319 : porosity *
320 (liquid_phase.property(MPL::PropertyType::density)
321 .template dValue<double>(
322 variables,
324 x_position, t, dt) /
325 rho_fr) +
326 (alpha - porosity) * (1.0 - alpha) /
327 _ip_data[ip].solid_material.getBulkModulus(
328 t, x_position, &C);
329
330 auto q = ip_data.darcy_velocity.head(DisplacementDim);
331 q.noalias() = -k_over_mu * (dNdx_p * p + rho_fr * gravity_vec);
332
333 laplace_p.noalias() +=
334 dNdx_p.transpose() * k_over_mu * dNdx_p * ip_w;
335 storage_p.noalias() += N_p.transpose() * S * N_p * ip_w;
336
337 rhs_p.noalias() +=
338 dNdx_p.transpose() * rho_fr * k_over_mu * gravity_vec * ip_w;
339 }
340 }
341
342 // displacement equation, pressure part
343 J_up.noalias() -= Kup;
344
345 // pressure equation, pressure part.
346 J_pp.noalias() += laplace_p + storage_p / dt;
347
348 // pressure equation, displacement part.
349 J_pu.noalias() += Kup.transpose() / dt;
350
351 // pressure equation
352 rhs_p.noalias() -= laplace_p * p + storage_p * (p - p_prev) / dt +
353 Kup.transpose() * (u - u_prev) / dt;
354
355 // displacement equation
356 rhs_u.noalias() -= -Kup * p;
357}
358
359template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
360 int DisplacementDim>
362 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
363 postTimestepConcreteWithVector(double const t, double const dt,
364 Eigen::VectorXd const& local_x)
365{
366 auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
368 if (_process_data.deactivate_matrix_in_flow)
369 {
371 }
372 auto u = local_x.segment(displacement_index, displacement_size);
373
375}
376
377template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
378 int DisplacementDim>
380 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
381 postTimestepConcreteWithBlockVectors(
382 double const t, double const dt,
383 Eigen::Ref<const Eigen::VectorXd> const& p,
384 Eigen::Ref<const Eigen::VectorXd> const& u)
385{
386 MPL::VariableArray variables;
387 MPL::VariableArray variables_prev;
389
390 auto const e_id = _element.getID();
391 x_position.setElementID(e_id);
392
394 KV sigma_avg = KV::Zero();
395 DisplacementDimVector velocity_avg;
396 velocity_avg.setZero();
397
398 unsigned const n_integration_points = _ip_data.size();
399
400 auto const& medium = _process_data.media_map.getMedium(_element.getID());
401 auto const& liquid_phase =
403 auto const T_ref =
405 .template value<double>(variables, x_position, t, dt);
406 variables.temperature = T_ref;
407 variables_prev.temperature = T_ref;
408
409 auto const B_dil_bar = getDilatationalBBarMatrix();
410
411 for (unsigned ip = 0; ip < n_integration_points; ip++)
412 {
413 auto& ip_data = _ip_data[ip];
414
415 auto const& eps_prev = ip_data.eps_prev;
416 auto const& sigma_eff_prev = ip_data.sigma_eff_prev;
417
418 auto& eps = ip_data.eps;
419 auto& sigma_eff = ip_data.sigma_eff;
420 auto& state = ip_data.material_state_variables;
421
422 auto const& N_u = ip_data.N_u;
423 auto const& N_p = ip_data.N_p;
424 auto const& dNdx_u = ip_data.dNdx_u;
425
426 variables.liquid_phase_pressure = N_p.dot(p);
427
428 x_position = {
429 std::nullopt, _element.getID(),
431 NumLib::interpolateCoordinates<ShapeFunctionDisplacement,
433 _element, N_u))};
434
435 auto const x_coord = x_position.getCoordinates().value()[0];
436
437 auto const B =
439 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
441 dNdx_u, N_u, B_dil_bar, x_coord, this->_is_axially_symmetric)
442 .eval();
443
444 eps.noalias() = B * u;
445
446 variables.mechanical_strain
448 eps);
449
450 variables_prev.stress
452 sigma_eff_prev);
453 variables_prev.mechanical_strain
455 eps_prev);
456
457 auto&& solution = _ip_data[ip].solid_material.integrateStress(
458 variables_prev, variables, t, x_position, dt, *state);
459
460 if (!solution)
461 {
462 OGS_FATAL("Computation of local constitutive relation failed.");
463 }
464
466 std::tie(sigma_eff, state, C) = std::move(*solution);
467
468 sigma_avg += ip_data.sigma_eff;
469
470 if (!_process_data.deactivate_matrix_in_flow) // Only for hydraulically
471 // active matrix
472 {
473 auto const rho_fr =
474 liquid_phase.property(MPL::PropertyType::density)
475 .template value<double>(variables, x_position, t, dt);
476 variables.density = rho_fr;
477 auto const mu =
478 liquid_phase.property(MPL::PropertyType::viscosity)
479 .template value<double>(variables, x_position, t, dt);
480
482 medium->property(MPL::PropertyType::permeability)
483 .value(variables, x_position, t, dt))
484 .eval();
485
486 GlobalDimMatrixType const k_over_mu = k / mu;
487
488 auto const& gravity_vec = _process_data.specific_body_force;
489 auto const& dNdx_p = ip_data.dNdx_p;
490
491 ip_data.darcy_velocity.head(DisplacementDim).noalias() =
492 -k_over_mu * (dNdx_p * p + rho_fr * gravity_vec);
493 velocity_avg.noalias() +=
494 ip_data.darcy_velocity.head(DisplacementDim);
495 }
496 }
497
498 sigma_avg /= n_integration_points;
499 velocity_avg /= n_integration_points;
500
501 Eigen::Map<KV>(
502 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
504
505 Eigen::Map<DisplacementDimVector>(
506 &(*_process_data.element_velocities)[e_id * DisplacementDim]) =
507 velocity_avg;
508
510 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
511 DisplacementDim>(_element, _is_axially_symmetric, p,
512 *_process_data.mesh_prop_nodal_p);
513}
514
515template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
516 int DisplacementDim>
518 ShapeFunctionDisplacement, ShapeFunctionPressure,
519 DisplacementDim>::setPressureOfInactiveNodes(double const t,
520 Eigen::Ref<Eigen::VectorXd> p)
521{
523 x_position.setElementID(_element.getID());
524 for (unsigned i = 0; i < pressure_size; i++)
525 {
526 // only inactive nodes
527 if (_process_data.p_element_status->isActiveNode(_element.getNode(i)))
528 {
529 continue;
530 }
531 x_position.setNodeID(getNodeIndex(_element, i));
532 auto const p0 = (*_process_data.p0)(t, x_position)[0];
533 p[i] = p0;
534 }
535}
536
537template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
538 int DisplacementDim>
539std::vector<double> const& HydroMechanicsLocalAssemblerMatrix<
540 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
541 getIntPtSigma(
542 const double /*t*/,
543 std::vector<GlobalVector*> const& /*x*/,
544 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
545 std::vector<double>& cache) const
546{
549}
550template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
551 int DisplacementDim>
552std::vector<double> const& HydroMechanicsLocalAssemblerMatrix<
553 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
554 getIntPtEpsilon(
555 const double /*t*/,
556 std::vector<GlobalVector*> const& /*x*/,
557 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
558 std::vector<double>& cache) const
559{
562}
563
564template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
565 int DisplacementDim>
566std::vector<double> const& HydroMechanicsLocalAssemblerMatrix<
567 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
568 getIntPtDarcyVelocity(
569 const double /*t*/,
570 std::vector<GlobalVector*> const& /*x*/,
571 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
572 std::vector<double>& cache) const
573{
574 unsigned const n_integration_points = _ip_data.size();
575
576 cache.clear();
577 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
578 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
579 cache, DisplacementDim, n_integration_points);
580
581 for (unsigned ip = 0; ip < n_integration_points; ip++)
582 {
583 cache_matrix.col(ip).noalias() = _ip_data[ip].darcy_velocity;
584 }
585
586 return cache;
587}
588} // namespace HydroMechanics
589} // namespace LIE
590} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
constexpr double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setNodeID(std::size_t node_id)
void setElementID(std::size_t element_id)
std::optional< MathLib::Point3d > const getCoordinates() const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
HydroMechanicsLocalAssemblerInterface(MeshLib::Element const &element, bool const is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method, std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
void postTimestepConcreteWithBlockVectors(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
void assembleBlockMatricesWithJacobian(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_prev, Eigen::Ref< const Eigen::VectorXd > const &u, Eigen::Ref< const Eigen::VectorXd > const &u_prev, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_u, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pu, Eigen::Ref< Eigen::MatrixXd > J_uu, Eigen::Ref< Eigen::MatrixXd > J_up)
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
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< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void interpolateToHigherOrderNodes(MeshLib::Element const &element, bool const is_axially_symmetric, Eigen::MatrixBase< EigenMatrixType > const &node_values, MeshLib::PropertyVector< double > &interpolated_values_global_vector)
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)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.