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