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
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#include <dune/common/exceptions.hh>
28
29// dumux environment
34
35namespace Dumux {
36
63template<class TypeTag> class MimeticPressure2P
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, nPhaseIdx = Indices::nPhaseIdx,
97 numPhases = getPropValue<TypeTag, Properties::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
107 using TraceType = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
109
112 using ScalarSolutionType = typename SolutionTypes::ScalarSolution;
113
116
117 using DimVector = Dune::FieldVector<Scalar, dim>;
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 default:
141 DUNE_THROW(Dune::NotImplemented, "Only saturation formulation sw and sn are implemented!");
142 }
143 if (sat > 1.0)
144 {
145 maxError = max(maxError, (sat - 1.0) / timeStep);
146 }
147 if (sat < 0.0)
148 {
149 maxError = max(maxError, (-sat) / timeStep);
150 }
151 }
152
153 lstiff_.setErrorInfo(maxError, timeStep);
154 A_.assemble(lstiff_, pressTrace_, f_);
155 return;
156 }
157
159 void solve();
160
161 void postprocess()
162 {
163 A_.calculatePressure(lstiff_, pressTrace_, problem_);
164
165 return;
166 }
167
168public:
170 void updateMaterialLaws();
171
180 void initialize(bool solveTwice = true)
181 {
182 const auto element = *problem_.gridView().template begin<0>();
183 FluidState fluidState;
184 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
185 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
186 fluidState.setTemperature(problem_.temperature(element));
187 fluidState.setSaturation(wPhaseIdx, 1.);
188 fluidState.setSaturation(nPhaseIdx, 0.);
189 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
190 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
191 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
192 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
193
195 A_.initialize();
196 pressTrace_.resize(problem_.gridView().size(1));
197 f_.resize(problem_.gridView().size(1));
198 lstiff_.initialize();
199 lstiff_.reset();
200
201 pressTrace_ = 0.0;
202 f_ = 0;
203
204 assemble(true);
205 solve();
206 postprocess();
207
208 return;
209 }
210
217 {
219 postprocess();
220 }
221
222 // TODO doc me!
223 void update()
224 {
225 lstiff_.reset();
226
227 assemble(false);
228 solve();
229 postprocess();
230
231 return;
232 }
233
240 template<class MultiWriter>
241 void addOutputVtkFields(MultiWriter &writer)
242 {
243 int size = problem_.gridView().size(0);
244 ScalarSolutionType *potential = writer.allocateManagedBuffer(size);
245 ScalarSolutionType *pressure = 0;
246 ScalarSolutionType *pressureSecond = 0;
247 ScalarSolutionType *potentialSecond = 0;
248 Dune::BlockVector < DimVector > *velocityWetting = 0;
249 Dune::BlockVector < DimVector > *velocityNonwetting = 0;
250
251 if (vtkOutputLevel_ > 0)
252 {
253 pressure = writer.allocateManagedBuffer(size);
254 pressureSecond = writer.allocateManagedBuffer(size);
255 potentialSecond = writer.allocateManagedBuffer(size);
256 velocityWetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
257 velocityNonwetting = writer.template allocateManagedBuffer<Scalar, dim>(size);
258 }
259
260
261 for (const auto& element : elements(problem_.gridView()))
262 {
263 int eIdxGlobal = problem_.variables().index(element);
264 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
265
266 if (pressureType == pw)
267 {
268 (*potential)[eIdxGlobal] = cellData.potential(wPhaseIdx);
269 }
270
271 if (pressureType == pn)
272 {
273 (*potential)[eIdxGlobal] = cellData.potential(nPhaseIdx);
274 }
275
276 if (vtkOutputLevel_ > 0)
277 {
278
279 if (pressureType == pw)
280 {
281 (*pressure)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
282 (*potentialSecond)[eIdxGlobal] = cellData.potential(nPhaseIdx);
283 (*pressureSecond)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
284 }
285
286 if (pressureType == pn)
287 {
288 (*pressure)[eIdxGlobal] = cellData.pressure(nPhaseIdx);
289 (*potentialSecond)[eIdxGlobal] = cellData.potential(wPhaseIdx);
290 (*pressureSecond)[eIdxGlobal] = cellData.pressure(wPhaseIdx);
291 }
292
293 const typename Element::Geometry& geometry = element.geometry();
294
295 // get corresponding reference element
296 const auto refElement = referenceElement(geometry);
297
298 const int numberOfFaces=refElement.size(1);
299 std::vector<Scalar> fluxW(numberOfFaces,0);
300 std::vector<Scalar> fluxNw(numberOfFaces,0);
301
302 // run through all intersections with neighbors and boundary
303 for (const auto& intersection : intersections(problem_.gridView(), element))
304 {
305 int isIndex = intersection.indexInInside();
306
307 fluxW[isIndex] += intersection.geometry().volume()
308 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(wPhaseIdx, isIndex));
309 fluxNw[isIndex] += intersection.geometry().volume()
310 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(nPhaseIdx, isIndex));
311 }
312
313 // calculate velocity on reference element as the Raviart-Thomas-0
314 // interpolant of the fluxes
315 Dune::FieldVector<Scalar, dim> refVelocity;
316 // simplices
317 if (refElement.type().isSimplex()) {
318 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
319 {
320 refVelocity[dimIdx] = -fluxW[dim - 1 - dimIdx];
321 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
322 {
323 refVelocity[dimIdx] += fluxW[fIdx]/(dim + 1);
324 }
325 }
326 }
327 // cubes
328 else if (refElement.type().isCube()){
329 for (int i = 0; i < dim; i++)
330 refVelocity[i] = 0.5 * (fluxW[2*i + 1] - fluxW[2*i]);
331 }
332 // 3D prism and pyramids
333 else {
334 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
335 }
336
337 const DimVector& localPos = refElement.position(0, 0);
338
339 // get the transposed Jacobian of the element mapping
340 const JacobianTransposed jacobianT = geometry.jacobianTransposed(localPos);
341
342 // calculate the element velocity by the Piola transformation
343 DimVector elementVelocity(0);
344 jacobianT.umtv(refVelocity, elementVelocity);
345 elementVelocity /= geometry.integrationElement(localPos);
346
347 (*velocityWetting)[eIdxGlobal] = elementVelocity;
348
349 // calculate velocity on reference element as the Raviart-Thomas-0
350 // interpolant of the fluxes
351 // simplices
352 if (refElement.type().isSimplex()) {
353 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
354 {
355 refVelocity[dimIdx] = -fluxNw[dim - 1 - dimIdx];
356 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
357 {
358 refVelocity[dimIdx] += fluxNw[fIdx]/(dim + 1);
359 }
360 }
361 }
362 // cubes
363 else if (refElement.type().isCube()){
364 for (int i = 0; i < dim; i++)
365 refVelocity[i] = 0.5 * (fluxNw[2*i + 1] - fluxNw[2*i]);
366 }
367 // 3D prism and pyramids
368 else {
369 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
370 }
371
372 // calculate the element velocity by the Piola transformation
373 elementVelocity = 0;
374 jacobianT.umtv(refVelocity, elementVelocity);
375 elementVelocity /= geometry.integrationElement(localPos);
376
377 (*velocityNonwetting)[eIdxGlobal] = elementVelocity;
378 }
379 }
380
381 if (pressureType == pw)
382 {
383 writer.attachCellData(*potential, "wetting potential");
384 }
385
386 if (pressureType == pn)
387 {
388 writer.attachCellData(*potential, "nonwetting potential");
389 }
390
391 if (vtkOutputLevel_ > 0)
392 {
393 if (pressureType == pw)
394 {
395 writer.attachCellData(*pressure, "wetting pressure");
396 writer.attachCellData(*pressureSecond, "nonwetting pressure");
397 writer.attachCellData(*potentialSecond, "nonwetting potential");
398 }
399
400 if (pressureType == pn)
401 {
402 writer.attachCellData(*pressure, "nonwetting pressure");
403 writer.attachCellData(*pressureSecond, "wetting pressure");
404 writer.attachCellData(*potentialSecond, "wetting potential");
405 }
406
407 writer.attachCellData(*velocityWetting, "wetting-velocity", dim);
408 writer.attachCellData(*velocityNonwetting, "nonwetting-velocity", dim);
409 }
410 }
411
420 void serializeEntity(std::ostream &outstream, const Element &element)
421 {
422 int numFaces = element.subEntities(1);
423 for (int i=0; i < numFaces; i++)
424 {
425 int fIdxGlobal = A_.faceMapper().subIndex(element, i, 1);
426 outstream << pressTrace_[fIdxGlobal][0];
427 }
428 }
429
436 void deserializeEntity(std::istream &instream, const Element &element)
437 {
438 int numFaces = element.subEntities(1);
439 for (int i=0; i < numFaces; i++)
440 {
441 int fIdxGlobal = A_.faceMapper().subIndex(element, i, 1);
442 instream >> pressTrace_[fIdxGlobal][0];
443 }
444 }
445
451 MimeticPressure2P(Problem& problem) :
452 problem_(problem),
453 A_(problem.gridView()), lstiff_(problem_, false, problem_.gridView())
454 {
455 if (pressureType != pw && pressureType != pn)
456 {
457 DUNE_THROW(Dune::NotImplemented, "Pressure type not supported!");
458 }
459 if (saturationType != sw && saturationType != sn)
460 {
461 DUNE_THROW(Dune::NotImplemented, "Saturation type not supported!");
462 }
463 if (getPropValue<TypeTag, Properties::EnableCompressibility>())
464 {
465 DUNE_THROW(Dune::NotImplemented, "Compressibility not supported!");
466 }
467
468 density_[wPhaseIdx] = 0.0;
469 density_[nPhaseIdx] = 0.0;
470 viscosity_[wPhaseIdx] = 0.0;
471 viscosity_[nPhaseIdx] = 0.0;
472
473 vtkOutputLevel_ = getParam<int>("Vtk.OutputLevel");
474 }
475
476private:
477 Problem& problem_;
478 TraceType pressTrace_; // vector of pressure traces
479 TraceType f_;
480 OperatorAssembler A_;
481 LocalStiffness lstiff_;
482
483 Scalar density_[numPhases];
484 Scalar viscosity_[numPhases];
485
486 int vtkOutputLevel_;
487};
488
490template<class TypeTag>
491void MimeticPressure2P<TypeTag>::solve()
492{
493 using Solver = GetPropType<TypeTag, Properties::LinearSolver>;
494
495 auto verboseLevelSolver = getParam<int>("LinearSolver.Verbosity", 0);
496
497 if (verboseLevelSolver)
498 std::cout << "MimeticPressure2P: solve for pressure" << std::endl;
499
500 auto solver = getSolver<Solver>(problem_);
501 solver.solve(*A_, pressTrace_, f_);
502
503 return;
504}
505
507template<class TypeTag>
509{
510 // iterate through leaf grid an evaluate c0 at cell center
511 for (const auto& element : elements(problem_.gridView()))
512 {
513 int eIdxGlobal = problem_.variables().index(element);
514
515 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
516
517 const Scalar satW = cellData.saturation(wPhaseIdx);
518
519 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
520
521 // initialize mobilities
522 const Scalar mobilityW = fluidMatrixInteraction.krw(satW) / viscosity_[wPhaseIdx];
523 const Scalar mobilityNw = fluidMatrixInteraction.krn(satW) / viscosity_[nPhaseIdx];
524
525 // initialize mobilities
526 cellData.setMobility(wPhaseIdx, mobilityW);
527 cellData.setMobility(nPhaseIdx, mobilityNw);
528
529 // initialize fractional flow functions
530 cellData.setFracFlowFunc(wPhaseIdx, mobilityW / (mobilityW + mobilityNw));
531 cellData.setFracFlowFunc(nPhaseIdx, mobilityNw / (mobilityW + mobilityNw));
532 }
533 return;
534}
535
536} // end namespace Dumux
537#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
const FaceMapper & faceMapper()
Definition: croperator.hh:199
void assemble(LocalStiffness &loc, Vector &u, Vector &f)
Assembles global stiffness matrix.
Definition: croperator.hh:227
void initialize()
Initialize the CR operator assembler.
Definition: croperator.hh:113
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:91
Mimetic method for the pressure equation.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:64
void updateVelocity()
Velocity update.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:216
void initialize(bool solveTwice=true)
Initializes the model.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:180
void addOutputVtkFields(MultiWriter &writer)
Write data file.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:241
void deserializeEntity(std::istream &instream, const Element &element)
General method for deserialization.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:436
void serializeEntity(std::ostream &outstream, const Element &element)
General method for serialization, output.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:420
void update()
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:223
void updateMaterialLaws()
Constitutive functions are initialized and stored in the variables object.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:508
MimeticPressure2P(Problem &problem)
Constructs a MimeticPressure2P object.
Definition: 2p/sequential/diffusion/mimetic/pressure.hh:451
Defines the basic properties required for a mimetic method.
Finite Volume Diffusion Model.