Loading [MathJax]/extensions/tex2jax.js
3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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#include <dune/common/exceptions.hh>
28
29// dumux environment
34
35namespace Dumux {
36
63template<class TypeTag> class MimeticPressure2PAdaptive
64{
68
70
72
75
76 enum
77 {
78 dim = GridView::dimension, dimWorld = GridView::dimensionworld
79 };
80 enum
81 {
82 pw = Indices::pressureW,
83 pn = Indices::pressureNw,
84 pGlobal = Indices::pressureGlobal,
85 Sw = Indices::saturationW,
86 Sn = Indices::saturationNw,
87 vw = Indices::velocityW,
88 vn = Indices::velocityNw,
90 pressureType = getPropValue<TypeTag, Properties::PressureFormulation>(),
92 saturationType = getPropValue<TypeTag, Properties::SaturationFormulation>(),
93 };
94 enum
95 {
96 wPhaseIdx = Indices::wPhaseIdx,
97 nPhaseIdx = Indices::nPhaseIdx,
98 numPhases = getPropValue<TypeTag, Properties::NumPhases>()
99 };
100
101 using Element = typename GridView::Traits::template Codim<0>::Entity;
102 using Grid = typename GridView::Grid;
103
104 using Geometry = typename Element::Geometry;
105 using JacobianTransposed = typename Geometry::JacobianTransposed ;
106
107 using DimVector = Dune::FieldVector<Scalar, dim>;
108
110 using TraceType = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
112
115 using ScalarSolutionType = typename SolutionTypes::ScalarSolution;
116
119
121 void initializeMatrix();
122
124 void assemble(bool first)
125 {
126 Scalar timeStep = problem_.timeManager().timeStepSize();
127 Scalar maxError = 0.0;
128 int size = problem_.gridView().size(0);
129 for (int i = 0; i < size; i++)
130 {
131 Scalar sat = 0;
132 using std::max;
133 switch (saturationType)
134 {
135 case Sw:
136 sat = problem_.variables().cellData(i).saturation(wPhaseIdx);
137 break;
138 case Sn:
139 sat = problem_.variables().cellData(i).saturation(nPhaseIdx);
140 break;
141 default:
142 DUNE_THROW(Dune::NotImplemented, "Only saturation formulation Sw and Sn are implemented!");
143 }
144 if (sat > 1.0)
145 {
146 maxError = max(maxError, (sat - 1.0) / timeStep);
147 }
148 if (sat < 0.0)
149 {
150 maxError = max(maxError, (-sat) / timeStep);
151 }
152 }
153
154 lstiff_.setErrorInfo(maxError, timeStep);
155 A_.assemble(lstiff_, pressTrace_, f_);
156 return;
157 }
158
160 void solve();
161
163 void postprocess()
164 {
165 A_.calculatePressure(lstiff_, pressTrace_, problem_);
166
167 return;
168 }
169
170public:
172 void updateMaterialLaws();
173
182 void initialize(bool solveTwice = true)
183 {
184 const auto element = *problem_.gridView().template begin<0>();
185 FluidState fluidState;
186 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
187 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
188 fluidState.setTemperature(problem_.temperature(element));
189 fluidState.setSaturation(wPhaseIdx, 1.);
190 fluidState.setSaturation(nPhaseIdx, 0.);
191 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
192 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
193 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
194 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
195
197 adapt();
198 lstiff_.initialize();
199 lstiff_.reset();
200
201 assemble(true);
202 solve();
203 postprocess();
204
205 return;
206 }
207
208 // TODO doc me!
209 void adapt()
210 {
211 A_.adapt();
212 pressTrace_.resize(A_.intersectionMapper().size());
213 f_.resize(A_.intersectionMapper().size());
214 pressTrace_ = 0.0;
215 f_ = 0;
216 lstiff_.adapt();
217 }
218
225 {
227 postprocess();
228 }
229
231 void update()
232 {
233 if (problem_.gridAdapt().wasAdapted())
234 {
235 adapt();
236 }
237
238 lstiff_.reset();
239 assemble(false);
240
241 solve();
242
243 postprocess();
244
245 return;
246 }
247
254 template<class MultiWriter>
255 void addOutputVtkFields(MultiWriter &writer)
256 {
257 int size = problem_.gridView().size(0);
258 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
259 ScalarSolutionType *pressure = 0;
260 ScalarSolutionType *pressureSecond = 0;
261 ScalarSolutionType *potentialSecond = 0;
262 Dune::BlockVector < DimVector > *velocityWetting = 0;
263 Dune::BlockVector < DimVector > *velocityNonwetting = 0;
264
265 if (vtkOutputLevel_ > 0)
266 {
267 pressure = writer.allocateManagedBuffer(size);
268 pressureSecond = writer.allocateManagedBuffer(size);
269 potentialSecond = writer.allocateManagedBuffer(size);
270 velocityWetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
271 velocityNonwetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
272 }
273
274
275 for (const auto& element : elements(problem_.gridView()))
276 {
277 int eIdxGlobal = problem_.variables().index(element);
278 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
279
280 if (pressureType == pw)
281 {
282 (*potential)[eIdxGlobal] = cellData.potential(wPhaseIdx);
283 }
284
285 if (pressureType == pn)
286 {
287 (*potential)[eIdxGlobal] = cellData.potential(nPhaseIdx);
288 }
289
290 if (vtkOutputLevel_ > 0)
291 {
292
293 if (pressureType == pw)
294 {
295 (*pressure)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
296 (*potentialSecond)[eIdxGlobal] = cellData.potential(nPhaseIdx);
297 (*pressureSecond)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
298 }
299
300 if (pressureType == pn)
301 {
302 (*pressure)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
303 (*potentialSecond)[eIdxGlobal] = cellData.potential(wPhaseIdx);
304 (*pressureSecond)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
305 }
306
307 const typename Element::Geometry& geometry = element.geometry();
308
309 // get corresponding reference element
310 const auto refElement = referenceElement(geometry);
311
312 const int numberOfFaces=refElement.size(1);
313 std::vector<Scalar> fluxW(numberOfFaces,0);
314 std::vector<Scalar> fluxNw(numberOfFaces,0);
315
316 // run through all intersections with neighbors and boundary
317 for (const auto& intersection : intersections(problem_.gridView(), element))
318 {
319 int isIndex = intersection.indexInInside();
320
321 fluxW[isIndex] += intersection.geometry().volume()
322 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
323 fluxNw[isIndex] += intersection.geometry().volume()
324 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
325 }
326
327 // calculate velocity on reference element as the Raviart-Thomas-0
328 // interpolant of the fluxes
329 Dune::FieldVector<Scalar, dim> refVelocity;
330 // simplices
331 if (refElement.type().isSimplex()) {
332 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
333 {
334 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
335 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
336 {
337 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
338 }
339 }
340 }
341 // cubes
342 else if (refElement.type().isCube()){
343 for (int i = 0; i < dim; i++)
344 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
345 }
346 // 3D prism and pyramids
347 else {
348 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
349 }
350
351 const DimVector& localPos = refElement.position(0, 0);
352
353 // get the transposed Jacobian of the element mapping
354 const JacobianTransposed jacobianT = geometry.jacobianTransposed(localPos);
355
356 // calculate the element velocity by the Piola transformation
357 DimVector elementVelocity(0);
358 jacobianT.umtv(refVelocity, elementVelocity);
359 elementVelocity /= geometry.integrationElement(localPos);
360
361 (*velocityWetting)[eIdxGlobal] = elementVelocity;
362
363 // calculate velocity on reference element as the Raviart-Thomas-0
364 // interpolant of the fluxes
365 // simplices
366 if (refElement.type().isSimplex()) {
367 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
368 {
369 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
370 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
371 {
372 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
373 }
374 }
375 }
376 // cubes
377 else if (refElement.type().isCube()){
378 for (int i = 0; i < dim; i++)
379 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
380 }
381 // 3D prism and pyramids
382 else {
383 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
384 }
385
386 // calculate the element velocity by the Piola transformation
387 elementVelocity = 0;
388 jacobianT.umtv(refVelocity, elementVelocity);
389 elementVelocity /= geometry.integrationElement(localPos);
390
391 (*velocityNonwetting)[eIdxGlobal] = elementVelocity;
392 }
393 }
394
395 if (pressureType == pw)
396 {
397 writer.attachCellData(*potential, "wetting potential");
398 }
399
400 if (pressureType == pn)
401 {
402 writer.attachCellData(*potential, "nonwetting potential");
403 }
404
405 if (vtkOutputLevel_ > 0)
406 {
407 if (pressureType == pw)
408 {
409 writer.attachCellData(*pressure, "wetting pressure");
410 writer.attachCellData(*pressureSecond, "nonwetting pressure");
411 writer.attachCellData(*potentialSecond, "nonwetting potential");
412 }
413
414 if (pressureType == pn)
415 {
416 writer.attachCellData(*pressure, "nonwetting pressure");
417 writer.attachCellData(*pressureSecond, "wetting pressure");
418 writer.attachCellData(*potentialSecond, "wetting potential");
419 }
420
421 writer.attachCellData(*velocityWetting, "wetting-velocity", dim);
422 writer.attachCellData(*velocityNonwetting, "nonwetting-velocity", dim);
423 }
424 }
425
434 void serializeEntity(std::ostream &outstream, const Element &element)
435 {
436 int numFaces = element.subEntities(1);
437 for (int i=0; i < numFaces; i++)
438 {
439 int isIdxGlobal = A_.intersectionMapper().subIndex(element, i);
440 outstream << pressTrace_[isIdxGlobal][0];
441 }
442 }
443
450 void deserializeEntity(std::istream &instream, const Element &element)
451 {
452 int numFaces = element.subEntities(1);
453 for (int i=0; i < numFaces; i++)
454 {
455 int isIdxGlobal = A_.intersectionMapper().subIndex(element, i);
456 instream >> pressTrace_[isIdxGlobal][0];
457 }
458 }
459
465 MimeticPressure2PAdaptive(Problem& problem) :
466 problem_(problem),
467 A_(problem.gridView()), lstiff_(problem_, false, problem_.gridView(), A_.intersectionMapper())
468 {
469 if (pressureType != pw && pressureType != pn)
470 {
471 DUNE_THROW(Dune::NotImplemented, "Pressure type not supported!");
472 }
473 if (saturationType != Sw && saturationType != Sn)
474 {
475 DUNE_THROW(Dune::NotImplemented, "Saturation type not supported!");
476 }
477 if (getPropValue<TypeTag, Properties::EnableCompressibility>())
478 {
479 DUNE_THROW(Dune::NotImplemented, "Compressibility not supported!");
480 }
481
482 density_[wPhaseIdx] = 0.0;
483 density_[nPhaseIdx] = 0.0;
484 viscosity_[wPhaseIdx] = 0.0;
485 viscosity_[nPhaseIdx] = 0.0;
486
487 vtkOutputLevel_ = getParam<int>("Vtk.OutputLevel");
488 }
489
490private:
491 Problem& problem_;
492 TraceType pressTrace_; // vector of pressure traces
493 TraceType f_;
494 OperatorAssembler A_;
495 LocalStiffness lstiff_;
496
497 Scalar density_[numPhases];
498 Scalar viscosity_[numPhases];
499
500 int vtkOutputLevel_;
501};
502
504template<class TypeTag>
505void MimeticPressure2PAdaptive<TypeTag>::solve()
506{
507 using Solver = GetPropType<TypeTag, Properties::LinearSolver>;
508
509 int verboseLevelSolver = getParam<int>("LinearSolver.Verbosity", 0);
510
511 if (verboseLevelSolver)
512 std::cout << "MimeticPressure2PAdaptive: solve for pressure" << std::endl;
513
514 auto solver = getSolver<Solver>(problem_);
515 solver.solve(*A_, pressTrace_, f_);
516
517// printmatrix(std::cout, *A_, "global stiffness matrix", "row", 11, 3);
518// printvector(std::cout, f_, "right hand side", "row", 200, 1, 3);
519// printvector(std::cout, pressTrace_, "pressure", "row", 200, 1, 3);
520 return;
521}
522
524template<class TypeTag>
526{
527 // iterate through leaf grid an evaluate c0 at cell center
528 for (const auto& element : elements(problem_.gridView()))
529 {
530 int eIdxGlobal = problem_.variables().index(element);
531
532 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
533
534 const Scalar satW = cellData.saturation(wPhaseIdx);
535
536 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
537
538 // initialize mobilities
539 const Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
540 const Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
541
542 // initialize mobilities
543 cellData.setMobility(wPhaseIdx, mobilityW);
544 cellData.setMobility(nPhaseIdx, mobilityNw);
545
546 //initialize fractional flow functions
547 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
548 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
549 }
550 return;
551}
552
553} // end namespace Dumux
554#endif
Local stiffness matrix for the diffusion equation discretized by mimetic FD.
An assembler for the Jacobian matrix based on mimetic FD.
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 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:275
unsigned int size() const
Definition: intersectionmapper.hh:331
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:90
Mimetic method for the pressure equation.
Definition: mimetic/pressureadaptive.hh:64
void adapt()
Definition: mimetic/pressureadaptive.hh:209
void updateVelocity()
Velocity update.
Definition: mimetic/pressureadaptive.hh:224
void initialize(bool solveTwice=true)
Initializes the model.
Definition: mimetic/pressureadaptive.hh:182
MimeticPressure2PAdaptive(Problem &problem)
Constructs a MimeticPressure2PAdaptive object.
Definition: mimetic/pressureadaptive.hh:465
void update()
updates the model
Definition: mimetic/pressureadaptive.hh:231
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: mimetic/pressureadaptive.hh:525
void serializeEntity(std::ostream &outstream, const Element &element)
General method for serialization, output.
Definition: mimetic/pressureadaptive.hh:434
void addOutputVtkFields(MultiWriter &writer)
Write data file.
Definition: mimetic/pressureadaptive.hh:255
void deserializeEntity(std::istream &instream, const Element &element)
General method for deserialization.
Definition: mimetic/pressureadaptive.hh:450
Defines the basic properties required for a mimetic method.
Finite Volume Diffusion Model.