38template <
typename BMatricesType,
typename ShapeMatrixType,
int DisplacementDim>
50 typename ShapeMatrixType::NodalRowVectorType
N;
51 typename ShapeMatrixType::GlobalDimNodalMatrixType
dNdx;
61 DisplacementDim>::MaterialStateVariables>
77 template <
typename DisplacementVectorType>
82 DisplacementVectorType
const& ,
83 double const degradation,
88 decltype(
sigma),
decltype(
D), DisplacementDim>(
100template <
typename ShapeMatrixType>
103 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>>
N;
106template <
typename ShapeFunction,
int DisplacementDim>
112 ShapeFunction::NPOINTS * DisplacementDim;
128 typename ShapeMatricesType::template VectorType<displacement_size>;
130 typename ShapeMatricesType::template VectorType<phasefield_size>;
146 bool const is_axially_symmetric,
153 unsigned const n_integration_points =
156 _ip_data.reserve(n_integration_points);
159 auto& solid_material =
165 auto const shape_matrices =
167 DisplacementDim>(e, is_axially_symmetric,
173 for (
unsigned ip = 0; ip < n_integration_points; ip++)
175 _ip_data.emplace_back(solid_material);
177 ip_data.integration_weight =
179 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
181 static const int kelvin_vector_size =
184 ip_data.eps_tensile.setZero(kelvin_vector_size);
185 ip_data.eps.setZero(kelvin_vector_size);
186 ip_data.eps_prev.resize(kelvin_vector_size);
187 ip_data.D.setZero(kelvin_vector_size, kelvin_vector_size);
189 ip_data.sigma_tensile.setZero(kelvin_vector_size);
190 ip_data.sigma_compressive.setZero(kelvin_vector_size);
191 ip_data.sigma.setZero(kelvin_vector_size);
192 ip_data.strain_energy_tensile = 0.0;
193 ip_data.elastic_energy = 0.0;
195 ip_data.N = shape_matrices[ip].N;
196 ip_data.dNdx = shape_matrices[ip].dNdx;
204 std::string_view
const name,
205 double const* values,
206 int const integration_order)
override;
209 std::vector<double>
const& ,
210 std::vector<double>
const& ,
211 std::vector<double>& ,
212 std::vector<double>& ,
213 std::vector<double>& )
override
216 "PhaseFieldLocalAssembler: assembly without jacobian is not "
221 double const t,
double const dt, Eigen::VectorXd
const& local_x,
222 Eigen::VectorXd
const& local_x_prev,
int const process_id,
223 std::vector<double>& local_b_data,
224 std::vector<double>& local_Jac_data)
override;
228 unsigned const n_integration_points =
231 for (
unsigned ip = 0; ip < n_integration_points; ip++)
249 double>::quiet_NaN() ,
253 ip_data.pushBackState();
258 Eigen::VectorXd
const& ,
259 double const ,
double const ,
262 unsigned const n_integration_points =
265 for (
unsigned ip = 0; ip < n_integration_points; ip++)
272 std::size_t mesh_item_id,
273 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
274 std::vector<GlobalVector*>
const& x,
double const t,
275 double& crack_volume)
override;
278 std::size_t mesh_item_id,
279 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
280 std::vector<GlobalVector*>
const& x,
double const t,
281 double& elastic_energy,
double& surface_energy,
282 double& pressure_work)
override;
285 const unsigned integration_point)
const override
290 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
296 std::vector<double>
getSigma()
const override;
298 std::vector<double>
getEpsilon()
const override;
303 std::vector<GlobalVector*>
const& ,
304 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
305 std::vector<double>& cache)
const override
313 std::vector<GlobalVector*>
const& ,
314 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
315 std::vector<double>& cache)
const override
323 std::vector<GlobalVector*>
const& ,
324 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
325 std::vector<double>& cache)
const override
333 std::vector<GlobalVector*>
const& ,
334 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
335 std::vector<double>& cache)
const override
343 std::vector<GlobalVector*>
const& ,
344 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
345 std::vector<double>& cache)
const override
352 double const t,
double const dt, Eigen::VectorXd
const& local_x,
353 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
356 double const t,
double const dt, Eigen::VectorXd
const& local_x,
357 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
361 std::vector<IpData, Eigen::aligned_allocator<IpData>>
_ip_data;
constexpr double getWeight() const
std::size_t getID() const
Returns the ID of the element.
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() const
void setElementID(std::size_t element_id)
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
void assembleWithJacobianForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
static constexpr int phasefield_size
static const int phase_process_id
typename ShapeMatricesType::template VectorType< phasefield_size > PhaseFieldVector
PhaseFieldProcessData< DisplacementDim > & _process_data
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
static constexpr int displacement_index
void assembleWithJacobianForDeformationEquations(double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
void assembleWithJacobianPhaseFieldEquations(double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
std::vector< double > const & getIntPtSigmaCompressive(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
MeshLib::Element const & _element
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
Returns number of read integration points.
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename BMatricesType::NodalForceVectorType NodalForceVectorType
NumLib::GenericIntegrationMethod const & _integration_method
static const int mechanics_process_id
static constexpr int displacement_size
PhaseFieldLocalAssembler(PhaseFieldLocalAssembler const &)=delete
void initializeConcrete() override
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtEpsilonTensile(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
static constexpr int phasefield_index
typename ShapeMatricesType::template VectorType< displacement_size > DeformationVector
typename ShapeMatricesType::template MatrixType< phasefield_size, phasefield_size > PhaseFieldMatrix
bool const _is_axially_symmetric
PhaseFieldLocalAssembler(MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, PhaseFieldProcessData< DisplacementDim > &process_data)
void computeCrackIntegral(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &crack_volume) override
PhaseFieldLocalAssembler(PhaseFieldLocalAssembler &&)=delete
void computeEnergy(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work) override
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
std::vector< double > getSigma() const override
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
std::vector< double > const & getIntPtSigmaTensile(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > getEpsilon() const override
typename ShapeMatricesType::NodalVectorType NodalVectorType
void calculateStress(T_VECTOR &sigma, T_VECTOR &sigma_tensile, T_VECTOR &sigma_compressive, T_VECTOR &eps_tensile, T_MATRIX &D, T_MATRIX &C_tensile, T_MATRIX &C_compressive, double &strain_energy_tensile, double &elastic_energy, double const degradation, T_VECTOR const &eps, EnergySplitModel const &energy_split_model, double const t, ParameterLib::SpatialPosition const &x, MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
BMatricesType::KelvinVectorType sigma_tensile
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material
double integration_weight
ShapeMatrixType::NodalRowVectorType N
BMatricesType::KelvinVectorType eps_prev
double strain_energy_tensile
BMatricesType::KelvinMatrixType C_compressive
BMatricesType::KelvinVectorType sigma_compressive
ShapeMatrixType::GlobalDimNodalMatrixType dNdx
BMatricesType::KelvinVectorType eps_tensile
BMatricesType::KelvinVectorType eps
BMatricesType::KelvinVectorType sigma
BMatricesType::KelvinMatrixType D
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
BMatricesType::KelvinMatrixType C_tensile
void updateConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x, double const, DisplacementVectorType const &, double const degradation, MaterialLib::Solids::Phasefield::EnergySplitModel const energy_split_model)
double history_variable_prev
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N