3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
2p/sequential/diffusion/mimetic/pressure.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_MIMETICPRESSURE2P_HH
25#define DUMUX_MIMETICPRESSURE2P_HH
26
27// dumux environment
32
33namespace Dumux {
34
61template<class TypeTag> class MimeticPressure2P
62{
63 using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
64 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
65 using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
66
67 using SpatialParams = typename GET_PROP_TYPE(TypeTag, SpatialParams);
68 using MaterialLaw = typename SpatialParams::MaterialLaw;
69
70 using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
71
72 using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
73 using FluidState = typename GET_PROP_TYPE(TypeTag, FluidState);
74
75 enum
76 {
77 dim = GridView::dimension, dimWorld = GridView::dimensionworld
78 };
79 enum
80 {
81 pw = Indices::pressureW,
82 pn = Indices::pressureNw,
83 pGlobal = Indices::pressureGlobal,
84 sw = Indices::saturationW,
85 sn = Indices::saturationNw,
86 vw = Indices::velocityW,
87 vn = Indices::velocityNw,
89 pressureType = GET_PROP_VALUE(TypeTag, PressureFormulation),
91 saturationType = GET_PROP_VALUE(TypeTag, SaturationFormulation),
92 };
93 enum
94 {
95 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
96 numPhases = GET_PROP_VALUE(TypeTag, NumPhases)
97 };
98
99 using Element = typename GridView::Traits::template Codim<0>::Entity;
100 using Grid = typename GridView::Grid;
101
102 using Geometry = typename Element::Geometry;
103 using JacobianTransposed = typename Geometry::JacobianTransposed ;
104
105 using LocalStiffness = typename GET_PROP_TYPE(TypeTag, LocalStiffness);
106 using TraceType = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
108
109 using CellData = typename GET_PROP_TYPE(TypeTag, CellData);
110 using SolutionTypes = typename GET_PROP(TypeTag, SolutionTypes);
111 using ScalarSolutionType = typename SolutionTypes::ScalarSolution;
112
113 using Matrix = typename GET_PROP_TYPE(TypeTag, PressureCoefficientMatrix);
114 using Vector = typename GET_PROP_TYPE(TypeTag, PressureRHSVector);
115
116 using DimVector = Dune::FieldVector<Scalar, dim>;
117
119 void initializeMatrix();
120
122 void assemble(bool first)
123 {
124 Scalar timeStep = problem_.timeManager().timeStepSize();
125 Scalar maxError = 0.0;
126 int size = problem_.gridView().size(0);
127 for (int i = 0; i < size; i++)
128 {
129 Scalar sat = 0;
130 using std::max;
131 switch (saturationType)
132 {
133 case sw:
134 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
135 break;
136 case sn:
137 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
138 break;
139 }
140 if (sat > 1.0)
141 {
142 maxError = max(maxError, (sat - 1.0) / timeStep);
143 }
144 if (sat < 0.0)
145 {
146 maxError = max(maxError, (-sat) / timeStep);
147 }
148 }
149
150 lstiff_.setErrorInfo(maxError, timeStep);
151 A_.assemble(lstiff_, pressTrace_, f_);
152 return;
153 }
154
156 void solve();
157
158 void postprocess()
159 {
160 A_.calculatePressure(lstiff_, pressTrace_, problem_);
161
162 return;
163 }
164
165public:
167 void updateMaterialLaws();
168
177 void initialize(bool solveTwice = true)
178 {
179 const auto element = *problem_.gridView().template begin<0>();
180 FluidState fluidState;
181 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
182 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
183 fluidState.setTemperature(problem_.temperature(element));
184 fluidState.setSaturation(wPhaseIdx, 1.);
185 fluidState.setSaturation(nPhaseIdx, 0.);
186 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
187 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
188 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
189 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
190
192 A_.initialize();
193 pressTrace_.resize(problem_.gridView().size(1));
194 f_.resize(problem_.gridView().size(1));
195 lstiff_.initialize();
196 lstiff_.reset();
197
198 pressTrace_ = 0.0;
199 f_ = 0;
200
201 assemble(true);
202 solve();
203 postprocess();
204
205 return;
206 }
207
214 {
216 postprocess();
217 }
218
219 // TODO doc me!
220 void update()
221 {
222 lstiff_.reset();
223
224 assemble(false);
225 solve();
226 postprocess();
227
228 return;
229 }
230
237 template<class MultiWriter>
238 void addOutputVtkFields(MultiWriter &writer)
239 {
240 int size = problem_.gridView().size(0);
241 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
242 ScalarSolutionType *pressure = 0;
243 ScalarSolutionType *pressureSecond = 0;
244 ScalarSolutionType *potentialSecond = 0;
245 Dune::BlockVector < DimVector > *velocityWetting = 0;
246 Dune::BlockVector < DimVector > *velocityNonwetting = 0;
247
248 if (vtkOutputLevel_ > 0)
249 {
250 pressure = writer.allocateManagedBuffer(size);
251 pressureSecond = writer.allocateManagedBuffer(size);
252 potentialSecond = writer.allocateManagedBuffer(size);
253 velocityWetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
254 velocityNonwetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
255 }
256
257
258 for (const auto& element : elements(problem_.gridView()))
259 {
260 int eIdxGlobal = problem_.variables().index(element);
261 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
262
263 if (pressureType == pw)
264 {
265 (*potential)[eIdxGlobal] = cellData.potential(wPhaseIdx);
266 }
267
268 if (pressureType == pn)
269 {
270 (*potential)[eIdxGlobal] = cellData.potential(nPhaseIdx);
271 }
272
273 if (vtkOutputLevel_ > 0)
274 {
275
276 if (pressureType == pw)
277 {
278 (*pressure)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
279 (*potentialSecond)[eIdxGlobal] = cellData.potential(nPhaseIdx);
280 (*pressureSecond)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
281 }
282
283 if (pressureType == pn)
284 {
285 (*pressure)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
286 (*potentialSecond)[eIdxGlobal] = cellData.potential(wPhaseIdx);
287 (*pressureSecond)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
288 }
289
290 const typename Element::Geometry& geometry = element.geometry();
291
292 // get corresponding reference element
293 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
294 const auto refElement = ReferenceElements::general(geometry.type());
295
296 const int numberOfFaces=refElement.size(1);
297 std::vector<Scalar> fluxW(numberOfFaces,0);
298 std::vector<Scalar> fluxNw(numberOfFaces,0);
299
300 // run through all intersections with neighbors and boundary
301 for (const auto& intersection : intersections(problem_.gridView(), element))
302 {
303 int isIndex = intersection.indexInInside();
304
305 fluxW[isIndex] += intersection.geometry().volume()
306 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
307 fluxNw[isIndex] += intersection.geometry().volume()
308 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
309 }
310
311 // calculate velocity on reference element as the Raviart-Thomas-0
312 // interpolant of the fluxes
313 Dune::FieldVector<Scalar, dim> refVelocity;
314 // simplices
315 if (refElement.type().isSimplex()) {
316 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
317 {
318 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
319 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
320 {
321 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
322 }
323 }
324 }
325 // cubes
326 else if (refElement.type().isCube()){
327 for (int i = 0; i < dim; i++)
328 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
329 }
330 // 3D prism and pyramids
331 else {
332 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
333 }
334
335 const DimVector& localPos = refElement.position(0, 0);
336
337 // get the transposed Jacobian of the element mapping
338 const JacobianTransposed jacobianT = geometry.jacobianTransposed(localPos);
339
340 // calculate the element velocity by the Piola transformation
341 DimVector elementVelocity(0);
342 jacobianT.umtv(refVelocity, elementVelocity);
343 elementVelocity /= geometry.integrationElement(localPos);
344
345 (*velocityWetting)[eIdxGlobal] = elementVelocity;
346
347 // calculate velocity on reference element as the Raviart-Thomas-0
348 // interpolant of the fluxes
349 // simplices
350 if (refElement.type().isSimplex()) {
351 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
352 {
353 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
354 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
355 {
356 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
357 }
358 }
359 }
360 // cubes
361 else if (refElement.type().isCube()){
362 for (int i = 0; i < dim; i++)
363 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
364 }
365 // 3D prism and pyramids
366 else {
367 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
368 }
369
370 // calculate the element velocity by the Piola transformation
371 elementVelocity = 0;
372 jacobianT.umtv(refVelocity, elementVelocity);
373 elementVelocity /= geometry.integrationElement(localPos);
374
375 (*velocityNonwetting)[eIdxGlobal] = elementVelocity;
376 }
377 }
378
379 if (pressureType == pw)
380 {
381 writer.attachCellData(*potential, "wetting potential");
382 }
383
384 if (pressureType == pn)
385 {
386 writer.attachCellData(*potential, "nonwetting potential");
387 }
388
389 if (vtkOutputLevel_ > 0)
390 {
391 if (pressureType == pw)
392 {
393 writer.attachCellData(*pressure, "wetting pressure");
394 writer.attachCellData(*pressureSecond, "nonwetting pressure");
395 writer.attachCellData(*potentialSecond, "nonwetting potential");
396 }
397
398 if (pressureType == pn)
399 {
400 writer.attachCellData(*pressure, "nonwetting pressure");
401 writer.attachCellData(*pressureSecond, "wetting pressure");
402 writer.attachCellData(*potentialSecond, "wetting potential");
403 }
404
405 writer.attachCellData(*velocityWetting, "wetting-velocity", dim);
406 writer.attachCellData(*velocityNonwetting, "non-wetting-velocity", dim);
407 }
408 }
409
418 void serializeEntity(std::ostream &outstream, const Element &element)
419 {
420 int numFaces = element.subEntities(1);
421 for (int i=0; i < numFaces; i++)
422 {
423 int fIdxGlobal = A_.faceMapper().subIndex(element, i, 1);
424 outstream << pressTrace_[fIdxGlobal][0];
425 }
426 }
427
434 void deserializeEntity(std::istream &instream, const Element &element)
435 {
436 int numFaces = element.subEntities(1);
437 for (int i=0; i < numFaces; i++)
438 {
439 int fIdxGlobal = A_.faceMapper().subIndex(element, i, 1);
440 instream >> pressTrace_[fIdxGlobal][0];
441 }
442 }
443
449 MimeticPressure2P(Problem& problem) :
450 problem_(problem),
451 A_(problem.gridView()), lstiff_(problem_, false, problem_.gridView())
452 {
453 if (pressureType != pw && pressureType != pn)
454 {
455 DUNE_THROW(Dune::NotImplemented, "Pressure type not supported!");
456 }
457 if (saturationType != sw && saturationType != sn)
458 {
459 DUNE_THROW(Dune::NotImplemented, "Saturation type not supported!");
460 }
462 {
463 DUNE_THROW(Dune::NotImplemented, "Compressibility not supported!");
464 }
465
466 density_[wPhaseIdx] = 0.0;
467 density_[nPhaseIdx] = 0.0;
468 viscosity_[wPhaseIdx] = 0.0;
469 viscosity_[nPhaseIdx] = 0.0;
470
471 vtkOutputLevel_ = getParam<int>("Vtk.OutputLevel");
472 }
473
474private:
475 Problem& problem_;
476 TraceType pressTrace_; // vector of pressure traces
477 TraceType f_;
478 OperatorAssembler A_;
479 LocalStiffness lstiff_;
480
481 Scalar density_[numPhases];
482 Scalar viscosity_[numPhases];
483
484 int vtkOutputLevel_;
485};
486
488template<class TypeTag>
489void MimeticPressure2P<TypeTag>::solve()
490{
491 using Solver = typename GET_PROP_TYPE(TypeTag, LinearSolver);
492
493 auto verboseLevelSolver = getParam<int>("LinearSolver.Verbosity");
494
495 if (verboseLevelSolver)
496 std::cout << "MimeticPressure2P: solve for pressure" << std::endl;
497
498 auto solver = getSolver<Solver>(problem_);
499 solver.solve(*A_, pressTrace_, f_);
500
501 return;
502}
503
505template<class TypeTag>
507{
508 // iterate through leaf grid an evaluate c0 at cell center
509 for (const auto& element : elements(problem_.gridView()))
510 {
511 int eIdxGlobal = problem_.variables().index(element);
512
513 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
514
515 Scalar satW = cellData.saturation(wPhaseIdx);
516
517 // initialize mobilities
518 Scalar mobilityW = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
519 / viscosity_[wPhaseIdx];
520 Scalar mobilityNw = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
521 / viscosity_[nPhaseIdx];
522
523 // initialize mobilities
524 cellData.setMobility(wPhaseIdx, mobilityW);
525 cellData.setMobility(nPhaseIdx, mobilityNw);
526
527 // initialize fractional flow functions
528 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
529 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
530 }
531 return;
532}
533
534} // end namespace Dumux
535#endif
#define GET_PROP_VALUE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:282
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
Local stiffness matrix for the diffusion equation discretized by mimetic FD.
An assembler for the Jacobian matrix based on mimetic FD.
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag PressureRHSVector
Type of the right hand side vector given to the linear solver.
Definition: sequential/pressureproperties.hh:60
Property tag EnableCompressibility
Returns whether compressibility is allowed.
Definition: porousmediumflow/2p/sequential/properties.hh:55
Property tag SaturationFormulation
The formulation of the saturation model.
Definition: porousmediumflow/2p/sequential/properties.hh:53
Property tag NumPhases
Number of phases in the system.
Definition: porousmediumflow/sequential/properties.hh:69
Property tag PressureCoefficientMatrix
Gives maximum number of intersections of an element and neighboring elements.
Definition: porousmediumflow/sequential/properties.hh:74
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
const FaceMapper & faceMapper()
Definition: croperator.hh:195
void assemble(LocalStiffness &loc, Vector &u, Vector &f)
Assembles global stiffness matrix.
Definition: croperator.hh:223
void initialize()
Initialize the CR operator assembler.
Definition: croperator.hh:112
Base class for local assemblers.
Definition: localstiffness.hh:60
Levelwise assembler.
Definition: operator.hh:45
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operator.hh:92
Mimetic method for the pressure equation.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:62
void updateVelocity()
Velocity update.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:213
void initialize(bool solveTwice=true)
Initializes the model.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:177
void addOutputVtkFields(MultiWriter &writer)
Write data file.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:238
void deserializeEntity(std::istream &instream, const Element &element)
General method for deserialization.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:434
void serializeEntity(std::ostream &outstream, const Element &element)
General method for serialization, output.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:418
void update()
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:220
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:506
MimeticPressure2P(Problem &problem)
Constructs a MimeticPressure2P object.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:449
Defines the basic properties required for a mimetic method.
Finite Volume Diffusion Model.