OGS
TwoPhaseFlowWithPPLocalAssembler-impl.h
Go to the documentation of this file.
1 
31 #pragma once
32 
33 #include "MaterialLib/MPL/Medium.h"
38 
39 namespace ProcessLib
40 {
41 namespace TwoPhaseFlowWithPP
42 {
43 namespace MPL = MaterialPropertyLib;
44 
45 template <typename ShapeFunction, typename IntegrationMethod, int GlobalDim>
46 void TwoPhaseFlowWithPPLocalAssembler<
47  ShapeFunction, IntegrationMethod,
48  GlobalDim>::assemble(double const t, double const dt,
49  std::vector<double> const& local_x,
50  std::vector<double> const& /*local_xdot*/,
51  std::vector<double>& local_M_data,
52  std::vector<double>& local_K_data,
53  std::vector<double>& local_b_data)
54 {
55  auto const local_matrix_size = local_x.size();
56 
57  assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
58 
59  auto local_M = MathLib::createZeroedMatrix<LocalMatrixType>(
60  local_M_data, local_matrix_size, local_matrix_size);
61  auto local_K = MathLib::createZeroedMatrix<LocalMatrixType>(
62  local_K_data, local_matrix_size, local_matrix_size);
63  auto local_b = MathLib::createZeroedVector<LocalVectorType>(
64  local_b_data, local_matrix_size);
65 
66  auto Mgp =
67  local_M.template block<nonwet_pressure_size, nonwet_pressure_size>(
68  nonwet_pressure_matrix_index, nonwet_pressure_matrix_index);
69  auto Mgpc = local_M.template block<nonwet_pressure_size, cap_pressure_size>(
70  nonwet_pressure_matrix_index, cap_pressure_matrix_index);
71 
72  auto Mlpc = local_M.template block<cap_pressure_size, cap_pressure_size>(
73  cap_pressure_matrix_index, cap_pressure_matrix_index);
74 
75  NodalMatrixType laplace_operator =
76  NodalMatrixType::Zero(ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
77 
78  auto Kgp =
79  local_K.template block<nonwet_pressure_size, nonwet_pressure_size>(
80  nonwet_pressure_matrix_index, nonwet_pressure_matrix_index);
81 
82  auto Klp = local_K.template block<cap_pressure_size, nonwet_pressure_size>(
83  cap_pressure_matrix_index, nonwet_pressure_matrix_index);
84 
85  auto Klpc = local_K.template block<cap_pressure_size, cap_pressure_size>(
86  cap_pressure_matrix_index, cap_pressure_matrix_index);
87 
88  auto Bg = local_b.template segment<nonwet_pressure_size>(
89  nonwet_pressure_matrix_index);
90 
91  auto Bl =
92  local_b.template segment<cap_pressure_size>(cap_pressure_matrix_index);
93 
94  auto const& medium = *_process_data.media_map->getMedium(_element.getID());
95  auto const& liquid_phase = medium.phase("AqueousLiquid");
96  auto const& gas_phase = medium.phase("Gas");
97 
98  unsigned const n_integration_points =
99  _integration_method.getNumberOfPoints();
100 
102  pos.setElementID(_element.getID());
103 
104  for (unsigned ip = 0; ip < n_integration_points; ip++)
105  {
106  pos.setIntegrationPoint(ip);
107  auto const& w = _ip_data[ip].integration_weight;
108  auto const& dNdx = _ip_data[ip].dNdx;
109  auto const& M = _ip_data[ip].massOperator;
110 
111  double pc_int_pt = 0.;
112  double pn_int_pt = 0.;
113  NumLib::shapeFunctionInterpolate(local_x, _ip_data[ip].N, pn_int_pt,
114  pc_int_pt);
115 
116  _pressure_wet[ip] = pn_int_pt - pc_int_pt;
117  MPL::VariableArray variables;
118 
119  variables[static_cast<int>(MPL::Variable::phase_pressure)] = pn_int_pt;
120  variables[static_cast<int>(MPL::Variable::capillary_pressure)] =
121  pc_int_pt;
122  variables[static_cast<int>(MPL::Variable::temperature)] =
123  _process_data.temperature(t, pos)[0];
124  variables[static_cast<int>(MPL::Variable::molar_mass)] =
125  gas_phase.property(MPL::PropertyType::molar_mass)
126  .template value<double>(variables, pos, t, dt);
127 
128  auto const rho_nonwet =
129  gas_phase.property(MPL::PropertyType::density)
130  .template value<double>(variables, pos, t, dt);
131 
132  auto const rho_wet = liquid_phase.property(MPL::PropertyType::density)
133  .template value<double>(variables, pos, t, dt);
134  auto& Sw = _saturation[ip];
135  Sw = medium.property(MPL::PropertyType::saturation)
136  .template value<double>(variables, pos, t, dt);
137 
138  auto const dSw_dpc =
139  medium.property(MPL::PropertyType::saturation)
140  .template dValue<double>(
141  variables, MPL::Variable::capillary_pressure, pos, t, dt);
142 
143  auto const porosity =
144  medium.property(MPL::PropertyType::porosity)
145  .template value<double>(variables, pos, t, dt);
146 
147  auto const drhononwet_dpn =
148  gas_phase.property(MPL::PropertyType::density)
149  .template dValue<double>(
150  variables, MPL::Variable::phase_pressure, pos, t, dt);
151 
152  auto const k_rel_wet =
154  .template value<double>(variables, pos, t, dt);
155  auto const k_rel_nonwet =
156  medium
157  .property(
159  .template value<double>(variables, pos, t, dt);
160 
161  auto const mu_nonwet =
162  gas_phase.property(MPL::PropertyType::viscosity)
163  .template value<double>(variables, pos, t, dt);
164 
165  auto const lambda_nonwet = k_rel_nonwet / mu_nonwet;
166 
167  auto const mu_wet = liquid_phase.property(MPL::PropertyType::viscosity)
168  .template value<double>(variables, pos, t, dt);
169 
170  auto const lambda_wet = k_rel_wet / mu_wet;
171 
172  auto const K = MPL::formEigenTensor<GlobalDim>(
173  medium.property(MPL::PropertyType::permeability)
174  .template value<double>(variables, pos, t, dt));
175 
176  Mgp.noalias() += porosity * (1 - Sw) * drhononwet_dpn * M;
177  Mgpc.noalias() += -porosity * rho_nonwet * dSw_dpc * M;
178  Mlpc.noalias() += porosity * dSw_dpc * rho_wet * M;
179 
180  laplace_operator.noalias() = dNdx.transpose() * K * dNdx * w;
181 
182  Kgp.noalias() += rho_nonwet * lambda_nonwet * laplace_operator;
183  Klp.noalias() += rho_wet * lambda_wet * laplace_operator;
184  Klpc.noalias() += -rho_wet * lambda_wet * laplace_operator;
185 
186  if (_process_data.has_gravity)
187  {
188  auto const& b = _process_data.specific_body_force;
189 
190  NodalVectorType gravity_operator = dNdx.transpose() * K * b * w;
191  Bg.noalias() +=
192  rho_nonwet * rho_nonwet * lambda_nonwet * gravity_operator;
193  Bl.noalias() += rho_wet * rho_wet * lambda_wet * gravity_operator;
194  } // end of has gravity
195  }
196  if (_process_data.has_mass_lumping)
197  {
198  for (unsigned row = 0; row < Mgpc.cols(); row++)
199  {
200  for (unsigned column = 0; column < Mgpc.cols(); column++)
201  {
202  if (row != column)
203  {
204  Mgpc(row, row) += Mgpc(row, column);
205  Mgpc(row, column) = 0.0;
206  Mgp(row, row) += Mgp(row, column);
207  Mgp(row, column) = 0.0;
208  Mlpc(row, row) += Mlpc(row, column);
209  Mlpc(row, column) = 0.0;
210  }
211  }
212  }
213  } // end of mass-lumping
214 }
215 
216 } // namespace TwoPhaseFlowWithPP
217 } // namespace ProcessLib
Definition of the PiecewiseLinearInterpolation class.
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
Definition: VariableType.h:108
@ relative_permeability_nonwetting_phase
Definition: PropertyType.h:77
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
Definition: Interpolation.h:79