3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
mimetic/pressureadaptive.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_MIMETICPRESSURE2PADAPTIVE_HH
25#define DUMUX_MIMETICPRESSURE2PADAPTIVE_HH
26
27// dumux environment
32
33namespace Dumux {
34
61template<class TypeTag> class MimeticPressure2PAdaptive
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,
96 nPhaseIdx = Indices::nPhaseIdx,
97 numPhases = GET_PROP_VALUE(TypeTag, NumPhases)
98 };
99
100 using Element = typename GridView::Traits::template Codim<0>::Entity;
101 using Grid = typename GridView::Grid;
102
103 using Geometry = typename Element::Geometry;
104 using JacobianTransposed = typename Geometry::JacobianTransposed ;
105
106 using DimVector = Dune::FieldVector<Scalar, dim>;
107
108 using LocalStiffness = typename GET_PROP_TYPE(TypeTag, LocalStiffness);
109 using TraceType = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
111
112 using CellData = typename GET_PROP_TYPE(TypeTag, CellData);
113 using SolutionTypes = typename GET_PROP(TypeTag, SolutionTypes);
114 using ScalarSolutionType = typename SolutionTypes::ScalarSolution;
115
116 using Matrix = typename GET_PROP_TYPE(TypeTag, PressureCoefficientMatrix);
117 using Vector = typename GET_PROP_TYPE(TypeTag, PressureRHSVector);
118
120 void initializeMatrix();
121
123 void assemble(bool first)
124 {
125 Scalar timeStep = problem_.timeManager().timeStepSize();
126 Scalar maxError = 0.0;
127 int size = problem_.gridView().size(0);
128 for (int i = 0; i < size; i++)
129 {
130 Scalar sat = 0;
131 using std::max;
132 switch (saturationType)
133 {
134 case Sw:
135 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
136 break;
137 case Sn:
138 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
139 break;
140 }
141 if (sat > 1.0)
142 {
143 maxError = max(maxError, (sat - 1.0) / timeStep);
144 }
145 if (sat < 0.0)
146 {
147 maxError = max(maxError, (-sat) / timeStep);
148 }
149 }
150
151 lstiff_.setErrorInfo(maxError, timeStep);
152 A_.assemble(lstiff_, pressTrace_, f_);
153 return;
154 }
155
157 void solve();
158
160 void postprocess()
161 {
162 A_.calculatePressure(lstiff_, pressTrace_, problem_);
163
164 return;
165 }
166
167public:
169 void updateMaterialLaws();
170
179 void initialize(bool solveTwice = true)
180 {
181 const auto element = *problem_.gridView().template begin<0>();
182 FluidState fluidState;
183 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
184 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
185 fluidState.setTemperature(problem_.temperature(element));
186 fluidState.setSaturation(wPhaseIdx, 1.);
187 fluidState.setSaturation(nPhaseIdx, 0.);
188 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
189 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
190 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
191 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
192
194 adapt();
195 lstiff_.initialize();
196 lstiff_.reset();
197
198 assemble(true);
199 solve();
200 postprocess();
201
202 return;
203 }
204
205 // TODO doc me!
206 void adapt()
207 {
208 A_.adapt();
209 pressTrace_.resize(A_.intersectionMapper().size());
210 f_.resize(A_.intersectionMapper().size());
211 pressTrace_ = 0.0;
212 f_ = 0;
213 lstiff_.adapt();
214 }
215
222 {
224 postprocess();
225 }
226
228 void update()
229 {
230 if (problem_.gridAdapt().wasAdapted())
231 {
232 adapt();
233 }
234
235 lstiff_.reset();
236 assemble(false);
237
238 solve();
239
240 postprocess();
241
242 return;
243 }
244
251 template<class MultiWriter>
252 void addOutputVtkFields(MultiWriter &writer)
253 {
254 int size = problem_.gridView().size(0);
255 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
256 ScalarSolutionType *pressure = 0;
257 ScalarSolutionType *pressureSecond = 0;
258 ScalarSolutionType *potentialSecond = 0;
259 Dune::BlockVector < DimVector > *velocityWetting = 0;
260 Dune::BlockVector < DimVector > *velocityNonwetting = 0;
261
262 if (vtkOutputLevel_ > 0)
263 {
264 pressure = writer.allocateManagedBuffer(size);
265 pressureSecond = writer.allocateManagedBuffer(size);
266 potentialSecond = writer.allocateManagedBuffer(size);
267 velocityWetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
268 velocityNonwetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
269 }
270
271
272 for (const auto& element : elements(problem_.gridView()))
273 {
274 int eIdxGlobal = problem_.variables().index(element);
275 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
276
277 if (pressureType == pw)
278 {
279 (*potential)[eIdxGlobal] = cellData.potential(wPhaseIdx);
280 }
281
282 if (pressureType == pn)
283 {
284 (*potential)[eIdxGlobal] = cellData.potential(nPhaseIdx);
285 }
286
287 if (vtkOutputLevel_ > 0)
288 {
289
290 if (pressureType == pw)
291 {
292 (*pressure)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
293 (*potentialSecond)[eIdxGlobal] = cellData.potential(nPhaseIdx);
294 (*pressureSecond)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
295 }
296
297 if (pressureType == pn)
298 {
299 (*pressure)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
300 (*potentialSecond)[eIdxGlobal] = cellData.potential(wPhaseIdx);
301 (*pressureSecond)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
302 }
303
304 const typename Element::Geometry& geometry = element.geometry();
305
306 // get corresponding reference element
307 using ReferenceElements = Dune::ReferenceElements<Scalar, dim>;
308 const auto refElement = ReferenceElements::general(geometry.type());
309
310 const int numberOfFaces=refElement.size(1);
311 std::vector<Scalar> fluxW(numberOfFaces,0);
312 std::vector<Scalar> fluxNw(numberOfFaces,0);
313
314 // run through all intersections with neighbors and boundary
315 for (const auto& intersection : intersections(problem_.gridView(), element))
316 {
317 int isIndex = intersection.indexInInside();
318
319 fluxW[isIndex] += intersection.geometry().volume()
320 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
321 fluxNw[isIndex] += intersection.geometry().volume()
322 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
323 }
324
325 // calculate velocity on reference element as the Raviart-Thomas-0
326 // interpolant of the fluxes
327 Dune::FieldVector<Scalar, dim> refVelocity;
328 // simplices
329 if (refElement.type().isSimplex()) {
330 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
331 {
332 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
333 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
334 {
335 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
336 }
337 }
338 }
339 // cubes
340 else if (refElement.type().isCube()){
341 for (int i = 0; i < dim; i++)
342 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
343 }
344 // 3D prism and pyramids
345 else {
346 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
347 }
348
349 const DimVector& localPos = refElement.position(0, 0);
350
351 // get the transposed Jacobian of the element mapping
352 const JacobianTransposed jacobianT = geometry.jacobianTransposed(localPos);
353
354 // calculate the element velocity by the Piola transformation
355 DimVector elementVelocity(0);
356 jacobianT.umtv(refVelocity, elementVelocity);
357 elementVelocity /= geometry.integrationElement(localPos);
358
359 (*velocityWetting)[eIdxGlobal] = elementVelocity;
360
361 // calculate velocity on reference element as the Raviart-Thomas-0
362 // interpolant of the fluxes
363 // simplices
364 if (refElement.type().isSimplex()) {
365 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
366 {
367 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
368 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
369 {
370 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
371 }
372 }
373 }
374 // cubes
375 else if (refElement.type().isCube()){
376 for (int i = 0; i < dim; i++)
377 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
378 }
379 // 3D prism and pyramids
380 else {
381 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
382 }
383
384 // calculate the element velocity by the Piola transformation
385 elementVelocity = 0;
386 jacobianT.umtv(refVelocity, elementVelocity);
387 elementVelocity /= geometry.integrationElement(localPos);
388
389 (*velocityNonwetting)[eIdxGlobal] = elementVelocity;
390 }
391 }
392
393 if (pressureType == pw)
394 {
395 writer.attachCellData(*potential, "wetting potential");
396 }
397
398 if (pressureType == pn)
399 {
400 writer.attachCellData(*potential, "nonwetting potential");
401 }
402
403 if (vtkOutputLevel_ > 0)
404 {
405 if (pressureType == pw)
406 {
407 writer.attachCellData(*pressure, "wetting pressure");
408 writer.attachCellData(*pressureSecond, "nonwetting pressure");
409 writer.attachCellData(*potentialSecond, "nonwetting potential");
410 }
411
412 if (pressureType == pn)
413 {
414 writer.attachCellData(*pressure, "nonwetting pressure");
415 writer.attachCellData(*pressureSecond, "wetting pressure");
416 writer.attachCellData(*potentialSecond, "wetting potential");
417 }
418
419 writer.attachCellData(*velocityWetting, "wetting-velocity", dim);
420 writer.attachCellData(*velocityNonwetting, "non-wetting-velocity", dim);
421 }
422 }
423
432 void serializeEntity(std::ostream &outstream, const Element &element)
433 {
434 int numFaces = element.subEntities(1);
435 for (int i=0; i < numFaces; i++)
436 {
437 int isIdxGlobal = A_.intersectionMapper().subIndex(element, i);
438 outstream << pressTrace_[isIdxGlobal][0];
439 }
440 }
441
448 void deserializeEntity(std::istream &instream, const Element &element)
449 {
450 int numFaces = element.subEntities(1);
451 for (int i=0; i < numFaces; i++)
452 {
453 int isIdxGlobal = A_.intersectionMapper().subIndex(element, i);
454 instream >> pressTrace_[isIdxGlobal][0];
455 }
456 }
457
463 MimeticPressure2PAdaptive(Problem& problem) :
464 problem_(problem),
465 A_(problem.gridView()), lstiff_(problem_, false, problem_.gridView(), A_.intersectionMapper())
466 {
467 if (pressureType != pw && pressureType != pn)
468 {
469 DUNE_THROW(Dune::NotImplemented, "Pressure type not supported!");
470 }
471 if (saturationType != Sw && saturationType != Sn)
472 {
473 DUNE_THROW(Dune::NotImplemented, "Saturation type not supported!");
474 }
476 {
477 DUNE_THROW(Dune::NotImplemented, "Compressibility not supported!");
478 }
479
480 density_[wPhaseIdx] = 0.0;
481 density_[nPhaseIdx] = 0.0;
482 viscosity_[wPhaseIdx] = 0.0;
483 viscosity_[nPhaseIdx] = 0.0;
484
485 vtkOutputLevel_ = getParam<int>("Vtk.OutputLevel");
486 }
487
488private:
489 Problem& problem_;
490 TraceType pressTrace_; // vector of pressure traces
491 TraceType f_;
492 OperatorAssembler A_;
493 LocalStiffness lstiff_;
494
495 Scalar density_[numPhases];
496 Scalar viscosity_[numPhases];
497
498 int vtkOutputLevel_;
499};
500
502template<class TypeTag>
503void MimeticPressure2PAdaptive<TypeTag>::solve()
504{
505 using Solver = typename GET_PROP_TYPE(TypeTag, LinearSolver);
506
507 int verboseLevelSolver = getParam<int>("LinearSolver.Verbosity");
508
509 if (verboseLevelSolver)
510 std::cout << "MimeticPressure2PAdaptive: solve for pressure" << std::endl;
511
512 auto solver = getSolver<Solver>(problem_);
513 solver.solve(*A_, pressTrace_, f_);
514
515// printmatrix(std::cout, *A_, "global stiffness matrix", "row", 11, 3);
516// printvector(std::cout, f_, "right hand side", "row", 200, 1, 3);
517// printvector(std::cout, pressTrace_, "pressure", "row", 200, 1, 3);
518 return;
519}
520
522template<class TypeTag>
524{
525 // iterate through leaf grid an evaluate c0 at cell center
526 for (const auto& element : elements(problem_.gridView()))
527 {
528 int eIdxGlobal = problem_.variables().index(element);
529
530 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
531
532 Scalar satW = cellData.saturation(wPhaseIdx);
533
534 // initialize mobilities
535 Scalar mobilityW = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satW)
536 / viscosity_[wPhaseIdx];
537 Scalar mobilityNw = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satW)
538 / viscosity_[nPhaseIdx];
539
540 // initialize mobilities
541 cellData.setMobility(wPhaseIdx, mobilityW);
542 cellData.setMobility(nPhaseIdx, mobilityNw);
543
544 //initialize fractional flow functions
545 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
546 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
547 }
548 return;
549}
550
551} // end namespace Dumux
552#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
int subIndex(int elemIdx, int fIdx)
Map interface fIdx'th interface of element index to array index.
Definition: intersectionmapper.hh:236
unsigned int size() const
Definition: intersectionmapper.hh:292
const IntersectionMapper & intersectionMapper()
Definition: croperatoradaptive.hh:130
void assemble(LocalStiffness &loc, Vector &u, Vector &f)
Assembles global stiffness matrix.
Definition: croperatoradaptive.hh:232
void adapt()
Definition: croperatoradaptive.hh:111
Base class for local assemblers.
Definition: localstiffness.hh:60
Levelwise assembler.
Definition: operatoradaptive.hh:45
void calculatePressure(LocalStiffness &loc, Vector &u, Problem &problem)
Definition: operatoradaptive.hh:91
Mimetic method for the pressure equation.
Definition: mimetic/pressureadaptive.hh:62
void adapt()
Definition: mimetic/pressureadaptive.hh:206
void updateVelocity()
Velocity update.
Definition: mimetic/pressureadaptive.hh:221
void initialize(bool solveTwice=true)
Initializes the model.
Definition: mimetic/pressureadaptive.hh:179
MimeticPressure2PAdaptive(Problem &problem)
Constructs a MimeticPressure2PAdaptive object.
Definition: mimetic/pressureadaptive.hh:463
void update()
updates the model
Definition: mimetic/pressureadaptive.hh:228
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: mimetic/pressureadaptive.hh:523
void serializeEntity(std::ostream &outstream, const Element &element)
General method for serialization, output.
Definition: mimetic/pressureadaptive.hh:432
void addOutputVtkFields(MultiWriter &writer)
Write data file.
Definition: mimetic/pressureadaptive.hh:252
void deserializeEntity(std::istream &instream, const Element &element)
General method for deserialization.
Definition: mimetic/pressureadaptive.hh:448
Defines the basic properties required for a mimetic method.
Finite Volume Diffusion Model.