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 43 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_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_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.
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
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 53 of file SteadyStateDiffusionFEM.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::SteadyStateDiffusion::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 48 of file SteadyStateDiffusionFEM.h.

◆ NodalMatrixType

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

Definition at line 51 of file SteadyStateDiffusionFEM.h.

◆ NodalVectorType

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

Definition at line 52 of file SteadyStateDiffusionFEM.h.

◆ ShapeMatrices

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

Definition at line 46 of file SteadyStateDiffusionFEM.h.

◆ ShapeMatricesType

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

Definition at line 45 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 58 of file SteadyStateDiffusionFEM.h.

References _element, _integration_method, _process_data, and _shape_matrices.

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 74 of file SteadyStateDiffusionFEM.h.

80 {
81 auto const local_matrix_size = local_x.size();
82 // This assertion is valid only if all nodal d.o.f. use the same shape
83 // matrices.
85
88
89 unsigned const n_integration_points =
90 _integration_method.getNumberOfPoints();
91
93 pos.setElementID(_element.getID());
94
95 auto const& medium =
96 *_process_data.media_map.getMedium(_element.getID());
98 vars.temperature =
99 medium
100 .property(
102 .template value<double>(vars, pos, t, dt);
103
104 for (unsigned ip = 0; ip < n_integration_points; ip++)
105 {
106 auto const& sm = _shape_matrices[ip];
107 auto const& wp = _integration_method.getWeightedPoint(ip);
108
109 pos = {std::nullopt, _element.getID(),
113 _element, sm.N))};
114
115 double p_int_pt = 0.0;
117 vars.liquid_phase_pressure = p_int_pt;
120 .value(vars, pos, t, dt));
121
122 local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
123 sm.integralMeasure * wp.getWeight();
124 }
125 }
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::diffusion, MaterialPropertyLib::formEigenTensor(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::NUM_NODAL_DOF, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), 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 129 of file SteadyStateDiffusionFEM.h.

132 {
133 // TODO (tf) Temporary value not used by current material models. Need
134 // extension of getFlux interface
136
137 // Eval shape matrices at given point
138 // Note: Axial symmetry is set to false here, because we only need dNdx
139 // here, which is not affected by axial symmetry.
140 auto const shape_matrices =
142 GlobalDim>(
143 _element, false /*is_axially_symmetric*/,
145
146 // fetch hydraulic conductivity
148 pos.setElementID(_element.getID());
149 auto const& medium =
150 *_process_data.media_map.getMedium(_element.getID());
151
153 vars.temperature =
154 medium
155 .property(
157 .template value<double>(vars, pos, t, dt);
158 double pressure = 0.0;
160 vars.liquid_phase_pressure = pressure;
161
164 .value(vars, pos, t, dt));
165
166 Eigen::Vector3d flux(0.0, 0.0, 0.0);
167 flux.head<GlobalDim>() =
168 -k * shape_matrices.dNdx *
170
171 return flux;
172 }

References _element, _process_data, NumLib::computeShapeMatrices(), MaterialPropertyLib::diffusion, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, 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 183 of file SteadyStateDiffusionFEM.h.

188 {
189 // TODO (tf) Temporary value not used by current material models. Need
190 // extension of secondary variable interface.
192
193 auto const n_integration_points =
194 _integration_method.getNumberOfPoints();
195
196 int const process_id = 0; // monolithic scheme
197 auto const indices =
199 assert(!indices.empty());
200 auto const local_x = x[process_id]->get(indices);
201 auto const local_x_vec =
204
205 cache.clear();
209
211 pos.setElementID(_element.getID());
212
213 auto const& medium =
214 *_process_data.media_map.getMedium(_element.getID());
215
217
218 double pressure = 0.0;
219 for (unsigned i = 0; i < n_integration_points; ++i)
220 {
222 pressure);
223 vars.liquid_phase_pressure = pressure;
224
225 pos = {std::nullopt, _element.getID(),
230
233 .value(vars, pos, t, dt));
234 // dimensions: (d x 1) = (d x n) * (n x 1)
235 cache_mat.col(i).noalias() =
236 -k * _shape_matrices[i].dNdx * local_x_vec;
237 }
238
239 return cache;
240 }
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::diffusion, MaterialPropertyLib::formEigenTensor(), NumLib::getIndices(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MathLib::toVector().

◆ 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 174 of file SteadyStateDiffusionFEM.h.

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

References _shape_matrices.

Member Data Documentation

◆ _element

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

◆ _integration_method

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

Definition at line 246 of file SteadyStateDiffusionFEM.h.

Referenced by LocalAssemblerData(), assemble(), and getIntPtDarcyVelocity().

◆ _process_data

template<typename ShapeFunction, int GlobalDim>
SteadyStateDiffusionData const& ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data
private

◆ _shape_matrices

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

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