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

Detailed Description

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

Definition at line 50 of file SteadyStateDiffusionFEM.h.

#include <SteadyStateDiffusionFEM.h>

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

Public Member Functions

 LocalAssemblerData (MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SteadyStateDiffusionData const &process_data)
 
void assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &, std::vector< double > &local_K_data, std::vector< double > &) override
 
Eigen::Vector3d getFlux (MathLib::Point3d const &p_local_coords, double const t, std::vector< double > const &local_x) const 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 & getIntPtDarcyVelocity (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 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 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< 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
 
SteadyStateDiffusionData 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::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 60 of file SteadyStateDiffusionFEM.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalAssemblerTraits
private

◆ NodalMatrixType

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

Definition at line 58 of file SteadyStateDiffusionFEM.h.

◆ NodalVectorType

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

Definition at line 59 of file SteadyStateDiffusionFEM.h.

◆ ShapeMatrices

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

Definition at line 53 of file SteadyStateDiffusionFEM.h.

◆ ShapeMatricesType

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

Definition at line 52 of file SteadyStateDiffusionFEM.h.

Constructor & Destructor Documentation

◆ LocalAssemblerData()

template<typename ShapeFunction , int GlobalDim>
ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::LocalAssemblerData ( MeshLib::Element const & element,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
SteadyStateDiffusionData const & process_data )
inline

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

Definition at line 65 of file SteadyStateDiffusionFEM.h.

71 : _element(element),
72 _process_data(process_data),
73 _integration_method(integration_method),
76 GlobalDim>(
77 element, is_axially_symmetric, _integration_method))
78 {
79 }
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)

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int GlobalDim>
void ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble ( double const t,
double const dt,
std::vector< double > const & local_x,
std::vector< double > const & ,
std::vector< double > & ,
std::vector< double > & local_K_data,
std::vector< double > &  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 81 of file SteadyStateDiffusionFEM.h.

87 {
88 auto const local_matrix_size = local_x.size();
89 // This assertion is valid only if all nodal d.o.f. use the same shape
90 // matrices.
91 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
92
93 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
94 local_K_data, local_matrix_size, local_matrix_size);
95
96 unsigned const n_integration_points =
98
101
102 auto const& medium =
105 vars.temperature =
106 medium
107 .property(
109 .template value<double>(vars, pos, t, dt);
110
111 for (unsigned ip = 0; ip < n_integration_points; ip++)
112 {
113 pos.setIntegrationPoint(ip);
114 auto const& sm = _shape_matrices[ip];
115 auto const& wp = _integration_method.getWeightedPoint(ip);
116
117 double p_int_pt = 0.0;
118 NumLib::shapeFunctionInterpolate(local_x, sm.N, p_int_pt);
119 vars.liquid_phase_pressure = p_int_pt;
120 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
122 .value(vars, pos, t, dt));
123
124 local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
125 sm.integralMeasure * wp.getWeight();
126 }
127 }
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)
void setIntegrationPoint(unsigned integration_point)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
MaterialPropertyLib::MaterialSpatialDistributionMap media_map

References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MaterialPropertyLib::diffusion, MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::GenericIntegrationMethod::getWeightedPoint(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, ProcessLib::SteadyStateDiffusion::NUM_NODAL_DOF, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.

◆ getFlux()

template<typename ShapeFunction , int GlobalDim>
Eigen::Vector3d ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getFlux ( MathLib::Point3d const & p_local_coords,
double const t,
std::vector< double > const & local_x ) const
inlineoverridevirtual

Computes the flux in the point p_local_coords that is given in local coordinates using the values from local_x.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 131 of file SteadyStateDiffusionFEM.h.

134 {
135 // TODO (tf) Temporary value not used by current material models. Need
136 // extension of getFlux interface
137 double const dt = std::numeric_limits<double>::quiet_NaN();
138
139 // Eval shape matrices at given point
140 // Note: Axial symmetry is set to false here, because we only need dNdx
141 // here, which is not affected by axial symmetry.
142 auto const shape_matrices =
144 GlobalDim>(
145 _element, false /*is_axially_symmetric*/,
146 std::array{p_local_coords})[0];
147
148 // fetch hydraulic conductivity
151 auto const& medium =
153
155 vars.temperature =
156 medium
157 .property(
159 .template value<double>(vars, pos, t, dt);
160 double pressure = 0.0;
161 NumLib::shapeFunctionInterpolate(local_x, shape_matrices.N, pressure);
162 vars.liquid_phase_pressure = pressure;
163
164 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
166 .value(vars, pos, t, dt));
167
168 Eigen::Vector3d flux(0.0, 0.0, 0.0);
169 flux.head<GlobalDim>() =
170 -k * shape_matrices.dNdx *
171 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
172
173 return flux;
174 }
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)

References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, NumLib::computeShapeMatrices(), MaterialPropertyLib::diffusion, MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.

◆ getIntPtDarcyVelocity()

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

Implements ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionLocalAssemblerInterface.

Definition at line 185 of file SteadyStateDiffusionFEM.h.

190 {
191 // TODO (tf) Temporary value not used by current material models. Need
192 // extension of secondary variable interface.
193 double const dt = std::numeric_limits<double>::quiet_NaN();
194
195 auto const n_integration_points =
197
198 int const process_id = 0; // monolithic scheme
199 auto const indices =
200 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
201 assert(!indices.empty());
202 auto const local_x = x[process_id]->get(indices);
203 auto const local_x_vec =
204 MathLib::toVector<Eigen::Matrix<double, ShapeFunction::NPOINTS, 1>>(
205 local_x, ShapeFunction::NPOINTS);
206
207 cache.clear();
208 auto cache_mat = MathLib::createZeroedMatrix<
209 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
210 cache, GlobalDim, n_integration_points);
211
214
215 auto const& medium =
217
219 vars.temperature =
220 medium
221 .property(
223 .template value<double>(vars, pos, t, dt);
224 double pressure = 0.0;
225 for (unsigned i = 0; i < n_integration_points; ++i)
226 {
227 pos.setIntegrationPoint(i);
229 pressure);
230 vars.liquid_phase_pressure = pressure;
231
232 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
234 .value(vars, pos, t, dt));
235 // dimensions: (d x 1) = (d x n) * (n x 1)
236 cache_mat.col(i).noalias() =
237 -k * _shape_matrices[i].dNdx * local_x_vec;
238 }
239
240 return cache;
241 }
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::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::diffusion, MeshLib::Element::getID(), NumLib::getIndices(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.

◆ getShapeMatrix()

template<typename ShapeFunction , int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::SteadyStateDiffusion::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 176 of file SteadyStateDiffusionFEM.h.

178 {
179 auto const& N = _shape_matrices[integration_point].N;
180
181 // assumes N is stored contiguously in memory
182 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
183 }

References ProcessLib::SteadyStateDiffusion::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: