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

Detailed Description

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

Definition at line 68 of file TESLocalAssembler.h.

#include <TESLocalAssembler.h>

Inheritance diagram for ProcessLib::TES::TESLocalAssembler< ShapeFunction_, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::TES::TESLocalAssembler< ShapeFunction_, 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, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, AssemblyParams const &asm_params)
 
void assemble (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) override
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
 
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
 
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< 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 LAT
 
using NodalMatrixType = typename LAT::LocalMatrix
 
using NodalVectorType = typename LAT::LocalVector
 

Private Attributes

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

Member Typedef Documentation

◆ LAT

template<typename ShapeFunction_ , int GlobalDim>
using ProcessLib::TES::TESLocalAssembler< ShapeFunction_, GlobalDim >::LAT
private
Initial value:
LocalAssemblerTraits<ShapeMatricesType, ShapeFunction::NPOINTS,
NODAL_DOF, GlobalDim>
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
const unsigned NODAL_DOF
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits

Definition at line 139 of file TESLocalAssembler.h.

◆ NodalMatrixType

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

Definition at line 144 of file TESLocalAssembler.h.

◆ NodalVectorType

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

Definition at line 145 of file TESLocalAssembler.h.

◆ ShapeFunction

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

Definition at line 71 of file TESLocalAssembler.h.

◆ ShapeMatrices

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

Definition at line 73 of file TESLocalAssembler.h.

◆ ShapeMatricesType

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

Definition at line 72 of file TESLocalAssembler.h.

Constructor & Destructor Documentation

◆ TESLocalAssembler()

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

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

125 : _element(e),
126 _integration_method(integration_method),
128 ShapeMatricesType, GlobalDim>(
129 e, is_axially_symmetric, _integration_method)),
130 _d(asm_params, _integration_method.getNumberOfPoints(), GlobalDim)
131{
132}
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
TESLocalAssemblerInner< LAT > _d
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::TES::TESLocalAssembler< ShapeFunction_, GlobalDim >::assemble ( 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 )
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
143 // matrices.
144 assert(local_matrix_size == ShapeFunction::NPOINTS * NODAL_DOF);
145
146 auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
147 local_M_data, local_matrix_size, local_matrix_size);
148 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
149 local_K_data, local_matrix_size, local_matrix_size);
150 auto local_b = MathLib::createZeroedVector<NodalVectorType>(
151 local_b_data, local_matrix_size);
152
153 unsigned const n_integration_points =
155
156 _d.preEachAssemble();
157
158 for (unsigned ip = 0; ip < n_integration_points; ip++)
159 {
160 auto const& sm = _shape_matrices[ip];
161 auto const& wp = _integration_method.getWeightedPoint(ip);
162 auto const weight = wp.getWeight();
163
164 _d.assembleIntegrationPoint(ip, local_x, sm, weight, local_M, local_K,
165 local_b);
166 }
167
168 if (_d.getAssemblyParameters().output_element_matrices)
169 {
170 std::puts("### Element: ?");
171
172 std::puts("---Velocity of water");
173 for (auto const& vs : _d.getData().velocity)
174 {
175 std::printf("| ");
176 for (auto v : vs)
177 {
178 std::printf("%23.16e ", v);
179 }
180 std::printf("|\n");
181 }
182
183 std::printf("\n---Mass matrix: \n");
184 ogs5OutMat(local_M);
185 std::printf("\n");
186
187 std::printf("---Laplacian + Advective + Content matrix: \n");
188 ogs5OutMat(local_K);
189 std::printf("\n");
190
191 std::printf("---RHS: \n");
192 ogs5OutVec(local_b);
193 std::printf("\n");
194 }
195}
static std::vector< double > getData(NcVar const &var, std::size_t const total_length, std::size_t const time_step, std::vector< std::size_t > const &length)
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const

References ProcessLib::TES::NODAL_DOF.

◆ checkBounds()

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

Implements ProcessLib::TES::TESLocalAssemblerInterface.

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

302{
303 return _d.getReactionAdaptor().checkBounds(local_x, local_x_prev_ts);
304}

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< ShapeFunction_, 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 257 of file TESLocalAssembler-impl.h.

262{
263 auto const n_integration_points = _integration_method.getNumberOfPoints();
264
265 constexpr int process_id = 0; // monolithic scheme
266 auto const indices =
267 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
268 assert(!indices.empty());
269 auto const local_x = x[process_id]->get(indices);
270 // local_x is ordered by component, local_x_mat is row major
271 auto const local_x_mat = MathLib::toMatrix<
272 Eigen::Matrix<double, NODAL_DOF, Eigen::Dynamic, Eigen::RowMajor>>(
273 local_x, NODAL_DOF, ShapeFunction_::NPOINTS);
274
275 cache.clear();
276 auto cache_mat = MathLib::createZeroedMatrix<
277 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
278 cache, GlobalDim, n_integration_points);
279
280 for (unsigned i = 0; i < n_integration_points; ++i)
281 {
282 double p, T, x;
284 x);
285 const double eta_GR = fluid_viscosity(p, T, x);
286
287 auto const& k = _d.getAssemblyParameters().solid_perm_tensor;
288 cache_mat.col(i).noalias() =
289 k *
290 (_shape_matrices[i].dNdx *
291 local_x_mat.row(COMPONENT_ID_PRESSURE).transpose()) /
292 -eta_GR;
293 }
294
295 return cache;
296}
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
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::detail::shapeFunctionInterpolate(), and MathLib::toMatrix().

◆ getIntPtLoading()

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

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

◆ getIntPtReactionDampingFactor()

template<typename ShapeFunction_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< ShapeFunction_, 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 230 of file TESLocalAssembler-impl.h.

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

◆ getIntPtReactionRate()

template<typename ShapeFunction_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< ShapeFunction_, 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 246 of file TESLocalAssembler-impl.h.

251{
252 return _d.getData().reaction_rate;
253}

◆ getIntPtSolidDensity()

template<typename ShapeFunction_ , int GlobalDim>
std::vector< double > const & ProcessLib::TES::TESLocalAssembler< ShapeFunction_, 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 199 of file TESLocalAssembler-impl.h.

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

◆ getShapeMatrix()

template<typename ShapeFunction_ , int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::TES::TESLocalAssembler< ShapeFunction_, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 89 of file TESLocalAssembler.h.

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

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

Member Data Documentation

◆ _d

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

Definition at line 142 of file TESLocalAssembler.h.

◆ _element

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

Definition at line 132 of file TESLocalAssembler.h.

◆ _integration_method

template<typename ShapeFunction_ , int GlobalDim>
NumLib::GenericIntegrationMethod const& ProcessLib::TES::TESLocalAssembler< ShapeFunction_, GlobalDim >::_integration_method
private

Definition at line 134 of file TESLocalAssembler.h.

◆ _shape_matrices

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

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