OGS
ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, typename IntegrationMethod, int DisplacementDim>
class ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >

Definition at line 111 of file PhaseFieldFEM.h.

#include <PhaseFieldFEM.h>

Inheritance diagram for ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesType = ShapeMatrixPolicyType< ShapeFunction, DisplacementDim >
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using BMatricesType = BMatrixPolicyType< ShapeFunction, DisplacementDim >
 
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType
 
using DeformationVector = typename ShapeMatricesType::template VectorType< displacement_size >
 
using PhaseFieldVector = typename ShapeMatricesType::template VectorType< phasefield_size >
 
using PhaseFieldMatrix = typename ShapeMatricesType::template MatrixType< phasefield_size, phasefield_size >
 
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
 
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
 
using IpData = IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >
 

Public Member Functions

 PhaseFieldLocalAssembler (PhaseFieldLocalAssembler const &)=delete
 
 PhaseFieldLocalAssembler (PhaseFieldLocalAssembler &&)=delete
 
 PhaseFieldLocalAssembler (MeshLib::Element const &e, std::size_t const, bool const is_axially_symmetric, unsigned const integration_order, 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 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
 
void initializeConcrete () override
 
void postTimestepConcrete (Eigen::VectorXd const &, double const, double 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
 
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
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point. More...
 
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
void setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, int const process_id)
 
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
 
virtual void preAssemble (double const, double const, std::vector< double > const &)
 
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double dxdot_dx, const double dx_dx, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
 
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
 
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
 
void postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
 
virtual std::vector< double > interpolateNodalValuesToIntegrationPoints (std::vector< double > const &)
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const
 Fits to staggered scheme. More...
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Private Member Functions

std::vector< double > const & getIntPtSigma (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilon (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
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)
 

Private Attributes

PhaseFieldProcessData< DisplacementDim > & _process_data
 
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
 
IntegrationMethod _integration_method
 
MeshLib::Element const & _element
 
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
 
bool const _is_axially_symmetric
 

Static Private Attributes

static constexpr int displacement_index = 0
 
static constexpr int displacement_size
 
static constexpr int phasefield_index
 
static constexpr int phasefield_size = ShapeFunction::NPOINTS
 
static const int phase_process_id = 1
 
static const int mechanics_process_id = 0
 

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 127 of file PhaseFieldFEM.h.

◆ DeformationVector

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::DeformationVector = typename ShapeMatricesType::template VectorType<displacement_size>

Definition at line 131 of file PhaseFieldFEM.h.

◆ IpData

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::IpData = IntegrationPointData<BMatricesType, ShapeMatricesType, DisplacementDim>

Definition at line 140 of file PhaseFieldFEM.h.

◆ NodalForceVectorType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType

Definition at line 129 of file PhaseFieldFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 138 of file PhaseFieldFEM.h.

◆ NodalVectorType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 139 of file PhaseFieldFEM.h.

◆ PhaseFieldMatrix

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::PhaseFieldMatrix = typename ShapeMatricesType::template MatrixType<phasefield_size, phasefield_size>

Definition at line 135 of file PhaseFieldFEM.h.

◆ PhaseFieldVector

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::PhaseFieldVector = typename ShapeMatricesType::template VectorType<phasefield_size>

Definition at line 133 of file PhaseFieldFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 126 of file PhaseFieldFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 122 of file PhaseFieldFEM.h.

Constructor & Destructor Documentation

◆ PhaseFieldLocalAssembler() [1/3]

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::PhaseFieldLocalAssembler ( PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim > const &  )
delete

◆ PhaseFieldLocalAssembler() [2/3]

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::PhaseFieldLocalAssembler ( PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim > &&  )
delete

◆ PhaseFieldLocalAssembler() [3/3]

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::PhaseFieldLocalAssembler ( MeshLib::Element const &  e,
std::size_t const  ,
bool const  is_axially_symmetric,
unsigned const  integration_order,
PhaseFieldProcessData< DisplacementDim > &  process_data 
)
inline

Definition at line 146 of file PhaseFieldFEM.h.

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  }
#define OGS_FATAL(...)
Definition: Error.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)
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
PhaseFieldProcessData< DisplacementDim > & _process_data
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)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23
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::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N
Definition: SecondaryData.h:28

References ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_element, ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_integration_method, ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data, ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_process_data, ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_secondary_data, MeshLib::Element::getID(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N, OGS_FATAL, MaterialLib::Solids::selectSolidConstitutiveRelation(), and ParameterLib::SpatialPosition::setElementID().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assemble ( double const  ,
double const  ,
std::vector< double > const &  ,
std::vector< double > const &  ,
std::vector< double > &  ,
std::vector< double > &  ,
std::vector< double > &   
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 213 of file PhaseFieldFEM.h.

219  {
220  OGS_FATAL(
221  "PhaseFieldLocalAssembler: assembly without jacobian is not "
222  "implemented.");
223  }

References OGS_FATAL.

◆ assembleWithJacobianForDeformationEquations()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobianForDeformationEquations ( double const  t,
double const  dt,
Eigen::VectorXd const &  local_x,
std::vector< double > &  local_b_data,
std::vector< double > &  local_Jac_data 
)
private

Definition at line 48 of file PhaseFieldFEM-impl.h.

52 {
53  using DeformationMatrix =
54  typename ShapeMatricesType::template MatrixType<displacement_size,
56  auto const d = local_x.template segment<phasefield_size>(phasefield_index);
57  auto const u =
58  local_x.template segment<displacement_size>(displacement_index);
59 
60  auto local_Jac = MathLib::createZeroedMatrix<DeformationMatrix>(
61  local_Jac_data, displacement_size, displacement_size);
62 
63  auto local_rhs = MathLib::createZeroedVector<DeformationVector>(
64  local_b_data, displacement_size);
65 
67  x_position.setElementID(_element.getID());
68 
69  auto local_pressure = 0.0;
70  if (_process_data.crack_pressure)
71  {
72  local_pressure = _process_data.unity_pressure;
73  }
74 
75  int const n_integration_points = _integration_method.getNumberOfPoints();
76  for (int ip = 0; ip < n_integration_points; ip++)
77  {
78  x_position.setIntegrationPoint(ip);
79  auto const& w = _ip_data[ip].integration_weight;
80  auto const& N = _ip_data[ip].N;
81  auto const& dNdx = _ip_data[ip].dNdx;
82 
83  auto const x_coord =
84  NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
85  _element, N);
86 
87  auto const& B =
88  LinearBMatrix::computeBMatrix<DisplacementDim,
89  ShapeFunction::NPOINTS,
91  dNdx, N, x_coord, _is_axially_symmetric);
92 
93  auto& eps = _ip_data[ip].eps;
94  eps.noalias() = B * u;
95  double const k = _process_data.residual_stiffness(t, x_position)[0];
96  double const d_ip = N.dot(d);
97  double const degradation = d_ip * d_ip * (1 - k) + k;
98  _ip_data[ip].updateConstitutiveRelation(t, x_position, dt, u,
99  degradation);
100 
101  auto& sigma = _ip_data[ip].sigma;
102  auto const& C_tensile = _ip_data[ip].C_tensile;
103  auto const& C_compressive = _ip_data[ip].C_compressive;
104 
105  typename ShapeMatricesType::template MatrixType<DisplacementDim,
107  N_u = ShapeMatricesType::template MatrixType<
108  DisplacementDim, displacement_size>::Zero(DisplacementDim,
110 
111  for (int i = 0; i < DisplacementDim; ++i)
112  {
113  N_u.template block<1, displacement_size / DisplacementDim>(
114  i, i * displacement_size / DisplacementDim)
115  .noalias() = N;
116  }
117 
118  auto const rho_sr = _process_data.solid_density(t, x_position)[0];
119  auto const& b = _process_data.specific_body_force;
120 
121  local_rhs.noalias() -=
122  (B.transpose() * sigma - N_u.transpose() * rho_sr * b -
123  local_pressure * N_u.transpose() * dNdx * d) *
124  w;
125 
126  local_Jac.noalias() +=
127  B.transpose() * (degradation * C_tensile + C_compressive) * B * w;
128  }
129 }
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
Definition: BMatrixPolicy.h:55
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
Definition: LinearBMatrix.h:42

References ProcessLib::LinearBMatrix::computeBMatrix(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

◆ assembleWithJacobianForStaggeredScheme()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::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 
)
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 23 of file PhaseFieldFEM-impl.h.

31 {
32  // For the equations with phase field.
33  if (process_id == phase_process_id)
34  {
35  assembleWithJacobianPhaseFieldEquations(t, dt, local_x, local_b_data,
36  local_Jac_data);
37  return;
38  }
39 
40  // For the equations with deformation
41  assembleWithJacobianForDeformationEquations(t, dt, local_x, local_b_data,
42  local_Jac_data);
43 }
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)

◆ assembleWithJacobianPhaseFieldEquations()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobianPhaseFieldEquations ( double const  t,
double const  dt,
Eigen::VectorXd const &  local_x,
std::vector< double > &  local_b_data,
std::vector< double > &  local_Jac_data 
)
private

Definition at line 134 of file PhaseFieldFEM-impl.h.

139 {
140  auto const d = local_x.template segment<phasefield_size>(phasefield_index);
141  auto const u =
142  local_x.template segment<displacement_size>(displacement_index);
143 
144  auto local_Jac = MathLib::createZeroedMatrix<PhaseFieldMatrix>(
145  local_Jac_data, phasefield_size, phasefield_size);
146  auto local_rhs = MathLib::createZeroedVector<PhaseFieldVector>(
147  local_b_data, phasefield_size);
148 
150  x_position.setElementID(_element.getID());
151 
152  auto local_pressure = 0.0;
153  if (_process_data.crack_pressure)
154  {
155  local_pressure = _process_data.unity_pressure;
156  }
157  else if (_process_data.hydro_crack)
158  {
159  local_pressure = _process_data.pressure;
160  }
161 
162  int const n_integration_points = _integration_method.getNumberOfPoints();
163  for (int ip = 0; ip < n_integration_points; ip++)
164  {
165  x_position.setIntegrationPoint(ip);
166  auto const& w = _ip_data[ip].integration_weight;
167  auto const& N = _ip_data[ip].N;
168  auto const& dNdx = _ip_data[ip].dNdx;
169 
170  double const gc = _process_data.crack_resistance(t, x_position)[0];
171  double const ls = _process_data.crack_length_scale(t, x_position)[0];
172 
173  double const k = _process_data.residual_stiffness(t, x_position)[0];
174  double const d_ip = N.dot(d);
175  double const degradation = d_ip * d_ip * (1 - k) + k;
176  auto const x_coord =
177  NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
178  _element, N);
179  auto const& B =
180  LinearBMatrix::computeBMatrix<DisplacementDim,
181  ShapeFunction::NPOINTS,
182  typename BMatricesType::BMatrixType>(
183  dNdx, N, x_coord, _is_axially_symmetric);
184 
185  auto& eps = _ip_data[ip].eps;
186  eps.noalias() = B * u;
187  _ip_data[ip].updateConstitutiveRelation(t, x_position, dt, u,
188  degradation);
189 
190  auto const& strain_energy_tensile = _ip_data[ip].strain_energy_tensile;
191 
192  auto& ip_data = _ip_data[ip];
193  ip_data.strain_energy_tensile = strain_energy_tensile;
194 
195  typename ShapeMatricesType::template MatrixType<DisplacementDim,
197  N_u = ShapeMatricesType::template MatrixType<
198  DisplacementDim, displacement_size>::Zero(DisplacementDim,
200 
201  for (int i = 0; i < DisplacementDim; ++i)
202  {
203  N_u.template block<1, displacement_size / DisplacementDim>(
204  i, i * displacement_size / DisplacementDim)
205  .noalias() = N;
206  }
207 
208  local_Jac.noalias() +=
209  (2 * N.transpose() * N * strain_energy_tensile) * w;
210 
211  local_rhs.noalias() -=
212  (N.transpose() * N * d * 2 * strain_energy_tensile -
213  local_pressure * dNdx.transpose() * N_u * u) *
214  w;
215 
216  switch (_process_data.phasefield_model)
217  {
219  {
220  local_Jac.noalias() +=
221  gc * 0.75 * dNdx.transpose() * dNdx * ls * w;
222 
223  local_rhs.noalias() -=
224  gc *
225  (-0.375 * N.transpose() / ls +
226  0.75 * dNdx.transpose() * dNdx * ls * d) *
227  w;
228  break;
229  }
231  {
232  local_Jac.noalias() +=
233  gc *
234  (N.transpose() * N / ls + dNdx.transpose() * dNdx * ls) * w;
235 
236  local_rhs.noalias() -=
237  gc *
238  ((N.transpose() * N / ls + dNdx.transpose() * dNdx * ls) *
239  d -
240  N.transpose() / ls) *
241  w;
242  break;
243  }
244  }
245  }
246 }

References ProcessLib::PhaseField::AT1, ProcessLib::PhaseField::AT2, ProcessLib::LinearBMatrix::computeBMatrix(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

◆ computeCrackIntegral()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::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 
)
overridevirtual

Implements ProcessLib::PhaseField::PhaseFieldLocalAssemblerInterface.

Definition at line 251 of file PhaseFieldFEM-impl.h.

258 {
259  assert(cpl_xs != nullptr);
260 
261  std::vector<std::vector<GlobalIndexType>> indices_of_processes;
262  indices_of_processes.reserve(dof_tables.size());
263  std::transform(dof_tables.begin(), dof_tables.end(),
264  std::back_inserter(indices_of_processes),
265  [&](NumLib::LocalToGlobalIndexMap const& dof_table)
266  { return NumLib::getIndices(mesh_item_id, dof_table); });
267 
268  auto local_coupled_xs =
269  getCoupledLocalSolutions(cpl_xs->coupled_xs, indices_of_processes);
270  assert(local_coupled_xs.size() == displacement_size + phasefield_size);
271 
272  auto const d = Eigen::Map<PhaseFieldVector const>(
273  &local_coupled_xs[phasefield_index], phasefield_size);
274  auto const u = Eigen::Map<DeformationVector const>(
275  &local_coupled_xs[displacement_index], displacement_size);
276 
278  x_position.setElementID(_element.getID());
279 
280  int const n_integration_points = _integration_method.getNumberOfPoints();
281  for (int ip = 0; ip < n_integration_points; ip++)
282  {
283  x_position.setIntegrationPoint(ip);
284  auto const& w = _ip_data[ip].integration_weight;
285  auto const& N = _ip_data[ip].N;
286  auto const& dNdx = _ip_data[ip].dNdx;
287 
288  typename ShapeMatricesType::template MatrixType<DisplacementDim,
290  N_u = ShapeMatricesType::template MatrixType<
291  DisplacementDim, displacement_size>::Zero(DisplacementDim,
293 
294  for (int i = 0; i < DisplacementDim; ++i)
295  {
296  N_u.template block<1, displacement_size / DisplacementDim>(
297  i, i * displacement_size / DisplacementDim)
298  .noalias() = N;
299  }
300 
301  crack_volume += (N_u * u).dot(dNdx * d) * w;
302  }
303 }
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType >> const &indices)

References ProcessLib::CoupledSolutionsForStaggeredScheme::coupled_xs, ProcessLib::getCoupledLocalSolutions(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

◆ computeEnergy()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::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 
)
overridevirtual

Implements ProcessLib::PhaseField::PhaseFieldLocalAssemblerInterface.

Definition at line 308 of file PhaseFieldFEM-impl.h.

316 {
317  assert(cpl_xs != nullptr);
318 
319  std::vector<std::vector<GlobalIndexType>> indices_of_processes;
320  indices_of_processes.reserve(dof_tables.size());
321  std::transform(dof_tables.begin(), dof_tables.end(),
322  std::back_inserter(indices_of_processes),
323  [&](NumLib::LocalToGlobalIndexMap const& dof_table)
324  { return NumLib::getIndices(mesh_item_id, dof_table); });
325 
326  auto const local_coupled_xs =
327  getCoupledLocalSolutions(cpl_xs->coupled_xs, indices_of_processes);
328  assert(local_coupled_xs.size() == displacement_size + phasefield_size);
329 
330  auto const d = Eigen::Map<PhaseFieldVector const>(
331  &local_coupled_xs[phasefield_index], phasefield_size);
332  auto const u = Eigen::Map<DeformationVector const>(
333  &local_coupled_xs[displacement_index], displacement_size);
334 
336  x_position.setElementID(_element.getID());
337 
338  double element_elastic_energy = 0.0;
339  double element_surface_energy = 0.0;
340  double element_pressure_work = 0.0;
341 
342  int const n_integration_points = _integration_method.getNumberOfPoints();
343  for (int ip = 0; ip < n_integration_points; ip++)
344  {
345  x_position.setIntegrationPoint(ip);
346  auto const& w = _ip_data[ip].integration_weight;
347  auto const& N = _ip_data[ip].N;
348  auto const& dNdx = _ip_data[ip].dNdx;
349  double const d_ip = N.dot(d);
350  auto pressure_ip = _process_data.pressure;
351  double const gc = _process_data.crack_resistance(t, x_position)[0];
352  double const ls = _process_data.crack_length_scale(t, x_position)[0];
353 
354  typename ShapeMatricesType::template MatrixType<DisplacementDim,
356  N_u = ShapeMatricesType::template MatrixType<
357  DisplacementDim, displacement_size>::Zero(DisplacementDim,
359 
360  for (int i = 0; i < DisplacementDim; ++i)
361  {
362  N_u.template block<1, displacement_size / DisplacementDim>(
363  i, i * displacement_size / DisplacementDim)
364  .noalias() = N;
365  }
366 
367  element_elastic_energy += _ip_data[ip].elastic_energy * w;
368 
369  switch (_process_data.phasefield_model)
370  {
372  {
373  element_surface_energy +=
374  gc * 0.375 *
375  ((1 - d_ip) / ls + (dNdx * d).dot((dNdx * d)) * ls) * w;
376 
377  break;
378  }
380  {
381  element_surface_energy += 0.5 * gc *
382  ((1 - d_ip) * (1 - d_ip) / ls +
383  (dNdx * d).dot((dNdx * d)) * ls) *
384  w;
385  break;
386  }
387  }
388 
389  if (_process_data.crack_pressure)
390  {
391  element_pressure_work += pressure_ip * (N_u * u).dot(dNdx * d) * w;
392  }
393  }
394 
395 #ifdef USE_PETSC
396  int const n_all_nodes = indices_of_processes[1].size();
397  int const n_regular_nodes = std::count_if(
398  begin(indices_of_processes[1]), end(indices_of_processes[1]),
399  [](GlobalIndexType const& index) { return index >= 0; });
400  if (n_all_nodes != n_regular_nodes)
401  {
402  element_elastic_energy *=
403  static_cast<double>(n_regular_nodes) / n_all_nodes;
404  element_surface_energy *=
405  static_cast<double>(n_regular_nodes) / n_all_nodes;
406  element_pressure_work *=
407  static_cast<double>(n_regular_nodes) / n_all_nodes;
408  }
409 #endif // USE_PETSC
410  elastic_energy += element_elastic_energy;
411  surface_energy += element_surface_energy;
412  pressure_work += element_pressure_work;
413 }
GlobalMatrix::IndexType GlobalIndexType

References ProcessLib::PhaseField::AT1, ProcessLib::PhaseField::AT2, ProcessLib::CoupledSolutionsForStaggeredScheme::coupled_xs, ProcessLib::getCoupledLocalSolutions(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

◆ getIntPtEpsilon()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector<double> const& ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilon ( const double  ,
std::vector< GlobalVector * > const &  ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  ,
std::vector< double > &  cache 
) const
inlineoverrideprivatevirtual

◆ getIntPtSigma()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector<double> const& ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigma ( const double  ,
std::vector< GlobalVector * > const &  ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  ,
std::vector< double > &  cache 
) const
inlineoverrideprivatevirtual

◆ getShapeMatrix()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
Eigen::Map<const Eigen::RowVectorXd> ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getShapeMatrix ( const unsigned  integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 273 of file PhaseFieldFEM.h.

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  }

References ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_secondary_data, and ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N.

◆ initializeConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 233 of file PhaseFieldFEM.h.

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  }

References ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_integration_method, and ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data.

◆ postTimestepConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::postTimestepConcrete ( Eigen::VectorXd const &  ,
double const  ,
double const   
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 244 of file PhaseFieldFEM.h.

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  }

References ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_integration_method, and ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data.

Member Data Documentation

◆ _element

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
MeshLib::Element const& ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_element
private

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
bool const ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_is_axially_symmetric
private

Definition at line 318 of file PhaseFieldFEM.h.

◆ _process_data

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
PhaseFieldProcessData<DisplacementDim>& ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_process_data
private

◆ _secondary_data

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
SecondaryData<typename ShapeMatrices::ShapeType> ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::_secondary_data
private

◆ displacement_index

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
constexpr int ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::displacement_index = 0
staticconstexprprivate

Definition at line 114 of file PhaseFieldFEM.h.

◆ displacement_size

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
constexpr int ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::displacement_size
staticconstexprprivate
Initial value:
=
ShapeFunction::NPOINTS * DisplacementDim

Definition at line 115 of file PhaseFieldFEM.h.

◆ mechanics_process_id

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
const int ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::mechanics_process_id = 0
staticprivate

Definition at line 321 of file PhaseFieldFEM.h.

◆ phase_process_id

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
const int ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::phase_process_id = 1
staticprivate

Definition at line 320 of file PhaseFieldFEM.h.

◆ phasefield_index

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
constexpr int ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::phasefield_index
staticconstexprprivate
Initial value:

Definition at line 117 of file PhaseFieldFEM.h.

◆ phasefield_size

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
constexpr int ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::phasefield_size = ShapeFunction::NPOINTS
staticconstexprprivate

Definition at line 119 of file PhaseFieldFEM.h.


The documentation for this class was generated from the following files: