Loading [MathJax]/extensions/tex2jax.js
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 31 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 std::optional< VectorSegmentgetVectorDeformationSegment () 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 49 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 47 of file StaggeredHTFEM.h.

◆ GlobalDimVectorType

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

Definition at line 46 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 36 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 39 of file StaggeredHTFEM.h.

◆ NodalMatrixType

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

Definition at line 43 of file StaggeredHTFEM.h.

◆ NodalRowVectorType

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

Definition at line 44 of file StaggeredHTFEM.h.

◆ NodalVectorType

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

Definition at line 42 of file StaggeredHTFEM.h.

◆ ShapeMatrices

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

Definition at line 34 of file StaggeredHTFEM.h.

◆ ShapeMatricesType

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

Definition at line 33 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 57 of file StaggeredHTFEM.h.

65 {
66 }
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:48

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 26 of file StaggeredHTFEM-impl.h.

31{
32 if (process_id == this->_process_data.heat_transport_process_id)
33 {
36 return;
37 }
38
41}
HTProcessData const & _process_data
Definition HTFEM.h:169
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 178 of file StaggeredHTFEM-impl.h.

184{
185 auto const local_p =
187 auto const local_T =
189
194
195 auto const& process_data = this->_process_data;
196 auto const& medium =
197 *process_data.media_map.getMedium(this->_element.getID());
198 auto const& liquid_phase = medium.phase("AqueousLiquid");
199
200 auto const& b =
202 .projected_specific_body_force_vectors[this->_element.getID()];
203
205
206 unsigned const n_integration_points =
207 this->_integration_method.getNumberOfPoints();
208
210 double average_velocity_norm = 0.0;
212
213 auto const& Ns =
214 process_data.shape_matrix_cache
216
217 for (unsigned ip(0); ip < n_integration_points; ip++)
218 {
219 auto const& ip_data = this->_ip_data[ip];
220 auto const& dNdx = ip_data.dNdx;
221 auto const& N = Ns[ip];
222 auto const& w = ip_data.integration_weight;
223
225 std::nullopt, this->_element.getID(),
228 this->_element, N))};
229
230 double p_at_xi = 0.;
232 double T_at_xi = 0.;
234
235 vars.temperature = T_at_xi;
236 vars.liquid_phase_pressure = p_at_xi;
237
238 vars.liquid_saturation = 1.0;
239
240 auto const porosity =
242 .template value<double>(vars, pos, t, dt);
243 vars.porosity = porosity;
244
245 // Use the fluid density model to compute the density
246 auto const fluid_density =
248 .template value<double>(vars, pos, t, dt);
249 vars.density = fluid_density;
252 .template value<double>(vars, pos, t, dt);
253
254 // Assemble mass matrix
255 local_M.noalias() += w *
259 N.transpose() * N;
260
261 // Assemble Laplace matrix
262 auto const viscosity =
264 .template value<double>(vars, pos, t, dt);
265
266 auto const intrinsic_permeability =
269 .value(vars, pos, t, dt));
270
274 process_data.has_gravity
276 (dNdx * local_p - fluid_density * b))
278
282 pos, t, dt);
283
284 local_K.noalias() +=
286
287 ip_flux_vector.emplace_back(velocity * fluid_density *
290 }
291
293 process_data.stabilizer, this->_ip_data,
294 process_data.shape_matrix_cache, ip_flux_vector,
295 average_velocity_norm / static_cast<double>(n_integration_points),
296 local_K);
297}
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:171
static const int temperature_index
Definition HTFEM.h:323
static const int temperature_size
Definition HTFEM.h:324
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:174
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:198
MeshLib::Element const & _element
Definition HTFEM.h:168
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Definition HTFEM.h:172
static const int pressure_index
Definition HTFEM.h:321
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
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 44 of file StaggeredHTFEM-impl.h.

48{
49 auto const local_p =
51 auto const local_T =
53
54 auto const local_T_prev =
56
63
64 auto const& process_data = this->_process_data;
65 auto const& medium =
66 *this->_process_data.media_map.getMedium(this->_element.getID());
67 auto const& liquid_phase = medium.phase("AqueousLiquid");
68 auto const& solid_phase = medium.phase("Solid");
69
70 auto const& b =
72 .projected_specific_body_force_vectors[this->_element.getID()];
73
75
76 unsigned const n_integration_points =
77 this->_integration_method.getNumberOfPoints();
78
79 auto const& Ns =
80 process_data.shape_matrix_cache
82
83 for (unsigned ip(0); ip < n_integration_points; ip++)
84 {
85 auto const& ip_data = this->_ip_data[ip];
86 auto const& dNdx = ip_data.dNdx;
87 auto const& N = Ns[ip];
88 auto const& w = ip_data.integration_weight;
89
91 std::nullopt, this->_element.getID(),
94 this->_element, N))};
95
96 double p_int_pt = 0.0;
97 double T_int_pt = 0.0;
100
101 vars.temperature = T_int_pt;
102 vars.liquid_phase_pressure = p_int_pt;
103
104 vars.liquid_saturation = 1.0;
105
106 auto const porosity =
108 .template value<double>(vars, pos, t, dt);
109 auto const fluid_density =
111 .template value<double>(vars, pos, t, dt);
112
113 vars.density = fluid_density;
114 const double dfluid_density_dp =
116 .template dValue<double>(
118 pos, t, dt);
119
120 // Use the viscosity model to compute the viscosity
121 auto const viscosity =
123 .template value<double>(vars, pos, t, dt);
124
125 // \todo the argument to getValue() has to be changed for non
126 // constant storage model
127 auto const specific_storage =
129 .template value<double>(vars, pos, t, dt);
130
131 auto const intrinsic_permeability =
134 .value(vars, pos, t, dt));
137
138 // matrix assembly
139 local_M.noalias() +=
140 w *
142 N.transpose() * N;
143
144 local_K.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
145
146 if (process_data.has_gravity)
147 {
148 local_b.noalias() +=
149 w * fluid_density * dNdx.transpose() * K_over_mu * b;
150 }
151
152 if (!process_data.has_fluid_thermal_expansion)
153 {
154 return;
155 }
156
157 // Add the thermal expansion term
158 {
159 auto const solid_thermal_expansion =
160 process_data.solid_thermal_expansion(t, pos)[0];
161 const double dfluid_density_dT =
164 .template dValue<double>(
166 t, dt);
167 double const Tdot_int_pt = (T_int_pt - local_T_prev.dot(N)) / dt;
168 auto const biot_constant = process_data.biot_constant(t, pos)[0];
169 const double eff_thermal_expansion =
172 local_b.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
173 }
174 }
175}
static const int pressure_size
Definition HTFEM.h:322
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 301 of file StaggeredHTFEM-impl.h.

306{
307 assert(x.size() == dof_table.size());
308 auto const n_processes = dof_table.size();
309
313 {
314 auto const indices =
316 assert(!indices.empty());
318 }
319 auto const local_xs =
321
323}
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
Definition HTFEM.h:239
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: