Loading [MathJax]/extensions/tex2jax.js
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>
79 double const t,
81 double const /*dt*/,
82 DisplacementVectorType const& /*u*/,
83 double const degradation,
85 energy_split_model)
86 {
88 decltype(sigma), decltype(D), DisplacementDim>(
91 eps, energy_split_model, t, x, solid_material);
92
95 }
96};
97
100template <typename ShapeMatrixType>
102{
103 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
104};
105
106template <typename ShapeFunction, int DisplacementDim>
108{
109private:
110 static constexpr int displacement_index = 0;
111 static constexpr int displacement_size =
112 ShapeFunction::NPOINTS * DisplacementDim;
113 static constexpr int phasefield_index =
115 static constexpr int phasefield_size = ShapeFunction::NPOINTS;
116
117public:
120 // Types for displacement.
121 // (Higher order elements = ShapeFunction).
124
126
128 typename ShapeMatricesType::template VectorType<displacement_size>;
130 typename ShapeMatricesType::template VectorType<phasefield_size>;
132 typename ShapeMatricesType::template MatrixType<phasefield_size,
136 using IpData =
138
141
143 MeshLib::Element const& e,
144 std::size_t const /*local_matrix_size*/,
145 NumLib::GenericIntegrationMethod const& integration_method,
146 bool const is_axially_symmetric,
148 : _process_data(process_data),
149 _integration_method(integration_method),
150 _element(e),
151 _is_axially_symmetric(is_axially_symmetric)
152 {
153 unsigned const n_integration_points =
155
156 _ip_data.reserve(n_integration_points);
157 _secondary_data.N.resize(n_integration_points);
158
159 auto& solid_material =
161 _process_data.solid_materials,
162 _process_data.material_ids,
163 e.getID());
164
165 auto const shape_matrices =
167 DisplacementDim>(e, is_axially_symmetric,
169
171 x_position.setElementID(_element.getID());
172
173 for (unsigned ip = 0; ip < n_integration_points; ip++)
174 {
175 _ip_data.emplace_back(solid_material);
176 auto& ip_data = _ip_data[ip];
177 ip_data.integration_weight =
179 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
180
181 static const int kelvin_vector_size =
183 DisplacementDim);
184 ip_data.eps_tensile.setZero(kelvin_vector_size);
185 ip_data.eps.setZero(kelvin_vector_size);
186 ip_data.eps_prev.resize(kelvin_vector_size);
187 ip_data.D.setZero(kelvin_vector_size, kelvin_vector_size);
188
189 ip_data.sigma_tensile.setZero(kelvin_vector_size);
190 ip_data.sigma_compressive.setZero(kelvin_vector_size);
191 ip_data.sigma.setZero(kelvin_vector_size);
192 ip_data.strain_energy_tensile = 0.0;
193 ip_data.elastic_energy = 0.0;
194
195 ip_data.N = shape_matrices[ip].N;
196 ip_data.dNdx = shape_matrices[ip].dNdx;
197
198 _secondary_data.N[ip] = shape_matrices[ip].N;
199 }
200 }
201
203 std::size_t setIPDataInitialConditions(
204 std::string_view const name,
205 double const* values,
206 int const integration_order) override;
207
208 void assemble(double const /*t*/, double const /*dt*/,
209 std::vector<double> const& /*local_x*/,
210 std::vector<double> const& /*local_x_prev*/,
211 std::vector<double>& /*local_M_data*/,
212 std::vector<double>& /*local_K_data*/,
213 std::vector<double>& /*local_rhs_data*/) override
214 {
215 OGS_FATAL(
216 "PhaseFieldLocalAssembler: assembly without jacobian is not "
217 "implemented.");
218 }
219
221 double const t, double const dt, Eigen::VectorXd const& local_x,
222 Eigen::VectorXd const& local_x_prev, int const process_id,
223 std::vector<double>& local_b_data,
224 std::vector<double>& local_Jac_data) override;
225
226 void initializeConcrete() override
227 {
228 unsigned const n_integration_points =
230
231 for (unsigned ip = 0; ip < n_integration_points; ip++)
232 {
233 auto& ip_data = _ip_data[ip];
234
235 ParameterLib::SpatialPosition const x_position{
236 std::nullopt, _element.getID(),
238 NumLib::interpolateCoordinates<ShapeFunction,
240 _element, ip_data.N))};
241
243 if (_process_data.initial_stress.value)
244 {
245 ip_data.sigma =
247 DisplacementDim>((*_process_data.initial_stress.value)(
248 std::numeric_limits<
249 double>::quiet_NaN() /* time independent */,
250 x_position));
251 }
252
253 ip_data.pushBackState();
254 }
255 }
256
257 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
258 Eigen::VectorXd const& /*local_x_prev*/,
259 double const /*t*/, double const /*dt*/,
260 int const /*process_id*/) override
261 {
262 unsigned const n_integration_points =
264
265 for (unsigned ip = 0; ip < n_integration_points; ip++)
266 {
267 _ip_data[ip].pushBackState();
268 }
269 }
270
272 std::size_t mesh_item_id,
273 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
274 std::vector<GlobalVector*> const& x, double const t,
275 double& crack_volume) override;
276
277 void computeEnergy(
278 std::size_t mesh_item_id,
279 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
280 std::vector<GlobalVector*> const& x, double const t,
281 double& elastic_energy, double& surface_energy,
282 double& pressure_work) override;
283
284 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
285 const unsigned integration_point) const override
286 {
287 auto const& N = _secondary_data.N[integration_point];
288
289 // assumes N is stored contiguously in memory
290 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
291 }
292
293 // TODO (naumov) This method is same as getIntPtSigma but for arguments and
294 // the ordering of the cache_mat.
295 // There should be only one.
296 std::vector<double> getSigma() const override;
297
298 std::vector<double> getEpsilon() const override;
299
300private:
301 std::vector<double> const& getIntPtSigma(
302 const double /*t*/,
303 std::vector<GlobalVector*> const& /*x*/,
304 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
305 std::vector<double>& cache) const override
306 {
308 _ip_data, &IpData::sigma, cache);
309 }
310
311 std::vector<double> const& getIntPtSigmaTensile(
312 const double /*t*/,
313 std::vector<GlobalVector*> const& /*x*/,
314 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
315 std::vector<double>& cache) const override
316 {
319 }
320
321 std::vector<double> const& getIntPtSigmaCompressive(
322 const double /*t*/,
323 std::vector<GlobalVector*> const& /*x*/,
324 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
325 std::vector<double>& cache) const override
326 {
329 }
330
331 std::vector<double> const& getIntPtEpsilon(
332 const double /*t*/,
333 std::vector<GlobalVector*> const& /*x*/,
334 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
335 std::vector<double>& cache) const override
336 {
338 _ip_data, &IpData::eps, cache);
339 }
340
341 std::vector<double> const& getIntPtEpsilonTensile(
342 const double /*t*/,
343 std::vector<GlobalVector*> const& /*x*/,
344 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
345 std::vector<double>& cache) const override
346 {
349 }
350
352 double const t, double const dt, Eigen::VectorXd const& local_x,
353 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
354
356 double const t, double const dt, Eigen::VectorXd const& local_x,
357 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
358
360
361 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
362
367
368 static const int phase_process_id = 1;
369 static const int mechanics_process_id = 0;
370};
371
372} // namespace PhaseField
373} // namespace ProcessLib
374
375#include "PhaseFieldFEM-impl.h"
#define OGS_FATAL(...)
Definition Error.h:26
constexpr 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
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
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
BMatricesType::KelvinMatrixType D
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
BMatricesType::KelvinMatrixType C_tensile
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