OGS
ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >

Definition at line 23 of file StaggeredHTFEM.h.

#include <StaggeredHTFEM.h>

Inheritance diagram for ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >:
[legend]

Public Member Functions

 StaggeredHTFEM (MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, HTProcessData const &process_data)
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) 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::HT::HTFEM< ShapeFunction, GlobalDim >
 HTFEM (MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HTProcessData const &process_data, const unsigned dof_per_node)
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
Eigen::Vector3d getFlux (MathLib::Point3d const &pnt_local_coords, double const t, std::vector< double > const &local_x) const override
Public Member Functions inherited from ProcessLib::HT::HTLocalAssemblerInterface
 HTLocalAssemblerInterface ()=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 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)
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 LocalMatrixType
using LocalVectorType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
using GlobalDimNodalMatrixType
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType

Private Member Functions

void assembleHydraulicEquation (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
void assembleHeatTransportEquation (double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_M_data, std::vector< double > &local_K_data)

Additional Inherited Members

Protected Member Functions inherited from ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >
double getHeatEnergyCoefficient (MaterialPropertyLib::VariableArray const &vars, const double porosity, const double fluid_density, const double specific_heat_capacity_fluid, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
GlobalDimMatrixType getThermalConductivityDispersivity (MaterialPropertyLib::VariableArray const &vars, const double fluid_density, const double specific_heat_capacity_fluid, const GlobalDimVectorType &velocity, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
std::vector< double > const & getIntPtDarcyVelocityLocal (const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
Protected Attributes inherited from ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >
MeshLib::Element const & _element
HTProcessData const & _process_data
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Static Protected Attributes inherited from ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >
static const int pressure_index = ShapeFunction::NPOINTS
static const int pressure_size = ShapeFunction::NPOINTS
static const int temperature_index = 0
static const int temperature_size = ShapeFunction::NPOINTS

Member Typedef Documentation

◆ GlobalDimMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 41 of file StaggeredHTFEM.h.

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::GlobalDimNodalMatrixType
private
Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 39 of file StaggeredHTFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 38 of file StaggeredHTFEM.h.

◆ LocalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::LocalMatrixType
private
Initial value:
typename ShapeMatricesType::template MatrixType<ShapeFunction::NPOINTS,
ShapeFunction::NPOINTS>

Definition at line 28 of file StaggeredHTFEM.h.

◆ LocalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::LocalVectorType
private
Initial value:
typename ShapeMatricesType::template VectorType<ShapeFunction::NPOINTS>

Definition at line 31 of file StaggeredHTFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
private

Definition at line 35 of file StaggeredHTFEM.h.

◆ NodalRowVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
private

Definition at line 36 of file StaggeredHTFEM.h.

◆ NodalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType
private

Definition at line 34 of file StaggeredHTFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 26 of file StaggeredHTFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 25 of file StaggeredHTFEM.h.

Constructor & Destructor Documentation

◆ StaggeredHTFEM()

template<typename ShapeFunction, int GlobalDim>
ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::StaggeredHTFEM ( MeshLib::Element const & element,
std::size_t const local_matrix_size,
NumLib::GenericIntegrationMethod const & integration_method,
bool is_axially_symmetric,
HTProcessData const & process_data )
inline

Definition at line 49 of file StaggeredHTFEM.h.

57 {
58 }
HTFEM(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HTProcessData const &process_data, const unsigned dof_per_node)
Definition HTFEM.h:41

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::HTFEM().

Member Function Documentation

◆ assembleForStaggeredScheme()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::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 )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 18 of file StaggeredHTFEM-impl.h.

23{
24 if (process_id == this->_process_data.heat_transport_process_id)
25 {
28 return;
29 }
30
33}
HTProcessData const & _process_data
Definition HTFEM.h:162
void assembleHydraulicEquation(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
void assembleHeatTransportEquation(double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_M_data, std::vector< double > &local_K_data)

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, assembleHeatTransportEquation(), and assembleHydraulicEquation().

◆ assembleHeatTransportEquation()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::assembleHeatTransportEquation ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
std::vector< double > & local_M_data,
std::vector< double > & local_K_data )
private

Definition at line 170 of file StaggeredHTFEM-impl.h.

176{
177 auto const local_p =
179 auto const local_T =
181
186
187 auto const& process_data = this->_process_data;
188 auto const& medium =
189 *process_data.media_map.getMedium(this->_element.getID());
190 auto const& liquid_phase = medium.phase("AqueousLiquid");
191
192 auto const& b =
194 .projected_specific_body_force_vectors[this->_element.getID()];
195
197
198 unsigned const n_integration_points =
199 this->_integration_method.getNumberOfPoints();
200
202 double average_velocity_norm = 0.0;
204
205 auto const& Ns =
206 process_data.shape_matrix_cache
208
209 for (unsigned ip(0); ip < n_integration_points; ip++)
210 {
211 auto const& ip_data = this->_ip_data[ip];
212 auto const& dNdx = ip_data.dNdx;
213 auto const& N = Ns[ip];
214 auto const& w = ip_data.integration_weight;
215
217 std::nullopt, this->_element.getID(),
220 this->_element, N))};
221
222 double p_at_xi = 0.;
224 double T_at_xi = 0.;
226
227 vars.temperature = T_at_xi;
228 vars.liquid_phase_pressure = p_at_xi;
229
230 vars.liquid_saturation = 1.0;
231
232 auto const porosity =
234 .template value<double>(vars, pos, t, dt);
235 vars.porosity = porosity;
236
237 // Use the fluid density model to compute the density
238 auto const fluid_density =
240 .template value<double>(vars, pos, t, dt);
241 vars.density = fluid_density;
244 .template value<double>(vars, pos, t, dt);
245
246 // Assemble mass matrix
247 local_M.noalias() += w *
251 N.transpose() * N;
252
253 // Assemble Laplace matrix
254 auto const viscosity =
256 .template value<double>(vars, pos, t, dt);
257
258 auto const intrinsic_permeability =
261 .value(vars, pos, t, dt));
262
266 process_data.has_gravity
268 (dNdx * local_p - fluid_density * b))
270
274 pos, t, dt);
275
276 local_K.noalias() +=
278
279 ip_flux_vector.emplace_back(velocity * fluid_density *
282 }
283
285 process_data.stabilizer, this->_ip_data,
286 process_data.shape_matrix_cache, ip_flux_vector,
287 average_velocity_norm / static_cast<double>(n_integration_points),
288 local_K);
289}
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:164
static const int temperature_index
Definition HTFEM.h:316
static const int temperature_size
Definition HTFEM.h:317
double getHeatEnergyCoefficient(MaterialPropertyLib::VariableArray const &vars, const double porosity, const double fluid_density, const double specific_heat_capacity_fluid, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
Definition HTFEM.h:167
GlobalDimMatrixType getThermalConductivityDispersivity(MaterialPropertyLib::VariableArray const &vars, const double fluid_density, const double specific_heat_capacity_fluid, const GlobalDimVectorType &velocity, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
Definition HTFEM.h:191
MeshLib::Element const & _element
Definition HTFEM.h:161
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Definition HTFEM.h:165
static const int pressure_index
Definition HTFEM.h:314
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
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 &)
void assembleAdvectionMatrix(IPData const &ip_data_vector, NumLib::ShapeMatrixCache const &shape_matrix_cache, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_element, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_ip_data, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, NumLib::detail::assembleAdvectionMatrix(), MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getHeatEnergyCoefficient(), ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getThermalConductivityDispersivity(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::VariableArray::porosity, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_index, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::VariableArray::temperature, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_index, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_size, and MaterialPropertyLib::viscosity.

Referenced by assembleForStaggeredScheme().

◆ assembleHydraulicEquation()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::assembleHydraulicEquation ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
std::vector< double > & local_M_data,
std::vector< double > & local_K_data,
std::vector< double > & local_b_data )
private

Definition at line 36 of file StaggeredHTFEM-impl.h.

40{
41 auto const local_p =
43 auto const local_T =
45
46 auto const local_T_prev =
48
55
56 auto const& process_data = this->_process_data;
57 auto const& medium =
58 *this->_process_data.media_map.getMedium(this->_element.getID());
59 auto const& liquid_phase = medium.phase("AqueousLiquid");
60 auto const& solid_phase = medium.phase("Solid");
61
62 auto const& b =
64 .projected_specific_body_force_vectors[this->_element.getID()];
65
67
68 unsigned const n_integration_points =
69 this->_integration_method.getNumberOfPoints();
70
71 auto const& Ns =
72 process_data.shape_matrix_cache
74
75 for (unsigned ip(0); ip < n_integration_points; ip++)
76 {
77 auto const& ip_data = this->_ip_data[ip];
78 auto const& dNdx = ip_data.dNdx;
79 auto const& N = Ns[ip];
80 auto const& w = ip_data.integration_weight;
81
83 std::nullopt, this->_element.getID(),
86 this->_element, N))};
87
88 double p_int_pt = 0.0;
89 double T_int_pt = 0.0;
92
93 vars.temperature = T_int_pt;
94 vars.liquid_phase_pressure = p_int_pt;
95
96 vars.liquid_saturation = 1.0;
97
98 auto const porosity =
100 .template value<double>(vars, pos, t, dt);
101 auto const fluid_density =
103 .template value<double>(vars, pos, t, dt);
104
105 vars.density = fluid_density;
106 const double dfluid_density_dp =
108 .template dValue<double>(
110 pos, t, dt);
111
112 // Use the viscosity model to compute the viscosity
113 auto const viscosity =
115 .template value<double>(vars, pos, t, dt);
116
117 // \todo the argument to getValue() has to be changed for non
118 // constant storage model
119 auto const specific_storage =
121 .template value<double>(vars, pos, t, dt);
122
123 auto const intrinsic_permeability =
126 .value(vars, pos, t, dt));
129
130 // matrix assembly
131 local_M.noalias() +=
132 w *
134 N.transpose() * N;
135
136 local_K.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
137
138 if (process_data.has_gravity)
139 {
140 local_b.noalias() +=
141 w * fluid_density * dNdx.transpose() * K_over_mu * b;
142 }
143
144 if (!process_data.has_fluid_thermal_expansion)
145 {
146 return;
147 }
148
149 // Add the thermal expansion term
150 {
151 auto const solid_thermal_expansion =
152 process_data.solid_thermal_expansion(t, pos)[0];
153 const double dfluid_density_dT =
156 .template dValue<double>(
158 t, dt);
159 double const Tdot_int_pt = (T_int_pt - local_T_prev.dot(N)) / dt;
160 auto const biot_constant = process_data.biot_constant(t, pos)[0];
161 const double eff_thermal_expansion =
164 local_b.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
165 }
166 }
167}
static const int pressure_size
Definition HTFEM.h:315
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_element, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_ip_data, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_process_data, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), NumLib::interpolateCoordinates(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_index, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::pressure_size, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::storage, MaterialPropertyLib::temperature, MaterialPropertyLib::VariableArray::temperature, ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::temperature_index, and MaterialPropertyLib::viscosity.

Referenced by assembleForStaggeredScheme().

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction, int GlobalDim>
std::vector< double > const & ProcessLib::HT::StaggeredHTFEM< 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::HT::HTLocalAssemblerInterface.

Definition at line 293 of file StaggeredHTFEM-impl.h.

298{
299 assert(x.size() == dof_table.size());
300 auto const n_processes = dof_table.size();
301
305 {
306 auto const indices =
308 assert(!indices.empty());
310 }
311 auto const local_xs =
313
315}
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
Definition HTFEM.h:232
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)

References ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::_element, ProcessLib::getCoupledLocalSolutions(), NumLib::getIndices(), and ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getIntPtDarcyVelocityLocal().


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