OGS
TwoPhaseFlowWithPPLocalAssembler-impl.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
12
13namespace ProcessLib
14{
15namespace TwoPhaseFlowWithPP
16{
17namespace MPL = MaterialPropertyLib;
18
39template <typename ShapeFunction, int GlobalDim>
41 setInitialConditionsConcrete(Eigen::VectorXd const local_x,
42 double const t,
43 int const /*process_id*/)
44{
45 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
46 unsigned const n_integration_points =
47 _integration_method.getNumberOfPoints();
48
49 auto const pc =
51 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
53 pos.setElementID(_element.getID());
54 for (unsigned ip = 0; ip < n_integration_points; ip++)
55 {
56 auto const& N = _ip_data[ip].N;
57 MPL::VariableArray variables;
58 variables.capillary_pressure = N.dot(pc);
59 _saturation[ip] = medium.property(MPL::PropertyType::saturation)
60 .template value<double>(variables, pos, t, dt);
61 }
62}
63
64template <typename ShapeFunction, int GlobalDim>
66 double const t, double const dt, std::vector<double> const& local_x,
67 std::vector<double> const& /*local_x_prev*/,
68 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
69 std::vector<double>& local_b_data)
70{
71 auto const local_matrix_size = local_x.size();
72
73 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
74
76 local_M_data, local_matrix_size, local_matrix_size);
78 local_K_data, local_matrix_size, local_matrix_size);
80 local_b_data, local_matrix_size);
81
82 auto Mgp =
83 local_M.template block<nonwet_pressure_size, nonwet_pressure_size>(
85 auto Mgpc = local_M.template block<nonwet_pressure_size, cap_pressure_size>(
87
88 auto Mlpc = local_M.template block<cap_pressure_size, cap_pressure_size>(
90
91 NodalMatrixType laplace_operator =
92 NodalMatrixType::Zero(ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
93
94 auto Kgp =
95 local_K.template block<nonwet_pressure_size, nonwet_pressure_size>(
97
98 auto Klp = local_K.template block<cap_pressure_size, nonwet_pressure_size>(
100
101 auto Klpc = local_K.template block<cap_pressure_size, cap_pressure_size>(
103
104 auto Bg = local_b.template segment<nonwet_pressure_size>(
106
107 auto Bl =
108 local_b.template segment<cap_pressure_size>(cap_pressure_matrix_index);
109
110 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
111 auto const& liquid_phase =
113 auto const& gas_phase = medium.phase(MaterialPropertyLib::PhaseName::Gas);
114
115 unsigned const n_integration_points =
116 _integration_method.getNumberOfPoints();
117
119 pos.setElementID(_element.getID());
120
121 for (unsigned ip = 0; ip < n_integration_points; ip++)
122 {
123 auto const& w = _ip_data[ip].integration_weight;
124 auto const& N = _ip_data[ip].N;
125 auto const& dNdx = _ip_data[ip].dNdx;
126
127 NodalMatrixType M = N.transpose() * N * w;
128
129 double pc_int_pt = 0.;
130 double pn_int_pt = 0.;
131 NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, pn_int_pt,
132 pc_int_pt);
133
134 _pressure_wet[ip] = pn_int_pt - pc_int_pt;
135 MPL::VariableArray variables;
136
137 variables.gas_phase_pressure = pn_int_pt;
138 variables.capillary_pressure = pc_int_pt;
139 variables.temperature = _process_data.temperature(t, pos)[0];
140 variables.molar_mass =
141 gas_phase.property(MPL::PropertyType::molar_mass)
142 .template value<double>(variables, pos, t, dt);
143
144 auto const rho_nonwet =
145 gas_phase.property(MPL::PropertyType::density)
146 .template value<double>(variables, pos, t, dt);
147
148 auto const rho_wet = liquid_phase.property(MPL::PropertyType::density)
149 .template value<double>(variables, pos, t, dt);
150 auto& Sw = _saturation[ip];
151 Sw = medium.property(MPL::PropertyType::saturation)
152 .template value<double>(variables, pos, t, dt);
153
154 variables.liquid_saturation = Sw;
155 auto const dSw_dpc =
156 medium.property(MPL::PropertyType::saturation)
157 .template dValue<double>(
158 variables, MPL::Variable::capillary_pressure, pos, t, dt);
159
160 auto const porosity =
161 medium.property(MPL::PropertyType::porosity)
162 .template value<double>(variables, pos, t, dt);
163
164 auto const drhononwet_dpn =
165 gas_phase.property(MPL::PropertyType::density)
166 .template dValue<double>(
167 variables, MPL::Variable::gas_phase_pressure, pos, t, dt);
168
169 auto const k_rel_wet =
171 .template value<double>(variables, pos, t, dt);
172 auto const k_rel_nonwet =
173 medium
174 .property(
176 .template value<double>(variables, pos, t, dt);
177
178 auto const mu_nonwet =
179 gas_phase.property(MPL::PropertyType::viscosity)
180 .template value<double>(variables, pos, t, dt);
181
182 auto const lambda_nonwet = k_rel_nonwet / mu_nonwet;
183
184 auto const mu_wet = liquid_phase.property(MPL::PropertyType::viscosity)
185 .template value<double>(variables, pos, t, dt);
186
187 auto const lambda_wet = k_rel_wet / mu_wet;
188
189 auto const K = MPL::formEigenTensor<GlobalDim>(
190 medium.property(MPL::PropertyType::permeability)
191 .template value<double>(variables, pos, t, dt));
192
193 Mgp.noalias() += porosity * (1 - Sw) * drhononwet_dpn * M;
194 Mgpc.noalias() += -porosity * rho_nonwet * dSw_dpc * M;
195 Mlpc.noalias() += porosity * dSw_dpc * rho_wet * M;
196
197 laplace_operator.noalias() = dNdx.transpose() * K * dNdx * w;
198
199 Kgp.noalias() += rho_nonwet * lambda_nonwet * laplace_operator;
200 Klp.noalias() += rho_wet * lambda_wet * laplace_operator;
201 Klpc.noalias() += -rho_wet * lambda_wet * laplace_operator;
202
203 if (_process_data.has_gravity)
204 {
205 auto const& b = _process_data.specific_body_force;
206
207 NodalVectorType gravity_operator = dNdx.transpose() * K * b * w;
208 Bg.noalias() +=
209 rho_nonwet * rho_nonwet * lambda_nonwet * gravity_operator;
210 Bl.noalias() += rho_wet * rho_wet * lambda_wet * gravity_operator;
211 } // end of has gravity
212 }
213 if (_process_data.has_mass_lumping)
214 {
215 for (unsigned row = 0; row < Mgpc.cols(); row++)
216 {
217 for (unsigned column = 0; column < Mgpc.cols(); column++)
218 {
219 if (row != column)
220 {
221 Mgpc(row, row) += Mgpc(row, column);
222 Mgpc(row, column) = 0.0;
223 Mgp(row, row) += Mgp(row, column);
224 Mgp(row, column) = 0.0;
225 Mlpc(row, row) += Mlpc(row, column);
226 Mlpc(row, column) = 0.0;
227 }
228 }
229 }
230 } // end of mass-lumping
231}
232
233} // namespace TwoPhaseFlowWithPP
234} // namespace ProcessLib
void setElementID(std::size_t element_id)
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
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) override
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
@ relative_permeability_nonwetting_phase
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)