OGS
ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int DisplacementDim>
class ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >

Definition at line 111 of file ThermoMechanicalPhaseFieldFEM.h.

#include <ThermoMechanicalPhaseFieldFEM.h>

Inheritance diagram for ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesType
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>
 
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType
 
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
 
using TemperatureVector
 
using DeformationVector
 
using PhaseFieldVector
 
using IpData
 

Public Member Functions

 ThermoMechanicalPhaseFieldLocalAssembler (ThermoMechanicalPhaseFieldLocalAssembler const &)=delete
 
 ThermoMechanicalPhaseFieldLocalAssembler (ThermoMechanicalPhaseFieldLocalAssembler &&)=delete
 
 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 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_M_data, std::vector< double > &local_K_data, 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
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
 
- 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_M_data, std::vector< double > &local_K_data, 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.
 
- 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 & getIntPtEpsilon (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtHeatFlux (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
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 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)
 
void assembleWithJacobianForPhaseFieldEquations (double const t, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 

Private Attributes

ThermoMechanicalPhaseFieldProcessData< 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
 
int const _mechanics_related_process_id
 ID of the processes that contains mechanical process.
 
int const _phase_field_process_id
 ID of phase field process.
 
int const _heat_conduction_process_id
 ID of heat conduction process.
 

Static Private Attributes

static constexpr int temperature_index = 0
 
static constexpr int temperature_size = ShapeFunction::NPOINTS
 
static constexpr int displacement_index
 
static constexpr int displacement_size
 
static constexpr int phasefield_index
 
static constexpr int phasefield_size = ShapeFunction::NPOINTS
 

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 131 of file ThermoMechanicalPhaseFieldFEM.h.

◆ DeformationVector

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::DeformationVector
Initial value:
typename ShapeMatricesType::template VectorType<displacement_size>

Definition at line 139 of file ThermoMechanicalPhaseFieldFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType

Definition at line 135 of file ThermoMechanicalPhaseFieldFEM.h.

◆ IpData

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::IpData
Initial value:
IntegrationPointData<BMatricesType, ShapeMatricesType, DisplacementDim>

Definition at line 143 of file ThermoMechanicalPhaseFieldFEM.h.

◆ NodalForceVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType

Definition at line 133 of file ThermoMechanicalPhaseFieldFEM.h.

◆ PhaseFieldVector

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::PhaseFieldVector
Initial value:
typename ShapeMatricesType::template VectorType<phasefield_size>

Definition at line 141 of file ThermoMechanicalPhaseFieldFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 130 of file ThermoMechanicalPhaseFieldFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType

◆ TemperatureVector

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::TemperatureVector
Initial value:
typename ShapeMatricesType::template VectorType<temperature_size>

Definition at line 137 of file ThermoMechanicalPhaseFieldFEM.h.

Constructor & Destructor Documentation

◆ ThermoMechanicalPhaseFieldLocalAssembler() [1/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicalPhaseFieldLocalAssembler ( ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim > const & )
delete

◆ ThermoMechanicalPhaseFieldLocalAssembler() [2/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicalPhaseFieldLocalAssembler ( ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim > && )
delete

◆ ThermoMechanicalPhaseFieldLocalAssembler() [3/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::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 )
inline

Definition at line 151 of file ThermoMechanicalPhaseFieldFEM.h.

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 }
#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)
int const _mechanics_related_process_id
ID of the processes that contains mechanical process.
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.
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< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

References ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_process_data, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N, OGS_FATAL, MaterialLib::Solids::selectSolidConstitutiveRelation(), and ParameterLib::SpatialPosition::setElementID().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< 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 225 of file ThermoMechanicalPhaseFieldFEM.h.

231 {
232 OGS_FATAL(
233 "ThermoMechanicalPhaseFieldLocalAssembler: assembly without "
234 "Jacobian is not implemented.");
235 }

References OGS_FATAL.

◆ assembleWithJacobianForDeformationEquations()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< 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 49 of file ThermoMechanicalPhaseFieldFEM-impl.h.

53{
54 assert(local_x.size() ==
56
57 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
58 auto const u =
59 local_x.template segment<displacement_size>(displacement_index);
60 auto const T =
61 local_x.template segment<temperature_size>(temperature_index);
62
63 auto local_Jac = MathLib::createZeroedMatrix<
64 typename ShapeMatricesType::template MatrixType<displacement_size,
66 local_Jac_data, displacement_size, displacement_size);
67
68 auto local_rhs = MathLib::createZeroedVector<
69 typename ShapeMatricesType::template VectorType<displacement_size>>(
70 local_b_data, displacement_size);
71
73 x_position.setElementID(_element.getID());
74
75 int const n_integration_points = _integration_method.getNumberOfPoints();
76 for (int ip = 0; ip < n_integration_points; ip++)
77 {
78 x_position.setIntegrationPoint(ip);
79 auto const& w = _ip_data[ip].integration_weight;
80 auto const& N = _ip_data[ip].N;
81 auto const& dNdx = _ip_data[ip].dNdx;
82
83 auto const x_coord =
84 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
85 _element, N);
86 auto const& B =
87 LinearBMatrix::computeBMatrix<DisplacementDim,
88 ShapeFunction::NPOINTS,
90 dNdx, N, x_coord, _is_axially_symmetric);
91
92 auto& eps = _ip_data[ip].eps;
93
94 auto const& D = _ip_data[ip].D;
95 auto const& sigma = _ip_data[ip].sigma;
96
97 auto const k = _process_data.residual_stiffness(t, x_position)[0];
98 auto rho_sr = _process_data.solid_density(t, x_position)[0];
99 auto const alpha = _process_data.linear_thermal_expansion_coefficient(
100 t, x_position)[0];
101 double const T0 = _process_data.reference_temperature;
102 auto const& b = _process_data.specific_body_force;
103
104 double const T_ip = N.dot(T);
105 double const delta_T = T_ip - T0;
106 // calculate real density
107 double const rho_s = rho_sr / (1 + 3 * alpha * delta_T);
108
109 double const d_ip = N.dot(d);
110 double const degradation = d_ip * d_ip * (1 - k) + k;
111
112 eps.noalias() = B * u;
113 _ip_data[ip].updateConstitutiveRelation(t, x_position, dt, u, alpha,
114 delta_T, degradation);
115
116 typename ShapeMatricesType::template MatrixType<DisplacementDim,
118 N_u = ShapeMatricesType::template MatrixType<
119 DisplacementDim, displacement_size>::Zero(DisplacementDim,
121
122 for (int i = 0; i < DisplacementDim; ++i)
123 {
124 N_u.template block<1, displacement_size / DisplacementDim>(
125 i, i * displacement_size / DisplacementDim)
126 .noalias() = N;
127 }
128
129 local_Jac.noalias() += B.transpose() * D * B * w;
130
131 local_rhs.noalias() -=
132 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
133 }
134}
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)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.

References ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), and ParameterLib::SpatialPosition::setElementID().

◆ assembleWithJacobianForHeatConductionEquations()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::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 )
private

Definition at line 137 of file ThermoMechanicalPhaseFieldFEM-impl.h.

142{
143 assert(local_x.size() ==
145
146 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
147 auto const T =
148 local_x.template segment<temperature_size>(temperature_index);
149 auto const T_prev =
150 local_x_prev.template segment<temperature_size>(temperature_index);
151 auto local_Jac = MathLib::createZeroedMatrix<
152 typename ShapeMatricesType::template MatrixType<temperature_size,
154 local_Jac_data, temperature_size, temperature_size);
155
156 auto local_rhs = MathLib::createZeroedVector<
157 typename ShapeMatricesType::template VectorType<temperature_size>>(
158 local_b_data, temperature_size);
159
161 x_position.setElementID(_element.getID());
162
163 int const n_integration_points = _integration_method.getNumberOfPoints();
164 for (int ip = 0; ip < n_integration_points; ip++)
165 {
166 x_position.setIntegrationPoint(ip);
167 auto const& w = _ip_data[ip].integration_weight;
168 auto const& N = _ip_data[ip].N;
169 auto const& dNdx = _ip_data[ip].dNdx;
170
171 auto& eps_m = _ip_data[ip].eps_m;
172 auto& heatflux = _ip_data[ip].heatflux;
173
174 auto rho_sr = _process_data.solid_density(t, x_position)[0];
175 auto const alpha = _process_data.linear_thermal_expansion_coefficient(
176 t, x_position)[0];
177 double const c = _process_data.specific_heat_capacity(t, x_position)[0];
178 auto const lambda =
179 _process_data.thermal_conductivity(t, x_position)[0];
180 auto const lambda_res =
181 _process_data.residual_thermal_conductivity(t, x_position)[0];
182 double const T0 = _process_data.reference_temperature;
183
184 double const d_ip = N.dot(d);
185 double const T_ip = N.dot(T);
186 double const T_dot_ip = (T_ip - N.dot(T_prev)) / dt;
187 double const delta_T = T_ip - T0;
188 // calculate real density
189 double const rho_s = rho_sr / (1 + 3 * alpha * delta_T);
190 // calculate effective thermal conductivity
191 auto const lambda_eff =
192 d_ip * d_ip * lambda + (1 - d_ip) * (1 - d_ip) * lambda_res;
193
194 using Invariants = MathLib::KelvinVector::Invariants<
196
197 double const eps_m_trace = Invariants::trace(eps_m);
198 if (eps_m_trace >= 0)
199 {
200 local_Jac.noalias() += (dNdx.transpose() * lambda_eff * dNdx +
201 N.transpose() * rho_s * c * N / dt) *
202 w;
203
204 local_rhs.noalias() -= (N.transpose() * rho_s * c * T_dot_ip +
205 dNdx.transpose() * lambda_eff * dNdx * T) *
206 w;
207
208 heatflux.noalias() = -(lambda_eff * dNdx * T) * w;
209 }
210 else
211 {
212 local_Jac.noalias() += (dNdx.transpose() * lambda * dNdx +
213 N.transpose() * rho_s * c * N / dt) *
214 w;
215
216 local_rhs.noalias() -= (N.transpose() * rho_s * c * T_dot_ip +
217 dNdx.transpose() * lambda * dNdx * T) *
218 w;
219
220 heatflux.noalias() = -(lambda * dNdx * T) * w;
221 }
222 }
223}

References MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MathLib::KelvinVector::kelvin_vector_dimensions(), and ParameterLib::SpatialPosition::setElementID().

◆ assembleWithJacobianForPhaseFieldEquations()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianForPhaseFieldEquations ( double const t,
Eigen::VectorXd const & local_x,
std::vector< double > & local_b_data,
std::vector< double > & local_Jac_data )
private

Definition at line 226 of file ThermoMechanicalPhaseFieldFEM-impl.h.

230{
231 assert(local_x.size() ==
233
234 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
235
236 auto local_Jac = MathLib::createZeroedMatrix<
237 typename ShapeMatricesType::template MatrixType<phasefield_size,
239 local_Jac_data, phasefield_size, phasefield_size);
240
241 auto local_rhs = MathLib::createZeroedVector<
242 typename ShapeMatricesType::template VectorType<phasefield_size>>(
243 local_b_data, phasefield_size);
244
246 x_position.setElementID(_element.getID());
247
248 int const n_integration_points = _integration_method.getNumberOfPoints();
249 for (int ip = 0; ip < n_integration_points; ip++)
250 {
251 x_position.setIntegrationPoint(ip);
252 auto const& w = _ip_data[ip].integration_weight;
253 auto const& N = _ip_data[ip].N;
254 auto const& dNdx = _ip_data[ip].dNdx;
255
256 auto const& strain_energy_tensile = _ip_data[ip].strain_energy_tensile;
257
258 auto const gc = _process_data.crack_resistance(t, x_position)[0];
259 auto const ls = _process_data.crack_length_scale(t, x_position)[0];
260
261 double const d_ip = N.dot(d);
262
263 local_Jac.noalias() += (dNdx.transpose() * gc * ls * dNdx +
264 N.transpose() * 2 * strain_energy_tensile * N +
265 N.transpose() * gc / ls * N) *
266 w;
267
268 local_rhs.noalias() -=
269 (dNdx.transpose() * gc * ls * dNdx * d +
270 N.transpose() * d_ip * 2 * strain_energy_tensile -
271 N.transpose() * gc / ls * (1 - d_ip)) *
272 w;
273 }
274}

References MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), and ParameterLib::SpatialPosition::setElementID().

◆ assembleWithJacobianForStaggeredScheme()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< 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_M_data,
std::vector< double > & local_K_data,
std::vector< double > & local_b_data,
std::vector< double > & local_Jac_data )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 21 of file ThermoMechanicalPhaseFieldFEM-impl.h.

28{
29 if (process_id == _phase_field_process_id)
30 {
31 assembleWithJacobianForPhaseFieldEquations(t, local_x, local_b_data,
32 local_Jac_data);
33 return;
34 }
35
36 if (process_id == _heat_conduction_process_id)
37 {
39 t, dt, local_x, local_x_prev, local_b_data, local_Jac_data);
40 return;
41 }
42
43 // For the equations with deformation
44 assembleWithJacobianForDeformationEquations(t, dt, local_x, local_b_data,
45 local_Jac_data);
46}
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 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)
void assembleWithJacobianForPhaseFieldEquations(double const t, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)

◆ getIntPtEpsilon()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsilon ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprivatevirtual

◆ getIntPtHeatFlux()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtHeatFlux ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverrideprivatevirtual

Implements ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssemblerInterface.

Definition at line 299 of file ThermoMechanicalPhaseFieldFEM.h.

304 {
306
307 auto const n_integration_points = _ip_data.size();
308
309 cache.clear();
310 auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
311 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
312 cache, DisplacementDim, n_integration_points);
313
314 for (unsigned ip = 0; ip < n_integration_points; ++ip)
315 {
316 auto const& heatflux = _ip_data[ip].heatflux;
317
318 for (typename KelvinVectorType::Index component = 0;
319 component < DisplacementDim;
320 ++component)
321 { // x, y, z components
322 cache_mat(component, ip) = heatflux[component];
323 }
324 }
325
326 return cache;
327 }
VectorType< _kelvin_vector_size > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType

References ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, and MathLib::createZeroedMatrix().

◆ getIntPtSigma()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtSigma ( 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::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 269 of file ThermoMechanicalPhaseFieldFEM.h.

271 {
272 auto const& N = _secondary_data.N[integration_point];
273
274 // assumes N is stored contiguously in memory
275 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
276 }

References ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N.

◆ initializeConcrete()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

◆ postTimestepConcrete()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete ( Eigen::VectorXd const & ,
Eigen::VectorXd const & ,
double const ,
double const ,
int const  )
inlineoverridevirtual

Member Data Documentation

◆ _element

◆ _heat_conduction_process_id

template<typename ShapeFunction , int DisplacementDim>
int const ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_heat_conduction_process_id
private

ID of heat conduction process.

Definition at line 358 of file ThermoMechanicalPhaseFieldFEM.h.

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

template<typename ShapeFunction , int DisplacementDim>
bool const ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_is_axially_symmetric
private

Definition at line 349 of file ThermoMechanicalPhaseFieldFEM.h.

◆ _mechanics_related_process_id

template<typename ShapeFunction , int DisplacementDim>
int const ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_mechanics_related_process_id
private

ID of the processes that contains mechanical process.

Definition at line 352 of file ThermoMechanicalPhaseFieldFEM.h.

◆ _phase_field_process_id

template<typename ShapeFunction , int DisplacementDim>
int const ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::_phase_field_process_id
private

ID of phase field process.

Definition at line 355 of file ThermoMechanicalPhaseFieldFEM.h.

◆ _process_data

◆ _secondary_data

◆ displacement_index

template<typename ShapeFunction , int DisplacementDim>
constexpr int ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::displacement_index
staticconstexprprivate
Initial value:

Definition at line 117 of file ThermoMechanicalPhaseFieldFEM.h.

◆ displacement_size

template<typename ShapeFunction , int DisplacementDim>
constexpr int ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::displacement_size
staticconstexprprivate
Initial value:
=
ShapeFunction::NPOINTS * DisplacementDim

Definition at line 119 of file ThermoMechanicalPhaseFieldFEM.h.

◆ phasefield_index

template<typename ShapeFunction , int DisplacementDim>
constexpr int ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::phasefield_index
staticconstexprprivate
Initial value:

Definition at line 121 of file ThermoMechanicalPhaseFieldFEM.h.

◆ phasefield_size

template<typename ShapeFunction , int DisplacementDim>
constexpr int ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::phasefield_size = ShapeFunction::NPOINTS
staticconstexprprivate

Definition at line 123 of file ThermoMechanicalPhaseFieldFEM.h.

◆ temperature_index

template<typename ShapeFunction , int DisplacementDim>
constexpr int ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::temperature_index = 0
staticconstexprprivate

Definition at line 115 of file ThermoMechanicalPhaseFieldFEM.h.

◆ temperature_size

template<typename ShapeFunction , int DisplacementDim>
constexpr int ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::temperature_size = ShapeFunction::NPOINTS
staticconstexprprivate

Definition at line 116 of file ThermoMechanicalPhaseFieldFEM.h.


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