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