3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
mimetic.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_MIMETIC2P_HH
25#define DUMUX_MIMETIC2P_HH
26
27#include<map>
28#include<iostream>
29#include<iomanip>
30#include<fstream>
31#include<vector>
32#include<sstream>
33
34#include<dune/common/exceptions.hh>
35#include<dune/grid/common/grid.hh>
36#include<dune/geometry/type.hh>
37#include<dune/geometry/quadraturerules.hh>
38
41#include "localstiffness.hh"
42
43#include <dune/common/dynvector.hh>
44
45namespace Dumux {
46
51template<class TypeTag>
53{
54 using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
55 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
56 using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
57
58 using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
59
60 // grid types
61 enum
62 {
63 dim = GridView::dimension
64 };
65 enum
66 {
67 wPhaseIdx = Indices::wPhaseIdx,
68 nPhaseIdx = Indices::nPhaseIdx,
69 pressureIdx = Indices::pressureIdx,
70 saturationIdx = Indices::saturationIdx,
71 pressureEqIdx = Indices::pressureEqIdx,
72 satEqIdx = Indices::satEqIdx,
73 numPhases = GET_PROP_VALUE(TypeTag, NumPhases)
74 };
75 enum
76 {
77 pw = Indices::pressureW,
78 pn = Indices::pressureNw,
79 pGlobal = Indices::pressureGlobal,
80 Sw = Indices::saturationW,
81 Sn = Indices::saturationNw,
82 vw = Indices::velocityW,
83 vn = Indices::velocityNw,
85 pressureType = GET_PROP_VALUE(TypeTag, PressureFormulation),
87 saturationType = GET_PROP_VALUE(TypeTag, SaturationFormulation)
88 };
89
90 using Grid = typename GridView::Grid;
91 using Element = typename GridView::Traits::template Codim<0>::Entity;
92
93 using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
94 using SolutionTypes = typename GET_PROP(TypeTag, SolutionTypes);
95 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
96 using CellData = typename GET_PROP_TYPE(TypeTag, CellData);
97 using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
98 using FluidState = typename GET_PROP_TYPE(TypeTag, FluidState);
99 using MaterialLaw = typename GET_PROP_TYPE(TypeTag, MaterialLaw);
100
101public:
102 // define the number of components of your system, this is used outside
103 // to allocate the correct size of (dense) blocks with a FieldMatrix
104 enum
105 {
106 m = 1
107 };
108 enum
109 {
110 size = 2 * dim
111 };
112
114 MimeticTwoPLocalStiffness(Problem& problem, bool levelBoundaryAsDirichlet,
115 const GridView& gridView, bool procBoundaryAsDirichlet = true) :
116 problem_(problem), gridView_(gridView), maxError_(0), timeStep_(1)
117 {
118 ErrorTermFactor_ = getParam<Scalar>("Impet.ErrorTermFactor");
119 ErrorTermLowerBound_ = getParam<Scalar>("Impet.ErrorTermLowerBound");
120 ErrorTermUpperBound_ = getParam<Scalar>("Impet.ErrorTermUpperBound");
121
122 density_[wPhaseIdx] = 0.0;
123 density_[nPhaseIdx] = 0.0;
124 viscosity_[wPhaseIdx] =0.0;
125 viscosity_[nPhaseIdx] = 0.0;
126 }
127
129 {
130 const auto element = *problem_.gridView().template begin<0>();
131 FluidState fluidState;
132 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
133 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
134 fluidState.setTemperature(problem_.temperature(element));
135 fluidState.setSaturation(wPhaseIdx, 1.);
136 fluidState.setSaturation(nPhaseIdx, 0.);
137 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
138 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
139 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
140 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
141
142 int size = gridView_.size(0);
143 rhs_.resize(size , 0.);
144 W_.resize(size);
145 }
146
147 void reset()
148 {
149 rhs_ = 0;
150 maxError_ = 0;
151 timeStep_ = 1;
152 }
153
154 void setErrorInfo(Scalar maxErr, Scalar dt)
155 {
156 maxError_ = maxErr;
157 timeStep_ = dt;
158 }
159
173 void assemble(const Element& element, int k = 1)
174 {
175 unsigned int numFaces = element.subEntities(1);
176 this->setcurrentsize(numFaces);
177
178 // clear assemble data
179 for (unsigned int i = 0; i < numFaces; i++)
180 {
181 this->b[i] = 0;
182 this->bctype[i].reset();
183 for (unsigned int j = 0; j < numFaces; j++)
184 this->A[i][j] = 0;
185 }
186
187 assembleV(element, k);
188 assembleBC(element, k);
189 }
190
191 // TODO/FIXME: this is only valid for linear problems where
192 // the local stiffness matrix is independend of the current
193 // solution. We need to implement this properly, but this
194 // should at least make the thing compile...
195 using VBlockType = Dune::FieldVector<Scalar, m>;
196 void assemble(const Element &cell, const Dune::BlockVector<VBlockType>& localSolution, int orderOfShapeFns = 1)
197 {
198 assemble(cell, orderOfShapeFns);
199 }
200
212 void assembleBoundaryCondition(const Element& element, int k = 1)
213 {
214 unsigned int numFaces = element.subEntities(1);
215 this->setcurrentsize(numFaces);
216
217 // clear assemble data
218 for (unsigned int i = 0; i < numFaces; i++)
219 {
220 this->b[i] = 0;
221 this->bctype[i].reset();
222 }
223
224 assembleBC(element, k);
225 }
226
227 // TODO doc me!
228 template<class Vector>
229 void completeRHS(const Element& element, Dune::FieldVector<int, 2*dim>& local2Global, Vector& f)
230 {
231 int eIdxGlobal = problem_.variables().index(element);
232 unsigned int numFaces = element.subEntities(1);
233
234 Dune::FieldVector<Scalar, 2 * dim> F(0.);
235 Scalar dInv = 0.;
236 computeReconstructionMatrices(element, W_[eIdxGlobal], F, dInv);
237
238 for (unsigned int i = 0; i < numFaces; i++)
239 {
240 if (!this->bc(i).isDirichlet(pressureEqIdx))
241 f[local2Global[i]][0] += (dInv * F[i] * rhs_[eIdxGlobal]);
242 }
243 }
244
245 // TODO doc me!
246 Scalar constructPressure(const Element& element, Dune::FieldVector<Scalar,2*dim>& pressTrace)
247 {
248 int eIdxGlobal = problem_.variables().index(element);
249 Scalar volume = element.geometry().volume();
250
251 PrimaryVariables sourceVec(0.0);
252 problem_.source(sourceVec, element);
253 Scalar qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
254 qmean += rhs_[eIdxGlobal];
255
256 Dune::FieldVector<Scalar, 2 * dim> F(0.);
257 Scalar dInv = 0.;
258 computeReconstructionMatrices(element, W_[eIdxGlobal], F, dInv);
259
260 return (dInv*(qmean + (F*pressTrace)));
261 }
262
263 // TODO doc me!
264 void constructVelocity(const Element& element, Dune::FieldVector<Scalar,2*dim>& vel,
265 Dune::FieldVector<Scalar,2*dim>& pressTrace, Scalar press)
266 {
267 int eIdxGlobal = problem_.variables().index(element);
268
269 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
270 for (const auto& intersection : intersections(gridView_, element))
271 {
272 faceVol[intersection.indexInInside()] = intersection.geometry().volume();
273 }
274
275 vel = 0;
276 for (int i = 0; i < 2*dim; i++)
277 for (int j = 0; j < 2*dim; j++)
278 vel[i] += W_[eIdxGlobal][i][j]*faceVol[j]*(press - pressTrace[j]);
279 }
280
281 // TODO doc me!
282 void constructVelocity(const Element& element, int fIdx, Scalar& vel,
283 Dune::FieldVector<Scalar,2*dim>& pressTrace, Scalar press)
284 {
285 int eIdxGlobal = problem_.variables().index(element);
286
287 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
288 for (const auto& intersection : intersections(gridView_, element))
289 {
290 faceVol[intersection.indexInInside()] = intersection.geometry().volume();
291 }
292
293 vel = 0;
294 for (int j = 0; j < 2*dim; j++)
295 vel += W_[eIdxGlobal][fIdx][j]*faceVol[j]*(press - pressTrace[j]);
296 }
297
298 // TODO doc me!
299 void computeReconstructionMatrices(const Element& element,
300 const Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& W,
301 Dune::FieldVector<Scalar, 2 * dim>& F, Scalar& dInv)
302 {
303 Dune::FieldVector<Scalar, 2 * dim> c(0);
304 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> Pi(0);
305
306 for (const auto& intersection : intersections(gridView_, element))
307 {
308 // local number of facet
309 int idx = intersection.indexInInside();
310
311 Scalar faceVol = intersection.geometry().volume();
312
313 // Corresponding to the element under consideration,
314 // calculate the part of the matrix C coupling velocities and element pressures.
315 // This is just a row vector of size numFaces.
316 // scale with volume
317 c[idx] = faceVol;
318 // Set up the element part of the matrix \Pi coupling velocities
319 // and pressure-traces. This is a diagonal matrix with entries given by faceVol.
320 Pi[idx][idx] = faceVol;
321 }
322
323 // Calculate the element part of the matrix D^{-1} = (c W c^T)^{-1} which is just a scalar value.
324 Dune::FieldVector<Scalar, 2 * dim> Wc(0);
325 W.umv(c, Wc);
326 dInv = 1.0 / (c * Wc);
327
328 // Calculate the element part of the matrix F = Pi W c^T which is a column vector.
329 F = 0;
330 Pi.umv(Wc, F);
331 }
332
333 void assembleElementMatrices(const Element& element, Dune::FieldVector<Scalar, 2 * dim>& faceVol,
334 Dune::FieldVector<Scalar, 2 * dim>& F,
335 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& Pi,
336 Scalar& dInv,
337 Scalar& qmean);
338
339
340 const Problem& problem() const
341 {
342 return problem_;
343 }
344
345private:
346 void assembleV(const Element& element, int k = 1);
347
348 void assembleBC(const Element& element, int k = 1);
349
350 // TODO doc me!
351 Scalar evaluateErrorTerm(CellData& cellData)
352 {
353 //error term for incompressible models to correct unphysical saturation over/undershoots due to saturation transport
354 // error reduction routine: volumetric error is damped and inserted to right hand side
355 Scalar sat = 0;
356 switch (saturationType)
357 {
358 case Sw:
359 sat = cellData.saturation(wPhaseIdx);
360 break;
361 case Sn:
362 sat = cellData.saturation(nPhaseIdx);
363 break;
364 }
365
366 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
367 if (sat < 0.0)
368 {
369 error = sat;
370 }
371 error /= timeStep_;
372
373 using std::abs;
374 Scalar errorAbs = abs(error);
375
376 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
377 && (!problem_.timeManager().willBeFinished()))
378 {
379 return ErrorTermFactor_ * error;
380 }
381 return 0.0;
382 }
383
384private:
385 Problem& problem_;
386 GridView gridView_;
387 Dune::DynamicVector<Scalar> rhs_;
388 std::vector<Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> > W_;
389 Scalar maxError_;
390 Scalar timeStep_;
391 Scalar ErrorTermFactor_; // Handling of error term: relaxation factor
392 Scalar ErrorTermLowerBound_; // Handling of error term: lower bound for error dampening
393 Scalar ErrorTermUpperBound_; //S Handling of error term: upper bound for error dampening
394
395 Scalar density_[numPhases];
396 Scalar viscosity_[numPhases];
397};
398
399// TODO doc me!
400template<class TypeTag>
401void MimeticTwoPLocalStiffness<TypeTag>::assembleV(const Element& element, int)
402{
403 unsigned int numFaces = element.subEntities(1);
404 this->setcurrentsize(numFaces);
405
406 int eIdxGlobal = problem_.variables().index(element);
407
408 // The notation is borrowed from Aarnes/Krogstadt/Lie 2006, Section 3.4.
409 // The matrix W developed here corresponds to one element-associated
410 // block of the matrix B^{-1} there.
411 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
412 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> Pi(0.);
413 Dune::FieldVector<Scalar, 2 * dim> F(0.);
414 Scalar dInv = 0.;
415 Scalar qmean = 0.;
416 this->assembleElementMatrices(element, faceVol, F, Pi, dInv, qmean);
417
418 // Calculate the element part of the matrix Pi W Pi^T.
419 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> PiWPiT(W_[eIdxGlobal]);
420 PiWPiT.rightmultiply(Pi);
421 PiWPiT.leftmultiply(Pi);
422
423 // Calculate the element part of the matrix F D^{-1} F^T.
424 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> FDinvFT(0);
425 for (unsigned int i = 0; i < numFaces; i++)
426 for (unsigned int j = 0; j < numFaces; j++)
427 FDinvFT[i][j] = dInv * F[i] * F[j];
428
429 // Calculate the element part of the matrix S = Pi W Pi^T - F D^{-1} F^T.
430 for (unsigned int i = 0; i < numFaces; i++)
431 for (unsigned int j = 0; j < numFaces; j++)
432 this->A[i][j] = PiWPiT[i][j] - FDinvFT[i][j];
433
434 // Calculate the source term F D^{-1} f
435 // NOT WORKING AT THE MOMENT
436 Scalar factor = dInv * qmean;
437 for (unsigned int i = 0; i < numFaces; i++)
438 this->b[i] = F[i] * factor;
439
440 // std::cout << "faceVol = " << faceVol << std::endl << "W = " << std::endl << W << std::endl
441 // << "c = " << c << std::endl << "Pi = " << std::endl << Pi << std::endl
442 // << "dinv = " << dinv << std::endl << "F = " << F << std::endl;
443 // std::cout << "dinvF = " << dinvF << ", q = " << qmean
444 // << ", b = " << this->b[0] << ", " << this->b[1] << ", " << this->b[2] << ", " << this->b[3] << std::endl;
445}
446
447// TODO doc me!
448template<class TypeTag>
450 Dune::FieldVector<Scalar, 2 * dim>& faceVol,
451 Dune::FieldVector<Scalar, 2 * dim>& F,
452 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& Pi,
453 Scalar& dInv,
454 Scalar& qmean)
455{
456 unsigned int numFaces = element.subEntities(1);
457 this->setcurrentsize(numFaces);
458
459 // get global coordinate of cell center
460 Dune::FieldVector<Scalar, dim> centerGlobal = element.geometry().center();
461
462 int eIdxGlobal = problem_.variables().index(element);
463
464 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
465
466 // cell volume
467 Scalar volume = element.geometry().volume();
468
469 // build the matrices R and ~N
470 Dune::FieldMatrix<Scalar, 2 * dim, dim> R(0), N(0);
471
472 // std::cout << "element " << elemId << ": center " << centerGlobal << std::endl;
473
474 // collect information needed for calculation of fluxes due to capillary-potential (pc + gravity!)
475 Scalar gravPotFace[2*dim];
476
477 Scalar gravPot = (problem_.bBoxMax() - centerGlobal) * problem_.gravity() * (density_[nPhaseIdx] - density_[wPhaseIdx]);
478
479 for (const auto& intersection : intersections(gridView_, element))
480 {
481 // local number of facet
482 int fIdx = intersection.indexInInside();
483
484 Dune::FieldVector<Scalar, dim> faceGlobal = intersection.geometry().center();
485 faceVol[fIdx] = intersection.geometry().volume();
486
487 // get normal vector
488 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
489
490 N[fIdx] = unitOuterNormal;
491
492 for (int k = 0; k < dim; k++)
493 // move origin to the center of gravity
494 R[fIdx][k] = faceVol[fIdx] * (faceGlobal[k] - centerGlobal[k]);
495
496 gravPotFace[fIdx] = (problem_.bBoxMax() - faceGlobal) * problem_.gravity() * (density_[nPhaseIdx] - density_[wPhaseIdx]);
497 }
498
499 // proceed along the lines of Algorithm 1 from
500 // Brezzi/Lipnikov/Simonicini M3AS 2005
501 // (1) orthonormalize columns of the matrix R
502 Scalar norm = R[0][0] * R[0][0];
503 using std::sqrt;
504 for (unsigned int i = 1; i < numFaces; i++)
505 norm += R[i][0] * R[i][0];
506 norm = sqrt(norm);
507 for (unsigned int i = 0; i < numFaces; i++)
508 R[i][0] /= norm;
509 Scalar weight = R[0][1] * R[0][0];
510 for (unsigned int i = 1; i < numFaces; i++)
511 weight += R[i][1] * R[i][0];
512 for (unsigned int i = 0; i < numFaces; i++)
513 R[i][1] -= weight * R[i][0];
514 norm = R[0][1] * R[0][1];
515 for (unsigned int i = 1; i < numFaces; i++)
516 norm += R[i][1] * R[i][1];
517 norm = sqrt(norm);
518 for (unsigned int i = 0; i < numFaces; i++)
519 R[i][1] /= norm;
520 if (dim == 3)
521 {
522 Scalar weight1 = R[0][2] * R[0][0];
523 Scalar weight2 = R[0][2] * R[0][1];
524 for (unsigned int i = 1; i < numFaces; i++)
525 {
526 weight1 += R[i][2] * R[i][0];
527 weight2 += R[i][2] * R[i][1];
528 }
529 for (unsigned int i = 0; i < numFaces; i++)
530 R[i][2] -= weight1 * R[i][0] + weight2 * R[i][1];
531 norm = R[0][2] * R[0][2];
532 for (unsigned int i = 1; i < numFaces; i++)
533 norm += R[i][2] * R[i][2];
534 using std::sqrt;
535 norm = sqrt(norm);
536 for (unsigned int i = 0; i < numFaces; i++)
537 R[i][2] /= norm;
538 }
539 // std::cout << "~R =\dim" << R;
540
541 // (2) Build the matrix ~D
542 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> D(0);
543 for (unsigned int s = 0; s < numFaces; s++)
544 {
545 Dune::FieldVector<Scalar, 2 * dim> es(0);
546 es[s] = 1;
547 for (unsigned int k = 0; k < numFaces; k++)
548 {
549 D[k][s] = es[k];
550 for (unsigned int i = 0; i < dim; i++)
551 {
552 D[k][s] -= R[s][i] * R[k][i];
553 }
554 }
555 }
556
557// std::cout<<"D = "<<D<<"\n";
558
559 // eval diffusion tensor, ASSUMING to be constant over each cell
560 Dune::FieldMatrix<Scalar, dim, dim> mobility(problem_.spatialParams().intrinsicPermeability(element));
561 mobility *= (cellData.mobility(wPhaseIdx) + cellData.mobility(nPhaseIdx));
562
563 Scalar traceMob = mobility[0][0];
564 for (int i = 1; i < dim; i++)
565 traceMob += mobility[i][i];
566 D *= 2.0 * traceMob / volume;
567 // std::cout << "u~D =\dim" << D;
568
569 // (3) Build the matrix W = Minv
570 Dune::FieldMatrix<Scalar, 2 * dim, dim> NK(N);
571 NK.rightmultiply(mobility);
572
573 for (unsigned int i = 0; i < numFaces; i++)
574 {
575 for (unsigned int j = 0; j < numFaces; j++)
576 {
577 W_[eIdxGlobal][i][j] = NK[i][0] * N[j][0];
578 for (unsigned int k = 1; k < dim; k++)
579 W_[eIdxGlobal][i][j] += NK[i][k] * N[j][k];
580 }
581 }
582
583 W_[eIdxGlobal] /= volume;
584 W_[eIdxGlobal] += D;
585
586 // std::cout << "W = \dim" << W;
587
588 // Now the notation is borrowed from Aarnes/Krogstadt/Lie 2006, Section 3.4.
589 // The matrix W developed so far corresponds to one element-associated
590 // block of the matrix B^{-1} there.
591
592 // Corresponding to the element under consideration,
593 // calculate the part of the matrix C coupling velocities and element pressures.
594 // This is just a row vector of size numFaces.
595 // scale with volume
596 Dune::FieldVector<Scalar, 2 * dim> c(0);
597 for (unsigned int i = 0; i < numFaces; i++)
598 c[i] = faceVol[i];
599
600 // Set up the element part of the matrix \Pi coupling velocities
601 // and pressure-traces. This is a diagonal matrix with entries given by faceVol.
602 for (unsigned int i = 0; i < numFaces; i++)
603 Pi[i][i] = faceVol[i];
604
605 // Calculate the element part of the matrix D^{-1} = (c W c^T)^{-1} which is just a scalar value.
606 Dune::FieldVector<Scalar, 2 * dim> Wc(0);
607 W_[eIdxGlobal].umv(c, Wc);
608 dInv = 1.0 / (c * Wc);
609
610 // Calculate the element part of the matrix F = Pi W c^T which is a column vector.
611 F = 0;
612 Pi.umv(Wc, F);
613 // std::cout << "Pi = \dim" << Pi << "c = " << c << ", F = " << F << std::endl;
614
615 // accumulate fluxes due to capillary potential (pc + gravity!)
616 for (const auto& intersection : intersections(gridView_, element))
617 {
618 int idx = intersection.indexInInside();
619
620 Scalar fracFlow = 0;
621
622 Scalar flux = 0;
623 for (int j = 0; j < 2 * dim; j++)
624 flux += W_[eIdxGlobal][idx][j] * faceVol[j] * (gravPot - gravPotFace[j]);
625
626 //it is enough to evaluate the capillary/gravity flux for neighbors -> not needed anyway at the boundaries!
627 if (intersection.neighbor())
628 {
629 int eIdxGlobalNeighbor = problem_.variables().index(intersection.outside());
630 if (flux >= 0.)
631 {
632 switch (pressureType)
633 {
634 case pw:
635 {
636 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
637 break;
638 }
639 case pn:
640 {
641 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
642 break;
643 }
644 }
645
646
647 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
648 rhs_[eIdxGlobalNeighbor] += (faceVol[idx] * fracFlow * flux);
649 }
650 }
651 else if (intersection.boundary())
652 {
653 BoundaryTypes bctype;
654 problem_.boundaryTypes(bctype, intersection);
655
656 if (bctype.isDirichlet(pressureEqIdx))
657 {
658 if (flux > 0. || !bctype.isDirichlet(satEqIdx))
659 {
660 switch (pressureType)
661 {
662 case pw:
663 {
664 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
665 break;
666 }
667 case pn:
668 {
669 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
670 break;
671 }
672 }
673 }
674 else if (flux < 0. && bctype.isDirichlet(satEqIdx))
675 {
676 PrimaryVariables boundValues(0.0);
677 problem_.dirichlet(boundValues, intersection);
678
679 Scalar krw = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element),
680 boundValues[saturationIdx]);
681 Scalar krn = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element),
682 boundValues[saturationIdx]);
683
684 switch (pressureType)
685 {
686 case pw:
687 {
688 fracFlow = krn / viscosity_[nPhaseIdx]
689 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
690 break;
691 }
692 case pn:
693 {
694 fracFlow = -krw / viscosity_[wPhaseIdx]
695 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
696 break;
697 }
698 }
699 }
700
701 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
702 }
703 }
704 }
705
706 PrimaryVariables sourceVec(0.0);
707 problem_.source(sourceVec, element);
708 qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
709
710 qmean += evaluateErrorTerm(cellData) * volume;
711}
712
713// TODO doc me!
714template<class TypeTag>
715void MimeticTwoPLocalStiffness<TypeTag>::assembleBC(const Element& element, int k)
716{
717 // evaluate boundary conditions via intersection iterator
718 for (const auto& intersection : intersections(gridView_, element))
719 {
720 int faceIndex = intersection.indexInInside();
721 if (intersection.boundary())
722 {
723 problem_.boundaryTypes(this->bctype[faceIndex], intersection);
724 PrimaryVariables boundValues(0.0);
725
726 if (this->bctype[faceIndex].isNeumann(pressureEqIdx))
727 {
728 problem_.neumann(boundValues, intersection);
729 Scalar J = (boundValues[wPhaseIdx]/density_[wPhaseIdx] + boundValues[nPhaseIdx]/density_[nPhaseIdx]);
730 this->b[faceIndex] -= J * intersection.geometry().volume();
731 }
732 else if (this->bctype[faceIndex].isDirichlet(pressureEqIdx))
733 {
734 problem_.dirichlet(boundValues, intersection);
735 if (pressureType == pw)
736 this->b[faceIndex] = boundValues[pressureIdx] +
737 (problem_.bBoxMax() - intersection.geometry().center()) * problem_.gravity() * density_[wPhaseIdx];
738 else if (pressureType == pn)
739 this->b[faceIndex] = boundValues[pressureIdx] +
740 (problem_.bBoxMax() - intersection.geometry().center()) * problem_.gravity() * density_[nPhaseIdx];
741 else
742 this->b[faceIndex] = boundValues[pressureIdx];
743 }
744 }
745 }
746
747}
748
749}
750#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
Definition of boundary condition types, extend if necessary.
Base class for assembling local stiffness matrices.
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 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
std::string mobility(int phaseIdx) noexcept
I/O name of mobility for multiphase systems.
Definition: name.hh:101
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Base class for local assemblers.
Definition: localstiffness.hh:60
void setcurrentsize(int s)
Sets the current size of the local stiffness matrix.
Definition: localstiffness.hh:233
const BoundaryTypes & bc(int i) const
Accesses boundary condition for each dof.
Definition: localstiffness.hh:223
Dune::Matrix< MBlockType > A
Definition: localstiffness.hh:250
std::vector< VBlockType > b
Definition: localstiffness.hh:251
std::vector< BoundaryTypes > bctype
Definition: localstiffness.hh:252
compute local stiffness matrix for conforming finite elements for the full 2-phase pressure equation
Definition: mimetic.hh:53
void reset()
Definition: mimetic.hh:147
void assemble(const Element &element, int k=1)
Assembles local stiffness matrix for given element and order.
Definition: mimetic.hh:173
MimeticTwoPLocalStiffness(Problem &problem, bool levelBoundaryAsDirichlet, const GridView &gridView, bool procBoundaryAsDirichlet=true)
Constructor.
Definition: mimetic.hh:114
@ m
Definition: mimetic.hh:106
void constructVelocity(const Element &element, int fIdx, Scalar &vel, Dune::FieldVector< Scalar, 2 *dim > &pressTrace, Scalar press)
Definition: mimetic.hh:282
void assembleBoundaryCondition(const Element &element, int k=1)
Assembles only boundary conditions for given element.
Definition: mimetic.hh:212
void constructVelocity(const Element &element, Dune::FieldVector< Scalar, 2 *dim > &vel, Dune::FieldVector< Scalar, 2 *dim > &pressTrace, Scalar press)
Definition: mimetic.hh:264
const Problem & problem() const
Definition: mimetic.hh:340
void completeRHS(const Element &element, Dune::FieldVector< int, 2 *dim > &local2Global, Vector &f)
Definition: mimetic.hh:229
void initialize()
Definition: mimetic.hh:128
void assemble(const Element &cell, const Dune::BlockVector< VBlockType > &localSolution, int orderOfShapeFns=1)
Definition: mimetic.hh:196
void setErrorInfo(Scalar maxErr, Scalar dt)
Definition: mimetic.hh:154
Scalar constructPressure(const Element &element, Dune::FieldVector< Scalar, 2 *dim > &pressTrace)
Definition: mimetic.hh:246
void assembleElementMatrices(const Element &element, Dune::FieldVector< Scalar, 2 *dim > &faceVol, Dune::FieldVector< Scalar, 2 *dim > &F, Dune::FieldMatrix< Scalar, 2 *dim, 2 *dim > &Pi, Scalar &dInv, Scalar &qmean)
Definition: mimetic.hh:449
void computeReconstructionMatrices(const Element &element, const Dune::FieldMatrix< Scalar, 2 *dim, 2 *dim > &W, Dune::FieldVector< Scalar, 2 *dim > &F, Scalar &dInv)
Definition: mimetic.hh:299
@ size
Definition: mimetic.hh:110
Dune::FieldVector< Scalar, m > VBlockType
Definition: mimetic.hh:195
Specifies the properties for immiscible 2p diffusion/pressure models.