3.2-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{
57
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 = getPropValue<TypeTag, Properties::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 = getPropValue<TypeTag, Properties::PressureFormulation>(),
87 saturationType = getPropValue<TypeTag, Properties::SaturationFormulation>()
88 };
89
90 using Grid = typename GridView::Grid;
91 using Element = typename GridView::Traits::template Codim<0>::Entity;
92
95 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
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 default:
365 DUNE_THROW(Dune::NotImplemented, "Only saturation formulation Sw and Sn are implemented!");
366 }
367
368 Scalar error = (sat > 1.0) ? sat - 1.0 : 0.0;
369 if (sat < 0.0)
370 {
371 error = sat;
372 }
373 error /= timeStep_;
374
375 using std::abs;
376 Scalar errorAbs = abs(error);
377
378 if ((errorAbs * timeStep_ > 1e-6) && (errorAbs > ErrorTermLowerBound_ * maxError_)
379 && (!problem_.timeManager().willBeFinished()))
380 {
381 return ErrorTermFactor_ * error;
382 }
383 return 0.0;
384 }
385
386private:
387 Problem& problem_;
388 GridView gridView_;
389 Dune::DynamicVector<Scalar> rhs_;
390 std::vector<Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> > W_;
391 Scalar maxError_;
392 Scalar timeStep_;
393 Scalar ErrorTermFactor_; // Handling of error term: relaxation factor
394 Scalar ErrorTermLowerBound_; // Handling of error term: lower bound for error dampening
395 Scalar ErrorTermUpperBound_; //S Handling of error term: upper bound for error dampening
396
397 Scalar density_[numPhases];
398 Scalar viscosity_[numPhases];
399};
400
401// TODO doc me!
402template<class TypeTag>
403void MimeticTwoPLocalStiffness<TypeTag>::assembleV(const Element& element, int)
404{
405 unsigned int numFaces = element.subEntities(1);
406 this->setcurrentsize(numFaces);
407
408 int eIdxGlobal = problem_.variables().index(element);
409
410 // The notation is borrowed from Aarnes/Krogstadt/Lie 2006, Section 3.4.
411 // The matrix W developed here corresponds to one element-associated
412 // block of the matrix B^{-1} there.
413 Dune::FieldVector<Scalar, 2 * dim> faceVol(0);
414 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> Pi(0.);
415 Dune::FieldVector<Scalar, 2 * dim> F(0.);
416 Scalar dInv = 0.;
417 Scalar qmean = 0.;
418 this->assembleElementMatrices(element, faceVol, F, Pi, dInv, qmean);
419
420 // Calculate the element part of the matrix Pi W Pi^T.
421 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> PiWPiT(W_[eIdxGlobal]);
422 PiWPiT.rightmultiply(Pi);
423 PiWPiT.leftmultiply(Pi);
424
425 // Calculate the element part of the matrix F D^{-1} F^T.
426 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> FDinvFT(0);
427 for (unsigned int i = 0; i < numFaces; i++)
428 for (unsigned int j = 0; j < numFaces; j++)
429 FDinvFT[i][j] = dInv * F[i] * F[j];
430
431 // Calculate the element part of the matrix S = Pi W Pi^T - F D^{-1} F^T.
432 for (unsigned int i = 0; i < numFaces; i++)
433 for (unsigned int j = 0; j < numFaces; j++)
434 this->A[i][j] = PiWPiT[i][j] - FDinvFT[i][j];
435
436 // Calculate the source term F D^{-1} f
437 // NOT WORKING AT THE MOMENT
438 Scalar factor = dInv * qmean;
439 for (unsigned int i = 0; i < numFaces; i++)
440 this->b[i] = F[i] * factor;
441
442 // std::cout << "faceVol = " << faceVol << std::endl << "W = " << std::endl << W << std::endl
443 // << "c = " << c << std::endl << "Pi = " << std::endl << Pi << std::endl
444 // << "dinv = " << dinv << std::endl << "F = " << F << std::endl;
445 // std::cout << "dinvF = " << dinvF << ", q = " << qmean
446 // << ", b = " << this->b[0] << ", " << this->b[1] << ", " << this->b[2] << ", " << this->b[3] << std::endl;
447}
448
449// TODO doc me!
450template<class TypeTag>
452 Dune::FieldVector<Scalar, 2 * dim>& faceVol,
453 Dune::FieldVector<Scalar, 2 * dim>& F,
454 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim>& Pi,
455 Scalar& dInv,
456 Scalar& qmean)
457{
458 unsigned int numFaces = element.subEntities(1);
459 this->setcurrentsize(numFaces);
460
461 // get global coordinate of cell center
462 Dune::FieldVector<Scalar, dim> centerGlobal = element.geometry().center();
463
464 int eIdxGlobal = problem_.variables().index(element);
465
466 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
467
468 // cell volume
469 Scalar volume = element.geometry().volume();
470
471 // build the matrices R and ~N
472 Dune::FieldMatrix<Scalar, 2 * dim, dim> R(0), N(0);
473
474 // std::cout << "element " << elemId << ": center " << centerGlobal << std::endl;
475
476 // collect information needed for calculation of fluxes due to capillary-potential (pc + gravity!)
477 Scalar gravPotFace[2*dim];
478
479 Scalar gravPot = (problem_.bBoxMax() - centerGlobal) * problem_.gravity() * (density_[nPhaseIdx] - density_[wPhaseIdx]);
480
481 for (const auto& intersection : intersections(gridView_, element))
482 {
483 // local number of facet
484 int fIdx = intersection.indexInInside();
485
486 Dune::FieldVector<Scalar, dim> faceGlobal = intersection.geometry().center();
487 faceVol[fIdx] = intersection.geometry().volume();
488
489 // get normal vector
490 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
491
492 N[fIdx] = unitOuterNormal;
493
494 for (int k = 0; k < dim; k++)
495 // move origin to the center of gravity
496 R[fIdx][k] = faceVol[fIdx] * (faceGlobal[k] - centerGlobal[k]);
497
498 gravPotFace[fIdx] = (problem_.bBoxMax() - faceGlobal) * problem_.gravity() * (density_[nPhaseIdx] - density_[wPhaseIdx]);
499 }
500
501 // proceed along the lines of Algorithm 1 from
502 // Brezzi/Lipnikov/Simonicini M3AS 2005
503 // (1) orthonormalize columns of the matrix R
504 Scalar norm = R[0][0] * R[0][0];
505 using std::sqrt;
506 for (unsigned int i = 1; i < numFaces; i++)
507 norm += R[i][0] * R[i][0];
508 norm = sqrt(norm);
509 for (unsigned int i = 0; i < numFaces; i++)
510 R[i][0] /= norm;
511 Scalar weight = R[0][1] * R[0][0];
512 for (unsigned int i = 1; i < numFaces; i++)
513 weight += R[i][1] * R[i][0];
514 for (unsigned int i = 0; i < numFaces; i++)
515 R[i][1] -= weight * R[i][0];
516 norm = R[0][1] * R[0][1];
517 for (unsigned int i = 1; i < numFaces; i++)
518 norm += R[i][1] * R[i][1];
519 norm = sqrt(norm);
520 for (unsigned int i = 0; i < numFaces; i++)
521 R[i][1] /= norm;
522 if (dim == 3)
523 {
524 Scalar weight1 = R[0][2] * R[0][0];
525 Scalar weight2 = R[0][2] * R[0][1];
526 for (unsigned int i = 1; i < numFaces; i++)
527 {
528 weight1 += R[i][2] * R[i][0];
529 weight2 += R[i][2] * R[i][1];
530 }
531 for (unsigned int i = 0; i < numFaces; i++)
532 R[i][2] -= weight1 * R[i][0] + weight2 * R[i][1];
533 norm = R[0][2] * R[0][2];
534 for (unsigned int i = 1; i < numFaces; i++)
535 norm += R[i][2] * R[i][2];
536 using std::sqrt;
537 norm = sqrt(norm);
538 for (unsigned int i = 0; i < numFaces; i++)
539 R[i][2] /= norm;
540 }
541 // std::cout << "~R =\dim" << R;
542
543 // (2) Build the matrix ~D
544 Dune::FieldMatrix<Scalar, 2 * dim, 2 * dim> D(0);
545 for (unsigned int s = 0; s < numFaces; s++)
546 {
547 Dune::FieldVector<Scalar, 2 * dim> es(0);
548 es[s] = 1;
549 for (unsigned int k = 0; k < numFaces; k++)
550 {
551 D[k][s] = es[k];
552 for (unsigned int i = 0; i < dim; i++)
553 {
554 D[k][s] -= R[s][i] * R[k][i];
555 }
556 }
557 }
558
559// std::cout<<"D = "<<D<<"\n";
560
561 // eval diffusion tensor, ASSUMING to be constant over each cell
562 Dune::FieldMatrix<Scalar, dim, dim> mobility(problem_.spatialParams().intrinsicPermeability(element));
563 mobility *= (cellData.mobility(wPhaseIdx) + cellData.mobility(nPhaseIdx));
564
565 Scalar traceMob = mobility[0][0];
566 for (int i = 1; i < dim; i++)
567 traceMob += mobility[i][i];
568 D *= 2.0 * traceMob / volume;
569 // std::cout << "u~D =\dim" << D;
570
571 // (3) Build the matrix W = Minv
572 Dune::FieldMatrix<Scalar, 2 * dim, dim> NK(N);
573 NK.rightmultiply(mobility);
574
575 for (unsigned int i = 0; i < numFaces; i++)
576 {
577 for (unsigned int j = 0; j < numFaces; j++)
578 {
579 W_[eIdxGlobal][i][j] = NK[i][0] * N[j][0];
580 for (unsigned int k = 1; k < dim; k++)
581 W_[eIdxGlobal][i][j] += NK[i][k] * N[j][k];
582 }
583 }
584
585 W_[eIdxGlobal] /= volume;
586 W_[eIdxGlobal] += D;
587
588 // std::cout << "W = \dim" << W;
589
590 // Now the notation is borrowed from Aarnes/Krogstadt/Lie 2006, Section 3.4.
591 // The matrix W developed so far corresponds to one element-associated
592 // block of the matrix B^{-1} there.
593
594 // Corresponding to the element under consideration,
595 // calculate the part of the matrix C coupling velocities and element pressures.
596 // This is just a row vector of size numFaces.
597 // scale with volume
598 Dune::FieldVector<Scalar, 2 * dim> c(0);
599 for (unsigned int i = 0; i < numFaces; i++)
600 c[i] = faceVol[i];
601
602 // Set up the element part of the matrix \Pi coupling velocities
603 // and pressure-traces. This is a diagonal matrix with entries given by faceVol.
604 for (unsigned int i = 0; i < numFaces; i++)
605 Pi[i][i] = faceVol[i];
606
607 // Calculate the element part of the matrix D^{-1} = (c W c^T)^{-1} which is just a scalar value.
608 Dune::FieldVector<Scalar, 2 * dim> Wc(0);
609 W_[eIdxGlobal].umv(c, Wc);
610 dInv = 1.0 / (c * Wc);
611
612 // Calculate the element part of the matrix F = Pi W c^T which is a column vector.
613 F = 0;
614 Pi.umv(Wc, F);
615 // std::cout << "Pi = \dim" << Pi << "c = " << c << ", F = " << F << std::endl;
616
617 // accumulate fluxes due to capillary potential (pc + gravity!)
618 for (const auto& intersection : intersections(gridView_, element))
619 {
620 int idx = intersection.indexInInside();
621
622 Scalar fracFlow = 0;
623
624 Scalar flux = 0;
625 for (int j = 0; j < 2 * dim; j++)
626 flux += W_[eIdxGlobal][idx][j] * faceVol[j] * (gravPot - gravPotFace[j]);
627
628 //it is enough to evaluate the capillary/gravity flux for neighbors -> not needed anyway at the boundaries!
629 if (intersection.neighbor())
630 {
631 int eIdxGlobalNeighbor = problem_.variables().index(intersection.outside());
632 if (flux >= 0.)
633 {
634 switch (pressureType)
635 {
636 case pw:
637 {
638 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
639 break;
640 }
641 case pn:
642 {
643 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
644 break;
645 }
646 default:
647 DUNE_THROW(Dune::NotImplemented, "Only pressure formulation pw and pn are implemented!");
648 }
649
650
651 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
652 rhs_[eIdxGlobalNeighbor] += (faceVol[idx] * fracFlow * flux);
653 }
654 }
655 else if (intersection.boundary())
656 {
657 BoundaryTypes bctype;
658 problem_.boundaryTypes(bctype, intersection);
659
660 if (bctype.isDirichlet(pressureEqIdx))
661 {
662 if (flux > 0. || !bctype.isDirichlet(satEqIdx))
663 {
664 switch (pressureType)
665 {
666 case pw:
667 {
668 fracFlow = cellData.fracFlowFunc(nPhaseIdx);
669 break;
670 }
671 case pn:
672 {
673 fracFlow = -cellData.fracFlowFunc(wPhaseIdx);
674 break;
675 }
676 default:
677 DUNE_THROW(Dune::NotImplemented, "Only pressure formulation pw and pn are implemented!");
678 }
679 }
680 else if (flux < 0. && bctype.isDirichlet(satEqIdx))
681 {
682 PrimaryVariables boundValues(0.0);
683 problem_.dirichlet(boundValues, intersection);
684
685 Scalar krw = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element),
686 boundValues[saturationIdx]);
687 Scalar krn = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element),
688 boundValues[saturationIdx]);
689
690 switch (pressureType)
691 {
692 case pw:
693 {
694 fracFlow = krn / viscosity_[nPhaseIdx]
695 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
696 break;
697 }
698 case pn:
699 {
700 fracFlow = -krw / viscosity_[wPhaseIdx]
701 / (krw / viscosity_[wPhaseIdx] + krn / viscosity_[nPhaseIdx]);
702 break;
703 }
704 default:
705 DUNE_THROW(Dune::NotImplemented, "Only pressure formulation pw and pn are implemented!");
706 }
707 }
708
709 rhs_[eIdxGlobal] -= (faceVol[idx] * fracFlow * flux);
710 }
711 }
712 }
713
714 PrimaryVariables sourceVec(0.0);
715 problem_.source(sourceVec, element);
716 qmean = volume * (sourceVec[wPhaseIdx]/density_[wPhaseIdx] + sourceVec[nPhaseIdx]/density_[nPhaseIdx]);
717
718 qmean += evaluateErrorTerm(cellData) * volume;
719}
720
721// TODO doc me!
722template<class TypeTag>
723void MimeticTwoPLocalStiffness<TypeTag>::assembleBC(const Element& element, int k)
724{
725 // evaluate boundary conditions via intersection iterator
726 for (const auto& intersection : intersections(gridView_, element))
727 {
728 int faceIndex = intersection.indexInInside();
729 if (intersection.boundary())
730 {
731 problem_.boundaryTypes(this->bctype[faceIndex], intersection);
732 PrimaryVariables boundValues(0.0);
733
734 if (this->bctype[faceIndex].isNeumann(pressureEqIdx))
735 {
736 problem_.neumann(boundValues, intersection);
737 Scalar J = (boundValues[wPhaseIdx]/density_[wPhaseIdx] + boundValues[nPhaseIdx]/density_[nPhaseIdx]);
738 this->b[faceIndex] -= J * intersection.geometry().volume();
739 }
740 else if (this->bctype[faceIndex].isDirichlet(pressureEqIdx))
741 {
742 problem_.dirichlet(boundValues, intersection);
743 if (pressureType == pw)
744 this->b[faceIndex] = boundValues[pressureIdx] +
745 (problem_.bBoxMax() - intersection.geometry().center()) * problem_.gravity() * density_[wPhaseIdx];
746 else if (pressureType == pn)
747 this->b[faceIndex] = boundValues[pressureIdx] +
748 (problem_.bBoxMax() - intersection.geometry().center()) * problem_.gravity() * density_[nPhaseIdx];
749 else
750 this->b[faceIndex] = boundValues[pressureIdx];
751 }
752 }
753 }
754
755}
756
757}
758#endif
Definition of boundary condition types, extend if necessary.
Base class for assembling local stiffness matrices.
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 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
@ size
Definition: mimetic.hh:110
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:451
@ m
Definition: mimetic.hh:106
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
Dune::FieldVector< Scalar, m > VBlockType
Definition: mimetic.hh:195
Specifies the properties for immiscible 2p diffusion/pressure models.