OGS
ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim > Class Template Referencefinal

Detailed Description

template<typename ShapeFunction_, typename IntegrationMethod_, int GlobalDim>
class ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >

Definition at line 68 of file TESLocalAssembler.h.

#include <TESLocalAssembler.h>

Inheritance diagram for ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >:
[legend]

Public Types

using ShapeFunction = ShapeFunction_
 
using ShapeMatricesType = ShapeMatrixPolicyType< ShapeFunction, GlobalDim >
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 

Public Member Functions

 TESLocalAssembler (MeshLib::Element const &e, std::size_t const local_matrix_size, bool is_axially_symmetric, unsigned const integration_order, AssemblyParams const &asm_params)
 
void assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point. More...
 
bool checkBounds (std::vector< double > const &local_x, std::vector< double > const &local_x_prev_ts) override
 
std::vector< double > const & getIntPtSolidDensity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtLoading (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtReactionDampingFactor (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtReactionRate (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
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::TES::TESLocalAssemblerInterface
 ~TESLocalAssemblerInterface () override=default
 
- 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< double > const &) 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 LAT = LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NODAL_DOF, GlobalDim >
 
using NodalMatrixType = typename LAT::LocalMatrix
 
using NodalVectorType = typename LAT::LocalVector
 

Private Attributes

MeshLib::Element const & _element
 
IntegrationMethod_ const _integration_method
 
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
 
TESLocalAssemblerInner< LAT_d
 

Member Typedef Documentation

◆ LAT

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
using ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::LAT = LocalAssemblerTraits<ShapeMatricesType, ShapeFunction::NPOINTS, NODAL_DOF, GlobalDim>
private

Definition at line 138 of file TESLocalAssembler.h.

◆ NodalMatrixType

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
using ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::NodalMatrixType = typename LAT::LocalMatrix
private

Definition at line 143 of file TESLocalAssembler.h.

◆ NodalVectorType

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
using ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::NodalVectorType = typename LAT::LocalVector
private

Definition at line 144 of file TESLocalAssembler.h.

◆ ShapeFunction

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
using ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::ShapeFunction = ShapeFunction_

Definition at line 71 of file TESLocalAssembler.h.

◆ ShapeMatrices

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
using ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 73 of file TESLocalAssembler.h.

◆ ShapeMatricesType

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
using ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>

Definition at line 72 of file TESLocalAssembler.h.

Constructor & Destructor Documentation

◆ TESLocalAssembler()

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::TESLocalAssembler ( MeshLib::Element const &  e,
std::size_t const  local_matrix_size,
bool  is_axially_symmetric,
unsigned const  integration_order,
AssemblyParams const &  asm_params 
)

Definition at line 119 of file TESLocalAssembler-impl.h.

125  : _element(e),
126  _integration_method(integration_order),
128  ShapeMatricesType, GlobalDim>(
129  e, is_axially_symmetric, _integration_method)),
130  _d(asm_params, _integration_method.getNumberOfPoints(), GlobalDim)
131 {
132 }
IntegrationMethod_ const _integration_method
TESLocalAssemblerInner< LAT > _d
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
MeshLib::Element const & _element
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::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::assemble ( double const  t,
double const  dt,
std::vector< double > const &  local_x,
std::vector< double > const &  local_xdot,
std::vector< double > &  local_M_data,
std::vector< double > &  local_K_data,
std::vector< double > &  local_b_data 
)
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 135 of file TESLocalAssembler-impl.h.

140 {
141  auto const local_matrix_size = local_x.size();
142  // This assertion is valid only if all nodal d.o.f. use the same shape matrices.
143  assert(local_matrix_size == ShapeFunction::NPOINTS * NODAL_DOF);
144 
145  auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
146  local_M_data, local_matrix_size, local_matrix_size);
147  auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
148  local_K_data, local_matrix_size, local_matrix_size);
149  auto local_b = MathLib::createZeroedVector<NodalVectorType>(local_b_data,
150  local_matrix_size);
151 
152  unsigned const n_integration_points =
153  _integration_method.getNumberOfPoints();
154 
155  _d.preEachAssemble();
156 
157  for (unsigned ip = 0; ip < n_integration_points; ip++)
158  {
159  auto const& sm = _shape_matrices[ip];
160  auto const& wp = _integration_method.getWeightedPoint(ip);
161  auto const weight = wp.getWeight();
162 
163  _d.assembleIntegrationPoint(ip, local_x, sm, weight, local_M, local_K,
164  local_b);
165  }
166 
167  if (_d.getAssemblyParameters().output_element_matrices)
168  {
169  std::puts("### Element: ?");
170 
171  std::puts("---Velocity of water");
172  for (auto const& vs : _d.getData().velocity)
173  {
174  std::printf("| ");
175  for (auto v : vs)
176  {
177  std::printf("%23.16e ", v);
178  }
179  std::printf("|\n");
180  }
181 
182  std::printf("\n---Mass matrix: \n");
183  ogs5OutMat(local_M);
184  std::printf("\n");
185 
186  std::printf("---Laplacian + Advective + Content matrix: \n");
187  ogs5OutMat(local_K);
188  std::printf("\n");
189 
190  std::printf("---RHS: \n");
191  ogs5OutVec(local_b);
192  std::printf("\n");
193  }
194 }
const unsigned NODAL_DOF

References ProcessLib::TES::NODAL_DOF, anonymous_namespace{TESLocalAssembler-impl.h}::ogs5OutMat(), and anonymous_namespace{TESLocalAssembler-impl.h}::ogs5OutVec().

◆ checkBounds()

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
bool ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::checkBounds ( std::vector< double > const &  local_x,
std::vector< double > const &  local_x_prev_ts 
)
overridevirtual

Implements ProcessLib::TES::TESLocalAssemblerInterface.

Definition at line 302 of file TESLocalAssembler-impl.h.

305 {
306  return _d.getReactionAdaptor().checkBounds(local_x, local_x_prev_ts);
307 }

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< 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
overridevirtual

Implements ProcessLib::TES::TESLocalAssemblerInterface.

Definition at line 260 of file TESLocalAssembler-impl.h.

266 {
267  auto const n_integration_points = _integration_method.getNumberOfPoints();
268 
269  constexpr int process_id = 0; // monolithic scheme
270  auto const indices =
271  NumLib::getIndices(_element.getID(), *dof_table[process_id]);
272  assert(!indices.empty());
273  auto const local_x = x[process_id]->get(indices);
274  // local_x is ordered by component, local_x_mat is row major
275  auto const local_x_mat = MathLib::toMatrix<
276  Eigen::Matrix<double, NODAL_DOF, Eigen::Dynamic, Eigen::RowMajor>>(
277  local_x, NODAL_DOF, ShapeFunction_::NPOINTS);
278 
279  cache.clear();
280  auto cache_mat = MathLib::createZeroedMatrix<
281  Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
282  cache, GlobalDim, n_integration_points);
283 
284  for (unsigned i = 0; i < n_integration_points; ++i)
285  {
286  double p, T, x;
288  x);
289  const double eta_GR = fluid_viscosity(p, T, x);
290 
291  auto const& k = _d.getAssemblyParameters().solid_perm_tensor;
292  cache_mat.col(i).noalias() =
293  k * (_shape_matrices[i].dNdx *
294  local_x_mat.row(COMPONENT_ID_PRESSURE).transpose()) /
295  -eta_GR;
296  }
297 
298  return cache;
299 }
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:32
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:72
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
Definition: Interpolation.h:79
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
double fluid_viscosity(const double p, const double T, const double x)
const unsigned COMPONENT_ID_PRESSURE

References ProcessLib::TES::COMPONENT_ID_PRESSURE, MathLib::createZeroedMatrix(), ProcessLib::TES::fluid_viscosity(), NumLib::getIndices(), ProcessLib::TES::NODAL_DOF, NumLib::shapeFunctionInterpolate(), and MathLib::toMatrix().

◆ getIntPtLoading()

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::getIntPtLoading ( const double  t,
std::vector< GlobalVector * > const &  x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  dof_table,
std::vector< double > &  cache 
) const
overridevirtual

Implements ProcessLib::TES::TESLocalAssemblerInterface.

Definition at line 210 of file TESLocalAssembler-impl.h.

216 {
217  auto const rho_SR = _d.getData().solid_density;
218  auto const rho_SR_dry = _d.getAssemblyParameters().rho_SR_dry;
219 
220  cache.clear();
221  cache.reserve(rho_SR.size());
222 
223  transform(cbegin(rho_SR), cend(rho_SR), back_inserter(cache),
224  [&](auto const rho) { return rho / rho_SR_dry - 1.0; });
225 
226  return cache;
227 }
@ rho
density. For some models, rho substitutes p

◆ getIntPtReactionDampingFactor()

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::getIntPtReactionDampingFactor ( const double  t,
std::vector< GlobalVector * > const &  x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  dof_table,
std::vector< double > &  cache 
) const
overridevirtual

Implements ProcessLib::TES::TESLocalAssemblerInterface.

Definition at line 231 of file TESLocalAssembler-impl.h.

236 {
237  auto const fac = _d.getData().reaction_adaptor->getReactionDampingFactor();
238  auto const num_integration_points = _d.getData().solid_density.size();
239 
240  cache.clear();
241  cache.resize(num_integration_points, fac);
242 
243  return cache;
244 }

◆ getIntPtReactionRate()

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::getIntPtReactionRate ( const double  t,
std::vector< GlobalVector * > const &  x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  dof_table,
std::vector< double > &  cache 
) const
overridevirtual

Implements ProcessLib::TES::TESLocalAssemblerInterface.

Definition at line 248 of file TESLocalAssembler-impl.h.

254 {
255  return _d.getData().reaction_rate;
256 }

◆ getIntPtSolidDensity()

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::getIntPtSolidDensity ( const double  t,
std::vector< GlobalVector * > const &  x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  dof_table,
std::vector< double > &  cache 
) const
overridevirtual

Implements ProcessLib::TES::TESLocalAssemblerInterface.

Definition at line 198 of file TESLocalAssembler-impl.h.

204 {
205  return _d.getData().solid_density;
206 }

◆ getShapeMatrix()

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
Eigen::Map<const Eigen::RowVectorXd> ProcessLib::TES::TESLocalAssembler< 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 88 of file TESLocalAssembler.h.

90  {
91  auto const& N = _shape_matrices[integration_point].N;
92 
93  // assumes N is stored contiguously in memory
94  return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
95  }

References ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::_shape_matrices.

Member Data Documentation

◆ _d

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
TESLocalAssemblerInner<LAT> ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::_d
private

Definition at line 141 of file TESLocalAssembler.h.

◆ _element

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
MeshLib::Element const& ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::_element
private

Definition at line 131 of file TESLocalAssembler.h.

◆ _integration_method

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
IntegrationMethod_ const ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::_integration_method
private

Definition at line 133 of file TESLocalAssembler.h.

◆ _shape_matrices

template<typename ShapeFunction_ , typename IntegrationMethod_ , int GlobalDim>
std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices> > ProcessLib::TES::TESLocalAssembler< ShapeFunction_, IntegrationMethod_, GlobalDim >::_shape_matrices
private

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