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