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