3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
operatoradaptive.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_MIMETICOPERATOR2PADAPTIVE_HH
25#define DUMUX_MIMETICOPERATOR2PADAPTIVE_HH
26
27#include "croperatoradaptive.hh"
30
31namespace Dumux {
32
43template<class TypeTag>
45{
47
51
52 enum
53 {
54 dim=GridView::dimension,
55 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 };
80
81 using FieldVector = Dune::FieldVector<Scalar, dimWorld>;
82public:
83
84 MimeticOperatorAssemblerTwoPAdaptive (const GridView& gridView)
85 : ParentType(gridView)
86 {}
87
88 // TODO doc me!
89 template<class Vector>
90 void calculatePressure(LocalStiffness& loc, Vector& u, Problem& problem)
91 {
92 Dune::DynamicVector<Scalar> velocityW(2*dim);
93 Dune::DynamicVector<Scalar> velocityNw(2*dim);
94 Dune::DynamicVector<Scalar> pressTraceW(2*dim);
95 Dune::DynamicVector<Scalar> pressTraceNw(2*dim);
96
97 const auto firstElement = *problem.gridView().template begin<0>();
98 FluidState fluidState;
99 fluidState.setPressure(wPhaseIdx, problem.referencePressure(firstElement));
100 fluidState.setPressure(nPhaseIdx, problem.referencePressure(firstElement));
101 fluidState.setTemperature(problem.temperature(firstElement));
102 fluidState.setSaturation(wPhaseIdx, 1.);
103 fluidState.setSaturation(nPhaseIdx, 0.);
104 Scalar densityDiff = FluidSystem::density(fluidState, nPhaseIdx) - FluidSystem::density(fluidState, wPhaseIdx);
105 Scalar viscosityW = FluidSystem::viscosity(fluidState, wPhaseIdx);
106 Scalar viscosityNw = FluidSystem::viscosity(fluidState, nPhaseIdx);
107
108 // reset velocity
109 for (int i = 0; i < problem.gridView().size(0); i++)
110 {
111 problem.variables().cellData(i).fluxData().resetVelocity();
112 }
113
114 // run over all level elements
115 for (const auto& element : elements(this->gridView_))
116 {
117 int eIdxGlobal = problem.variables().index(element);
118
119 unsigned int numFaces = this->intersectionMapper_.size(eIdxGlobal);
120
121 // get local to global id map and pressure traces
122 velocityW.resize(numFaces);
123 velocityNw.resize(numFaces);
124 pressTraceW.resize(numFaces);
125 pressTraceNw.resize(numFaces);
126
127 CellData& cellData = problem.variables().cellData(eIdxGlobal);
128 FieldVector globalPos = element.geometry().center();
129
130 int intersectionIdx = -1;
131 // get local to global id map and pressure traces
132 for (const auto& intersection : intersections(problem.gridView(), element))
133 {
134 ++intersectionIdx;
135
136 int fIdxGlobal = this->intersectionMapper_.subIndex(element, intersectionIdx);
137
138 Scalar pcPotFace = (problem.bBoxMax() - intersection.geometry().center()) * problem.gravity() * densityDiff;
139
140 switch (pressureType)
141 {
142 case pw:
143 {
144 pressTraceW[intersectionIdx] = u[fIdxGlobal];
145 pressTraceNw[intersectionIdx] = u[fIdxGlobal] + pcPotFace;
146 break;
147 }
148 case pn:
149 {
150 pressTraceNw[intersectionIdx] = u[fIdxGlobal];
151 pressTraceW[intersectionIdx] = u[fIdxGlobal] - pcPotFace;
152
153 break;
154 }
155 }
156 }
157
158 switch (pressureType)
159 {
160 case pw:
161 {
162 Scalar potW = loc.constructPressure(element, pressTraceW);
163 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
164 Scalar potNw = potW + gravPot;
165
166 cellData.setPotential(wPhaseIdx, potW);
167 cellData.setPotential(nPhaseIdx, potNw);
168
169 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, wPhaseIdx);
170
171 cellData.setPressure(wPhaseIdx, potW - gravPot);
172
173 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, nPhaseIdx);
174
175 cellData.setPressure(nPhaseIdx, potNw - gravPot);
176
177 break;
178 }
179 case pn:
180 {
181 Scalar potNw = loc.constructPressure(element, pressTraceNw);
182 Scalar gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * densityDiff;
183 Scalar potW = potNw - gravPot;
184
185 cellData.setPotential(nPhaseIdx, potNw);
186 cellData.setPotential(wPhaseIdx, potW);
187
188 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, wPhaseIdx);
189
190 cellData.setPressure(wPhaseIdx, potW - gravPot);
191
192 gravPot = (problem.bBoxMax() - globalPos) * problem.gravity() * FluidSystem::density(fluidState, nPhaseIdx);
193
194 cellData.setPressure(nPhaseIdx, potNw - gravPot);
195
196 break;
197 }
198 }
199
200 // velocity reconstruction: !!! The velocity which is not reconstructed from the primary
201 // pressure variable can be slightly wrong and not conservative!!!!
202 // -> Should not be used for transport!!
203 loc.constructVelocity(element, velocityW, pressTraceW, cellData.potential(wPhaseIdx));
204 loc.constructVelocity(element, velocityNw, pressTraceNw, cellData.potential(nPhaseIdx));
205
206 intersectionIdx = -1;
207 for (const auto& intersection : intersections(problem.gridView(), element))
208 {
209 ++intersectionIdx;
210 int idxInInside = intersection.indexInInside();
211
212 cellData.fluxData().addUpwindPotential(wPhaseIdx, idxInInside, velocityW[intersectionIdx]);
213 cellData.fluxData().addUpwindPotential(nPhaseIdx, idxInInside, velocityNw[intersectionIdx]);
214
215 Scalar mobilityW = 0;
216 Scalar mobilityNw = 0;
217
218 if (intersection.neighbor())
219 {
220 int neighborIdx = problem.variables().index(intersection.outside());
221
222 CellData& cellDataNeighbor = problem.variables().cellData(neighborIdx);
223
224 mobilityW =
225 (velocityW[intersectionIdx] >= 0.) ? cellData.mobility(wPhaseIdx) :
226 cellDataNeighbor.mobility(wPhaseIdx);
227 mobilityNw =
228 (velocityNw[intersectionIdx] >= 0.) ? cellData.mobility(nPhaseIdx) :
229 cellDataNeighbor.mobility(nPhaseIdx);
230
231 if (velocityW[intersectionIdx] >= 0.)
232 {
233 FieldVector velocity(intersection.centerUnitOuterNormal());
234 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[intersectionIdx];
235 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
236 cellDataNeighbor.fluxData().addVelocity(wPhaseIdx, intersection.indexInOutside(), velocity);
237 }
238 if (velocityNw[intersectionIdx] >= 0.)
239 {
240 FieldVector velocity(intersection.centerUnitOuterNormal());
241 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[intersectionIdx];
242 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
243 cellDataNeighbor.fluxData().addVelocity(nPhaseIdx, intersection.indexInOutside(), velocity);
244 }
245
246 cellData.fluxData().setVelocityMarker(idxInInside);
247 }
248 else
249 {
250 BoundaryTypes bctype;
251 problem.boundaryTypes(bctype, intersection);
252 if (bctype.isDirichlet(satEqIdx))
253 {
254 PrimaryVariables boundValues(0.0);
255 problem.dirichlet(boundValues, intersection);
256
257 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
258 if (velocityW[intersectionIdx] >= 0.)
259 {
260 mobilityW = cellData.mobility(wPhaseIdx);
261 }
262 else
263 {
264 mobilityW = fluidMatrixInteraction.krw(boundValues[saturationIdx]) / viscosityW;
265 }
266
267 if (velocityNw[intersectionIdx] >= 0.)
268 {
269 mobilityNw = cellData.mobility(nPhaseIdx);
270 }
271 else
272 {
273 mobilityNw = fluidMatrixInteraction.krn(boundValues[saturationIdx]) / viscosityNw;
274 }
275 }
276 else
277 {
278 mobilityW = cellData.mobility(wPhaseIdx);
279 mobilityNw = cellData.mobility(nPhaseIdx);
280 }
281
282 FieldVector velocity(intersection.centerUnitOuterNormal());
283 velocity *= mobilityW/(mobilityW+mobilityNw) * velocityW[intersectionIdx];
284 cellData.fluxData().addVelocity(wPhaseIdx, idxInInside, velocity);
285
286
287 velocity = intersection.centerUnitOuterNormal();
288 velocity *= mobilityNw/(mobilityW+mobilityNw) * velocityNw[intersectionIdx];
289 cellData.fluxData().addVelocity(nPhaseIdx, idxInInside, velocity);
290 cellData.fluxData().setVelocityMarker(idxInInside);
291 }
292 }
293 }
294 }
295};
296} // end namespace Dumux
297#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
int subIndex(int elemIdx, int fIdx)
Map interface fIdx'th interface of element index to array index.
Definition: intersectionmapper.hh:275
unsigned int size() const
Definition: intersectionmapper.hh:331
Extends CROperatorBase by a generic methods to assemble global stiffness matrix from local stiffness ...
Definition: croperatoradaptive.hh:69
const GridView gridView_
Definition: croperatoradaptive.hh:163
IntersectionMapper intersectionMapper_
Definition: croperatoradaptive.hh:165
Levelwise assembler.
Definition: operatoradaptive.hh:45
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operatoradaptive.hh:90
MimeticOperatorAssemblerTwoPAdaptive(const GridView &gridView)
Definition: operatoradaptive.hh:84
Specifies the properties for immiscible 2p diffusion/pressure models.
Defines the basic properties required for a mimetic method.