32namespace HeatConduction
43 std::vector<GlobalVector*>
const& x,
44 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
45 std::vector<double>& cache)
const = 0;
48template <
typename ShapeFunction,
int GlobalDim>
66 std::size_t
const local_matrix_size,
68 bool is_axially_symmetric,
80 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
81 (void)local_matrix_size;
85 std::vector<double>
const& local_x,
86 std::vector<double>
const& ,
87 std::vector<double>& local_M_data,
88 std::vector<double>& local_K_data,
89 std::vector<double>& )
override
91 auto const local_matrix_size = local_x.size();
94 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
97 local_M_data, local_matrix_size, local_matrix_size);
99 local_K_data, local_matrix_size, local_matrix_size);
101 unsigned const n_integration_points =
108 for (
unsigned ip = 0; ip < n_integration_points; ip++)
122 double T_int_pt = 0.0;
130 .value(vars, pos, t, dt));
131 auto const specific_heat_capacity =
134 specific_heat_capacity)
135 .
template value<double>(vars, pos, t, dt);
138 .template value<double>(vars, pos, t, dt);
140 local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
141 wp.getWeight() * sm.integralMeasure;
142 local_M.noalias() += sm.N.transpose() * density *
143 specific_heat_capacity * sm.N * sm.detJ *
144 wp.getWeight() * sm.integralMeasure;
148 local_M = local_M.colwise().sum().eval().asDiagonal();
153 std::vector<double>
const& local_x,
154 std::vector<double>
const& local_x_prev,
155 std::vector<double>& local_rhs_data,
156 std::vector<double>& local_Jac_data)
override
158 auto const local_matrix_size = local_x.size();
161 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
163 auto x = Eigen::Map<NodalVectorType const>(local_x.data(),
166 auto x_prev = Eigen::Map<NodalVectorType const>(local_x_prev.data(),
170 local_Jac_data, local_matrix_size, local_matrix_size);
172 local_rhs_data, local_matrix_size);
175 NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
177 NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
179 unsigned const n_integration_points =
189 for (
unsigned ip = 0; ip < n_integration_points; ip++)
199 double T_int_pt = 0.0;
207 .value(vars, pos, t, dt));
208 auto const specific_heat_capacity =
211 specific_heat_capacity)
212 .
template value<double>(vars, pos, t, dt);
215 .template value<double>(vars, pos, t, dt);
217 laplace.noalias() += sm.dNdx.transpose() * k * sm.dNdx * w;
219 sm.N.transpose() * density * specific_heat_capacity * sm.N * w;
223 storage = storage.colwise().sum().eval().asDiagonal();
226 local_Jac.noalias() += laplace + storage / dt;
227 local_rhs.noalias() -= laplace * x + storage * (x - x_prev) / dt;
231 const unsigned integration_point)
const override
236 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
241 std::vector<GlobalVector*>
const& x,
242 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
243 std::vector<double>& cache)
const override
245 int const process_id = 0;
248 assert(!indices.empty());
249 auto const& local_x = x[process_id]->get(indices);
251 auto const T_nodal_values = Eigen::Map<const NodalVectorType>(
252 local_x.data(), ShapeFunction::NPOINTS);
254 unsigned const n_integration_points =
264 double const dt = std::numeric_limits<double>::quiet_NaN();
267 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
268 cache, GlobalDim, n_integration_points);
270 for (
unsigned ip = 0; ip < n_integration_points; ip++)
282 .value(vars, pos, t, dt));
285 cache_mat.col(ip).noalias() = -k * sm.dNdx * T_nodal_values;
296 std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>>
Medium * getMedium(std::size_t element_id)
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)
void setIntegrationPoint(unsigned integration_point)
virtual std::vector< double > const & getIntPtHeatFlux(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
HeatConductionProcessData const & _process_data
typename LocalAssemblerTraits::LocalVector NodalVectorType
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
std::vector< double > const & getIntPtHeatFlux(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
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_rhs_data, std::vector< double > &local_Jac_data) override
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &) override
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
LocalAssemblerData(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, HeatConductionProcessData const &process_data)
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
MeshLib::Element const & _element
NumLib::GenericIntegrationMethod const & _integration_method
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
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)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
const unsigned NUM_NODAL_DOF
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
VectorType< GlobalDim > GlobalDimVectorType
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
bool const mass_lumping
If set mass lumping will be applied to the equation.
Vector< NNodes *NodalDOF > LocalVector
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix