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

Detailed Description

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

Definition at line 43 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 47 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 68 of file WellboreCompensateNeumannBoundaryConditionLocalAssembler.h.

73 {
75 _local_rhs.setZero();
76
77 unsigned const n_integration_points =
78 Base::_integration_method.getNumberOfPoints();
79
80 auto const indices_current_variable =
83 mesh_item_id, *_data.dof_table_boundary_pressure);
85 mesh_item_id, *_data.dof_table_boundary_velocity);
87 mesh_item_id, *_data.dof_table_boundary_enthalpy);
88
95
96 auto const& medium = *_data.media_map.getMedium(_element.getID());
97 auto const& liquid_phase = medium.phase("AqueousLiquid");
98 auto const& gas_phase = medium.phase("Gas");
99
101 pos.setElementID(_element.getID());
102
104
105 for (unsigned ip = 0; ip < n_integration_points; ip++)
106 {
108 auto const& N = n_and_weight.N;
109 auto const& w = n_and_weight.weight;
110
111 double pressure_int_pt = 0.0;
112 double velocity_int_pt = 0.0;
113 double enthalpy_int_pt = 0.0;
114
121
122 vars.liquid_phase_pressure = pressure_int_pt;
123 vars.enthalpy = enthalpy_int_pt;
124
125 double liquid_water_density =
127 .property(
129 .template value<double>(vars, pos, 0, 0);
130
131 double const vapour_water_density =
133 .property(
135 .template value<double>(vars, pos, 0, 0);
136
137 double const h_sat_liq_w =
139 .property(
141 .template value<double>(vars, pos, 0, 0);
142
143 double const h_sat_vap_w =
145 .property(
147 .template value<double>(vars, pos, 0, 0);
148
149 double const dryness =
152
153 double const T_int_pt =
154 (dryness == 0)
156 .property(
158 .template value<double>(vars, pos, 0, 0)
159 : gas_phase
162 .template value<double>(vars, pos, 0, 0);
163
164 vars.temperature = T_int_pt;
165
166 // For the calculation of the void fraction of vapour,
167 // see Rohuani, Z., and E. Axelsson. "Calculation of volume void
168 // fraction in a subcooled and quality region." International
169 // Journal of Heat and Mass Transfer 17 (1970): 383-393.
170
171 // Profile parameter of drift flux
172 double const C_0 = 1 + 0.12 * (1 - dryness);
173
174 // For the surface tension calculation, see
175 // Cooper, J. R., and R. B. Dooley. "IAPWS release on surface
176 // tension of ordinary water substance." International Association
177 // for the Properties of Water and Steam (1994).
178 double const sigma_gl = 0.2358 *
179 std::pow((1 - T_int_pt / 647.096), 1.256) *
180 (1 - 0.625 * (1 - T_int_pt / 647.096));
181 // drift flux velocity
182 double const u_gu =
183 1.18 * (1 - dryness) *
184 std::pow((9.81) * sigma_gl *
186 0.25) /
188
189 // solving void fraction of vapour using local Newton
190 // iteration.
191 double alpha = 0;
192 if (dryness != 0)
193 {
194 // Local Newton solver
195 using LocalJacobianMatrix =
200
202
204 {
207 (1 - alpha) * liquid_water_density) *
211 (1 - alpha) * liquid_water_density) *
213 alpha * C_0 * (1 - dryness) *
216 (1 - alpha) * liquid_water_density) *
220 };
221
223 {
224 jacobian(0) =
228 C_0 * (1 - dryness) * vapour_water_density) *
230 (1 - 2 * alpha) * liquid_water_density) *
233 };
234
235 auto const update_solution =
236 [&](LocalUnknownVector const& increment)
237 {
238 // increment solution vectors
239 alpha += increment[0];
240 };
241
242 const int maximum_iterations(20);
243 const double residuum_tolerance(1.e-10);
244 const double increment_tolerance(0);
245
251
252 auto const success_iterations = newton_solver.solve(J_loc);
253
255 {
256 WARN(
257 "Attention! Steam void fraction has not been correctly "
258 "calculated!");
259 }
260 }
261
262 if (alpha == 0)
263 {
267 .template value<double>(vars, pos, 0, 0);
268 }
269
270 double const mix_density = vapour_water_density * alpha +
272
273 // slip parameter between two phases
274 double const gamma =
276 mix_density / (1 - alpha) /
278 (1 - alpha * C_0) * liquid_water_density),
279 2) *
280 std::pow((C_0 - 1) * velocity_int_pt + u_gu, 2);
281
282 double const neumann_ip_values =
283 _data.pressure_coefficient * mix_density * velocity_int_pt +
284 _data.velocity_coefficient *
286 _data.enthalpy_coefficient * mix_density * velocity_int_pt *
288 _local_rhs.noalias() += N.transpose() * neumann_ip_values * w;
289 }
290
292 }
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:42
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

template<typename ShapeFunction, int GlobalDim>
MeshLib::Element const& ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_element
private

◆ _local_matrix_size

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

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