OGS
ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim > Class Template Referencefinal

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >

Definition at line 42 of file WellboreCompensateNeumannBoundaryConditionLocalAssembler.h.

#include <WellboreCompensateNeumannBoundaryConditionLocalAssembler.h>

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

Public Member Functions

 WellboreCompensateNeumannBoundaryConditionLocalAssembler (MeshLib::Element const &e, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, WellboreCompensateNeumannBoundaryConditionData const &data)
void assemble (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, double const, std::vector< GlobalVector * > const &x, int const process_id, GlobalMatrix *, GlobalVector &b, GlobalMatrix *) override
Public Member Functions inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >
 GenericNaturalBoundaryConditionLocalAssembler (MeshLib::Element const &e, bool is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method)
Public Member Functions inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssemblerInterface
virtual ~GenericNaturalBoundaryConditionLocalAssemblerInterface ()=default

Private Types

using Base
using NodalVectorType = typename Base::NodalVectorType
using NodalMatrixType = typename Base::NodalMatrixType

Private Attributes

MeshLib::Element const & _element
WellboreCompensateNeumannBoundaryConditionData const & _data
Base::NodalVectorType _local_matrix_size

Additional Inherited Members

Protected Types inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >
using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
Protected Attributes inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< NAndWeight, Eigen::aligned_allocator< NAndWeight > > const _ns_and_weights
MeshLib::Element const & _element

Member Typedef Documentation

◆ Base

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::Base
private
Initial value:
GenericNaturalBoundaryConditionLocalAssembler(MeshLib::Element const &e, bool is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method)

Definition at line 46 of file WellboreCompensateNeumannBoundaryConditionLocalAssembler.h.

◆ NodalMatrixType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::NodalMatrixType = typename Base::NodalMatrixType
private

◆ NodalVectorType

template<typename ShapeFunction, int GlobalDim>
using ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::NodalVectorType = typename Base::NodalVectorType
private

Constructor & Destructor Documentation

◆ WellboreCompensateNeumannBoundaryConditionLocalAssembler()

template<typename ShapeFunction, int GlobalDim>
ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::WellboreCompensateNeumannBoundaryConditionLocalAssembler ( MeshLib::Element const & e,
std::size_t const local_matrix_size,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
WellboreCompensateNeumannBoundaryConditionData const & data )
inline

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int GlobalDim>
void ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::assemble ( std::size_t const mesh_item_id,
NumLib::LocalToGlobalIndexMap const & dof_table_boundary,
double const ,
std::vector< GlobalVector * > const & x,
int const process_id,
GlobalMatrix * ,
GlobalVector & b,
GlobalMatrix *  )
inlineoverridevirtual

Implements ProcessLib::GenericNaturalBoundaryConditionLocalAssemblerInterface.

Definition at line 67 of file WellboreCompensateNeumannBoundaryConditionLocalAssembler.h.

72 {
74 _local_rhs.setZero();
75
76 unsigned const n_integration_points =
77 Base::_integration_method.getNumberOfPoints();
78
79 auto const indices_current_variable =
82 mesh_item_id, *_data.dof_table_boundary_pressure);
84 mesh_item_id, *_data.dof_table_boundary_velocity);
86 mesh_item_id, *_data.dof_table_boundary_enthalpy);
87
94
95 auto const& medium = *_data.media_map.getMedium(_element.getID());
96 auto const& liquid_phase = medium.phase("AqueousLiquid");
97 auto const& gas_phase = medium.phase("Gas");
98
100 pos.setElementID(_element.getID());
101
103
104 for (unsigned ip = 0; ip < n_integration_points; ip++)
105 {
107 auto const& N = n_and_weight.N;
108 auto const& w = n_and_weight.weight;
109
110 double pressure_int_pt = 0.0;
111 double velocity_int_pt = 0.0;
112 double enthalpy_int_pt = 0.0;
113
120
121 vars.liquid_phase_pressure = pressure_int_pt;
122 vars.enthalpy = enthalpy_int_pt;
123
124 double liquid_water_density =
126 .property(
128 .template value<double>(vars, pos, 0, 0);
129
130 double const vapour_water_density =
132 .property(
134 .template value<double>(vars, pos, 0, 0);
135
136 double const h_sat_liq_w =
138 .property(
140 .template value<double>(vars, pos, 0, 0);
141
142 double const h_sat_vap_w =
144 .property(
146 .template value<double>(vars, pos, 0, 0);
147
148 double const dryness =
151
152 double const T_int_pt =
153 (dryness == 0)
155 .property(
157 .template value<double>(vars, pos, 0, 0)
158 : gas_phase
161 .template value<double>(vars, pos, 0, 0);
162
163 vars.temperature = T_int_pt;
164
165 // For the calculation of the void fraction of vapour,
166 // see Rohuani, Z., and E. Axelsson. "Calculation of volume void
167 // fraction in a subcooled and quality region." International
168 // Journal of Heat and Mass Transfer 17 (1970): 383-393.
169
170 // Profile parameter of drift flux
171 double const C_0 = 1 + 0.12 * (1 - dryness);
172
173 // For the surface tension calculation, see
174 // Cooper, J. R., and R. B. Dooley. "IAPWS release on surface
175 // tension of ordinary water substance." International Association
176 // for the Properties of Water and Steam (1994).
177 double const sigma_gl = 0.2358 *
178 std::pow((1 - T_int_pt / 647.096), 1.256) *
179 (1 - 0.625 * (1 - T_int_pt / 647.096));
180 // drift flux velocity
181 double const u_gu =
182 1.18 * (1 - dryness) *
183 std::pow((9.81) * sigma_gl *
185 0.25) /
187
188 // solving void fraction of vapour using local Newton
189 // iteration.
190 double alpha = 0;
191 if (dryness != 0)
192 {
193 // Local Newton solver
194 using LocalJacobianMatrix =
199
201
203 {
206 (1 - alpha) * liquid_water_density) *
210 (1 - alpha) * liquid_water_density) *
212 alpha * C_0 * (1 - dryness) *
215 (1 - alpha) * liquid_water_density) *
219 };
220
222 {
223 jacobian(0) =
227 C_0 * (1 - dryness) * vapour_water_density) *
229 (1 - 2 * alpha) * liquid_water_density) *
232 };
233
234 auto const update_solution =
235 [&](LocalUnknownVector const& increment)
236 {
237 // increment solution vectors
238 alpha += increment[0];
239 };
240
241 const int maximum_iterations(20);
242 const double residuum_tolerance(1.e-10);
243 const double increment_tolerance(0);
244
250
251 auto const success_iterations = newton_solver.solve(J_loc);
252
254 {
255 WARN(
256 "Attention! Steam void fraction has not been correctly "
257 "calculated!");
258 }
259 }
260
261 if (alpha == 0)
262 {
266 .template value<double>(vars, pos, 0, 0);
267 }
268
269 double const mix_density = vapour_water_density * alpha +
271
272 // slip parameter between two phases
273 double const gamma =
275 mix_density / (1 - alpha) /
277 (1 - alpha * C_0) * liquid_water_density),
278 2) *
279 std::pow((C_0 - 1) * velocity_int_pt + u_gu, 2);
280
281 double const neumann_ip_values =
282 _data.coefficients.pressure * mix_density * velocity_int_pt +
283 _data.coefficients.velocity *
285 _data.coefficients.enthalpy * mix_density * velocity_int_pt *
287 _local_rhs.noalias() += N.transpose() * neumann_ip_values * w;
288 }
289
291 }
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
std::vector< NAndWeight, Eigen::aligned_allocator< NAndWeight > > const _ns_and_weights
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References _data, _element, ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_integration_method, _local_matrix_size, ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_ns_and_weights, MathLib::EigenVector::add(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::enthalpy, NumLib::getIndices(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::saturation_density, MaterialPropertyLib::saturation_enthalpy, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), NumLib::NewtonRaphson< LinearSolver, JacobianMatrixUpdate, ResidualUpdate, SolutionUpdate >::solve(), MaterialPropertyLib::temperature, MaterialPropertyLib::VariableArray::temperature, and WARN().

Member Data Documentation

◆ _data

◆ _element

◆ _local_matrix_size


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