38template <
typename BMatricesType,
typename ShapeMatrixType,
int DisplacementDim>
50 typename ShapeMatrixType::NodalRowVectorType
N;
51 typename ShapeMatrixType::GlobalDimNodalMatrixType
dNdx;
61 DisplacementDim>::MaterialStateVariables>
77 template <
typename DisplacementVectorType>
81 DisplacementVectorType
const& ,
82 double const degradation,
85 auto linear_elastic_mp =
88 .getMaterialProperties();
90 auto const bulk_modulus = linear_elastic_mp.bulk_modulus(t, x);
91 auto const mu = linear_elastic_mp.mu(t, x);
93 switch (energy_split_model)
99 MaterialLib::Solids::Phasefield::
100 calculateIsotropicDegradedStress<DisplacementDim>(
101 degradation, bulk_modulus, mu,
eps);
108 MaterialLib::Solids::Phasefield::
109 calculateVolDevDegradedStress<DisplacementDim>(
110 degradation, bulk_modulus, mu,
eps);
118 calculateIsotropicDegradedStressWithRankineEnergy<
119 DisplacementDim>(degradation, bulk_modulus, mu,
eps);
124 double temperature = 0.;
126 C_ortho =
static_cast<
129 .getElasticTensor(t, x, temperature);
134 calculateOrthoVolDevDegradedStress<DisplacementDim>(
135 degradation,
eps, C_ortho);
140 double temperature = 0.;
142 C_ortho =
static_cast<
145 .getElasticTensor(t, x, temperature);
150 calculateOrthoMasonryDegradedStress<DisplacementDim>(
151 degradation,
eps, C_ortho);
163template <
typename ShapeMatrixType>
166 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>>
N;
169template <
typename ShapeFunction,
int DisplacementDim>
175 ShapeFunction::NPOINTS * DisplacementDim;
191 typename ShapeMatricesType::template VectorType<displacement_size>;
193 typename ShapeMatricesType::template VectorType<phasefield_size>;
209 bool const is_axially_symmetric,
216 unsigned const n_integration_points =
219 _ip_data.reserve(n_integration_points);
222 auto& solid_material =
228 auto const shape_matrices =
230 DisplacementDim>(e, is_axially_symmetric,
236 for (
unsigned ip = 0; ip < n_integration_points; ip++)
238 _ip_data.emplace_back(solid_material);
240 ip_data.integration_weight =
242 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
244 static const int kelvin_vector_size =
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);
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;
258 ip_data.N = shape_matrices[ip].N;
259 ip_data.dNdx = shape_matrices[ip].dNdx;
266 std::vector<double>
const& ,
267 std::vector<double>
const& ,
268 std::vector<double>& ,
269 std::vector<double>& ,
270 std::vector<double>& )
override
273 "PhaseFieldLocalAssembler: assembly without jacobian is not "
278 double const t,
double const dt, Eigen::VectorXd
const& local_x,
279 Eigen::VectorXd
const& local_x_prev,
int const process_id,
280 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
281 std::vector<double>& local_b_data,
282 std::vector<double>& local_Jac_data)
override;
286 unsigned const n_integration_points =
289 for (
unsigned ip = 0; ip < n_integration_points; ip++)
296 Eigen::VectorXd
const& ,
297 double const ,
double const ,
301 unsigned const n_integration_points =
304 for (
unsigned ip = 0; ip < n_integration_points; ip++)
311 std::size_t mesh_item_id,
313 std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
const&
315 std::vector<GlobalVector*>
const& x,
double const t,
316 double& crack_volume)
override;
319 std::vector<std::reference_wrapper<
321 std::vector<GlobalVector*>
const& x,
double const t,
322 double& elastic_energy,
double& surface_energy,
323 double& pressure_work)
override;
326 const unsigned integration_point)
const override
331 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
337 std::vector<GlobalVector*>
const& ,
338 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
339 std::vector<double>& cache)
const override
341 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
347 std::vector<GlobalVector*>
const& ,
348 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
349 std::vector<double>& cache)
const override
351 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
357 std::vector<GlobalVector*>
const& ,
358 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
359 std::vector<double>& cache)
const override
361 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
367 std::vector<GlobalVector*>
const& ,
368 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
369 std::vector<double>& cache)
const override
371 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
377 std::vector<GlobalVector*>
const& ,
378 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
379 std::vector<double>& cache)
const override
381 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
386 double const t,
double const dt, Eigen::VectorXd
const& local_x,
387 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
390 double const t,
double const dt, Eigen::VectorXd
const& local_x,
391 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
395 std::vector<IpData, Eigen::aligned_allocator<IpData>>
_ip_data;
virtual std::size_t getID() const final
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
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
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_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename BMatricesType::NodalForceVectorType NodalForceVectorType
void computeCrackIntegral(std::size_t mesh_item_id, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap > > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &crack_volume) override
NumLib::GenericIntegrationMethod const & _integration_method
typename ShapeMatricesType::template VectorType< displacement_size > DeformationVector
static const int mechanics_process_id
static constexpr int displacement_size
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
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 MatrixType< phasefield_size, phasefield_size > PhaseFieldMatrix
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, bool const, int const) override
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)
PhaseFieldLocalAssembler(PhaseFieldLocalAssembler &&)=delete
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void computeEnergy(std::size_t mesh_item_id, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap > > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work) override
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
std::vector< double > const & getIntPtSigmaTensile(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
typename ShapeMatricesType::NodalVectorType NodalVectorType
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, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
VectorType< ShapeFunction::NPOINTS > NodalVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N
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
void updateConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x, double const, DisplacementVectorType const &, double const degradation, EnergySplitModel const energy_split_model)
BMatricesType::KelvinMatrixType D
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
BMatricesType::KelvinMatrixType C_tensile
double history_variable_prev
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N