30template <
typename ShapeFunction,
int DisplacementDim>
36 bool const is_axially_symmetric,
43 unsigned const n_integration_points =
46 _ip_data.reserve(n_integration_points);
49 auto const shape_matrices =
51 DisplacementDim>(e, is_axially_symmetric,
57 for (
unsigned ip = 0; ip < n_integration_points; ip++)
59 _ip_data.emplace_back(solid_material);
61 auto const& sm = shape_matrices[ip];
63 ip_data.dNdx_u = sm.dNdx;
64 ip_data.integration_weight =
66 sm.integralMeasure * sm.detJ;
69 static const int kelvin_vector_size =
71 ip_data._sigma.setZero(kelvin_vector_size);
72 ip_data._eps.setZero(kelvin_vector_size);
75 ip_data._sigma_prev.resize(kelvin_vector_size);
76 ip_data._eps_prev.resize(kelvin_vector_size);
78 ip_data._C.resize(kelvin_vector_size, kelvin_vector_size);
84template <
typename ShapeFunction,
int DisplacementDim>
87 std::vector<double>
const& local_x,
88 std::vector<double>
const& ,
89 std::vector<double>& local_b_data,
90 std::vector<double>& local_Jac_data)
92 assert(
_element.getDimension() == DisplacementDim);
94 auto const local_matrix_size = local_x.size();
97 local_Jac_data, local_matrix_size, local_matrix_size);
100 local_b_data, local_matrix_size);
102 unsigned const n_integration_points =
110 for (
unsigned ip = 0; ip < n_integration_points; ip++)
112 auto const& w =
_ip_data[ip].integration_weight;
115 auto const& dNdx =
_ip_data[ip].dNdx_u;
118 std::nullopt, this->
_element.getID(),
130 auto const& eps_prev =
_ip_data[ip]._eps_prev;
131 auto const& sigma_prev =
_ip_data[ip]._sigma_prev;
135 auto& state =
_ip_data[ip]._material_state_variables;
138 B * Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
139 local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
152 auto&& solution =
_ip_data[ip]._solid_material.integrateStress(
153 variables_prev, variables, t, x_position, dt, *state);
157 OGS_FATAL(
"Computation of local constitutive relation failed.");
161 std::tie(sigma, state, C) = std::move(*solution);
163 local_b.noalias() -= B.transpose() * sigma * w;
164 local_Jac.noalias() += B.transpose() * C * B * w;
168template <
typename ShapeFunction,
int DisplacementDim>
171 double const , Eigen::VectorXd
const& )
175 KV sigma_avg = KV::Zero();
178 unsigned const n_integration_points =
180 for (
unsigned ip = 0; ip < n_integration_points; ip++)
184 sigma_avg /= n_integration_points;
187 &(*
_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
191template <
typename ShapeFunction,
int DisplacementDim>
192std::vector<double>
const&
196 std::vector<GlobalVector*>
const& ,
197 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
198 std::vector<double>& cache)
const
203template <
typename ShapeFunction,
int DisplacementDim>
204std::vector<double>
const&
208 std::vector<GlobalVector*>
const& ,
209 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
210 std::vector<double>& cache)
const
KelvinVector mechanical_strain
std::size_t getID() const
Returns the ID of the element.
std::optional< MathLib::Point3d > const getCoordinates() const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
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)