OGS
ThermoMechanicalPhaseFieldFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
30
31namespace ProcessLib
32{
33namespace ThermoMechanicalPhaseField
34{
35template <typename BMatricesType, typename ShapeMatrixType, int DisplacementDim>
37{
46
47 typename ShapeMatrixType::NodalRowVectorType N;
48 typename ShapeMatrixType::GlobalDimNodalMatrixType dNdx;
49
50 typename BMatricesType::KelvinVectorType eps, eps_prev;
51 typename BMatricesType::KelvinVectorType eps_m;
52
53 typename BMatricesType::KelvinVectorType sigma_tensile, sigma_compressive,
56 typename ShapeMatrixType::GlobalDimVectorType heatflux;
57
59 std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
60 DisplacementDim>::MaterialStateVariables>
62
63 typename BMatricesType::KelvinMatrixType D, C_tensile, C_compressive;
65
67 {
68 eps_prev = eps;
69 material_state_variables->pushBackState();
70 }
71
74
75 template <typename DisplacementVectorType>
76 void updateConstitutiveRelation(double const t,
78 double const /*dt*/,
79 DisplacementVectorType const& /*u*/,
80 double const alpha,
81 double const delta_T,
82 double const degradation)
83 {
84 eps_m.noalias() = eps - alpha * delta_T * Invariants::identity2;
85
86 auto linear_elastic_mp =
88 DisplacementDim> const&>(solid_material)
89 .getMaterialProperties();
90
91 auto const bulk_modulus = linear_elastic_mp.bulk_modulus(t, x);
92 auto const mu = linear_elastic_mp.mu(t, x);
93
97 DisplacementDim>(degradation, bulk_modulus, mu, eps_m);
98 }
100};
101
104template <typename ShapeMatrixType>
106{
107 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
108};
109
110template <typename ShapeFunction, int DisplacementDim>
113{
114private:
115 static constexpr int temperature_index = 0;
116 static constexpr int temperature_size = ShapeFunction::NPOINTS;
117 static constexpr int displacement_index =
119 static constexpr int displacement_size =
120 ShapeFunction::NPOINTS * DisplacementDim;
121 static constexpr int phasefield_index =
123 static constexpr int phasefield_size = ShapeFunction::NPOINTS;
124
125public:
128 // Types for displacement.
129 // (Higher order elements = ShapeFunction).
132
134
136
138 typename ShapeMatricesType::template VectorType<temperature_size>;
140 typename ShapeMatricesType::template VectorType<displacement_size>;
142 typename ShapeMatricesType::template VectorType<phasefield_size>;
143 using IpData =
145
150
152 MeshLib::Element const& e,
153 std::size_t const /*local_matrix_size*/,
154 NumLib::GenericIntegrationMethod const& integration_method,
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)
160 : _process_data(process_data),
161 _integration_method(integration_method),
162 _element(e),
163 _is_axially_symmetric(is_axially_symmetric),
164 _mechanics_related_process_id(mechanics_related_process_id),
165 _phase_field_process_id(phase_field_process_id),
166 _heat_conduction_process_id(heat_conduction_process_id)
167 {
168 unsigned const n_integration_points =
170
171 _ip_data.reserve(n_integration_points);
172 _secondary_data.N.resize(n_integration_points);
173
174 auto& solid_material =
176 _process_data.solid_materials,
177 _process_data.material_ids,
178 e.getID());
180 DisplacementDim> const*>(&solid_material))
181 {
182 OGS_FATAL(
183 "For phasefield process only linear elastic solid material "
184 "support is implemented.");
185 }
186
187 auto const shape_matrices =
189 DisplacementDim>(e, is_axially_symmetric,
191
193 x_position.setElementID(_element.getID());
194
195 for (unsigned ip = 0; ip < n_integration_points; ip++)
196 {
197 _ip_data.emplace_back(solid_material);
198 auto& ip_data = _ip_data[ip];
199 ip_data.integration_weight =
201 shape_matrices[ip].integralMeasure * shape_matrices[ip].detJ;
202
203 static const int kelvin_vector_size =
205 DisplacementDim);
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);
209
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;
217
218 ip_data.N = shape_matrices[ip].N;
219 ip_data.dNdx = shape_matrices[ip].dNdx;
220
221 _secondary_data.N[ip] = shape_matrices[ip].N;
222 }
223 }
224
225 void assemble(double const /*t*/, double const /*dt*/,
226 std::vector<double> const& /*local_x*/,
227 std::vector<double> const& /*local_x_prev*/,
228 std::vector<double>& /*local_M_data*/,
229 std::vector<double>& /*local_K_data*/,
230 std::vector<double>& /*local_rhs_data*/) override
231 {
232 OGS_FATAL(
233 "ThermoMechanicalPhaseFieldLocalAssembler: assembly without "
234 "Jacobian is not implemented.");
235 }
236
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;
242
243 void initializeConcrete() override
244 {
245 unsigned const n_integration_points =
247
248 for (unsigned ip = 0; ip < n_integration_points; ip++)
249 {
250 _ip_data[ip].pushBackState();
251 }
252 }
253
254 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
255 Eigen::VectorXd const& /*local_x_prev*/,
256 double const /*t*/, double const /*dt*/,
257 int const /*process_id*/) override
258 {
259 unsigned const n_integration_points =
261
262 for (unsigned ip = 0; ip < n_integration_points; ip++)
263 {
264 _ip_data[ip].pushBackState();
265 }
266 }
267
268 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
269 const unsigned integration_point) const override
270 {
271 auto const& N = _secondary_data.N[integration_point];
272
273 // assumes N is stored contiguously in memory
274 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
275 }
276
277private:
278 std::vector<double> const& getIntPtSigma(
279 const double /*t*/,
280 std::vector<GlobalVector*> const& /*x*/,
281 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
282 std::vector<double>& cache) const override
283 {
284 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
285 _ip_data, &IpData::sigma, cache);
286 }
287
288 std::vector<double> const& getIntPtEpsilon(
289 const double /*t*/,
290 std::vector<GlobalVector*> const& /*x*/,
291 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
292 std::vector<double>& cache) const override
293 {
294 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
295 _ip_data, &IpData::eps, cache);
296 }
297
298 std::vector<double> const& getIntPtHeatFlux(
299 const double /*t*/,
300 std::vector<GlobalVector*> const& /*x*/,
301 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
302 std::vector<double>& cache) const override
303 {
304 using KelvinVectorType = typename BMatricesType::KelvinVectorType;
305
306 auto const n_integration_points = _ip_data.size();
307
308 cache.clear();
309 auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
310 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
311 cache, DisplacementDim, n_integration_points);
312
313 for (unsigned ip = 0; ip < n_integration_points; ++ip)
314 {
315 auto const& heatflux = _ip_data[ip].heatflux;
316
317 for (typename KelvinVectorType::Index component = 0;
318 component < DisplacementDim;
319 ++component)
320 { // x, y, z components
321 cache_mat(component, ip) = heatflux[component];
322 }
323 }
324
325 return cache;
326 }
327
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);
331
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);
336
338 double const t, Eigen::VectorXd const& local_x,
339 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data);
340
342
343 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
344
349
352
355
358};
359
360} // namespace ThermoMechanicalPhaseField
361} // namespace ProcessLib
362
#define OGS_FATAL(...)
Definition Error.h:26
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)
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
VectorType< _kelvin_vector_size > KelvinVectorType
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
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)
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
std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void assembleWithJacobianForPhaseFieldEquations(double const t, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
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)
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.
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N
void updateConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x, double const, DisplacementVectorType const &, double const alpha, double const delta_T, double const degradation)
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
IntegrationPointData(MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N