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