OGS
PhaseFieldFEM.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 <memory>
7#include <vector>
8
26
27namespace ProcessLib
28{
29namespace PhaseField
30{
31template <typename BMatricesType, typename ShapeMatrixType, int DisplacementDim>
33{
42
43 typename ShapeMatrixType::NodalRowVectorType N;
44 typename ShapeMatrixType::GlobalDimNodalMatrixType dNdx;
45
46 typename BMatricesType::KelvinVectorType eps, eps_prev, eps_tensile;
47
48 typename BMatricesType::KelvinVectorType sigma_tensile, sigma_compressive,
51
53 std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
54 DisplacementDim>::MaterialStateVariables>
56
57 typename BMatricesType::KelvinMatrixType D, C_tensile, C_compressive;
58
61
69
70 template <typename DisplacementVectorType>
72 double const t,
74 double const /*dt*/,
75 DisplacementVectorType const& /*u*/,
76 double const degradation,
78 energy_split_model)
79 {
81 decltype(sigma), decltype(D), DisplacementDim>(
84 eps, energy_split_model, t, x, solid_material);
85
88 }
89};
90
93template <typename ShapeMatrixType>
95{
96 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
97};
98
99template <typename ShapeFunction, int DisplacementDim>
101{
102private:
103 static constexpr int displacement_index = 0;
104 static constexpr int displacement_size =
105 ShapeFunction::NPOINTS * DisplacementDim;
106 static constexpr int phasefield_index =
108 static constexpr int phasefield_size = ShapeFunction::NPOINTS;
109
110public:
113 // Types for displacement.
114 // (Higher order elements = ShapeFunction).
117
119
121 typename ShapeMatricesType::template VectorType<displacement_size>;
123 typename ShapeMatricesType::template VectorType<phasefield_size>;
125 typename ShapeMatricesType::template MatrixType<phasefield_size,
129 using IpData =
131
134
136 MeshLib::Element const& e,
137 std::size_t const /*local_matrix_size*/,
138 NumLib::GenericIntegrationMethod const& integration_method,
139 bool const is_axially_symmetric,
141 : _process_data(process_data),
142 _integration_method(integration_method),
143 _element(e),
144 _is_axially_symmetric(is_axially_symmetric)
145 {
146 unsigned const n_integration_points =
147 _integration_method.getNumberOfPoints();
148
149 _ip_data.reserve(n_integration_points);
150 _secondary_data.N.resize(n_integration_points);
151
152 auto& solid_material =
154 _process_data.solid_materials,
155 _process_data.material_ids,
156 e.getID());
157
158 auto const shape_matrices =
160 DisplacementDim>(e, is_axially_symmetric,
162
164 x_position.setElementID(_element.getID());
165
166 for (unsigned ip = 0; ip < n_integration_points; ip++)
167 {
168 _ip_data.emplace_back(solid_material);
169 auto& ip_data = _ip_data[ip];
170 ip_data.integration_weight =
171 _integration_method.getWeightedPoint(ip).getWeight() *
172 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
173
174 static const int kelvin_vector_size =
176 DisplacementDim);
177 ip_data.eps_tensile.setZero(kelvin_vector_size);
178 ip_data.eps.setZero(kelvin_vector_size);
179 ip_data.eps_prev.resize(kelvin_vector_size);
180 ip_data.D.setZero(kelvin_vector_size, kelvin_vector_size);
181
182 ip_data.sigma_tensile.setZero(kelvin_vector_size);
183 ip_data.sigma_compressive.setZero(kelvin_vector_size);
184 ip_data.sigma.setZero(kelvin_vector_size);
185 ip_data.strain_energy_tensile = 0.0;
186 ip_data.elastic_energy = 0.0;
187
188 ip_data.N = shape_matrices[ip].N;
189 ip_data.dNdx = shape_matrices[ip].dNdx;
190
191 _secondary_data.N[ip] = shape_matrices[ip].N;
192 }
193 }
194
196 std::size_t setIPDataInitialConditions(
197 std::string_view const name,
198 double const* values,
199 int const integration_order) override;
200
201 void assemble(double const /*t*/, double const /*dt*/,
202 std::vector<double> const& /*local_x*/,
203 std::vector<double> const& /*local_x_prev*/,
204 std::vector<double>& /*local_M_data*/,
205 std::vector<double>& /*local_K_data*/,
206 std::vector<double>& /*local_rhs_data*/) override
207 {
208 OGS_FATAL(
209 "PhaseFieldLocalAssembler: assembly without jacobian is not "
210 "implemented.");
211 }
212
214 double const t, double const dt, Eigen::VectorXd const& local_x,
215 Eigen::VectorXd const& local_x_prev, int const process_id,
216 std::vector<double>& local_b_data,
217 std::vector<double>& local_Jac_data) override;
218
219 void initializeConcrete() override
220 {
221 unsigned const n_integration_points =
222 _integration_method.getNumberOfPoints();
223
224 for (unsigned ip = 0; ip < n_integration_points; ip++)
225 {
226 auto& ip_data = _ip_data[ip];
227
228 ParameterLib::SpatialPosition const x_position{
229 std::nullopt, _element.getID(),
233 _element, ip_data.N))};
234
236 if (_process_data.initial_stress.value)
237 {
238 ip_data.sigma =
240 DisplacementDim>((*_process_data.initial_stress.value)(
241 std::numeric_limits<
242 double>::quiet_NaN() /* time independent */,
243 x_position));
244 }
245
246 ip_data.pushBackState();
247 }
248 }
249
250 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
251 Eigen::VectorXd const& /*local_x_prev*/,
252 double const /*t*/, double const /*dt*/,
253 int const /*process_id*/) override
254 {
255 unsigned const n_integration_points =
256 _integration_method.getNumberOfPoints();
257
258 for (unsigned ip = 0; ip < n_integration_points; ip++)
259 {
260 _ip_data[ip].pushBackState();
261 }
262 }
263
265 std::size_t mesh_item_id,
266 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
267 std::vector<GlobalVector*> const& x, double const t,
268 double& crack_volume) override;
269
270 void computeEnergy(
271 std::size_t mesh_item_id,
272 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
273 std::vector<GlobalVector*> const& x, double const t,
274 double& elastic_energy, double& surface_energy,
275 double& pressure_work) override;
276
277 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
278 const unsigned integration_point) const override
279 {
280 auto const& N = _secondary_data.N[integration_point];
281
282 // assumes N is stored contiguously in memory
283 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
284 }
285
286 // TODO (naumov) This method is same as getIntPtSigma but for arguments and
287 // the ordering of the cache_mat.
288 // There should be only one.
289 std::vector<double> getSigma() const override;
290
291 std::vector<double> getEpsilon() const override;
292
293private:
294 std::vector<double> const& getIntPtSigma(
295 const double /*t*/,
296 std::vector<GlobalVector*> const& /*x*/,
297 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
298 std::vector<double>& cache) const override
299 {
301 _ip_data, &IpData::sigma, cache);
302 }
303
304 std::vector<double> const& getIntPtSigmaTensile(
305 const double /*t*/,
306 std::vector<GlobalVector*> const& /*x*/,
307 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
308 std::vector<double>& cache) const override
309 {
312 }
313
314 std::vector<double> const& getIntPtSigmaCompressive(
315 const double /*t*/,
316 std::vector<GlobalVector*> const& /*x*/,
317 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
318 std::vector<double>& cache) const override
319 {
322 }
323
324 std::vector<double> const& getIntPtEpsilon(
325 const double /*t*/,
326 std::vector<GlobalVector*> const& /*x*/,
327 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
328 std::vector<double>& cache) const override
329 {
331 _ip_data, &IpData::eps, cache);
332 }
333
334 std::vector<double> const& getIntPtEpsilonTensile(
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 override
339 {
342 }
343
345 double const t, double const dt, Eigen::VectorXd const& local_x,
346 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
347
349 double const t, double const dt, Eigen::VectorXd const& local_x,
350 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
351
353
354 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
355
360
361 static const int phase_process_id = 1;
362 static const int mechanics_process_id = 0;
363};
364
365} // namespace PhaseField
366} // namespace ProcessLib
367
368#include "PhaseFieldFEM-impl.h"
#define OGS_FATAL(...)
Definition Error.h:19
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
void setElementID(std::size_t element_id)
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
void assembleWithJacobianForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
typename ShapeMatricesType::template VectorType< phasefield_size > PhaseFieldVector
PhaseFieldProcessData< DisplacementDim > & _process_data
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
void assembleWithJacobianForDeformationEquations(double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
void assembleWithJacobianPhaseFieldEquations(double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
std::vector< double > const & getIntPtSigmaCompressive(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
Returns number of read integration points.
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename BMatricesType::NodalForceVectorType NodalForceVectorType
NumLib::GenericIntegrationMethod const & _integration_method
PhaseFieldLocalAssembler(PhaseFieldLocalAssembler const &)=delete
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtEpsilonTensile(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
typename ShapeMatricesType::template VectorType< displacement_size > DeformationVector
typename ShapeMatricesType::template MatrixType< phasefield_size, phasefield_size > PhaseFieldMatrix
BMatrixPolicyType< ShapeFunction, DisplacementDim > BMatricesType
PhaseFieldLocalAssembler(MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, PhaseFieldProcessData< DisplacementDim > &process_data)
void computeCrackIntegral(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &crack_volume) override
PhaseFieldLocalAssembler(PhaseFieldLocalAssembler &&)=delete
void computeEnergy(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work) override
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim > IpData
std::vector< double > getSigma() const override
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
std::vector< double > const & getIntPtSigmaTensile(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > getEpsilon() const override
typename ShapeMatricesType::NodalVectorType NodalVectorType
void calculateStress(T_VECTOR &sigma, T_VECTOR &sigma_tensile, T_VECTOR &sigma_compressive, T_VECTOR &eps_tensile, T_MATRIX &D, T_MATRIX &C_tensile, T_MATRIX &C_compressive, double &strain_energy_tensile, double &elastic_energy, double const degradation, T_VECTOR const &eps, EnergySplitModel const &energy_split_model, double const t, ParameterLib::SpatialPosition const &x, MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
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)
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
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
void updateConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x, double const, DisplacementVectorType const &, double const degradation, MaterialLib::Solids::Phasefield::EnergySplitModel const energy_split_model)
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N