OGS
PhaseFieldFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
33
34namespace ProcessLib
35{
36namespace PhaseField
37{
38template <typename BMatricesType, typename ShapeMatrixType, int DisplacementDim>
40{
49
50 typename ShapeMatrixType::NodalRowVectorType N;
51 typename ShapeMatrixType::GlobalDimNodalMatrixType dNdx;
52
53 typename BMatricesType::KelvinVectorType eps, eps_prev, eps_tensile;
54
55 typename BMatricesType::KelvinVectorType sigma_tensile, sigma_compressive,
58
60 std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
61 DisplacementDim>::MaterialStateVariables>
63
64 typename BMatricesType::KelvinMatrixType D, C_tensile, C_compressive;
65
68
76
77 template <typename DisplacementVectorType>
78 void updateConstitutiveRelation(double const t,
80 double const /*dt*/,
81 DisplacementVectorType const& /*u*/,
82 double const degradation,
83 EnergySplitModel const energy_split_model)
84 {
85 auto linear_elastic_mp =
87 DisplacementDim> const&>(solid_material)
88 .getMaterialProperties();
89
90 auto const bulk_modulus = linear_elastic_mp.bulk_modulus(t, x);
91 auto const mu = linear_elastic_mp.mu(t, x);
92
93 switch (energy_split_model)
94 {
96 {
101 degradation, bulk_modulus, mu, eps);
102 break;
103 }
105 {
110 degradation, bulk_modulus, mu, eps);
111 break;
112 }
114 {
117 C_compressive) = MaterialLib::Solids::Phasefield::
118 calculateIsotropicDegradedStressWithRankineEnergy<
119 DisplacementDim>(degradation, bulk_modulus, mu, eps);
120 break;
121 }
123 {
124 double temperature = 0.;
126 C_ortho = static_cast<
128 DisplacementDim> const&>(solid_material)
129 .getElasticTensor(t, x, temperature);
130
135 degradation, eps, C_ortho);
136 break;
137 }
139 {
140 double temperature = 0.;
142 C_ortho = static_cast<
144 DisplacementDim> const&>(solid_material)
145 .getElasticTensor(t, x, temperature);
146
147 std::tie(eps_tensile, sigma, sigma_tensile, D,
151 degradation, eps, C_ortho);
152 break;
153 }
154 }
155
158 }
159};
160
163template <typename ShapeMatrixType>
165{
166 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
167};
168
169template <typename ShapeFunction, int DisplacementDim>
171{
172private:
173 static constexpr int displacement_index = 0;
174 static constexpr int displacement_size =
175 ShapeFunction::NPOINTS * DisplacementDim;
176 static constexpr int phasefield_index =
178 static constexpr int phasefield_size = ShapeFunction::NPOINTS;
179
180public:
183 // Types for displacement.
184 // (Higher order elements = ShapeFunction).
187
189
191 typename ShapeMatricesType::template VectorType<displacement_size>;
193 typename ShapeMatricesType::template VectorType<phasefield_size>;
195 typename ShapeMatricesType::template MatrixType<phasefield_size,
199 using IpData =
201
204
206 MeshLib::Element const& e,
207 std::size_t const /*local_matrix_size*/,
208 NumLib::GenericIntegrationMethod const& integration_method,
209 bool const is_axially_symmetric,
211 : _process_data(process_data),
212 _integration_method(integration_method),
213 _element(e),
214 _is_axially_symmetric(is_axially_symmetric)
215 {
216 unsigned const n_integration_points =
218
219 _ip_data.reserve(n_integration_points);
220 _secondary_data.N.resize(n_integration_points);
221
222 auto& solid_material =
224 _process_data.solid_materials,
225 _process_data.material_ids,
226 e.getID());
227
228 auto const shape_matrices =
230 DisplacementDim>(e, is_axially_symmetric,
232
234 x_position.setElementID(_element.getID());
235
236 for (unsigned ip = 0; ip < n_integration_points; ip++)
237 {
238 _ip_data.emplace_back(solid_material);
239 auto& ip_data = _ip_data[ip];
240 ip_data.integration_weight =
242 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
243
244 static const int kelvin_vector_size =
246 DisplacementDim);
247 ip_data.eps_tensile.setZero(kelvin_vector_size);
248 ip_data.eps.setZero(kelvin_vector_size);
249 ip_data.eps_prev.resize(kelvin_vector_size);
250 ip_data.D.setZero(kelvin_vector_size, kelvin_vector_size);
251
252 ip_data.sigma_tensile.setZero(kelvin_vector_size);
253 ip_data.sigma_compressive.setZero(kelvin_vector_size);
254 ip_data.sigma.setZero(kelvin_vector_size);
255 ip_data.strain_energy_tensile = 0.0;
256 ip_data.elastic_energy = 0.0;
257
258 ip_data.N = shape_matrices[ip].N;
259 ip_data.dNdx = shape_matrices[ip].dNdx;
260
261 _secondary_data.N[ip] = shape_matrices[ip].N;
262 }
263 }
264
266 std::size_t setIPDataInitialConditions(
267 std::string_view const name,
268 double const* values,
269 int const integration_order) override;
270
271 void assemble(double const /*t*/, double const /*dt*/,
272 std::vector<double> const& /*local_x*/,
273 std::vector<double> const& /*local_x_prev*/,
274 std::vector<double>& /*local_M_data*/,
275 std::vector<double>& /*local_K_data*/,
276 std::vector<double>& /*local_rhs_data*/) override
277 {
278 OGS_FATAL(
279 "PhaseFieldLocalAssembler: assembly without jacobian is not "
280 "implemented.");
281 }
282
284 double const t, double const dt, Eigen::VectorXd const& local_x,
285 Eigen::VectorXd const& local_x_prev, int const process_id,
286 std::vector<double>& local_b_data,
287 std::vector<double>& local_Jac_data) override;
288
289 void initializeConcrete() override
290 {
291 unsigned const n_integration_points =
293
294 for (unsigned ip = 0; ip < n_integration_points; ip++)
295 {
296 auto& ip_data = _ip_data[ip];
297
298 ParameterLib::SpatialPosition const x_position{
299 std::nullopt, _element.getID(), ip,
301 NumLib::interpolateCoordinates<ShapeFunction,
303 _element, ip_data.N))};
304
306 if (_process_data.initial_stress.value)
307 {
308 ip_data.sigma =
310 DisplacementDim>((*_process_data.initial_stress.value)(
311 std::numeric_limits<
312 double>::quiet_NaN() /* time independent */,
313 x_position));
314 }
315
316 ip_data.pushBackState();
317 }
318 }
319
320 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
321 Eigen::VectorXd const& /*local_x_prev*/,
322 double const /*t*/, double const /*dt*/,
323 int const /*process_id*/) override
324 {
325 unsigned const n_integration_points =
327
328 for (unsigned ip = 0; ip < n_integration_points; ip++)
329 {
330 _ip_data[ip].pushBackState();
331 }
332 }
333
335 std::size_t mesh_item_id,
336 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
337 std::vector<GlobalVector*> const& x, double const t,
338 double& crack_volume) override;
339
340 void computeEnergy(
341 std::size_t mesh_item_id,
342 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
343 std::vector<GlobalVector*> const& x, double const t,
344 double& elastic_energy, double& surface_energy,
345 double& pressure_work) override;
346
347 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
348 const unsigned integration_point) const override
349 {
350 auto const& N = _secondary_data.N[integration_point];
351
352 // assumes N is stored contiguously in memory
353 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
354 }
355
356 // TODO (naumov) This method is same as getIntPtSigma but for arguments and
357 // the ordering of the cache_mat.
358 // There should be only one.
359 std::vector<double> getSigma() const override;
360
361 std::vector<double> getEpsilon() const override;
362
363private:
364 std::vector<double> const& getIntPtSigma(
365 const double /*t*/,
366 std::vector<GlobalVector*> const& /*x*/,
367 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
368 std::vector<double>& cache) const override
369 {
371 _ip_data, &IpData::sigma, cache);
372 }
373
374 std::vector<double> const& getIntPtSigmaTensile(
375 const double /*t*/,
376 std::vector<GlobalVector*> const& /*x*/,
377 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
378 std::vector<double>& cache) const override
379 {
382 }
383
384 std::vector<double> const& getIntPtSigmaCompressive(
385 const double /*t*/,
386 std::vector<GlobalVector*> const& /*x*/,
387 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
388 std::vector<double>& cache) const override
389 {
392 }
393
394 std::vector<double> const& getIntPtEpsilon(
395 const double /*t*/,
396 std::vector<GlobalVector*> const& /*x*/,
397 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
398 std::vector<double>& cache) const override
399 {
401 _ip_data, &IpData::eps, cache);
402 }
403
404 std::vector<double> const& getIntPtEpsilonTensile(
405 const double /*t*/,
406 std::vector<GlobalVector*> const& /*x*/,
407 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
408 std::vector<double>& cache) const override
409 {
412 }
413
415 double const t, double const dt, Eigen::VectorXd const& local_x,
416 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
417
419 double const t, double const dt, Eigen::VectorXd const& local_x,
420 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
421
423
424 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
425
430
431 static const int phase_process_id = 1;
432 static const int mechanics_process_id = 0;
433};
434
435} // namespace PhaseField
436} // namespace ProcessLib
437
438#include "PhaseFieldFEM-impl.h"
#define OGS_FATAL(...)
Definition Error.h:26
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)
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
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
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
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateOrthoVolDevDegradedStress(double const degradation, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_ortho)
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateIsotropicDegradedStress(double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateVolDevDegradedStress(double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateOrthoMasonryDegradedStress(double const degradation, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_ortho)
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)
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)
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
BMatricesType::KelvinVectorType sigma_tensile
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material
ShapeMatrixType::NodalRowVectorType N
BMatricesType::KelvinVectorType eps_prev
BMatricesType::KelvinMatrixType C_compressive
BMatricesType::KelvinVectorType sigma_compressive
ShapeMatrixType::GlobalDimNodalMatrixType dNdx
BMatricesType::KelvinVectorType eps_tensile
BMatricesType::KelvinVectorType eps
BMatricesType::KelvinVectorType sigma
void updateConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x, double const, DisplacementVectorType const &, double const degradation, EnergySplitModel const energy_split_model)
BMatricesType::KelvinMatrixType D
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
BMatricesType::KelvinMatrixType C_tensile
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N