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 100 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::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.
virtual int getNumberOfVectorElementsForDeformation () const
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 116 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 120 of file PhaseFieldFEM.h.

◆ IpData

◆ NodalForceVectorType

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

Definition at line 118 of file PhaseFieldFEM.h.

◆ NodalMatrixType

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

Definition at line 127 of file PhaseFieldFEM.h.

◆ NodalVectorType

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

Definition at line 128 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 124 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 122 of file PhaseFieldFEM.h.

◆ ShapeMatrices

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

Definition at line 115 of file PhaseFieldFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 111 of file PhaseFieldFEM.h.

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 135 of file PhaseFieldFEM.h.

143 _element(e),
145 {
146 unsigned const n_integration_points =
147 _integration_method.getNumberOfPoints();
148
151
152 auto& solid_material =
154 _process_data.solid_materials,
155 _process_data.material_ids,
156 e.getID());
157
158 auto const shape_matrices =
162
164 x_position.setElementID(_element.getID());
165
166 for (unsigned ip = 0; ip < n_integration_points; ip++)
167 {
168 _ip_data.emplace_back(solid_material);
169 auto& ip_data = _ip_data[ip];
170 ip_data.integration_weight =
171 _integration_method.getWeightedPoint(ip).getWeight() *
172 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
173
174 static const int kelvin_vector_size =
177 ip_data.eps_tensile.setZero(kelvin_vector_size);
178 ip_data.eps.setZero(kelvin_vector_size);
179 ip_data.eps_prev.resize(kelvin_vector_size);
181
182 ip_data.sigma_tensile.setZero(kelvin_vector_size);
183 ip_data.sigma_compressive.setZero(kelvin_vector_size);
184 ip_data.sigma.setZero(kelvin_vector_size);
185 ip_data.strain_energy_tensile = 0.0;
186 ip_data.elastic_energy = 0.0;
187
189 ip_data.dNdx = shape_matrices[ip].dNdx;
190
192 }
193 }
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.

References _element, _integration_method, _ip_data, _is_axially_symmetric, _process_data, _secondary_data, MeshLib::Element::getID(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), 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 201 of file PhaseFieldFEM.h.

207 {
208 OGS_FATAL(
209 "PhaseFieldLocalAssembler: assembly without jacobian is not "
210 "implemented.");
211 }
#define OGS_FATAL(...)
Definition Error.h:19

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 36 of file PhaseFieldFEM-impl.h.

40{
41 using DeformationMatrix =
45 auto const u =
47
50
53
54 auto local_pressure = 0.0;
55 if (_process_data.pressurized_crack)
56 {
57 local_pressure = _process_data.unity_pressure;
58 }
59
60 int const n_integration_points = _integration_method.getNumberOfPoints();
61 for (int ip = 0; ip < n_integration_points; ip++)
62 {
63 auto const& w = _ip_data[ip].integration_weight;
64 auto const& N = _ip_data[ip].N;
65 auto const& dNdx = _ip_data[ip].dNdx;
66
68 std::nullopt, this->_element.getID(),
71 this->_element, N))};
72
73 auto const x_coord =
74 x_position.getCoordinates().value()[0]; // r for axisymmetry
75
76 auto const& B =
81
82 auto& eps = _ip_data[ip].eps;
83 eps.noalias() = B * u;
84 double const k = _process_data.residual_stiffness(t, x_position)[0];
85 double const ls = _process_data.crack_length_scale(t, x_position)[0];
86 double const d_ip = N.dot(d);
87 double const degradation =
88 _process_data.degradation_derivative->degradation(d_ip, k, ls);
89 _ip_data[ip].updateConstitutiveRelation(
91 _process_data.energy_split_model);
92
93 auto& sigma = _ip_data[ip].sigma;
94 auto const& D = _ip_data[ip].D;
95
101
102 for (int i = 0; i < DisplacementDim; ++i)
103 {
106 .noalias() = N;
107 }
108
109 auto const rho_sr = _process_data.solid_density(t, x_position)[0];
110 auto const& b = _process_data.specific_body_force;
111
112 local_rhs.noalias() -=
113 (B.transpose() * sigma - N_u.transpose() * rho_sr * b -
114 local_pressure * N_u.transpose() * dNdx * d) *
115 w;
116 local_Jac.noalias() += B.transpose() * D * B * w;
117 }
118}
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)

References _element, _integration_method, _ip_data, _is_axially_symmetric, _process_data, ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), displacement_index, displacement_size, ParameterLib::SpatialPosition::getCoordinates(), NumLib::interpolateCoordinates(), and phasefield_index.

Referenced by assembleWithJacobianForStaggeredScheme().

◆ 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 16 of file PhaseFieldFEM-impl.h.

21{
22 // For the equations with phase field.
24 {
27 return;
28 }
29
30 // For the equations with deformation
33}
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)

References assembleWithJacobianForDeformationEquations(), assembleWithJacobianPhaseFieldEquations(), and phase_process_id.

◆ 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 121 of file PhaseFieldFEM-impl.h.

126{
127 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
128 auto const u =
130
135
136 auto local_pressure = 0.0;
137 if (_process_data.static_pressurized_crack)
138 {
139 local_pressure = _process_data.unity_pressure;
140 }
141 else if (_process_data.propagating_pressurized_crack)
142 {
143 local_pressure = _process_data.pressure;
144 }
145
146 int const n_integration_points = _integration_method.getNumberOfPoints();
147 for (int ip = 0; ip < n_integration_points; ip++)
148 {
149 auto const& w = _ip_data[ip].integration_weight;
150 auto const& N = _ip_data[ip].N;
151 auto const& dNdx = _ip_data[ip].dNdx;
152
153 double const d_ip = N.dot(d);
154
156 std::nullopt, this->_element.getID(),
159 this->_element, N))};
160
161 double const k = _process_data.residual_stiffness(t, x_position)[0];
162 double const ls = _process_data.crack_length_scale(t, x_position)[0];
163 double const gc = _process_data.crack_resistance(t, x_position)[0];
164 double const degradation =
165 _process_data.degradation_derivative->degradation(d_ip, k, ls);
166 double const degradation_df1 =
167 _process_data.degradation_derivative->degradationDf1(d_ip, k, ls);
168 double const degradation_df2 =
169 _process_data.degradation_derivative->degradationDf2(d_ip, k, ls);
170
171 auto const x_coord =
172 x_position.getCoordinates().value()[0]; // r for axisymmetry
173 auto const& B =
178
179 auto& eps = _ip_data[ip].eps;
180 eps.noalias() = B * u;
181 _ip_data[ip].updateConstitutiveRelation(
183 _process_data.energy_split_model);
184
185 auto const& strain_energy_tensile = _ip_data[ip].strain_energy_tensile;
186
187 auto& ip_data = _ip_data[ip];
188 ip_data.strain_energy_tensile = strain_energy_tensile;
189
195
196 for (int i = 0; i < DisplacementDim; ++i)
197 {
200 .noalias() = N;
201 }
202
203 local_Jac.noalias() +=
204 (N.transpose() * N * degradation_df2 * strain_energy_tensile) * w;
205
206 local_rhs.noalias() -=
207 (N.transpose() * degradation_df1 * strain_energy_tensile -
208 local_pressure * dNdx.transpose() * N_u * u) *
209 w;
210
212 decltype(dNdx), decltype(N), decltype(w), decltype(d),
213 decltype(local_Jac), decltype(local_rhs)>(
214 dNdx, N, w, d, local_Jac, local_rhs, gc, ls,
215 _process_data.phasefield_model);
216 }
217}

References _element, _integration_method, _ip_data, _is_axially_symmetric, _process_data, ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), displacement_index, displacement_size, ParameterLib::SpatialPosition::getCoordinates(), NumLib::interpolateCoordinates(), phasefield_index, and phasefield_size.

Referenced by assembleWithJacobianForStaggeredScheme().

◆ 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 220 of file PhaseFieldFEM-impl.h.

226{
228 indices_of_processes.reserve(dof_tables.size());
229 std::transform(dof_tables.begin(), dof_tables.end(),
231 [&](auto const dof_table)
232 { return NumLib::getIndices(mesh_item_id, *dof_table); });
233
236
241
242 int const n_integration_points = _integration_method.getNumberOfPoints();
243 for (int ip = 0; ip < n_integration_points; ip++)
244 {
245 auto const& w = _ip_data[ip].integration_weight;
246 auto const& N = _ip_data[ip].N;
247 auto const& dNdx = _ip_data[ip].dNdx;
248
254
255 for (int i = 0; i < DisplacementDim; ++i)
256 {
259 .noalias() = N;
260 }
261
262 crack_volume += (N_u * u).dot(dNdx * d) * w;
263 }
264}
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)

References _integration_method, _ip_data, displacement_index, displacement_size, ProcessLib::getCoupledLocalSolutions(), phasefield_index, and phasefield_size.

◆ 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 267 of file PhaseFieldFEM-impl.h.

272{
274 indices_of_processes.reserve(dof_tables.size());
275 std::transform(dof_tables.begin(), dof_tables.end(),
277 [&](auto const dof_table)
278 { return NumLib::getIndices(mesh_item_id, *dof_table); });
279
280 auto const local_coupled_xs =
283
288
289 double element_elastic_energy = 0.0;
290 double element_surface_energy = 0.0;
291 double element_pressure_work = 0.0;
292
293 int const n_integration_points = _integration_method.getNumberOfPoints();
294 for (int ip = 0; ip < n_integration_points; ip++)
295 {
296 auto const& w = _ip_data[ip].integration_weight;
297 auto const& N = _ip_data[ip].N;
298 auto const& dNdx = _ip_data[ip].dNdx;
299 auto pressure_ip = _process_data.pressure;
300
302 std::nullopt, this->_element.getID(),
305 this->_element, N))};
306
312
313 for (int i = 0; i < DisplacementDim; ++i)
314 {
317 .noalias() = N;
318 }
319
320 element_elastic_energy += _ip_data[ip].elastic_energy * w;
321
322 double const gc = _process_data.crack_resistance(t, x_position)[0];
323 double const ls = _process_data.crack_length_scale(t, x_position)[0];
324 double const d_ip = N.dot(d);
325 switch (_process_data.phasefield_model)
326 {
328 {
330 gc * 0.375 *
331 ((1 - d_ip) / ls + (dNdx * d).dot((dNdx * d)) * ls) * w;
332
333 break;
334 }
336 {
337 element_surface_energy += 0.5 * gc *
338 ((1 - d_ip) * (1 - d_ip) / ls +
339 (dNdx * d).dot((dNdx * d)) * ls) *
340 w;
341 break;
342 }
344 {
347 ((1 - d_ip * d_ip) / ls + (dNdx * d).dot((dNdx * d)) * ls) *
348 w;
349 break;
350 }
351 }
352
353 if (_process_data.pressurized_crack)
354 {
356 }
357 }
358
359#ifdef USE_PETSC
360 int const n_all_nodes = indices_of_processes[1].size();
361 int const n_regular_nodes = std::count_if(
363 [](GlobalIndexType const& index) { return index >= 0; });
365 {
367 static_cast<double>(n_regular_nodes) / n_all_nodes;
369 static_cast<double>(n_regular_nodes) / n_all_nodes;
371 static_cast<double>(n_regular_nodes) / n_all_nodes;
372 }
373#endif // USE_PETSC
377}

References _element, _integration_method, _ip_data, _process_data, MaterialLib::Solids::Phasefield::AT1, MaterialLib::Solids::Phasefield::AT2, MaterialLib::Solids::Phasefield::COHESIVE, displacement_index, displacement_size, ProcessLib::getCoupledLocalSolutions(), NumLib::interpolateCoordinates(), phasefield_index, and phasefield_size.

◆ getEpsilon()

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

◆ 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 277 of file PhaseFieldFEM.h.

279 {
280 auto const& N = _secondary_data.N[integration_point];
281
282 // assumes N is stored contiguously in memory
283 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
284 }

References _secondary_data.

◆ 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 219 of file PhaseFieldFEM.h.

220 {
221 unsigned const n_integration_points =
222 _integration_method.getNumberOfPoints();
223
224 for (unsigned ip = 0; ip < n_integration_points; ip++)
225 {
226 auto& ip_data = _ip_data[ip];
227
229 std::nullopt, _element.getID(),
233 _element, ip_data.N))};
234
236 if (_process_data.initial_stress.value)
237 {
238 ip_data.sigma =
240 DisplacementDim>((*_process_data.initial_stress.value)(
242 double>::quiet_NaN() /* time independent */,
243 x_position));
244 }
245
246 ip_data.pushBackState();
247 }
248 }

References _element, _integration_method, _ip_data, _process_data, 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 250 of file PhaseFieldFEM.h.

254 {
255 unsigned const n_integration_points =
256 _integration_method.getNumberOfPoints();
257
258 for (unsigned ip = 0; ip < n_integration_points; ip++)
259 {
260 _ip_data[ip].pushBackState();
261 }
262 }

References _integration_method, and _ip_data.

◆ 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 382 of file PhaseFieldFEM-impl.h.

385{
386 if (integration_order !=
387 static_cast<int>(_integration_method.getIntegrationOrder()))
388 {
389 OGS_FATAL(
390 "Setting integration point initial conditions; The integration "
391 "order of the local assembler for element {:d} is different from "
392 "the integration order in the initial condition.",
393 _element.getID());
394 }
395
396 if (name == "sigma")
397 {
398 if (_process_data.initial_stress.value != nullptr)
399 {
400 OGS_FATAL(
401 "Setting initial conditions for stress from integration "
402 "point data and from a parameter '{:s}' is not possible "
403 "simultaneously.",
404 _process_data.initial_stress.value->name);
405 }
406
409 }
410
411 return 0;
412}
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)

References _element, _integration_method, _ip_data, _process_data, OGS_FATAL, ProcessLib::setIntegrationPointKelvinVectorData(), and ProcessLib::PhaseField::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::sigma.

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

◆ _process_data

◆ _secondary_data

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

Definition at line 358 of file PhaseFieldFEM.h.

Referenced by PhaseFieldLocalAssembler(), and getShapeMatrix().

◆ displacement_index

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

◆ displacement_size

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

Definition at line 104 of file PhaseFieldFEM.h.

Referenced by assembleWithJacobianForDeformationEquations(), assembleWithJacobianPhaseFieldEquations(), computeCrackIntegral(), and computeEnergy().

◆ mechanics_process_id

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

Definition at line 362 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 361 of file PhaseFieldFEM.h.

Referenced by assembleWithJacobianForStaggeredScheme().

◆ phasefield_index

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

◆ phasefield_size

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

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