OGS
ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >

Definition at line 42 of file HeatConductionFEM.h.

#include <HeatConductionFEM.h>

Inheritance diagram for ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >:
[legend]

Public Member Functions

 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
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 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.
virtual int getNumberOfVectorElementsForDeformation () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Private Types

using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
using LocalAssemblerTraits
using NodalMatrixType = typename LocalAssemblerTraits::LocalMatrix
using NodalVectorType = typename LocalAssemblerTraits::LocalVector
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType

Private Attributes

MeshLib::Element const & _element
HeatConductionProcessData const & _process_data
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices

Member Typedef Documentation

◆ GlobalDimVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 52 of file HeatConductionFEM.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalAssemblerTraits
private
Initial value:
ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits

Definition at line 47 of file HeatConductionFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::NodalMatrixType = typename LocalAssemblerTraits::LocalMatrix
private

Definition at line 50 of file HeatConductionFEM.h.

◆ NodalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::NodalVectorType = typename LocalAssemblerTraits::LocalVector
private

Definition at line 51 of file HeatConductionFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 45 of file HeatConductionFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 44 of file HeatConductionFEM.h.

Constructor & Destructor Documentation

◆ LocalAssemblerData()

template<typename ShapeFunction, int GlobalDim>
ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::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 )
inline

The thermal_conductivity factor is directly integrated into the local element matrix.

Definition at line 57 of file HeatConductionFEM.h.

68 GlobalDim>(
70 {
71 // This assertion is valid only if all nodal d.o.f. use the same shape
72 // matrices.
75 }
HeatConductionProcessData const & _process_data
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
NumLib::GenericIntegrationMethod const & _integration_method

References _element, _integration_method, _process_data, _shape_matrices, and ProcessLib::HeatConduction::NUM_NODAL_DOF.

Member Function Documentation

◆ assemble()

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 77 of file HeatConductionFEM.h.

83 {
84 auto const local_matrix_size = local_x.size();
85 // This assertion is valid only if all nodal d.o.f. use the same shape
86 // matrices.
88
93
94 unsigned const n_integration_points =
95 _integration_method.getNumberOfPoints();
96
97 auto const& medium =
98 *_process_data.media_map.getMedium(_element.getID());
100
101 for (unsigned ip = 0; ip < n_integration_points; ip++)
102 {
103 auto const& sm = _shape_matrices[ip];
105 std::nullopt, _element.getID(),
109 sm.N))};
110
111 auto const& wp = _integration_method.getWeightedPoint(ip);
112
113 // get the local temperature and put it in the variable array for
114 // access in MPL
115 double T_int_pt = 0.0;
117 vars.temperature = T_int_pt;
118
120 medium
121 .property(
123 .value(vars, pos, t, dt));
124 auto const specific_heat_capacity =
125 medium
128 .template value<double>(vars, pos, t, dt);
129 auto const density =
131 .template value<double>(vars, pos, t, dt);
132
133 local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
134 wp.getWeight() * sm.integralMeasure;
135 local_M.noalias() += sm.N.transpose() * density *
136 specific_heat_capacity * sm.N * sm.detJ *
137 wp.getWeight() * sm.integralMeasure;
138 }
139 if (_process_data.mass_lumping)
140 {
141 local_M = local_M.colwise().sum().eval().asDiagonal();
142 }
143 }
constexpr 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 &)

References _element, _integration_method, _process_data, _shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), NumLib::interpolateCoordinates(), ProcessLib::HeatConduction::NUM_NODAL_DOF, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.

◆ assembleWithJacobian()

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 145 of file HeatConductionFEM.h.

150 {
151 auto const local_matrix_size = local_x.size();
152 // This assertion is valid only if all nodal d.o.f. use the same shape
153 // matrices.
155
158
161
166
171
172 unsigned const n_integration_points =
173 _integration_method.getNumberOfPoints();
174
175 auto const& medium =
176 *_process_data.media_map.getMedium(_element.getID());
178
179 for (unsigned ip = 0; ip < n_integration_points; ip++)
180 {
181 auto const& sm = _shape_matrices[ip];
183 std::nullopt, _element.getID(),
187 sm.N))};
188
189 double const w =
190 _integration_method.getWeightedPoint(ip).getWeight() * sm.detJ *
191 sm.integralMeasure;
192
193 // get the local temperature and put it in the variable array for
194 // access in MPL
195 double T_int_pt = 0.0;
197 vars.temperature = T_int_pt;
198
200 medium
201 .property(
203 .value(vars, pos, t, dt));
204 auto const specific_heat_capacity =
205 medium
208 .template value<double>(vars, pos, t, dt);
209 auto const density =
211 .template value<double>(vars, pos, t, dt);
212
213 laplace.noalias() += sm.dNdx.transpose() * k * sm.dNdx * w;
214 storage.noalias() +=
215 sm.N.transpose() * density * specific_heat_capacity * sm.N * w;
216 }
217 if (_process_data.mass_lumping)
218 {
219 storage = storage.colwise().sum().eval().asDiagonal();
220 }
221
222 local_Jac.noalias() += laplace + storage / dt;
223 local_rhs.noalias() -= laplace * x + storage * (x - x_prev) / dt;
224 }
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)

References _element, _integration_method, _process_data, _shape_matrices, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), NumLib::interpolateCoordinates(), ProcessLib::HeatConduction::NUM_NODAL_DOF, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.

◆ getIntPtHeatFlux()

template<typename ShapeFunction, int GlobalDim>
std::vector< double > const & ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtHeatFlux ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface.

Definition at line 235 of file HeatConductionFEM.h.

240 {
241 int const process_id = 0; // monolithic case.
242 auto const indices =
244 assert(!indices.empty());
245 auto const& local_x = x[process_id]->get(indices);
246
249
250 unsigned const n_integration_points =
251 _integration_method.getNumberOfPoints();
252
253 auto const& medium =
254 *_process_data.media_map.getMedium(_element.getID());
256
258 cache.clear();
262
263 for (unsigned ip = 0; ip < n_integration_points; ip++)
264 {
265 auto const& sm = _shape_matrices[ip];
266
268 std::nullopt, _element.getID(),
272 sm.N))};
273
274 // get the local temperature and put it in the variable array for
275 // access in MPL
276 vars.temperature = sm.N.dot(T_nodal_values);
277
279 medium
280 .property(
282 .value(vars, pos, t, dt));
283
284 // heat flux only computed for output.
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 _element, _integration_method, _process_data, _shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::formEigenTensor(), NumLib::getIndices(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.

◆ getShapeMatrix()

template<typename ShapeFunction, int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 226 of file HeatConductionFEM.h.

228 {
229 auto const& N = _shape_matrices[integration_point].N;
230
231 // assumes N is stored contiguously in memory
232 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
233 }

References _shape_matrices.

Member Data Documentation

◆ _element

template<typename ShapeFunction, int GlobalDim>
MeshLib::Element const& ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_element
private

◆ _integration_method

template<typename ShapeFunction, int GlobalDim>
NumLib::GenericIntegrationMethod const& ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method
private

◆ _process_data

template<typename ShapeFunction, int GlobalDim>
HeatConductionProcessData const& ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data
private

◆ _shape_matrices

template<typename ShapeFunction, int GlobalDim>
std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices> > ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices
private

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