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

Detailed Description

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

Definition at line 170 of file PhaseFieldFEM.h.

#include <PhaseFieldFEM.h>

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

Public Types

using ShapeMatricesType
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>
 
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType
 
using DeformationVector
 
using PhaseFieldVector
 
using PhaseFieldMatrix
 
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
 
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
 
using IpData
 

Public Member Functions

 PhaseFieldLocalAssembler (PhaseFieldLocalAssembler const &)=delete
 
 PhaseFieldLocalAssembler (PhaseFieldLocalAssembler &&)=delete
 
 PhaseFieldLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, PhaseFieldProcessData< DisplacementDim > &process_data)
 
std::size_t setIPDataInitialConditions (std::string_view const name, double const *values, int const integration_order) override
 Returns number of read integration points.
 
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_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
 
void initializeConcrete () override
 
void postTimestepConcrete (Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
 
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
 
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
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
 
std::vector< double > getSigma () const override
 
std::vector< double > getEpsilon () const override
 
- Public Member Functions inherited from ProcessLib::PhaseField::PhaseFieldLocalAssemblerInterface
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, 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_x_prev, 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_x_prev, 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_prev, 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, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
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.
 
- 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 & getIntPtSigmaTensile (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtSigmaCompressive (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
 
std::vector< double > const & getIntPtEpsilonTensile (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
 
NumLib::GenericIntegrationMethod const & _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 , int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 186 of file PhaseFieldFEM.h.

◆ DeformationVector

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

Definition at line 190 of file PhaseFieldFEM.h.

◆ IpData

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

Definition at line 199 of file PhaseFieldFEM.h.

◆ NodalForceVectorType

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

Definition at line 188 of file PhaseFieldFEM.h.

◆ NodalMatrixType

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

Definition at line 197 of file PhaseFieldFEM.h.

◆ NodalVectorType

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

Definition at line 198 of file PhaseFieldFEM.h.

◆ PhaseFieldMatrix

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

Definition at line 194 of file PhaseFieldFEM.h.

◆ PhaseFieldVector

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

Definition at line 192 of file PhaseFieldFEM.h.

◆ ShapeMatrices

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

Definition at line 185 of file PhaseFieldFEM.h.

◆ ShapeMatricesType

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

Constructor & Destructor Documentation

◆ PhaseFieldLocalAssembler() [1/3]

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

◆ PhaseFieldLocalAssembler() [2/3]

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

◆ PhaseFieldLocalAssembler() [3/3]

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

Definition at line 205 of file PhaseFieldFEM.h.

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 }
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)
PhaseFieldProcessData< DisplacementDim > & _process_data
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
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.
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

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

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, 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 271 of file PhaseFieldFEM.h.

277 {
278 OGS_FATAL(
279 "PhaseFieldLocalAssembler: assembly without jacobian is not "
280 "implemented.");
281 }
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ assembleWithJacobianForDeformationEquations()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, 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 43 of file PhaseFieldFEM-impl.h.

47{
48 using DeformationMatrix =
49 typename ShapeMatricesType::template MatrixType<displacement_size,
51 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
52 auto const u =
53 local_x.template segment<displacement_size>(displacement_index);
54
56 local_Jac_data, displacement_size, displacement_size);
57
59 local_b_data, displacement_size);
60
62 x_position.setElementID(_element.getID());
63
64 auto local_pressure = 0.0;
65 if (_process_data.pressurized_crack)
66 {
67 local_pressure = _process_data.unity_pressure;
68 }
69
70 int const n_integration_points = _integration_method.getNumberOfPoints();
71 for (int ip = 0; ip < n_integration_points; ip++)
72 {
73 x_position.setIntegrationPoint(ip);
74 auto const& w = _ip_data[ip].integration_weight;
75 auto const& N = _ip_data[ip].N;
76 auto const& dNdx = _ip_data[ip].dNdx;
77
78 auto const x_coord =
80 _element, N);
81
82 auto const& B =
83 LinearBMatrix::computeBMatrix<DisplacementDim,
84 ShapeFunction::NPOINTS,
86 dNdx, N, x_coord, _is_axially_symmetric);
87
88 auto& eps = _ip_data[ip].eps;
89 eps.noalias() = B * u;
90 double const k = _process_data.residual_stiffness(t, x_position)[0];
91 double const ls = _process_data.crack_length_scale(t, x_position)[0];
92 double const d_ip = N.dot(d);
93 double const degradation =
94 _process_data.degradation_derivative->degradation(d_ip, k, ls);
95 _ip_data[ip].updateConstitutiveRelation(
96 t, x_position, dt, u, degradation,
97 _process_data.energy_split_model);
98
99 auto& sigma = _ip_data[ip].sigma;
100 auto const& D = _ip_data[ip].D;
101
102 typename ShapeMatricesType::template MatrixType<DisplacementDim,
104 N_u = ShapeMatricesType::template MatrixType<
105 DisplacementDim, displacement_size>::Zero(DisplacementDim,
107
108 for (int i = 0; i < DisplacementDim; ++i)
109 {
110 N_u.template block<1, displacement_size / DisplacementDim>(
111 i, i * displacement_size / DisplacementDim)
112 .noalias() = N;
113 }
114
115 auto const rho_sr = _process_data.solid_density(t, x_position)[0];
116 auto const& b = _process_data.specific_body_force;
117
118 local_rhs.noalias() -=
119 (B.transpose() * sigma - N_u.transpose() * rho_sr * b -
120 local_pressure * N_u.transpose() * dNdx * d) *
121 w;
122 local_Jac.noalias() += B.transpose() * D * B * w;
123 }
124}
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
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.

References ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), NumLib::interpolateXCoordinate(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

◆ assembleWithJacobianForStaggeredScheme()

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

Reimplemented from ProcessLib::LocalAssemblerInterface.

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

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

◆ assembleWithJacobianPhaseFieldEquations()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, 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 127 of file PhaseFieldFEM-impl.h.

132{
133 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
134 auto const u =
135 local_x.template segment<displacement_size>(displacement_index);
136
138 local_Jac_data, phasefield_size, phasefield_size);
140 local_b_data, phasefield_size);
141
143 x_position.setElementID(_element.getID());
144
145 auto local_pressure = 0.0;
146 if (_process_data.static_pressurized_crack)
147 {
148 local_pressure = _process_data.unity_pressure;
149 }
150 else if (_process_data.propagating_pressurized_crack)
151 {
152 local_pressure = _process_data.pressure;
153 }
154
155 double const k = _process_data.residual_stiffness(t, x_position)[0];
156 double const ls = _process_data.crack_length_scale(t, x_position)[0];
157 double const gc = _process_data.crack_resistance(t, x_position)[0];
158
159 int const n_integration_points = _integration_method.getNumberOfPoints();
160 for (int ip = 0; ip < n_integration_points; ip++)
161 {
162 x_position.setIntegrationPoint(ip);
163 auto const& w = _ip_data[ip].integration_weight;
164 auto const& N = _ip_data[ip].N;
165 auto const& dNdx = _ip_data[ip].dNdx;
166
167 double const d_ip = N.dot(d);
168 double const degradation =
169 _process_data.degradation_derivative->degradation(d_ip, k, ls);
170 double const degradation_df1 =
171 _process_data.degradation_derivative->degradation_df1(d_ip, ls);
172 double const degradation_df2 =
173 _process_data.degradation_derivative->degradation_df2(d_ip, ls);
174 auto const x_coord =
176 _element, N);
177 auto const& B =
178 LinearBMatrix::computeBMatrix<DisplacementDim,
179 ShapeFunction::NPOINTS,
181 dNdx, N, x_coord, _is_axially_symmetric);
182
183 auto& eps = _ip_data[ip].eps;
184 eps.noalias() = B * u;
185 _ip_data[ip].updateConstitutiveRelation(
186 t, x_position, dt, u, degradation,
187 _process_data.energy_split_model);
188
189 auto const& strain_energy_tensile = _ip_data[ip].strain_energy_tensile;
190
191 auto& ip_data = _ip_data[ip];
192 ip_data.strain_energy_tensile = strain_energy_tensile;
193
194 typename ShapeMatricesType::template MatrixType<DisplacementDim,
196 N_u = ShapeMatricesType::template MatrixType<
197 DisplacementDim, displacement_size>::Zero(DisplacementDim,
199
200 for (int i = 0; i < DisplacementDim; ++i)
201 {
202 N_u.template block<1, displacement_size / DisplacementDim>(
203 i, i * displacement_size / DisplacementDim)
204 .noalias() = N;
205 }
206
207 local_Jac.noalias() +=
208 (N.transpose() * N * degradation_df2 * strain_energy_tensile) * w;
209
210 local_rhs.noalias() -=
211 (N.transpose() * degradation_df1 * strain_energy_tensile -
212 local_pressure * dNdx.transpose() * N_u * u) *
213 w;
214
215 switch (_process_data.phasefield_model)
216 {
218 {
219 auto const local_Jac_AT1 =
220 (gc * 0.75 * dNdx.transpose() * dNdx * ls * w).eval();
221 local_Jac.noalias() += local_Jac_AT1;
222
223 local_rhs.noalias() -=
224 gc * (-0.375 * N.transpose() / ls) * w + local_Jac_AT1 * d;
225 break;
226 }
228 {
229 auto const local_Jac_AT2 =
230 (gc *
231 (N.transpose() * N / ls + dNdx.transpose() * dNdx * ls) *
232 w)
233 .eval();
234 local_Jac.noalias() += local_Jac_AT2;
235
236 local_rhs.noalias() -=
237 local_Jac_AT2 * d - gc * (N.transpose() / ls) * w;
238 break;
239 }
241 {
242 auto const local_Jac_COHESIVE =
243 (2.0 / std::numbers::pi * gc *
244 (-N.transpose() * N / ls + dNdx.transpose() * dNdx * ls) *
245 w)
246 .eval();
247
248 local_Jac.noalias() += local_Jac_COHESIVE;
249
250 local_rhs.noalias() -= local_Jac_COHESIVE * d;
251 break;
252 }
253 }
254 }
255}

References ProcessLib::PhaseField::AT1, ProcessLib::PhaseField::AT2, ProcessLib::PhaseField::COHESIVE, ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), NumLib::interpolateXCoordinate(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

◆ computeCrackIntegral()

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

Implements ProcessLib::PhaseField::PhaseFieldLocalAssemblerInterface.

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

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

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

◆ computeEnergy()

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

Implements ProcessLib::PhaseField::PhaseFieldLocalAssemblerInterface.

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

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

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

◆ getEpsilon()

template<typename ShapeFunctionDisplacement , int DisplacementDim>
std::vector< double > ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunctionDisplacement, DisplacementDim >::getEpsilon ( ) const
overridevirtual

Implements ProcessLib::PhaseField::PhaseFieldLocalAssemblerInterface.

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

465{
466 auto const kelvin_vector_size =
468 unsigned const n_integration_points =
470
471 std::vector<double> ip_epsilon_values;
472 auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
473 double, Eigen::Dynamic, kelvin_vector_size, Eigen::RowMajor>>(
474 ip_epsilon_values, n_integration_points, kelvin_vector_size);
475
476 for (unsigned ip = 0; ip < n_integration_points; ++ip)
477 {
478 auto const& eps = _ip_data[ip].eps;
479 cache_mat.row(ip) =
481 }
482
483 return ip_epsilon_values;
484}
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)

References MathLib::createZeroedMatrix(), MathLib::KelvinVector::kelvin_vector_dimensions(), and MathLib::KelvinVector::kelvinVectorToSymmetricTensor().

◆ getIntPtEpsilon()

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

◆ getIntPtEpsilonTensile()

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

◆ getIntPtSigma()

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

◆ getIntPtSigmaCompressive()

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

◆ getIntPtSigmaTensile()

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

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 347 of file PhaseFieldFEM.h.

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 }

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

◆ getSigma()

template<typename ShapeFunctionDisplacement , int DisplacementDim>
std::vector< double > ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunctionDisplacement, DisplacementDim >::getSigma ( ) const
overridevirtual

◆ initializeConcrete()

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

Set initial stress from parameter.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 289 of file PhaseFieldFEM.h.

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 }
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)

References ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_process_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::interpolateCoordinates(), and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ postTimestepConcrete()

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

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 320 of file PhaseFieldFEM.h.

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 }

References ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, and NumLib::GenericIntegrationMethod::getNumberOfPoints().

◆ setIPDataInitialConditions()

template<typename ShapeFunctionDisplacement , int DisplacementDim>
std::size_t ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunctionDisplacement, DisplacementDim >::setIPDataInitialConditions ( std::string_view const name,
double const * values,
int const integration_order )
overridevirtual

Returns number of read integration points.

Implements ProcessLib::PhaseField::PhaseFieldLocalAssemblerInterface.

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

425{
426 if (integration_order !=
427 static_cast<int>(_integration_method.getIntegrationOrder()))
428 {
429 OGS_FATAL(
430 "Setting integration point initial conditions; The integration "
431 "order of the local assembler for element {:d} is different from "
432 "the integration order in the initial condition.",
433 _element.getID());
434 }
435
436 if (name == "sigma")
437 {
438 if (_process_data.initial_stress.value != nullptr)
439 {
440 OGS_FATAL(
441 "Setting initial conditions for stress from integration "
442 "point data and from a parameter '{:s}' is not possible "
443 "simultaneously.",
444 _process_data.initial_stress.value->name);
445 }
446
448 values, _ip_data, &IpData::sigma);
449 }
450
451 return 0;
452}
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)

References OGS_FATAL, and ProcessLib::setIntegrationPointKelvinVectorData().

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

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

Definition at line 429 of file PhaseFieldFEM.h.

◆ _process_data

◆ _secondary_data

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

◆ displacement_index

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

Definition at line 173 of file PhaseFieldFEM.h.

◆ displacement_size

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

Definition at line 174 of file PhaseFieldFEM.h.

◆ mechanics_process_id

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

Definition at line 432 of file PhaseFieldFEM.h.

◆ phase_process_id

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

Definition at line 431 of file PhaseFieldFEM.h.

◆ phasefield_index

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

Definition at line 176 of file PhaseFieldFEM.h.

◆ phasefield_size

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

Definition at line 178 of file PhaseFieldFEM.h.


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