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

Detailed Description

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

Definition at line 49 of file SteadyStateDiffusionFEM.h.

#include <SteadyStateDiffusionFEM.h>

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

Public Member Functions

 LocalAssemblerData (MeshLib::Element const &element, std::size_t const, bool is_axially_symmetric, unsigned const integration_order, 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. More...
 
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
 
void setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, 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_xdot, 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_xdot, const double dxdot_dx, const double dx_dx, 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_xdot, const double dxdot_dx, const double dx_dx, 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_dot, 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, double const t, double const dt)
 
void postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
 
virtual std::vector< double > interpolateNodalValuesToIntegrationPoints (std::vector< double > const &)
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const
 Fits to staggered scheme. More...
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Private Types

using ShapeMatricesType = ShapeMatrixPolicyType< ShapeFunction, GlobalDim >
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using LocalAssemblerTraits = ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim >
 
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
 
IntegrationMethod const _integration_method
 
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
 

Member Typedef Documentation

◆ GlobalDimVectorType

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

Definition at line 59 of file SteadyStateDiffusionFEM.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalAssemblerTraits = ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>
private

Definition at line 54 of file SteadyStateDiffusionFEM.h.

◆ NodalMatrixType

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

Definition at line 57 of file SteadyStateDiffusionFEM.h.

◆ NodalVectorType

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

Definition at line 58 of file SteadyStateDiffusionFEM.h.

◆ ShapeMatrices

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

Definition at line 52 of file SteadyStateDiffusionFEM.h.

◆ ShapeMatricesType

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

Definition at line 51 of file SteadyStateDiffusionFEM.h.

Constructor & Destructor Documentation

◆ LocalAssemblerData()

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

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

Definition at line 64 of file SteadyStateDiffusionFEM.h.

69  : _element(element),
70  _process_data(process_data),
71  _integration_method(integration_order),
74  GlobalDim>(
75  element, is_axially_symmetric, _integration_method))
76  {
77  }
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
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 , typename IntegrationMethod , int GlobalDim>
void ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, 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 79 of file SteadyStateDiffusionFEM.h.

85  {
86  auto const local_matrix_size = local_x.size();
87  // This assertion is valid only if all nodal d.o.f. use the same shape
88  // matrices.
89  assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
90 
91  auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
92  local_K_data, local_matrix_size, local_matrix_size);
93 
94  unsigned const n_integration_points =
95  _integration_method.getNumberOfPoints();
96 
99 
100  auto const& medium =
101  *_process_data.media_map->getMedium(_element.getID());
103  vars[static_cast<int>(MaterialPropertyLib::Variable::temperature)] =
104  medium
105  .property(
107  .template value<double>(vars, pos, t, dt);
108 
109  for (unsigned ip = 0; ip < n_integration_points; ip++)
110  {
111  pos.setIntegrationPoint(ip);
112  auto const& sm = _shape_matrices[ip];
113  auto const& wp = _integration_method.getWeightedPoint(ip);
114 
115  double p_int_pt = 0.0;
116  NumLib::shapeFunctionInterpolate(local_x, sm.N, p_int_pt);
117  vars[static_cast<int>(
119  auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
121  .value(vars, pos, t, dt));
122 
123  local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
124  sm.integralMeasure * wp.getWeight();
125  }
126  }
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
Definition: VariableType.h:108
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
Definition: Interpolation.h:79
std::unique_ptr< MaterialPropertyLib::MaterialSpatialDistributionMap > media_map

References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_shape_matrices, MaterialPropertyLib::diffusion, MeshLib::Element::getID(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, ProcessLib::SteadyStateDiffusion::NUM_NODAL_DOF, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::shapeFunctionInterpolate(), and MaterialPropertyLib::temperature.

◆ getFlux()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
Eigen::Vector3d ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, 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 130 of file SteadyStateDiffusionFEM.h.

133  {
134  // TODO (tf) Temporary value not used by current material models. Need
135  // extension of getFlux interface
136  double const dt = std::numeric_limits<double>::quiet_NaN();
137 
138  // Eval shape matrices at given point
139  // Note: Axial symmetry is set to false here, because we only need dNdx
140  // here, which is not affected by axial symmetry.
141  auto const shape_matrices =
143  GlobalDim>(
144  _element, false /*is_axially_symmetric*/,
145  std::array{p_local_coords})[0];
146 
147  // fetch hydraulic conductivity
149  pos.setElementID(_element.getID());
150  auto const& medium =
151  *_process_data.media_map->getMedium(_element.getID());
152 
154  vars[static_cast<int>(MaterialPropertyLib::Variable::temperature)] =
155  medium
156  .property(
158  .template value<double>(vars, pos, t, dt);
159  double pressure = 0.0;
160  NumLib::shapeFunctionInterpolate(local_x, shape_matrices.N, pressure);
161  vars[static_cast<int>(MaterialPropertyLib::Variable::phase_pressure)] =
162  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, IntegrationMethod, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, NumLib::computeShapeMatrices(), MaterialPropertyLib::diffusion, MeshLib::Element::getID(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::shapeFunctionInterpolate(), and MaterialPropertyLib::temperature.

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
std::vector<double> const& ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, 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 =
196  _integration_method.getNumberOfPoints();
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 
213  pos.setElementID(_element.getID());
214 
215  auto const& medium =
216  *_process_data.media_map->getMedium(_element.getID());
217 
219  vars[static_cast<int>(MaterialPropertyLib::Variable::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[static_cast<int>(
232 
233  auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
235  .value(vars, pos, t, dt));
236  // dimensions: (d x 1) = (d x n) * (n x 1)
237  cache_mat.col(i).noalias() =
238  -k * _shape_matrices[i].dNdx * local_x_vec;
239  }
240 
241  return cache;
242  }
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:32
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::diffusion, MeshLib::Element::getID(), NumLib::getIndices(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::shapeFunctionInterpolate(), and MaterialPropertyLib::temperature.

◆ getShapeMatrix()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
Eigen::Map<const Eigen::RowVectorXd> ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, 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, IntegrationMethod, GlobalDim >::_shape_matrices.

Member Data Documentation

◆ _element

◆ _integration_method

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
IntegrationMethod const ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method
private

◆ _process_data

◆ _shape_matrices


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