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.
 
- 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.

62 : HTFEM<ShapeFunction, GlobalDim>(element, local_matrix_size,
63 integration_method,
64 is_axially_symmetric, process_data, 1)
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

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 {
34 assembleHeatTransportEquation(t, dt, local_x, local_M_data,
35 local_K_data);
36 return;
37 }
38
39 assembleHydraulicEquation(t, dt, local_x, local_x_prev, local_M_data,
40 local_K_data, local_b_data);
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)

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

183{
184 auto const local_p =
185 local_x.template segment<pressure_size>(pressure_index);
186 auto const local_T =
187 local_x.template segment<temperature_size>(temperature_index);
188
190 local_M_data, temperature_size, temperature_size);
192 local_K_data, temperature_size, temperature_size);
193
195 pos.setElementID(this->_element.getID());
196
197 auto const& process_data = this->_process_data;
198 auto const& medium =
199 *process_data.media_map.getMedium(this->_element.getID());
200 auto const& liquid_phase = medium.phase("AqueousLiquid");
201
202 auto const& b =
203 process_data
204 .projected_specific_body_force_vectors[this->_element.getID()];
205
207
208 unsigned const n_integration_points =
210
211 std::vector<GlobalDimVectorType> ip_flux_vector;
212 double average_velocity_norm = 0.0;
213 ip_flux_vector.reserve(n_integration_points);
214
215 auto const& Ns =
216 process_data.shape_matrix_cache
217 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
218
219 for (unsigned ip(0); ip < n_integration_points; ip++)
220 {
221 pos.setIntegrationPoint(ip);
222
223 auto const& ip_data = this->_ip_data[ip];
224 auto const& dNdx = ip_data.dNdx;
225 auto const& N = Ns[ip];
226 auto const& w = ip_data.integration_weight;
227
228 double p_at_xi = 0.;
229 NumLib::shapeFunctionInterpolate(local_p, N, p_at_xi);
230 double T_at_xi = 0.;
231 NumLib::shapeFunctionInterpolate(local_T, N, T_at_xi);
232
233 vars.temperature = T_at_xi;
234 vars.liquid_phase_pressure = p_at_xi;
235
236 vars.liquid_saturation = 1.0;
237
238 auto const porosity =
240 .template value<double>(vars, pos, t, dt);
241 vars.porosity = porosity;
242
243 // Use the fluid density model to compute the density
244 auto const fluid_density =
245 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
246 .template value<double>(vars, pos, t, dt);
247 vars.density = fluid_density;
248 auto const specific_heat_capacity_fluid =
250 .template value<double>(vars, pos, t, dt);
251
252 // Assemble mass matrix
253 local_M.noalias() += w *
255 vars, porosity, fluid_density,
256 specific_heat_capacity_fluid, pos, t, dt) *
257 N.transpose() * N;
258
259 // Assemble Laplace matrix
260 auto const viscosity =
262 .template value<double>(vars, pos, t, dt);
263
264 auto const intrinsic_permeability =
267 .value(vars, pos, t, dt));
268
269 GlobalDimMatrixType const K_over_mu =
270 intrinsic_permeability / viscosity;
271 GlobalDimVectorType const velocity =
272 process_data.has_gravity
273 ? GlobalDimVectorType(-K_over_mu *
274 (dNdx * local_p - fluid_density * b))
275 : GlobalDimVectorType(-K_over_mu * dNdx * local_p);
276
277 GlobalDimMatrixType const thermal_conductivity_dispersivity =
279 vars, fluid_density, specific_heat_capacity_fluid, velocity,
280 pos, t, dt);
281
282 local_K.noalias() +=
283 w * dNdx.transpose() * thermal_conductivity_dispersivity * dNdx;
284
285 ip_flux_vector.emplace_back(velocity * fluid_density *
286 specific_heat_capacity_fluid);
287 average_velocity_norm += velocity.norm();
288 }
289
291 process_data.stabilizer, this->_ip_data,
292 process_data.shape_matrix_cache, ip_flux_vector,
293 average_velocity_norm / static_cast<double>(n_integration_points),
294 local_K);
295}
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:171
static const int temperature_index
Definition HTFEM.h:325
static const int temperature_size
Definition HTFEM.h:326
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:323
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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)
double fluid_density(const double p, const double T, const double x)
MaterialPropertyLib::MaterialSpatialDistributionMap media_map

References NumLib::detail::assembleAdvectionMatrix(), MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::VariableArray::porosity, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

◆ 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 =
50 local_x.template segment<pressure_size>(pressure_index);
51 auto const local_T =
52 local_x.template segment<temperature_size>(temperature_index);
53
54 auto const local_T_prev =
55 local_x_prev.template segment<temperature_size>(temperature_index);
56
58 local_M_data, pressure_size, pressure_size);
60 local_K_data, pressure_size, pressure_size);
61 auto local_b = MathLib::createZeroedVector<LocalVectorType>(local_b_data,
63
65 pos.setElementID(this->_element.getID());
66
67 auto const& process_data = this->_process_data;
68 auto const& medium =
70 auto const& liquid_phase = medium.phase("AqueousLiquid");
71 auto const& solid_phase = medium.phase("Solid");
72
73 auto const& b =
74 process_data
75 .projected_specific_body_force_vectors[this->_element.getID()];
76
78
79 unsigned const n_integration_points =
81
82 auto const& Ns =
83 process_data.shape_matrix_cache
84 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
85
86 for (unsigned ip(0); ip < n_integration_points; ip++)
87 {
88 pos.setIntegrationPoint(ip);
89
90 auto const& ip_data = this->_ip_data[ip];
91 auto const& dNdx = ip_data.dNdx;
92 auto const& N = Ns[ip];
93 auto const& w = ip_data.integration_weight;
94
95 double p_int_pt = 0.0;
96 double T_int_pt = 0.0;
97 NumLib::shapeFunctionInterpolate(local_p, N, p_int_pt);
98 NumLib::shapeFunctionInterpolate(local_T, N, T_int_pt);
99
100 vars.temperature = T_int_pt;
101 vars.liquid_phase_pressure = p_int_pt;
102
103 vars.liquid_saturation = 1.0;
104
105 auto const porosity =
107 .template value<double>(vars, pos, t, dt);
108 auto const fluid_density =
109 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
110 .template value<double>(vars, pos, t, dt);
111
112 vars.density = fluid_density;
113 const double dfluid_density_dp =
114 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
115 .template dValue<double>(
117 pos, t, dt);
118
119 // Use the viscosity model to compute the viscosity
120 auto const viscosity =
122 .template value<double>(vars, pos, t, dt);
123
124 // \todo the argument to getValue() has to be changed for non
125 // constant storage model
126 auto const specific_storage =
128 .template value<double>(vars, pos, t, dt);
129
130 auto const intrinsic_permeability =
133 .value(vars, pos, t, dt));
134 GlobalDimMatrixType const K_over_mu =
135 intrinsic_permeability / viscosity;
136
137 // matrix assembly
138 local_M.noalias() +=
139 w *
140 (porosity * dfluid_density_dp / fluid_density + specific_storage) *
141 N.transpose() * N;
142
143 local_K.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
144
145 if (process_data.has_gravity)
146 {
147 local_b.noalias() +=
148 w * fluid_density * dNdx.transpose() * K_over_mu * b;
149 }
150
151 if (!process_data.has_fluid_thermal_expansion)
152 {
153 return;
154 }
155
156 // Add the thermal expansion term
157 {
158 auto const solid_thermal_expansion =
159 process_data.solid_thermal_expansion(t, pos)[0];
160 const double dfluid_density_dT =
161 liquid_phase
163 .template dValue<double>(
165 t, dt);
166 double const Tdot_int_pt = (T_int_pt - local_T_prev.dot(N)) / dt;
167 auto const biot_constant = process_data.biot_constant(t, pos)[0];
168 const double eff_thermal_expansion =
169 3.0 * (biot_constant - porosity) * solid_thermal_expansion -
170 porosity * dfluid_density_dT / fluid_density;
171 local_b.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
172 }
173 }
174}
static const int pressure_size
Definition HTFEM.h:324
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)

References MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::storage, MaterialPropertyLib::temperature, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

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

304{
305 assert(x.size() == dof_table.size());
306 auto const n_processes = dof_table.size();
307
308 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes;
309 indices_of_all_coupled_processes.reserve(n_processes);
310 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
311 {
312 auto const indices =
313 NumLib::getIndices(this->_element.getID(), *dof_table[process_id]);
314 assert(!indices.empty());
315 indices_of_all_coupled_processes.push_back(indices);
316 }
317 auto const local_xs =
318 getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
319
320 return this->getIntPtDarcyVelocityLocal(t, local_xs, cache);
321}
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::getCoupledLocalSolutions(), and NumLib::getIndices().


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