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