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 49 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 > &, std::vector< double > &, 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_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 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 59 of file HeatConductionFEM.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalAssemblerTraits
private
Initial value:

Definition at line 54 of file HeatConductionFEM.h.

◆ NodalMatrixType

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

Definition at line 57 of file HeatConductionFEM.h.

◆ NodalVectorType

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

Definition at line 58 of file HeatConductionFEM.h.

◆ ShapeMatrices

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

Definition at line 52 of file HeatConductionFEM.h.

◆ ShapeMatricesType

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

Definition at line 51 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 64 of file HeatConductionFEM.h.

70 : _element(element),
71 _process_data(process_data),
72 _integration_method(integration_method),
75 GlobalDim>(
76 element, is_axially_symmetric, _integration_method))
77 {
78 // This assertion is valid only if all nodal d.o.f. use the same shape
79 // matrices.
80 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
81 (void)local_matrix_size;
82 }
HeatConductionProcessData const & _process_data
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
NumLib::GenericIntegrationMethod const & _integration_method
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)

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

90 {
91 auto const local_matrix_size = local_x.size();
92 // This assertion is valid only if all nodal d.o.f. use the same shape
93 // matrices.
94 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
95
96 auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
97 local_M_data, local_matrix_size, local_matrix_size);
98 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
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 {
110 auto const& sm = _shape_matrices[ip];
112 std::nullopt, _element.getID(), ip,
114 NumLib::interpolateCoordinates<ShapeFunction,
116 sm.N))};
117
118 auto const& wp = _integration_method.getWeightedPoint(ip);
119
120 // get the local temperature and put it in the variable array for
121 // access in MPL
122 double T_int_pt = 0.0;
123 NumLib::shapeFunctionInterpolate(local_x, sm.N, T_int_pt);
124 vars.temperature = T_int_pt;
125
126 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
127 medium
128 .property(
130 .value(vars, pos, t, dt));
131 auto const specific_heat_capacity =
132 medium
135 .template value<double>(vars, pos, t, dt);
136 auto const density =
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 *
143 specific_heat_capacity * sm.N * sm.detJ *
144 wp.getWeight() * sm.integralMeasure;
145 }
147 {
148 local_M = local_M.colwise().sum().eval().asDiagonal();
149 }
150 }
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
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, MaterialPropertyLib::density, 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.

◆ 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 > & ,
std::vector< double > & ,
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.

159 {
160 auto const local_matrix_size = local_x.size();
161 // This assertion is valid only if all nodal d.o.f. use the same shape
162 // matrices.
163 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
164
165 auto x = Eigen::Map<NodalVectorType const>(local_x.data(),
166 local_matrix_size);
167
168 auto x_prev = Eigen::Map<NodalVectorType const>(local_x_prev.data(),
169 local_matrix_size);
170
171 auto local_Jac = MathLib::createZeroedMatrix<NodalMatrixType>(
172 local_Jac_data, local_matrix_size, local_matrix_size);
173 auto local_rhs = MathLib::createZeroedVector<NodalVectorType>(
174 local_rhs_data, local_matrix_size);
175
176 NodalMatrixType laplace =
177 NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
179 NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
180
181 unsigned const n_integration_points =
183
186
187 auto const& medium =
190
191 for (unsigned ip = 0; ip < n_integration_points; ip++)
192 {
193 pos.setIntegrationPoint(ip);
194 auto const& sm = _shape_matrices[ip];
195 double const w =
197 sm.integralMeasure;
198
199 // get the local temperature and put it in the variable array for
200 // access in MPL
201 double T_int_pt = 0.0;
202 NumLib::shapeFunctionInterpolate(local_x, sm.N, T_int_pt);
203 vars.temperature = T_int_pt;
204
205 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
206 medium
207 .property(
209 .value(vars, pos, t, dt));
210 auto const specific_heat_capacity =
211 medium
214 .template value<double>(vars, pos, t, dt);
215 auto const density =
217 .template value<double>(vars, pos, t, dt);
218
219 laplace.noalias() += sm.dNdx.transpose() * k * sm.dNdx * w;
220 storage.noalias() +=
221 sm.N.transpose() * density * specific_heat_capacity * sm.N * w;
222 }
224 {
225 storage = storage.colwise().sum().eval().asDiagonal();
226 }
227
228 local_Jac.noalias() += laplace + storage / dt;
229 local_rhs.noalias() -= laplace * x + storage * (x - x_prev) / dt;
230 }
double getWeight() const
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType

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, MaterialPropertyLib::density, 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.

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

246 {
247 int const process_id = 0; // monolithic case.
248 auto const indices =
249 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
250 assert(!indices.empty());
251 auto const& local_x = x[process_id]->get(indices);
252
253 auto const T_nodal_values = Eigen::Map<const NodalVectorType>(
254 local_x.data(), ShapeFunction::NPOINTS);
255
256 unsigned const n_integration_points =
258
261
262 auto const& medium =
265
266 double const dt = std::numeric_limits<double>::quiet_NaN();
267 cache.clear();
268 auto cache_mat = MathLib::createZeroedMatrix<
269 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
270 cache, GlobalDim, n_integration_points);
271
272 for (unsigned ip = 0; ip < n_integration_points; ip++)
273 {
274 pos.setIntegrationPoint(ip);
275 auto const& sm = _shape_matrices[ip];
276 // get the local temperature and put it in the variable array for
277 // access in MPL
278 vars.temperature = sm.N.dot(T_nodal_values);
279
280 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
281 medium
282 .property(
284 .value(vars, pos, t, dt));
285
286 // heat flux only computed for output.
287 cache_mat.col(ip).noalias() = -k * sm.dNdx * T_nodal_values;
288 }
289
290 return cache;
291 }
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
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(), 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.

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

234 {
235 auto const& N = _shape_matrices[integration_point].N;
236
237 // assumes N is stored contiguously in memory
238 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
239 }

References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices.

Member Data Documentation

◆ _element

◆ _integration_method

◆ _process_data

◆ _shape_matrices


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