template<typename ShapeFunction, int GlobalDim>
class ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >
Definition at line 49 of file HeatConductionFEM.h.
|
| 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) |
|
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 |
|
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 |
|
Eigen::Map< const Eigen::RowVectorXd > | getShapeMatrix (const unsigned integration_point) const override |
| Provides the shape matrix at the given integration point.
|
|
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 |
|
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 | 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) |
|
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.
|
|
template<typename ShapeFunction , int GlobalDim>
void ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::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 > & | ) |
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 84 of file HeatConductionFEM.h.
90 {
91 auto const local_matrix_size = local_x.size();
92
93
94 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
95
97 local_M_data, local_matrix_size, local_matrix_size);
99 local_K_data, local_matrix_size, local_matrix_size);
100
101 unsigned const n_integration_points =
103
104 auto const& medium =
107
108 for (unsigned ip = 0; ip < n_integration_points; ip++)
109 {
116 sm.N))};
117
119
120
121
122 double T_int_pt = 0.0;
125
127 medium
128 .property(
130 .value(vars, pos, t, dt));
132 medium
134 specific_heat_capacity)
135 .template value<double>(vars, pos, t, dt);
138 .template value<double>(vars, pos, t, dt);
139
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 *
144 wp.getWeight() * sm.integralMeasure;
145 }
147 {
148 local_M = local_M.colwise().sum().eval().asDiagonal();
149 }
150 }
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
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
bool const mass_lumping
If set mass lumping will be applied to the equation.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::interpolateCoordinates(), ProcessLib::HeatConduction::HeatConductionProcessData::mass_lumping, ProcessLib::HeatConduction::HeatConductionProcessData::media_map, ProcessLib::HeatConduction::NUM_NODAL_DOF, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.
template<typename ShapeFunction , int GlobalDim>
void ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::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 ) |
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 152 of file HeatConductionFEM.h.
157 {
158 auto const local_matrix_size = local_x.size();
159
160
161 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
162
163 auto x = Eigen::Map<NodalVectorType const>(local_x.data(),
164 local_matrix_size);
165
166 auto x_prev = Eigen::Map<NodalVectorType const>(local_x_prev.data(),
167 local_matrix_size);
168
170 local_Jac_data, local_matrix_size, local_matrix_size);
172 local_rhs_data, local_matrix_size);
173
175 NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
177 NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
178
179 unsigned const n_integration_points =
181
184
185 auto const& medium =
188
189 for (unsigned ip = 0; ip < n_integration_points; ip++)
190 {
193 double const w =
195 sm.integralMeasure;
196
197
198
199 double T_int_pt = 0.0;
202
204 medium
205 .property(
207 .value(vars, pos, t, dt));
209 medium
211 specific_heat_capacity)
212 .template value<double>(vars, pos, t, dt);
215 .template value<double>(vars, pos, t, dt);
216
217 laplace.noalias() += sm.dNdx.transpose() * k * sm.dNdx * w;
220 }
222 {
224 }
225
226 local_Jac.noalias() += laplace +
storage / dt;
227 local_rhs.noalias() -= laplace * x +
storage * (x - x_prev) / dt;
228 }
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), ProcessLib::HeatConduction::HeatConductionProcessData::mass_lumping, ProcessLib::HeatConduction::HeatConductionProcessData::media_map, ProcessLib::HeatConduction::NUM_NODAL_DOF, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.
template<typename ShapeFunction , int GlobalDim>
Implements ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface.
Definition at line 239 of file HeatConductionFEM.h.
244 {
245 int const process_id = 0;
246 auto const indices =
248 assert(!indices.empty());
249 auto const& local_x = x[process_id]->get(indices);
250
251 auto const T_nodal_values = Eigen::Map<const NodalVectorType>(
252 local_x.data(), ShapeFunction::NPOINTS);
253
254 unsigned const n_integration_points =
256
259
260 auto const& medium =
263
264 double const dt = std::numeric_limits<double>::quiet_NaN();
265 cache.clear();
267 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
268 cache, GlobalDim, n_integration_points);
269
270 for (unsigned ip = 0; ip < n_integration_points; ip++)
271 {
274
275
277
279 medium
280 .property(
282 .value(vars, pos, t, dt));
283
284
285 cache_mat.col(ip).noalias() = -k * sm.dNdx * T_nodal_values;
286 }
287
288 return cache;
289 }
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), NumLib::getIndices(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::HeatConduction::HeatConductionProcessData::media_map, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.