OGS
PhaseFieldFEM.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <memory>
14 #include <vector>
15 
25 #include "PhaseFieldProcessData.h"
29 
30 namespace ProcessLib
31 {
32 namespace PhaseField
33 {
34 template <typename BMatricesType, typename ShapeMatrixType, int DisplacementDim>
36 {
42  solid_material.createMaterialStateVariables())
43  {
44  }
45 
46  typename ShapeMatrixType::NodalRowVectorType N;
47  typename ShapeMatrixType::GlobalDimNodalMatrixType dNdx;
48 
50 
54 
56  std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
57  DisplacementDim>::MaterialStateVariables>
59 
63 
65  {
68  eps_prev = eps;
69  material_state_variables->pushBackState();
70  }
71 
72  template <typename DisplacementVectorType>
73  void updateConstitutiveRelation(double const t,
75  double const /*dt*/,
76  DisplacementVectorType const& /*u*/,
77  double const degradation)
78  {
79  auto linear_elastic_mp =
81  DisplacementDim> const&>(solid_material)
82  .getMaterialProperties();
83 
84  auto const bulk_modulus = linear_elastic_mp.bulk_modulus(t, x);
85  auto const mu = linear_elastic_mp.mu(t, x);
86 
90  DisplacementDim>(degradation, bulk_modulus, mu, eps);
91 
94  }
96 
97  static constexpr int kelvin_vector_size =
99 };
100 
103 template <typename ShapeMatrixType>
105 {
106  std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
107 };
108 
109 template <typename ShapeFunction, typename IntegrationMethod,
110  int DisplacementDim>
112 {
113 private:
114  static constexpr int displacement_index = 0;
115  static constexpr int displacement_size =
116  ShapeFunction::NPOINTS * DisplacementDim;
117  static constexpr int phasefield_index =
119  static constexpr int phasefield_size = ShapeFunction::NPOINTS;
120 
121 public:
124  // Types for displacement.
125  // (Higher order elements = ShapeFunction).
128 
130 
132  typename ShapeMatricesType::template VectorType<displacement_size>;
134  typename ShapeMatricesType::template VectorType<phasefield_size>;
136  typename ShapeMatricesType::template MatrixType<phasefield_size,
140  using IpData =
142 
145 
147  MeshLib::Element const& e,
148  std::size_t const /*local_matrix_size*/,
149  bool const is_axially_symmetric,
150  unsigned const integration_order,
152  : _process_data(process_data),
153  _integration_method(integration_order),
154  _element(e),
155  _is_axially_symmetric(is_axially_symmetric)
156  {
157  unsigned const n_integration_points =
158  _integration_method.getNumberOfPoints();
159 
160  _ip_data.reserve(n_integration_points);
161  _secondary_data.N.resize(n_integration_points);
162 
163  auto& solid_material =
165  _process_data.solid_materials,
166  _process_data.material_ids,
167  e.getID());
169  DisplacementDim> const*>(&solid_material))
170  {
171  OGS_FATAL(
172  "For phasefield process only linear elastic solid material "
173  "support is implemented.");
174  }
175 
176  auto const shape_matrices =
178  DisplacementDim>(e, is_axially_symmetric,
180 
182  x_position.setElementID(_element.getID());
183 
184  for (unsigned ip = 0; ip < n_integration_points; ip++)
185  {
186  _ip_data.emplace_back(solid_material);
187  auto& ip_data = _ip_data[ip];
188  ip_data.integration_weight =
189  _integration_method.getWeightedPoint(ip).getWeight() *
190  shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
191 
192  static const int kelvin_vector_size =
194  DisplacementDim);
195  ip_data.eps.setZero(kelvin_vector_size);
196  ip_data.eps_prev.resize(kelvin_vector_size);
197  ip_data.C_tensile.setZero(kelvin_vector_size, kelvin_vector_size);
198  ip_data.C_compressive.setZero(kelvin_vector_size,
199  kelvin_vector_size);
200  ip_data.sigma_tensile.setZero(kelvin_vector_size);
201  ip_data.sigma_compressive.setZero(kelvin_vector_size);
202  ip_data.sigma.setZero(kelvin_vector_size);
203  ip_data.strain_energy_tensile = 0.0;
204  ip_data.elastic_energy = 0.0;
205 
206  ip_data.N = shape_matrices[ip].N;
207  ip_data.dNdx = shape_matrices[ip].dNdx;
208 
209  _secondary_data.N[ip] = shape_matrices[ip].N;
210  }
211  }
212 
213  void assemble(double const /*t*/, double const /*dt*/,
214  std::vector<double> const& /*local_x*/,
215  std::vector<double> const& /*local_xdot*/,
216  std::vector<double>& /*local_M_data*/,
217  std::vector<double>& /*local_K_data*/,
218  std::vector<double>& /*local_rhs_data*/) override
219  {
220  OGS_FATAL(
221  "PhaseFieldLocalAssembler: assembly without jacobian is not "
222  "implemented.");
223  }
224 
226  double const t, double const dt, Eigen::VectorXd const& local_x,
227  Eigen::VectorXd const& local_xdot, const double dxdot_dx,
228  const double dx_dx, int const process_id,
229  std::vector<double>& local_M_data, std::vector<double>& local_K_data,
230  std::vector<double>& local_b_data,
231  std::vector<double>& local_Jac_data) override;
232 
233  void initializeConcrete() override
234  {
235  unsigned const n_integration_points =
236  _integration_method.getNumberOfPoints();
237 
238  for (unsigned ip = 0; ip < n_integration_points; ip++)
239  {
240  _ip_data[ip].pushBackState();
241  }
242  }
243 
244  void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
245  double const /*t*/, double const /*dt*/) override
246  {
247  unsigned const n_integration_points =
248  _integration_method.getNumberOfPoints();
249 
250  for (unsigned ip = 0; ip < n_integration_points; ip++)
251  {
252  _ip_data[ip].pushBackState();
253  }
254  }
255 
257  std::size_t mesh_item_id,
258  std::vector<
259  std::reference_wrapper<NumLib::LocalToGlobalIndexMap>> const&
260  dof_tables,
261  GlobalVector const& x, double const t, double& crack_volume,
262  CoupledSolutionsForStaggeredScheme const* const cpl_xs) override;
263 
264  void computeEnergy(
265  std::size_t mesh_item_id,
266  std::vector<
267  std::reference_wrapper<NumLib::LocalToGlobalIndexMap>> const&
268  dof_tables,
269  GlobalVector const& x, double const t, double& elastic_energy,
270  double& surface_energy, double& pressure_work,
271  CoupledSolutionsForStaggeredScheme const* const cpl_xs) override;
272 
273  Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
274  const unsigned integration_point) const override
275  {
276  auto const& N = _secondary_data.N[integration_point];
277 
278  // assumes N is stored contiguously in memory
279  return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
280  }
281 
282 private:
283  std::vector<double> const& getIntPtSigma(
284  const double /*t*/,
285  std::vector<GlobalVector*> const& /*x*/,
286  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
287  std::vector<double>& cache) const override
288  {
289  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
290  _ip_data, &IpData::sigma, cache);
291  }
292 
293  std::vector<double> const& getIntPtEpsilon(
294  const double /*t*/,
295  std::vector<GlobalVector*> const& /*x*/,
296  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
297  std::vector<double>& cache) const override
298  {
299  return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
300  _ip_data, &IpData::eps, cache);
301  }
302 
304  double const t, double const dt, Eigen::VectorXd const& local_x,
305  std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
306 
308  double const t, double const dt, Eigen::VectorXd const& local_x,
309  std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
310 
312 
313  std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
314 
315  IntegrationMethod _integration_method;
319 
320  static const int phase_process_id = 1;
321  static const int mechanics_process_id = 0;
322 };
323 
324 } // namespace PhaseField
325 } // namespace ProcessLib
326 
327 #include "PhaseFieldFEM-impl.h"
#define OGS_FATAL(...)
Definition: Error.h:26
Global vector based on Eigen vector.
Definition: EigenVector.h:26
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
void setElementID(std::size_t element_id)
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
Definition: BMatrixPolicy.h:44
typename BMatricesType::NodalForceVectorType NodalForceVectorType
std::vector< double > const & getIntPtSigma(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
std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void computeCrackIntegral(std::size_t mesh_item_id, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, GlobalVector const &x, double const t, double &crack_volume, CoupledSolutionsForStaggeredScheme const *const cpl_xs) override
PhaseFieldLocalAssembler(MeshLib::Element const &e, std::size_t const, bool const is_axially_symmetric, unsigned const integration_order, PhaseFieldProcessData< DisplacementDim > &process_data)
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
typename ShapeMatricesType::template VectorType< phasefield_size > PhaseFieldVector
typename ShapeMatricesType::template MatrixType< phasefield_size, phasefield_size > PhaseFieldMatrix
PhaseFieldLocalAssembler(PhaseFieldLocalAssembler const &)=delete
void postTimestepConcrete(Eigen::VectorXd const &, double const, double const) override
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
typename ShapeMatricesType::NodalVectorType NodalVectorType
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
PhaseFieldLocalAssembler(PhaseFieldLocalAssembler &&)=delete
void assembleWithJacobianForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
void computeEnergy(std::size_t mesh_item_id, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, GlobalVector const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work, CoupledSolutionsForStaggeredScheme const *const cpl_xs) override
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
PhaseFieldProcessData< DisplacementDim > & _process_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)
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)
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double > calculateDegradedStress(double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
MechanicsBase< DisplacementDim > & selectSolidConstitutiveRelation(std::map< int, std::unique_ptr< MechanicsBase< DisplacementDim >>> const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Definition: KelvinVector.h:56
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)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
VectorType< ShapeFunction::NPOINTS > NodalVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N
Definition: SecondaryData.h:28
BMatricesType::KelvinVectorType sigma_tensile
Definition: PhaseFieldFEM.h:51
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material
Definition: PhaseFieldFEM.h:55
ShapeMatrixType::NodalRowVectorType N
Definition: PhaseFieldFEM.h:46
BMatricesType::KelvinVectorType eps_prev
Definition: PhaseFieldFEM.h:49
void updateConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x, double const, DisplacementVectorType const &, double const degradation)
Definition: PhaseFieldFEM.h:73
BMatricesType::KelvinMatrixType C_compressive
Definition: PhaseFieldFEM.h:60
BMatricesType::KelvinVectorType sigma_compressive
Definition: PhaseFieldFEM.h:51
ShapeMatrixType::GlobalDimNodalMatrixType dNdx
Definition: PhaseFieldFEM.h:47
BMatricesType::KelvinVectorType eps
Definition: PhaseFieldFEM.h:49
BMatricesType::KelvinVectorType sigma
Definition: PhaseFieldFEM.h:52
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
Definition: PhaseFieldFEM.h:37
BMatricesType::KelvinMatrixType C_tensile
Definition: PhaseFieldFEM.h:60
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
Definition: PhaseFieldFEM.h:58
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N