3.3.0
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
32
33namespace Dumux {
34
45template<class TypeTag>
47{
49
53 enum
54 {
55 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
56 };
58
59 using Element = typename GridView::template Codim<0>::Entity;
60
66 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
67
69
70 enum
71 {
72 pw = Indices::pressureW,
73 pn = Indices::pressureNw,
74 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
75 wPhaseIdx = Indices::wPhaseIdx,
76 nPhaseIdx = Indices::nPhaseIdx,
77 saturationIdx = Indices::saturationIdx,
78 satEqIdx = Indices::satEqIdx,
79 pressureEqIdx = Indices::pressureEqIdx
80 };
81
82 using FieldVector = Dune::FieldVector<Scalar, dimWorld>;
83
84public:
85
86 MimeticOperatorAssemblerTwoP(const GridView& gridView) :
87 ParentType(gridView)
88 {
89 }
90
91 // TODO doc me!
92 template<class Vector>
93 void calculatePressure(LocalStiffness& loc, Vector& u, Problem& problem)
94 {
95 Dune::FieldVector<Scalar, 2 * dim> velocityW(0);
96 Dune::FieldVector<Scalar, 2 * dim> velocityNw(0);
97 Dune::FieldVector<Scalar, 2 * dim> pressTrace(0);
98 Dune::FieldVector<Scalar, 2 * dim> gravPotTrace(0);
99
100 const auto firstElement = *this->gridView_.template begin<0>();
101 FluidState fluidState;
102 fluidState.setPressure(wPhaseIdx, problem.referencePressure(firstElement));
103 fluidState.setPressure(nPhaseIdx, problem.referencePressure(firstElement));
104 fluidState.setTemperature(problem.temperature(firstElement));
105 fluidState.setSaturation(wPhaseIdx, 1.);
106 fluidState.setSaturation(nPhaseIdx, 0.);
107 Scalar densityDiff = FluidSystem::density(fluidState, nPhaseIdx) - FluidSystem::density(fluidState, wPhaseIdx);
108 Scalar viscosityW = FluidSystem::viscosity(fluidState, wPhaseIdx);
109 Scalar viscosityNw = FluidSystem::viscosity(fluidState, nPhaseIdx);
110
111 //reset velocity
112 for (int i = 0; i < problem.gridView().size(0); i++)
113 {
114 problem.variables().cellData(i).fluxData().resetVelocity();
115 }
116
117 // run over all level elements
118 for (const auto& element : elements(this->gridView_))
119 {
120 int eIdxGlobal = problem.variables().index(element);
121
122 CellData& cellData = problem.variables().cellData(eIdxGlobal);
123 FieldVector globalPos = element.geometry().center();
124
125 // get local to global id map and pressure traces
126 for (const auto& intersection : intersections(problem.gridView(), element))
127 {
128 int indexInInside = intersection.indexInInside();
129
130 int fIdxGlobal = this->faceMapper_.subIndex(element, indexInInside, 1);
131
132 pressTrace[indexInInside] = u[fIdxGlobal];
133
134 gravPotTrace[indexInInside] = (problem.bBoxMax() - intersection.geometry().center()) * problem.gravity() * densityDiff;
135 }
136
137 switch (pressureType)
138 {
139 case pw:
140 {
141 Scalar potW = loc.constructPressure(element, pressTrace);
142 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
143 Scalar potNw = potW + gravPot;
144
145 cellData.setPotential(wPhaseIdx, potW);
146 cellData.setPotential(nPhaseIdx, potNw);
147
148 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, wPhaseIdx);
149
150 cellData.setPressure(wPhaseIdx, potW - gravPot);
151
152 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, nPhaseIdx);
153
154 cellData.setPressure(nPhaseIdx, potNw - gravPot);
155
156 break;
157 }
158 case pn:
159 {
160 Scalar potNw = loc.constructPressure(element, pressTrace);
161 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
162 Scalar potW = potNw - gravPot;
163
164 cellData.setPotential(nPhaseIdx, potNw);
165 cellData.setPotential(wPhaseIdx, potW);
166
167 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, wPhaseIdx);
168
169 cellData.setPressure(wPhaseIdx, potW - gravPot);
170
171 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, nPhaseIdx);
172
173 cellData.setPressure(nPhaseIdx, potNw - gravPot);
174
175 break;
176 }
177 }
178
179 // velocity reconstruction: !!! The velocity which is not reconstructed from the primary
180 // pressure variable can be slightly wrong and not conservative!!!!
181 // -> Should not be used for transport!!
182 switch (pressureType)
183 {
184 case pw:
185 {
186 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(wPhaseIdx));
187 pressTrace += gravPotTrace;
188 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(nPhaseIdx));
189
190 break;
191 }
192 case pn:
193 {
194 loc.constructVelocity(element, velocityW, pressTrace, cellData.potential(nPhaseIdx));
195 pressTrace -= gravPotTrace;
196 loc.constructVelocity(element, velocityNw, pressTrace, cellData.potential(wPhaseIdx));
197
198 break;
199 }
200 }
201
202 for (const auto& intersection : intersections(problem.gridView(), element))
203 {
204 int idxInInside = intersection.indexInInside();
205
206 cellData.fluxData().setUpwindPotential(wPhaseIdx, idxInInside, velocityW[idxInInside]);
207 cellData.fluxData().setUpwindPotential(nPhaseIdx, idxInInside, velocityNw[idxInInside]);
208
209 Scalar mobilityW = 0;
210 Scalar mobilityNw = 0;
211
212 if (intersection.neighbor())
213 {
214 int neighborIdx = problem.variables().index(intersection.outside());
215
216 CellData& cellDataNeighbor = problem.variables().cellData(neighborIdx);
217
218 mobilityW =
219 (velocityW[idxInInside] >= 0.) ? cellData.mobility(wPhaseIdx) :
220 cellDataNeighbor.mobility(wPhaseIdx);
221 mobilityNw =
222 (velocityNw[idxInInside] >= 0.) ? cellData.mobility(nPhaseIdx) :
223 cellDataNeighbor.mobility(nPhaseIdx);
224
225 if (velocityW[idxInInside] >= 0.)
226 {
227 FieldVector velocity(intersection.centerUnitOuterNormal());
228 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[idxInInside];
229 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
230 cellDataNeighbor.fluxData().addVelocity(wPhaseIdx, intersection.indexInOutside(), velocity);
231 }
232 if (velocityNw[idxInInside] >= 0.)
233 {
234 FieldVector velocity(intersection.centerUnitOuterNormal());
235 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[idxInInside];
236 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
237 cellDataNeighbor.fluxData().addVelocity(nPhaseIdx, intersection.indexInOutside(), velocity);
238 }
239
240 cellData.fluxData().setVelocityMarker(idxInInside);
241 }
242 else
243 {
244 BoundaryTypes bctype;
245 problem.boundaryTypes(bctype, intersection);
246 if (bctype.isDirichlet(satEqIdx))
247 {
248 PrimaryVariables boundValues(0.0);
249 problem.dirichlet(boundValues, intersection);
250
251 // old material law interface is deprecated: Replace this by
252 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
253 // after the release of 3.3, when the deprecated interface is no longer supported
254 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem.spatialParams(), element);
255
256 if (velocityW[idxInInside] >= 0.)
257 {
258 mobilityW = cellData.mobility(wPhaseIdx);
259 }
260 else
261 {
262 mobilityW = fluidMatrixInteraction.krw(boundValues[saturationIdx]) / viscosityW;
263 }
264
265 if (velocityNw[idxInInside] >= 0.)
266 {
267 mobilityNw = cellData.mobility(nPhaseIdx);
268 }
269 else
270 {
271 mobilityNw = fluidMatrixInteraction.krn(boundValues[saturationIdx]) / viscosityNw;
272 }
273 }
274 else
275 {
276 mobilityW = cellData.mobility(wPhaseIdx);
277 mobilityNw = cellData.mobility(nPhaseIdx);
278 }
279
280 FieldVector velocity(intersection.centerUnitOuterNormal());
281 velocity *= mobilityW / (mobilityW + mobilityNw) * velocityW[idxInInside];
282 cellData.fluxData().setVelocity(wPhaseIdx, idxInInside, velocity);
283
284 velocity = 0;
285 velocity = intersection.centerUnitOuterNormal();
286 velocity *= mobilityNw / (mobilityW + mobilityNw) * velocityNw[idxInInside];
287 cellData.fluxData().setVelocity(nPhaseIdx, idxInInside, velocity);
288 cellData.fluxData().setVelocityMarker(idxInInside);
289 }
290 }
291 }
292 }
293};
294} // end namespace Dumux
295#endif
Helpers for deprecation.
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 (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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:66
const GridView gridView_
Definition: croperator.hh:327
FaceMapper faceMapper_
Definition: croperator.hh:328
Levelwise assembler.
Definition: operator.hh:47
MimeticOperatorAssemblerTwoP(const GridView &gridView)
Definition: operator.hh:86
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operator.hh:93
Specifies the properties for immiscible 2p diffusion/pressure models.
Defines the basic properties required for a mimetic method.