33namespace ThermoMechanicalPhaseField
35template <
typename BMatricesType,
typename ShapeMatrixType,
int DisplacementDim>
47 typename ShapeMatrixType::NodalRowVectorType
N;
48 typename ShapeMatrixType::GlobalDimNodalMatrixType
dNdx;
51 typename BMatricesType::KelvinVectorType
eps_m;
56 typename ShapeMatrixType::GlobalDimVectorType
heatflux;
60 DisplacementDim>::MaterialStateVariables>
75 template <
typename DisplacementVectorType>
79 DisplacementVectorType
const& ,
82 double const degradation)
86 auto linear_elastic_mp =
89 .getMaterialProperties();
91 auto const bulk_modulus = linear_elastic_mp.bulk_modulus(t, x);
92 auto const mu = linear_elastic_mp.mu(t, x);
97 DisplacementDim>(degradation, bulk_modulus, mu,
eps_m);
104template <
typename ShapeMatrixType>
107 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>>
N;
110template <
typename ShapeFunction,
int DisplacementDim>
120 ShapeFunction::NPOINTS * DisplacementDim;
138 typename ShapeMatricesType::template VectorType<temperature_size>;
140 typename ShapeMatricesType::template VectorType<displacement_size>;
142 typename ShapeMatricesType::template VectorType<phasefield_size>;
155 bool const is_axially_symmetric,
157 int const mechanics_related_process_id,
158 int const phase_field_process_id,
159 int const heat_conduction_process_id)
168 unsigned const n_integration_points =
171 _ip_data.reserve(n_integration_points);
174 auto& solid_material =
180 DisplacementDim
> const*>(&solid_material))
183 "For phasefield process only linear elastic solid material "
184 "support is implemented.");
187 auto const shape_matrices =
189 DisplacementDim>(e, is_axially_symmetric,
195 for (
unsigned ip = 0; ip < n_integration_points; ip++)
197 _ip_data.emplace_back(solid_material);
199 ip_data.integration_weight =
201 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
203 static const int kelvin_vector_size =
206 ip_data.eps.setZero(kelvin_vector_size);
207 ip_data.eps_prev.resize(kelvin_vector_size);
208 ip_data.eps_m.setZero(kelvin_vector_size);
210 ip_data.D.setZero(kelvin_vector_size, kelvin_vector_size);
211 ip_data.sigma_tensile.setZero(kelvin_vector_size);
212 ip_data.sigma_compressive.setZero(kelvin_vector_size);
213 ip_data.heatflux.setZero(DisplacementDim);
214 ip_data.sigma.setZero(kelvin_vector_size);
215 ip_data.strain_energy_tensile = 0.0;
216 ip_data.elastic_energy = 0.0;
218 ip_data.N = shape_matrices[ip].N;
219 ip_data.dNdx = shape_matrices[ip].dNdx;
226 std::vector<double>
const& ,
227 std::vector<double>
const& ,
228 std::vector<double>& ,
229 std::vector<double>& ,
230 std::vector<double>& )
override
233 "ThermoMechanicalPhaseFieldLocalAssembler: assembly without "
234 "Jacobian is not implemented.");
238 double const t,
double const dt, Eigen::VectorXd
const& local_x,
239 Eigen::VectorXd
const& local_x_prev,
int const process_id,
240 std::vector<double>& local_b_data,
241 std::vector<double>& local_Jac_data)
override;
245 unsigned const n_integration_points =
248 for (
unsigned ip = 0; ip < n_integration_points; ip++)
255 Eigen::VectorXd
const& ,
256 double const ,
double const ,
259 unsigned const n_integration_points =
262 for (
unsigned ip = 0; ip < n_integration_points; ip++)
269 const unsigned integration_point)
const override
274 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
280 std::vector<GlobalVector*>
const& ,
281 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
282 std::vector<double>& cache)
const override
290 std::vector<GlobalVector*>
const& ,
291 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
292 std::vector<double>& cache)
const override
300 std::vector<GlobalVector*>
const& ,
301 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
302 std::vector<double>& cache)
const override
306 auto const n_integration_points =
_ip_data.size();
310 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
311 cache, DisplacementDim, n_integration_points);
313 for (
unsigned ip = 0; ip < n_integration_points; ++ip)
315 auto const& heatflux =
_ip_data[ip].heatflux;
317 for (
typename KelvinVectorType::Index component = 0;
318 component < DisplacementDim;
321 cache_mat(component, ip) = heatflux[component];
329 double const t,
double const dt, Eigen::VectorXd
const& local_x,
330 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
333 double const t,
double const dt, Eigen::VectorXd
const& local_x,
334 Eigen::VectorXd
const& local_x_prev, std::vector<double>& local_b_data,
335 std::vector<double>& local_Jac_data);
338 double const t, Eigen::VectorXd
const& local_x,
339 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
343 std::vector<IpData, Eigen::aligned_allocator<IpData>>
_ip_data;
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.
VectorType< _kelvin_vector_size > KelvinVectorType
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
typename ShapeMatricesType::template VectorType< displacement_size > DeformationVector
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)
typename ShapeMatricesType::template VectorType< temperature_size > TemperatureVector
static constexpr int displacement_index
ThermoMechanicalPhaseFieldLocalAssembler(ThermoMechanicalPhaseFieldLocalAssembler &&)=delete
std::vector< double > const & getIntPtHeatFlux(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void assembleWithJacobianForHeatConductionEquations(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
int const _mechanics_related_process_id
ID of the processes that contains mechanical process.
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
ThermoMechanicalPhaseFieldLocalAssembler(MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoMechanicalPhaseFieldProcessData< DisplacementDim > &process_data, int const mechanics_related_process_id, int const phase_field_process_id, int const heat_conduction_process_id)
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
bool const _is_axially_symmetric
int const _heat_conduction_process_id
ID of heat conduction process.
int const _phase_field_process_id
ID of phase field process.
static constexpr int phasefield_index
static constexpr int displacement_size
void initializeConcrete() 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
ThermoMechanicalPhaseFieldLocalAssembler(ThermoMechanicalPhaseFieldLocalAssembler const &)=delete
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatricesType::template VectorType< phasefield_size > PhaseFieldVector
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
static constexpr int phasefield_size
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
static constexpr int temperature_size
typename BMatricesType::NodalForceVectorType NodalForceVectorType
void assembleWithJacobianForPhaseFieldEquations(double const t, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
ThermoMechanicalPhaseFieldProcessData< DisplacementDim > & _process_data
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
MeshLib::Element const & _element
NumLib::GenericIntegrationMethod const & _integration_method
static constexpr int temperature_index
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateVolDevDegradedStress(double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
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::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
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< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
VectorType< GlobalDim > GlobalDimVectorType
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.
BMatricesType::KelvinVectorType eps_prev
void updateConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x, double const, DisplacementVectorType const &, double const alpha, double const delta_T, double const degradation)
ShapeMatrixType::GlobalDimVectorType heatflux
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
BMatricesType::KelvinVectorType eps
BMatricesType::KelvinVectorType eps_m
BMatricesType::KelvinVectorType sigma_tensile
BMatricesType::KelvinMatrixType C_tensile
ShapeMatrixType::NodalRowVectorType N
double strain_energy_tensile
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
BMatricesType::KelvinVectorType sigma
BMatricesType::KelvinMatrixType D
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
BMatricesType::KelvinMatrixType C_compressive
BMatricesType::KelvinVectorType sigma_compressive
ShapeMatrixType::GlobalDimNodalMatrixType dNdx
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material
double integration_weight
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N