3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
operator.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18*****************************************************************************/
24#ifndef DUMUX_MIMETICOPERATOR2P_HH
25#define DUMUX_MIMETICOPERATOR2P_HH
26
27#include "croperator.hh"
30
31namespace Dumux {
32
43template<class TypeTag>
45{
47
51 enum
52 {
53 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
54 };
56
57 using Element = typename GridView::template Codim<0>::Entity;
58
64 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
65
67
68 enum
69 {
70 pw = Indices::pressureW,
71 pn = Indices::pressureNw,
72 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
73 wPhaseIdx = Indices::wPhaseIdx,
74 nPhaseIdx = Indices::nPhaseIdx,
75 saturationIdx = Indices::saturationIdx,
76 satEqIdx = Indices::satEqIdx,
77 pressureEqIdx = Indices::pressureEqIdx
78 };
79
80 using FieldVector = Dune::FieldVector<Scalar, dimWorld>;
81
82public:
83
84 MimeticOperatorAssemblerTwoP(const GridView& gridView) :
85 ParentType(gridView)
86 {
87 }
88
89 // TODO doc me!
90 template<class Vector>
91 void calculatePressure(LocalStiffness& loc, Vector& u, Problem& problem)
92 {
93 Dune::FieldVector<Scalar, 2 * dim> velocityW(0);
94 Dune::FieldVector<Scalar, 2 * dim> velocityNw(0);
95 Dune::FieldVector<Scalar, 2 * dim> pressTrace(0);
96 Dune::FieldVector<Scalar, 2 * dim> gravPotTrace(0);
97
98 const auto firstElement = *this->gridView_.template begin<0>();
99 FluidState fluidState;
100 fluidState.setPressure(wPhaseIdx, problem.referencePressure(firstElement));
101 fluidState.setPressure(nPhaseIdx, problem.referencePressure(firstElement));
102 fluidState.setTemperature(problem.temperature(firstElement));
103 fluidState.setSaturation(wPhaseIdx, 1.);
104 fluidState.setSaturation(nPhaseIdx, 0.);
105 Scalar densityDiff = FluidSystem::density(fluidState, nPhaseIdx) - FluidSystem::density(fluidState, wPhaseIdx);
106 Scalar viscosityW = FluidSystem::viscosity(fluidState, wPhaseIdx);
107 Scalar viscosityNw = FluidSystem::viscosity(fluidState, nPhaseIdx);
108
109 //reset velocity
110 for (int i = 0; i < problem.gridView().size(0); i++)
111 {
112 problem.variables().cellData(i).fluxData().resetVelocity();
113 }
114
115 // run over all level elements
116 for (const auto& element : elements(this->gridView_))
117 {
118 int eIdxGlobal = problem.variables().index(element);
119
120 CellData& cellData = problem.variables().cellData(eIdxGlobal);
121 FieldVector globalPos = element.geometry().center();
122
123 // get local to global id map and pressure traces
124 for (const auto& intersection : intersections(problem.gridView(), element))
125 {
126 int indexInInside = intersection.indexInInside();
127
128 int fIdxGlobal = this->faceMapper_.subIndex(element, indexInInside, 1);
129
130 pressTrace[indexInInside] = u[fIdxGlobal];
131
132 gravPotTrace[indexInInside] = (problem.bBoxMax() - intersection.geometry().center()) * problem.gravity() * densityDiff;
133 }
134
135 switch (pressureType)
136 {
137 case pw:
138 {
139 Scalar potW = loc.constructPressure(element, pressTrace);
140 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
141 Scalar potNw = potW + gravPot;
142
143 cellData.setPotential(wPhaseIdx, potW);
144 cellData.setPotential(nPhaseIdx, potNw);
145
146 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, wPhaseIdx);
147
148 cellData.setPressure(wPhaseIdx, potW - gravPot);
149
150 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, nPhaseIdx);
151
152 cellData.setPressure(nPhaseIdx, potNw - gravPot);
153
154 break;
155 }
156 case pn:
157 {
158 Scalar potNw = loc.constructPressure(element, pressTrace);
159 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
160 Scalar potW = potNw - gravPot;
161
162 cellData.setPotential(nPhaseIdx, potNw);
163 cellData.setPotential(wPhaseIdx, potW);
164
165 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, wPhaseIdx);
166
167 cellData.setPressure(wPhaseIdx, potW - gravPot);
168
169 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, nPhaseIdx);
170
171 cellData.setPressure(nPhaseIdx, potNw - gravPot);
172
173 break;
174 }
175 }
176
177 // velocity reconstruction: !!! The velocity which is not reconstructed from the primary
178 // pressure variable can be slightly wrong and not conservative!!!!
179 // -> Should not be used for transport!!
180 switch (pressureType)
181 {
182 case pw:
183 {
184 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(wPhaseIdx));
185 pressTrace += gravPotTrace;
186 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(nPhaseIdx));
187
188 break;
189 }
190 case pn:
191 {
192 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(nPhaseIdx));
193 pressTrace -= gravPotTrace;
194 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(wPhaseIdx));
195
196 break;
197 }
198 }
199
200 for (const auto& intersection : intersections(problem.gridView(), element))
201 {
202 int idxInInside = intersection.indexInInside();
203
204 cellData.fluxData().setUpwindPotential(wPhaseIdx, idxInInside, velocityW[idxInInside]);
205 cellData.fluxData().setUpwindPotential(nPhaseIdx, idxInInside, velocityNw[idxInInside]);
206
207 Scalar mobilityW = 0;
208 Scalar mobilityNw = 0;
209
210 if (intersection.neighbor())
211 {
212 int neighborIdx = problem.variables().index(intersection.outside());
213
214 CellData& cellDataNeighbor = problem.variables().cellData(neighborIdx);
215
216 mobilityW =
217 (velocityW[idxInInside] >= 0.) ? cellData.mobility(wPhaseIdx) :
218 cellDataNeighbor.mobility(wPhaseIdx);
219 mobilityNw =
220 (velocityNw[idxInInside] >= 0.) ? cellData.mobility(nPhaseIdx) :
221 cellDataNeighbor.mobility(nPhaseIdx);
222
223 if (velocityW[idxInInside] >= 0.)
224 {
225 FieldVector velocity(intersection.centerUnitOuterNormal());
226 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[idxInInside];
227 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
228 cellDataNeighbor.fluxData().addVelocity(wPhaseIdx, intersection.indexInOutside(), velocity);
229 }
230 if (velocityNw[idxInInside] >= 0.)
231 {
232 FieldVector velocity(intersection.centerUnitOuterNormal());
233 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[idxInInside];
234 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
235 cellDataNeighbor.fluxData().addVelocity(nPhaseIdx, intersection.indexInOutside(), velocity);
236 }
237
238 cellData.fluxData().setVelocityMarker(idxInInside);
239 }
240 else
241 {
242 BoundaryTypes bctype;
243 problem.boundaryTypes(bctype, intersection);
244 if (bctype.isDirichlet(satEqIdx))
245 {
246 PrimaryVariables boundValues(0.0);
247 problem.dirichlet(boundValues, intersection);
248
249 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
250 if (velocityW[idxInInside] >= 0.)
251 {
252 mobilityW = cellData.mobility(wPhaseIdx);
253 }
254 else
255 {
256 mobilityW = fluidMatrixInteraction.krw(boundValues[saturationIdx]) / viscosityW;
257 }
258
259 if (velocityNw[idxInInside] >= 0.)
260 {
261 mobilityNw = cellData.mobility(nPhaseIdx);
262 }
263 else
264 {
265 mobilityNw = fluidMatrixInteraction.krn(boundValues[saturationIdx]) / viscosityNw;
266 }
267 }
268 else
269 {
270 mobilityW = cellData.mobility(wPhaseIdx);
271 mobilityNw = cellData.mobility(nPhaseIdx);
272 }
273
274 FieldVector velocity(intersection.centerUnitOuterNormal());
275 velocity *= mobilityW / (mobilityW + mobilityNw) * velocityW[idxInInside];
276 cellData.fluxData().setVelocity(wPhaseIdx, idxInInside, velocity);
277
278 velocity = 0;
279 velocity = intersection.centerUnitOuterNormal();
280 velocity *= mobilityNw / (mobilityW + mobilityNw) * velocityNw[idxInInside];
281 cellData.fluxData().setVelocity(nPhaseIdx, idxInInside, velocity);
282 cellData.fluxData().setVelocityMarker(idxInInside);
283 }
284 }
285 }
286 }
287};
288} // end namespace Dumux
289#endif
Defines a class for Crozieux-Raviart piecewise linear finite element functions.
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property
Definition: propertysystem.hh:141
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Extends CROperatorBase by a generic methods to assemble global stiffness matrix from local stiffness ...
Definition: croperator.hh:67
const GridView gridView_
Definition: croperator.hh:331
FaceMapper faceMapper_
Definition: croperator.hh:332
Levelwise assembler.
Definition: operator.hh:45
MimeticOperatorAssemblerTwoP(const GridView &gridView)
Definition: operator.hh:84
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operator.hh:91
Specifies the properties for immiscible 2p diffusion/pressure models.
Defines the basic properties required for a mimetic method.