3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fv2dpressureadaptive.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 *****************************************************************************/
25#ifndef DUMUX_FV2DPRESSURE2P2C_ADAPTIVE_HH
26#define DUMUX_FV2DPRESSURE2P2C_ADAPTIVE_HH
27
28// dune environent:
29#include <dune/istl/bvector.hh>
30#include <dune/istl/operators.hh>
31#include <dune/istl/solvers.hh>
32#include <dune/istl/preconditioners.hh>
33
34// dumux environment
37// include 2p mpfa pressure model
39
40#include <dumux/common/math.hh>
42
43namespace Dumux {
44
76template<class TypeTag> class FV2dPressure2P2CAdaptive
77: public FVPressure2P2C<TypeTag>
78{
79 //the model implementation
82
86
88
90 enum
91 {
92 dim = GridView::dimension, dimWorld = GridView::dimensionworld
93 };
94 enum
95 {
96 pw = Indices::pressureW,
97 pn = Indices::pressureN
98 };
99 enum
100 {
101 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
102 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx
103 };
104 enum
105 {
106 rhs = BaseType::rhs, matrix = BaseType::matrix,
107 };
108
109 // using declarations to abbreviate several dune classes...
110 using Intersection = typename GridView::Intersection;
111 using IntersectionIterator = typename GridView::IntersectionIterator;
112
113 // convenience shortcuts for Vectors/Matrices
114 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
115 using TransmissivityMatrix = Dune::FieldVector<Scalar,dim+1>;
116 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
117 using PhaseVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumPhases>()>;
118
119 // the typenames used for the stiffness matrix and solution vector
121
123protected:
124 Problem& problem()
125 {
126 return this->problem_;
127 }
128 const Problem& problem() const
129 {
130 return this->problem_;
131 }
132
133public:
134 //initializes the matrix to store the system of equations
135 void initializeMatrix();
136 //function which assembles the system of equations to be solved
137 void assemble(bool first);
138
139 void getMpfaFlux(const IntersectionIterator&, const CellData&);
140
141 // mpfa transmissibilities
142 int computeTransmissibilities(const IntersectionIterator&,
143 IntersectionIterator&,
144 TransmissivityMatrix&,
145 TransmissivityMatrix&,
146 GlobalPosition&,
147 int&);
148
151 {
152 int gridSize = problem().gridView().size(0);
153 this->pressure().resize(gridSize);
154
155 for(int i=0; i< gridSize; i++)
156 {
157 this->pressure()[i]
158 = problem().variables().cellData(i).pressure(this->pressureType);
159 }
160 }
161
163
168 {
169 enableVolumeIntegral = getParam<bool>("Impet.EnableVolumeIntegral");
170 enableMPFA = getParam<bool>("GridAdapt.EnableMultiPointFluxApproximation");
171
172 if(getParam<int>("GridAdapt.MaxInteractionVolumes") != 1)
174 else
175 enableSecondHalfEdge = false;
176
177 // prepare rotation matrix: R_ is initialized as zero matrix
178 R_=0.;
179 // evaluate matrix R
180 if (dim == 2)
181 {
182 R_[0][1] = 1;
183 R_[1][0] = -1;
184 }
185 }
186
187private:
188 Problem& problem_;
190 Implementation &asImp_()
191 { return *static_cast<Implementation *>(this);}
192
194 const Implementation &asImp_() const
195 { return *static_cast<const Implementation *>(this);}
196
197protected:
198 // mpfa transmissibilities from mpfal2pfa
199 int transmissibilityAdapter_(const IntersectionIterator&,
200 const IntersectionIterator&,
201 IntersectionIterator&,
202 GlobalPosition&,
203 TransmissivityMatrix&);
204
206 DimMatrix R_;
212};
213
214
216
222template<class TypeTag>
224{
225 int gridSize_ = problem().gridView().size(0);
226 // update RHS vector, matrix
227 this->A_.setSize (gridSize_,gridSize_); //
228 this->f_.resize(gridSize_);
229
230 // determine matrix row sizes
231 for (const auto& element : elements(problem().gridView()))
232 {
233 // cell index
234 int globalIdxI = problem().variables().index(element);
235 CellData& cellDataI = problem().variables().cellData(globalIdxI);
236
237 // initialize row size
238 int rowSize = 1;
239
240 // set perimeter to zero
241 cellDataI.perimeter() = 0;
242
243 // prepare storage for all found 3rd cells
244 std::vector<int> foundAdditionals;
245
246 int numberOfIntersections = 0;
247 // run through all intersections with neighbors
248 const auto isEndIt = problem().gridView().iend(element);
249 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
250 {
251 const auto& intersection = *isIt;
252
253 cellDataI.perimeter() += intersection.geometry().volume();
254 numberOfIntersections++;
255 if (intersection.neighbor())
256 {
257 rowSize++;
258
259 // if mpfa is used, more entries might be needed if both halfedges are regarded
260 if (enableMPFA && (enableSecondHalfEdge && intersection.outside().level() != element.level()))
261 {
262 GlobalPosition globalPos3(0.);
263 int globalIdx3=-1;
264 TransmissivityMatrix T(0.);
265 auto additionalIsIt = isIt;
266 TransmissivityMatrix additionalT(0.);
267 // compute Transmissibilities: also examines subcontrolvolume information
268 int halfedgesStored
269 = problem().variables().getMpfaData(intersection, additionalIsIt, T, additionalT, globalPos3, globalIdx3);
270 if (halfedgesStored == 0)
271 halfedgesStored = problem().pressureModel().computeTransmissibilities(isIt,additionalIsIt, T,additionalT,
272 globalPos3, globalIdx3 );
273
274 if(halfedgesStored == 2)
275 {
276 bool increaseRowSize = true;
277 //check if additional cell is ordinary neighbor of eIt
278 for (const auto& intersection2 : intersections(problem_.gridView(), element))
279 {
280 if(!intersection2.neighbor())
281 continue;
282 if(additionalIsIt->outside() == intersection2.outside() )
283 increaseRowSize = false;
284 }
285 //also check if additional cell was already used for another interaction triangle
286 for (unsigned int i = 0; i < foundAdditionals.size(); i++)
287 if(foundAdditionals[i] == problem().variables().index(additionalIsIt->outside()))
288 increaseRowSize = false;
289
290 if (increaseRowSize)
291 {
292 rowSize++;
293 foundAdditionals.push_back(problem().variables().index(additionalIsIt->outside()));
294 }
295 }
296 }
297 }
298 }
299
300 cellDataI.fluxData().resize(numberOfIntersections);
301 this->A_.setrowsize(globalIdxI, rowSize);
302 }
303 this->A_.endrowsizes();
304
305 // determine position of matrix entries
306 for (const auto& element : elements(problem().gridView()))
307 {
308 // cell index
309 int globalIdxI = problem().variables().index(element);
310
311 // add diagonal index
312 this->A_.addindex(globalIdxI, globalIdxI);
313
314 // run through all intersections with neighbors
315 const auto isEndIt = problem().gridView().iend(element);
316 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
317 {
318 const auto& intersection = *isIt;
319
320 if (intersection.neighbor())
321 {
322 // access neighbor
323 int globalIdxJ = problem().variables().index(intersection.outside());
324
325 // add off diagonal index
326 this->A_.addindex(globalIdxI, globalIdxJ);
327
328 // if mpfa is used, more entries might be needed if both halfedges are regarded
329 if (enableMPFA && (enableSecondHalfEdge && intersection.outside().level() != element.level()))
330 {
331 GlobalPosition globalPos3(0.);
332 int globalIdx3=-1;
333 TransmissivityMatrix T(0.);
334 auto additionalIsIt = isIt;
335 TransmissivityMatrix additionalT(0.);
336 // compute Transmissibilities: also examines subcontrolvolume information
337 int halfedgesStored
338 = problem().variables().getMpfaData(intersection, additionalIsIt, T, additionalT, globalPos3, globalIdx3);
339 if (halfedgesStored == 0)
340 halfedgesStored = problem().pressureModel().computeTransmissibilities(isIt,additionalIsIt, T,additionalT,
341 globalPos3, globalIdx3 );
342 // add off diagonal index if 2 half-edges regarded
343 if(halfedgesStored == 2)
344 this->A_.addindex(globalIdxI, problem().variables().index(additionalIsIt->outside()));
345 }
346 }
347 }
348 }
349 this->A_.endindices();
350
351 return;
352}
353
355
363template<class TypeTag>
365{
366 if(first)
367 {
368 BaseType::assemble(true);
369 return;
370 }
371
372 initializeMatrix();
373 // initialization: set matrix A_ to zero
374 this->A_ = 0;
375 this->f_ = 0;
376
377 for (const auto& element : elements(problem().gridView()))
378 {
379 // get the global index of the cell
380 int globalIdxI = problem().variables().index(element);
381
382 // assemble interior element contributions
383 if (element.partitionType() == Dune::InteriorEntity)
384 {
385 // get the cell data
386 CellData& cellDataI = problem().variables().cellData(globalIdxI);
387
388 Dune::FieldVector<Scalar, 2> entries(0.);
389
390 /***** source term ***********/
391 problem().pressureModel().getSource(entries, element, cellDataI, first);
392 this->f_[globalIdxI] += entries[rhs];
393
394 /***** flux term ***********/
395 // iterate over all faces of the cell
396 auto isEndIt = problem().gridView().iend(element);
397 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
398 {
399 const auto& intersection = *isIt;
400
401 /************* handle interior face *****************/
402 if (intersection.neighbor())
403 {
404 auto elementNeighbor = intersection.outside();
405
406 int globalIdxJ = problem().variables().index(elementNeighbor);
407
408 //check for hanging nodes
409 //take a hanging node never from the element with smaller level!
410 bool haveSameLevel = (element.level() == elementNeighbor.level());
411 // calculate only from one side, but add matrix entries for both sides
412 // the last condition is needed to properly assemble in the presence
413 // of ghost elements
414 if (getPropValue<TypeTag, Properties::VisitFacesOnlyOnce>()
415 && (globalIdxI > globalIdxJ) && haveSameLevel
416 && elementNeighbor.partitionType() == Dune::InteriorEntity)
417 continue;
418
419 entries = 0;
420 //check for hanging nodes
421 if(!haveSameLevel && enableMPFA)
422 {
423 problem().pressureModel().getMpfaFlux(isIt, cellDataI);
424 }
425 else
426 {
427 problem().pressureModel().getFlux(entries, intersection, cellDataI, first);
428
429 //set right hand side
430 this->f_[globalIdxI] -= entries[rhs];
431
432 // set diagonal entry
433 this->A_[globalIdxI][globalIdxI] += entries[matrix];
434
435 // set off-diagonal entry
436 this->A_[globalIdxI][globalIdxJ] -= entries[matrix];
437
438 // The second condition is needed to not spoil the ghost element entries
439 if (getPropValue<TypeTag, Properties::VisitFacesOnlyOnce>()
440 && elementNeighbor.partitionType() == Dune::InteriorEntity)
441 {
442 this->f_[globalIdxJ] += entries[rhs];
443 this->A_[globalIdxJ][globalIdxJ] += entries[matrix];
444 this->A_[globalIdxJ][globalIdxI] -= entries[matrix];
445 }
446 }
447
448 } // end neighbor
449
450 /************* boundary face ************************/
451 else
452 {
453 entries = 0;
454 problem().pressureModel().getFluxOnBoundary(entries, intersection, cellDataI, first);
455
456 //set right hand side
457 this->f_[globalIdxI] += entries[rhs];
458 // set diagonal entry
459 this->A_[globalIdxI][globalIdxI] += entries[matrix];
460 }
461 } //end interfaces loop
462// printmatrix(std::cout, this->A_, "global stiffness matrix", "row", 11, 3);
463
464 /***** storage term ***********/
465 entries = 0;
466 problem().pressureModel().getStorage(entries, element, cellDataI, first);
467 this->f_[globalIdxI] += entries[rhs];
468 // set diagonal entry
469 this->A_[globalIdxI][globalIdxI] += entries[matrix];
470 }
471 // assemble overlap and ghost element contributions
472 else
473 {
474 this->A_[globalIdxI] = 0.0;
475 this->A_[globalIdxI][globalIdxI] = 1.0;
476 this->f_[globalIdxI] = this->pressure()[globalIdxI];
477 }
478
479 } // end grid traversal
480// printmatrix(std::cout, this->A_, "global stiffness matrix after assempling", "row", 11,3);
481// printvector(std::cout, this->f_, "right hand side", "row", 10);
482 return;
483}
484
486
506template<class TypeTag>
507void FV2dPressure2P2CAdaptive<TypeTag>::getMpfaFlux(const IntersectionIterator& intersectionIterator,
508 const CellData& cellDataI)
509{
510 // acess Cell I
511 auto elementI = intersectionIterator->inside();
512 int globalIdxI = problem().variables().index(elementI);
513
514 // get global coordinate of cell center
515 const GlobalPosition& globalPos = elementI.geometry().center();
516
517 // cell volume & perimeter, assume linear map here
518 Scalar volume = elementI.geometry().volume();
519 Scalar perimeter = cellDataI.perimeter();
520
521 const GlobalPosition& gravity_ = problem().gravity();
522
523 // get absolute permeability
524 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
525
526 // access neighbor
527 auto neighbor = intersectionIterator->outside();
528 int globalIdxJ = problem().variables().index(neighbor);
529 CellData& cellDataJ = problem().variables().cellData(globalIdxJ);
530
531 // gemotry info of neighbor
532 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
533
534 // distance vector between barycenters
535 GlobalPosition distVec = globalPosNeighbor - globalPos;
536
537 // compute distance between cell centers
538 Scalar dist = distVec.two_norm();
539
540 GlobalPosition unitDistVec(distVec);
541 unitDistVec /= dist;
542
543 DimMatrix permeabilityJ
544 = problem().spatialParams().intrinsicPermeability(neighbor);
545
546 // compute vectorized permeabilities
547 DimMatrix meanPermeability(0);
548 harmonicMeanMatrix(meanPermeability, permeabilityI, permeabilityJ);
549
550 Dune::FieldVector<Scalar, dim> permeability(0);
551 meanPermeability.mv(unitDistVec, permeability);
552
553 // get average density for gravity flux
554 Scalar rhoMeanW = 0.5 * (cellDataI.density(wPhaseIdx) + cellDataJ.density(wPhaseIdx));
555 Scalar rhoMeanNW = 0.5 * (cellDataI.density(nPhaseIdx) + cellDataJ.density(nPhaseIdx));
556
557 // reset potential gradients
558 Scalar potentialW = 0;
559 Scalar potentialNW = 0;
560
561 // determine volume derivatives in neighbor
562 if (!cellDataJ.hasVolumeDerivatives())
563 this->volumeDerivatives(globalPosNeighbor, neighbor);
564
565 Scalar dv_dC1 = (cellDataJ.dv(wPhaseIdx)
566 + cellDataI.dv(wPhaseIdx)) / 2; // dV/dm1= dv/dC^1
567 Scalar dv_dC2 = (cellDataJ.dv(nPhaseIdx)
568 + cellDataI.dv(nPhaseIdx)) / 2;
569
570 Scalar graddv_dC1 = (cellDataJ.dv(wPhaseIdx)
571 - cellDataI.dv(wPhaseIdx)) / dist;
572 Scalar graddv_dC2 = (cellDataJ.dv(nPhaseIdx)
573 - cellDataI.dv(nPhaseIdx)) / dist;
574
575
576// potentialW = problem().variables().potentialWetting(globalIdxI, isIndex);
577// potentialNW = problem().variables().potentialNonwetting(globalIdxI, isIndex);
578//
579// densityW = (potentialW > 0.) ? densityWI : densityWJ;
580// densityNW = (potentialNW > 0.) ? densityNWI : densityNWJ;
581//
582// densityW = (potentialW == 0.) ? rhoMeanW : densityW;
583// densityNW = (potentialNW == 0.) ? rhoMeanNW : densityNW;
584 //jochen: central weighting for gravity term
585 Scalar densityW = rhoMeanW;
586 Scalar densityNW = rhoMeanNW;
587
588 // Prepare MPFA
590 GlobalPosition globalPos3(0.);
591 int globalIdx3=-1;
592 TransmissivityMatrix T(0.);
593 // prepare second half-edge
594 auto additionalIsIt = intersectionIterator;
595 TransmissivityMatrix additionalT(0.);
596
597 int halfedgesStored = problem().variables().getMpfaData(*intersectionIterator,
598 additionalIsIt, T, additionalT,
599 globalPos3, globalIdx3 );
600 if (halfedgesStored == 0)
601 halfedgesStored = problem().pressureModel().computeTransmissibilities(intersectionIterator,
602 additionalIsIt, T, additionalT,
603 globalPos3, globalIdx3 );
604
605 if(!halfedgesStored)
606 Dune::dgrave << "something went wrong getting mpfa data on cell " << globalIdxI << std::endl;
607
608 // shortcurts mpfa case
609 CellData& cellData3 = problem().variables().cellData(globalIdx3);
610 Scalar temp1 = globalPos * gravity_;
611 Scalar temp2 = globalPosNeighbor * gravity_;
612 Scalar temp3 = globalPos3 * gravity_;
613
614 potentialW = (cellDataI.pressure(wPhaseIdx)-temp1*densityW) * T[2]
615 + (cellDataJ.pressure(wPhaseIdx)-temp2*densityW) * T[0]
616 + (cellData3.pressure(wPhaseIdx)-temp3*densityW) * T[1];
617 potentialNW = (cellDataI.pressure(nPhaseIdx)-temp1*densityNW) * T[2]
618 + (cellDataJ.pressure(nPhaseIdx)-temp2*densityNW) * T[0]
619 + (cellData3.pressure(nPhaseIdx)-temp3*densityNW) * T[1];
620 // regard second half edge, if there is one
621 if(halfedgesStored == 2)
622 {
623 int AdditionalIdx = problem().variables().index(additionalIsIt->outside());
624 CellData& cellDataAdditional = problem().variables().cellData(AdditionalIdx);
625 potentialW += (cellDataI.pressure(wPhaseIdx)-temp1*densityW) * additionalT[2]
626 + (cellDataJ.pressure(wPhaseIdx)-temp2*densityW) * additionalT[0]
627 + (cellDataAdditional.pressure(wPhaseIdx)
628 -(additionalIsIt->outside().geometry().center()*gravity_)
629 *densityW) * additionalT[1];
630 potentialNW += (cellDataI.pressure(nPhaseIdx)-temp1*densityNW) * additionalT[2]
631 + (cellDataJ.pressure(nPhaseIdx)-temp2*densityNW) * additionalT[0]
632 + (cellDataAdditional.pressure(nPhaseIdx)
633 -(additionalIsIt->outside().geometry().center()*gravity_)
634 *densityNW) * additionalT[1];
635 }
636
637
638 // initialize convenience shortcuts
639 Scalar lambdaW(0.), lambdaN(0.);
640 Scalar dV_w(0.), dV_n(0.); // dV_a = \sum_k \rho_a * dv/dC^k * X^k_a
641 Scalar gV_w(0.), gV_n(0.); // multipaper eq(3.3) line 3 analogon dV_w
642
643
644 //do the upwinding of the mobility depending on the phase potentials
645 if (potentialW >= 0.)
646 {
647 dV_w = (dv_dC1 * cellDataI.massFraction(wPhaseIdx, wCompIdx)
648 + dv_dC2 * cellDataI.massFraction(wPhaseIdx, nCompIdx));
649 lambdaW = cellDataI.mobility(wPhaseIdx);
650 gV_w = (graddv_dC1 * cellDataI.massFraction(wPhaseIdx, wCompIdx)
651 + graddv_dC2 * cellDataI.massFraction(wPhaseIdx, nCompIdx));
652 dV_w *= cellDataI.density(wPhaseIdx);
653 gV_w *= cellDataI.density(wPhaseIdx);
654 }
655 else
656 {
657 dV_w = (dv_dC1 * cellDataJ.massFraction(wPhaseIdx, wCompIdx)
658 + dv_dC2 * cellDataJ.massFraction(wPhaseIdx, nCompIdx));
659 lambdaW = cellDataJ.mobility(wPhaseIdx);
660 gV_w = (graddv_dC1 * cellDataJ.massFraction(wPhaseIdx, wCompIdx)
661 + graddv_dC2 * cellDataJ.massFraction(wPhaseIdx, nCompIdx));
662 dV_w *= cellDataJ.density(wPhaseIdx);
663 gV_w *= cellDataJ.density(wPhaseIdx);
664 }
665 if (potentialNW >= 0.)
666 {
667 dV_n = (dv_dC1 * cellDataI.massFraction(nPhaseIdx, wCompIdx)
668 + dv_dC2 * cellDataI.massFraction(nPhaseIdx, nCompIdx));
669 lambdaN = cellDataI.mobility(nPhaseIdx);
670 gV_n = (graddv_dC1 * cellDataI.massFraction(nPhaseIdx, wCompIdx)
671 + graddv_dC2 * cellDataI.massFraction(nPhaseIdx, nCompIdx));
672 dV_n *= cellDataI.density(nPhaseIdx);
673 gV_n *= cellDataI.density(nPhaseIdx);
674 }
675 else
676 {
677 dV_n = (dv_dC1 * cellDataJ.massFraction(nPhaseIdx, wCompIdx)
678 + dv_dC2 * cellDataJ.massFraction(nPhaseIdx, nCompIdx));
679 lambdaN = cellDataJ.mobility(nPhaseIdx);
680 gV_n = (graddv_dC1 * cellDataJ.massFraction(nPhaseIdx, wCompIdx)
681 + graddv_dC2 * cellDataJ.massFraction(nPhaseIdx, nCompIdx));
682 dV_n *= cellDataJ.density(nPhaseIdx);
683 gV_n *= cellDataJ.density(nPhaseIdx);
684 }
685
687 /* extend T with other matrix entries and assemble to A_ */
688 this->A_[globalIdxI][globalIdxJ] += (lambdaW * dV_w + lambdaN * dV_n) * T[0];
689 this->A_[globalIdxI][globalIdx3] += (lambdaW * dV_w + lambdaN * dV_n) * T[1];
690 this->A_[globalIdxI][globalIdxI] += (lambdaW * dV_w + lambdaN * dV_n) * T[2];
691
692 // add gravity to RHS vector
693 this->f_[globalIdxI] += (densityW * lambdaW * dV_w + densityNW * lambdaN * dV_n) * temp2 * T[0];
694 this->f_[globalIdxI] += (densityW * lambdaW * dV_w + densityNW * lambdaN * dV_n) * temp3 * T[1];
695 this->f_[globalIdxI] += (densityW * lambdaW * dV_w + densityNW * lambdaN * dV_n) * temp1 * T[2];
696
697 // weithing accounts for the fraction of the subcontrol volume
698 Scalar weightingFactor = volume / perimeter; // transforms flux through area A -> V * A/perimeter
699 if(enableVolumeIntegral) // switch off volume integral for mpfa case
700 {
701 // correct for area integral
702 this->A_[globalIdxI][globalIdxJ] -= weightingFactor * (lambdaW * gV_w + lambdaN * gV_n)*T[0];
703 this->A_[globalIdxI][globalIdx3] -= weightingFactor * (lambdaW * gV_w + lambdaN * gV_n)*T[1];
704 this->A_[globalIdxI][globalIdxI] -= weightingFactor * (lambdaW * gV_w + lambdaN * gV_n)*T[2];
705
706 // add gravity to RHS vector
707 this->f_[globalIdxI] -= weightingFactor * (densityW * lambdaW * gV_w + densityNW * lambdaN * gV_n) * temp2 * T[0];
708 this->f_[globalIdxI] -= weightingFactor * (densityW * lambdaW * gV_w + densityNW * lambdaN * gV_n) * temp3 * T[1];
709 this->f_[globalIdxI] -= weightingFactor * (densityW * lambdaW * gV_w + densityNW * lambdaN * gV_n) * temp1 * T[2];
710 }
711
712 // capillary pressure flux
713 Scalar pcGradient = cellDataI.capillaryPressure() * T[2]
714 + cellDataJ.capillaryPressure() * T[0]
715 + cellData3.capillaryPressure() * T[1];
716
717 if (this->pressureType == pw)
718 pcGradient *= + lambdaN * dV_n - enableVolumeIntegral * weightingFactor * lambdaN * gV_n;
719 else if (this->pressureType == pn)
720 pcGradient *= - lambdaW * dV_w + enableVolumeIntegral * weightingFactor * lambdaW * gV_w;
721
722 this->f_[globalIdxI] += pcGradient;
723
724 // include second half-edge
725 if(halfedgesStored == 2)
726 {
727 int AdditionalIdx = problem().variables().index(additionalIsIt->outside());
728 CellData& cellDataAdditional = problem().variables().cellData(AdditionalIdx);
729
730 /* extend T with other matrix entries and assemble to A_ */
731 this->A_[globalIdxI][globalIdxJ] += (lambdaW * dV_w + lambdaN * dV_n) * additionalT[0];
732 this->A_[globalIdxI][AdditionalIdx] += (lambdaW * dV_w + lambdaN * dV_n) * additionalT[1];
733 this->A_[globalIdxI][globalIdxI] += (lambdaW * dV_w + lambdaN * dV_n) * additionalT[2];
734
735 // add gravity to RHS vector
736 this->f_[globalIdxI] += (densityW * lambdaW * dV_w + densityNW * lambdaN * dV_n) * temp2 * additionalT[0];
737 this->f_[globalIdxI] += (densityW * lambdaW * dV_w + densityNW * lambdaN * dV_n)
738 * (additionalIsIt->outside().geometry().center()*gravity_) * additionalT[1];
739 this->f_[globalIdxI] += (densityW * lambdaW * dV_w + densityNW * lambdaN * dV_n) * temp1 * additionalT[2];
740
741 if(enableVolumeIntegral) // switch off volume integral for mpfa case
742 {
743 // correct for area integral
744 this->A_[globalIdxI][globalIdxJ] -= weightingFactor * (lambdaW * gV_w + lambdaN * gV_n)*additionalT[0];
745 this->A_[globalIdxI][AdditionalIdx] -= weightingFactor * (lambdaW * gV_w + lambdaN * gV_n)*additionalT[1];
746 this->A_[globalIdxI][globalIdxI] -= weightingFactor * (lambdaW * gV_w + lambdaN * gV_n)*additionalT[2];
747
748 // add gravity to RHS vector
749 this->f_[globalIdxI] -= weightingFactor * (densityW * lambdaW * gV_w + densityNW * lambdaN * gV_n) * temp2 * additionalT[0];
750 this->f_[globalIdxI] -= weightingFactor * (densityW * lambdaW * gV_w + densityNW * lambdaN * gV_n)
751 * (additionalIsIt->outside().geometry().center()*gravity_) * additionalT[1];
752 this->f_[globalIdxI] -= weightingFactor * (densityW * lambdaW * gV_w + densityNW * lambdaN * gV_n) * temp1 * additionalT[2];
753 }
754
755 // capillary pressure flux
756 pcGradient = cellDataI.capillaryPressure() * additionalT[2]
757 + cellDataJ.capillaryPressure() * additionalT[0]
758 + cellDataAdditional.capillaryPressure() * additionalT[1];
759 if (this->pressureType == pw)
760 pcGradient *= + lambdaN * dV_n - enableVolumeIntegral * weightingFactor * lambdaN * gV_n;
761 else if (this->pressureType == pn)
762 pcGradient *= - lambdaW * dV_w + enableVolumeIntegral * weightingFactor * lambdaW * gV_w;
763
764 this->f_[globalIdxI] += pcGradient;
765 }
766}
767
769
804template <class TypeTag>
806 IntersectionIterator& additionalIntersectionIt,
807 TransmissivityMatrix& T,
808 TransmissivityMatrix& additionalT,
809 GlobalPosition& globalPos3,
810 int& globalIdx3)
811{
812 const auto& intersection = *isIt;
813
814 // get geometry information of cellI = cell1, cellJ = cell2
815 auto element = intersection.inside();
816 auto neighbor = intersection.outside();
817 GlobalPosition globalPos1 = element.geometry().center();
818 GlobalPosition globalPos2 = neighbor.geometry().center();
819 DimMatrix K1(problem().spatialParams().intrinsicPermeability(element));
820 DimMatrix K2(problem().spatialParams().intrinsicPermeability(neighbor));
821
823 // geometry and Data of face IJ in nomenclature of mpfa
824 GlobalPosition globalPosFace12 = intersection.geometry().center();
825 GlobalPosition integrationOuterNormaln12 = intersection.centerUnitOuterNormal();
826 integrationOuterNormaln12 *= intersection.geometry().volume() / 2.0; // TODO: 2.0 only in 2D
827
828
829 // nextIs points to next intersection
830 auto nextIs = isIt;
831 ++nextIs;
832 if (nextIs == problem().gridView().iend(element))
833 nextIs = problem().gridView().ibegin(element);
834
835 // get last intersection : --intersection does not exist
836 // paceingIt loops one IS bevore prevIs
837 auto prevIs = problem().gridView().ibegin(element);
838 auto paceingIt = prevIs;
839 for (++paceingIt; paceingIt != problem().gridView().iend(element); ++paceingIt)
840 {
841 if (!paceingIt->neighbor()) // continue if no neighbor found
842 ++prevIs; // we investigate next paceingIt -> prevIs is also increased
843 else if (paceingIt->outside() == intersection.outside()) // we already found prevIs
844 break;
845 else if (paceingIt == problem().gridView().iend(element))
846 prevIs = paceingIt; // this could only happen if isIt is begin, so prevIs has to be last.
847 else
848 ++prevIs; // we investigate next paceingIt -> prevIs is also increased
849 }
850
852 auto face13 = isIt; // as long as face13 == intersection, it is still not found!
853 auto face23 = isIt; // as long as face23 == intersection, it is still not found!
854
855 // store other intersection for the other interaction region for the other half-edge
856 auto isEndIt = problem().gridView().iend(neighbor);
857 for (auto isIt23 = problem().gridView().ibegin(neighbor); isIt23 != isEndIt; ++isIt23)
858 {
859 const auto& intersection23 = *isIt23;
860
861 // stop search if found
862 if(face13->outside() != intersection.outside())
863 break;
864
865 if(!intersection23.neighbor())
866 continue;
867
868 // either prevIs or nextIs is face13, it is that with common interface with cell2
869 // investigate if prevIs points to cell 3
870 if (prevIs->neighbor())
871 {
872 if (prevIs->outside() == intersection23.outside())
873 {
874 face23 = isIt23;
875 face13 = prevIs;
876 additionalIntersectionIt = nextIs;
877 }
878 }
879 // investigate if nextIs points to cell 3
880 if (nextIs->neighbor())
881 {
882 if (nextIs->outside() == intersection23.outside())
883 {
884 face23 = isIt23;
885 face13 = nextIs;
886 additionalIntersectionIt = prevIs;
887 }
888 }
889 }
890 if(face13->outside() == intersection.outside()) //means isIt13 not found yet
891 Dune::dgrave << "is 13 not found!!!" << std::endl;
892
893 // get information of cell3
894 globalPos3 = face13->outside().geometry().center();
895 globalIdx3 = problem().variables().index(face13->outside());
896 // get absolute permeability of neighbor cell 3
897 DimMatrix K3(problem().spatialParams().intrinsicPermeability(face13->outside()));
898
899
900 // get the intersection node /bar^{x_3} between 'isIt' and 'isIt13', denoted as 'corner123'
901 // initialization of corner123
902 GlobalPosition corner123(0.);
903
904 // get the global coordinate of corner123
905 for (int i = 0; i < intersection.geometry().corners(); ++i)
906 {
907 for (int j = 0; j < face13->geometry().corners(); ++j)
908 {
909 if (face13->geometry().corner(j) == intersection.geometry().corner(i))
910 {
911 corner123 = face13->geometry().corner(j);
912
913 // stop outer (i) and inner (j) loop
914 i = intersection.geometry().corners();
915 break;
916 }
917 }
918 }
919
921 // center of face in global coordinates, i.e., the midpoint of edge 'isIt24'
922 GlobalPosition globalPosFace23 = face23->geometry().center();
923
924 // get face volume
925 Scalar face23vol = face23->geometry().volume();
926
927 // get outer normal vector scaled with half volume of face 'isIt24'
928 Dune::FieldVector<Scalar,dimWorld> integrationOuterNormaln23 = face23->centerUnitOuterNormal();
929 integrationOuterNormaln23 *= face23vol / 2.0; //TODO: 2.0 only for 2D
930
931 // compute normal vectors nu1-nu7 in triangle R for first half edge
932 GlobalPosition nu1(0);
933 R_.umv(globalPosFace12-globalPos2 ,nu1); //globalPos2 = globalPosNeighbor
934
935 GlobalPosition nu2(0);
936 R_.umv(globalPos2-globalPosFace23, nu2);
937
938 GlobalPosition nu3(0);
939 R_.umv(globalPosFace23-globalPos3, nu3);
940
941 GlobalPosition nu4(0);
942 R_.umv(globalPos3-corner123, nu4);
943
944 GlobalPosition nu5(0);
945 R_.umv(corner123-globalPos1, nu5); //globalPos1 = globalPos
946
947 GlobalPosition nu6(0);
948 R_.umv(globalPos1-globalPosFace12, nu6);
949
950 GlobalPosition nu7(0);
951 R_.umv(corner123-globalPos2, nu7);
952
953 // compute T, i.e., the area of quadrilateral made by normal vectors 'nu'
954 GlobalPosition Rnu2(0);
955 R_.umv(nu2, Rnu2);
956 Scalar T1 = nu1 * Rnu2;
957
958 GlobalPosition Rnu4(0);
959 R_.umv(nu4, Rnu4);
960 Scalar T2 = nu3 * Rnu4;
961
962 GlobalPosition Rnu6(0);
963 R_.umv(nu6, Rnu6);
964 Scalar T3 = nu5 * Rnu6;
965
966 // compute components needed for flux calculation, denoted as 'omega' and 'chi'
967 GlobalPosition K2nu1(0);
968 K2.umv(nu1, K2nu1); // permeabilityJ = K2 = perm of cell 2 around x_1
969 GlobalPosition K2nu2(0);
970 K2.umv(nu2, K2nu2);
971 GlobalPosition K3nu3(0);
972 K3.umv(nu3, K3nu3);
973 GlobalPosition K3nu4(0);
974 K3.umv(nu4, K3nu4);
975 GlobalPosition K1nu5(0);
976 K1.umv(nu5, K1nu5); // permeabilityI = K1 = perm of cell 1 around x_3
977 GlobalPosition K1nu6(0);
978 K1.umv(nu6, K1nu6);
979 GlobalPosition Rnu1(0);
980 R_.umv(nu1, Rnu1);
981
982 double omega111 = /*lambda2 */ (integrationOuterNormaln23 * K2nu1)/T1;
983 double omega112 = /*lambda2 */ (integrationOuterNormaln23 * K2nu2)/T1;
984 double omega211 = /*lambda2 */ (integrationOuterNormaln12 * K2nu1)/T1;
985 double omega212 = /*lambda2 */ (integrationOuterNormaln12 * K2nu2)/T1;
986 double omega123 = /*lambda3 */ (integrationOuterNormaln23 * K3nu3)/T2;
987 double omega124 = /*lambda3 */ (integrationOuterNormaln23 * K3nu4)/T2;
988 double omega235 = /*lambda1 */ (integrationOuterNormaln12 * K1nu5)/T3;
989 double omega236 = /*lambda1 */ (integrationOuterNormaln12 * K1nu6)/T3;
990 double chi711 = (nu7 * Rnu1)/T1;
991 double chi712 = (nu7 * Rnu2)/T1;
992
994 // compute transmissibility matrix T = CA^{-1}B+D
995 DimMatrix C(0), A(0);
996 Dune::FieldMatrix<Scalar,dim,2*dim-dim+1> D(0), B(0);
997
998 // evaluate matrix C, D, A, B
999 C[0][0] = -omega111;
1000 C[0][1] = -omega112;
1001 C[1][0] = -omega211;
1002 C[1][1] = -omega212;
1003
1004 D[0][0] = omega111 + omega112;
1005 D[1][0] = omega211 + omega212;
1006
1007 A[0][0] = omega111 - omega124 - omega123*chi711;
1008 A[0][1] = omega112 - omega123*chi712;
1009 A[1][0] = omega211 - omega236*chi711;
1010 A[1][1] = omega212 - omega235 - omega236*chi712;
1011
1012 B[0][0] = omega111 + omega112 + omega123*(1.0 - chi711 - chi712);
1013 B[0][1] = -omega123 - omega124;
1014 B[1][0] = omega211 + omega212 + omega236*(1.0 - chi711 - chi712);
1015 B[1][2] = -omega235 - omega236;
1016
1017 // compute T
1018 A.invert();
1019 D += B.leftmultiply(C.rightmultiply(A));
1020 T = D[1];
1021 if(!enableSecondHalfEdge )//or abs(intersection.centerUnitOuterNormal()[0])<0.5) // [0]<0.5 => switch off vertical 2hes
1022 {
1023 T *= 2;
1024 // set your map entry
1025 problem().variables().storeMpfaData(intersection, T, globalPos3, globalIdx3);
1026 return 1; // indicates that only 1 halfedge was regarded
1027 }
1028 else
1029 {
1030 // derive additional T for second half edge
1031
1032 // get the intersection node /bar^{x_3} between 'isIt' and 'aditionalIntersection', denoted as 'corner1245'
1033 // initialization of corner1245
1034 GlobalPosition corner1245(0.);
1035
1036 // get the global coordinate of corner1245, which is not connected with face13
1037 auto tempIntersection = face13;
1038 bool corner1245found = false;
1039 // ensure iterator increases over local end
1040 if (tempIntersection == problem().gridView().iend(element))
1041 tempIntersection = problem().gridView().ibegin(element);
1042 while (!corner1245found)
1043 {
1044 ++tempIntersection;
1045
1046 // ensure iterator increases over local end
1047 if (tempIntersection == problem().gridView().iend(element))
1048 tempIntersection = problem().gridView().ibegin(element);
1049 // enshure we do not arrive at isIt
1050 if (tempIntersection == isIt)
1051 continue;
1052
1053 // loop over both corners of is
1054 for (int i = 0; i < intersection.geometry().corners(); ++i)
1055 {
1056 // test if a corner of additionalIntersectionIt also lies on is
1057 for (int j = 0; j < tempIntersection->geometry().corners(); ++j)
1058 {
1059 if (tempIntersection->geometry().corner(j) == intersection.geometry().corner(i))
1060 {
1061 corner1245 = tempIntersection->geometry().corner(j);
1062 additionalIntersectionIt = tempIntersection;
1063
1064 // stop outer (i) and inner (j) loop
1065 i = intersection.geometry().corners();
1066 // stop also Intersection loop
1067 corner1245found = true;
1068 break;
1069 }
1070 }
1071 }
1072 }
1073 if (!additionalIntersectionIt->neighbor())// or (additionalIntersectionIt->outside().level() == intersection.inside().level()))
1074 {
1075 T *= 2;
1076 // set your map entry
1077 problem().variables().storeMpfaData(intersection, T, globalPos3, globalIdx3);
1078 return 1;
1079 }
1080
1081 // use Markus mpfa-l implementation
1082 this->transmissibilityAdapter_(isIt, face23, additionalIntersectionIt,
1083 corner1245, additionalT);
1084
1085 // store transmissivity data for second half edge
1086 problem().variables().storeMpfaData(intersection, additionalIntersectionIt, T,additionalT, globalPos3, globalIdx3);
1087 // return that information about two interaction volumes have been stored
1088 return 2;
1089 }
1090}
1091
1093
1119template<class TypeTag>
1121 const IntersectionIterator& face23_2p2cnaming,
1122 IntersectionIterator& additionalIntersectionIt,
1123 GlobalPosition& corner1234,
1124 TransmissivityMatrix& additionalT)
1125{
1126 const auto& intersection = *isIt;
1127
1128 /****** find all 4 faces *********/
1129 // store additionalIntersection, which is face23 in 2p mpfa-l naming scheme
1130 auto isIt23 = additionalIntersectionIt;
1131 auto isIt14 = isIt; //has to be initialized with something
1132
1133 // search for 4th intersection that connects to corner1245, which is needed
1134 // for markus implementation of the mpfa
1135
1136 // reuse corner1245found bool
1137 bool face14found = false;
1138 // repeat search procedure with isIt and face23 (in 2p2c naming)
1139 auto tempIntersection = face23_2p2cnaming;
1140
1141 while (!face14found)
1142 {
1143 ++tempIntersection;
1144
1145 // ensure iterator increases over local end of neighbor J
1146 if (tempIntersection== problem().gridView().iend(intersection.outside()))
1147 tempIntersection = problem().gridView().ibegin(intersection.outside());
1148
1149 if(!tempIntersection->neighbor())
1150 continue;
1151
1152 // enshure we hace not arrived at isIt (but seen from other side)
1153 if (tempIntersection->outside() == intersection.inside())
1154 continue; // restart loop to recheck after next increase for local end
1155
1156 // loop over both corners of is
1157 for (int i = 0; i < intersection.geometry().corners(); ++i)
1158 {
1159 // test if a corner of additionalIntersectionIt also lies on is
1160 for (int j = 0; j < tempIntersection->geometry().corners(); ++j)
1161 {
1162 if (tempIntersection->geometry().corner(j) == intersection.geometry().corner(i))
1163 {
1164// // test if this tempIntersection is better than additionalIntersectionIt
1165// if (tempIntersection->outside().level() > additionalIntersectionIt->outside().level())
1166// {
1167// // select this as element for second half edge
1168// additionalIntersectionIt = tempIntersection;
1169// corner1245 = additionalIntersectionIt->geometry().corner(j);
1170// }
1171 isIt14 = tempIntersection;
1172 // stop outer (i) and inner (j) loop
1173 i = intersection.geometry().corners();
1174 // stop also Intersection loop
1175 face14found = true;
1176 break;
1177 }
1178 }
1179 }
1180 }
1181 /**** end find 4 faces **/
1182
1183 // create Interaction Volume object
1184 FVMPFALInteractionVolume<TypeTag> interactionVolume(problem().gridView().grid());
1185
1186 interactionVolume.setCenterPosition(corner1234);
1187
1188 //*************** store pointer 1
1189 interactionVolume.setSubVolumeElement(intersection.outside(), 0);
1190 interactionVolume.setIndexOnElement(intersection.indexInOutside(), 0, 0);
1191 interactionVolume.setIndexOnElement(isIt14->indexInInside(), 0, 1);
1192
1193 // center of face in global coordinates, i.e., the midpoint of edge 'isIt12'
1194 const GlobalPosition& globalPosFace12 = intersection.geometry().center();
1195
1196 // get face volume
1197 Scalar faceVol12 = intersection.geometry().volume() / 2.0;
1198
1199 // get outer normal vector scaled with half volume of face 'isIt12'
1200 Dune::FieldVector<Scalar, dimWorld> unitOuterNormal12 = intersection.centerUnitOuterNormal();
1201 unitOuterNormal12 *=-1;
1202
1203 // center of face in global coordinates, i.e., the midpoint of edge 'isIt14'
1204 const GlobalPosition& globalPosFace41 = isIt14->geometry().center();
1205
1206 // get face volume
1207 Scalar faceVol41 = isIt14->geometry().volume() / 2.0;
1208
1209 // get outer normal vector scaled with half volume of face 'isIt14': for numbering of n see Aavatsmark, Eigestad
1210 Dune::FieldVector<Scalar, dimWorld> unitOuterNormal14 = isIt14->centerUnitOuterNormal();
1211
1212 interactionVolume.setNormal(unitOuterNormal12, 0, 0);
1213 interactionVolume.setNormal(unitOuterNormal14, 0, 1);
1214 //get the normals of from cell 2 and 4
1215 unitOuterNormal14 *= -1;
1216 unitOuterNormal12 *= -1;
1217 interactionVolume.setFaceArea(faceVol12, 0, 0);
1218 interactionVolume.setFaceArea(faceVol41, 0, 1);
1219 interactionVolume.setFacePosition(globalPosFace12, 0, 0);
1220 interactionVolume.setFacePosition(globalPosFace41, 0, 1);
1221
1222 // access neighbor cell 2 of 'isIt12'
1223 auto element2 = intersection.inside();
1224
1225 //**************** store pointer 2
1226 interactionVolume.setSubVolumeElement(element2, 1);
1227// interactionVolume.setIndexOnElement(intersection.indexInInside(), 1, 1);
1228 interactionVolume.setNormal(unitOuterNormal12, 1, 1);
1229 interactionVolume.setFaceArea(faceVol12, 1, 1);
1230 interactionVolume.setFacePosition(globalPosFace12, 1, 1);
1231
1232
1233 //**************** data for cell 4
1234 auto element4 = isIt14->outside();
1235
1236 //store pointer 4
1237 interactionVolume.setSubVolumeElement(element4, 3);
1238// interactionVolume.setIndexOnElement(globalIdx4, 3, 0);
1239
1240 interactionVolume.setNormal(unitOuterNormal14, 3, 0);
1241 interactionVolume.setFaceArea(faceVol41, 3, 0);
1242 interactionVolume.setFacePosition(globalPosFace41, 3, 0);
1243
1244 //**************** data for cell 3
1245 const auto& intersection23 = *isIt23;
1246 auto element3 = intersection23.outside();
1247 //store pointer 3
1248 interactionVolume.setSubVolumeElement(element3, 2);
1249
1250 GlobalPosition globalPosFace23 = intersection23.geometry().center();
1251 Scalar faceVol23 = intersection23.geometry().volume() / 2.0;
1252 // get outer normal vector scaled with half volume of face : for numbering of n see Aavatsmark, Eigestad
1253 GlobalPosition unitOuterNormal23 = intersection23.centerUnitOuterNormal();
1254
1255 interactionVolume.setNormal(unitOuterNormal23, 1, 0);
1256 unitOuterNormal23 *= -1;
1257 interactionVolume.setNormal(unitOuterNormal23, 2, 1);
1258 interactionVolume.setFaceArea(faceVol23, 1, 0);
1259 interactionVolume.setFaceArea(faceVol23, 2, 1);
1260 Scalar dummy=NAN;
1261 interactionVolume.setFaceArea(dummy, 2, 0);
1262 interactionVolume.setFaceArea(dummy, 3, 1);
1263 interactionVolume.setFacePosition(globalPosFace23, 1, 0);
1264 interactionVolume.setFacePosition(globalPosFace23, 2, 1);
1265
1266 Dune::FieldVector<Scalar, dim> unity(1.);
1267 std::vector<Dune::FieldVector<Scalar, dim> > lambda(4, unity);
1268
1269 Dune::FieldMatrix<Scalar,dim,2*dim-dim+1> T;
1270 int triangleType = transmissibilityCalculator_.calculateTransmissibility(
1271 T, interactionVolume, lambda,
1272 0, 1, 2, 3);
1273
1274 // 3.decide which triangle (which transmissibility coefficients) to use
1275 if (triangleType == TransmissibilityCalculator::rightTriangle)
1276 {
1277 additionalIntersectionIt = isIt23;
1278 // Translate flux from 2p mpfa-l local indexing to
1279 // 2p2c naming with i and j (reverse direction)
1280 // and adjust for the fact that it is right Triangle
1281 additionalT[0] = -T[1][2];
1282 additionalT[1] = -T[1][1];
1283 additionalT[2] = -T[1][0];
1284 }
1285 else if (triangleType == TransmissibilityCalculator::leftTriangle)
1286 {
1287 additionalIntersectionIt = isIt14;
1288 // Translate flux from 2p mpfa-l local indexing to
1289 // 2p2c naming with i and j (reverse direction)
1290 // and adjust for the fact that it is left Triangle
1291 additionalT[0] = -T[1][0];
1292 additionalT[1] = -T[1][1];
1293 additionalT[2] = -T[1][2];
1294
1295 }
1296 else
1297 DUNE_THROW(Dune::MathError, "No transmissivity for second half edge found!");
1298
1299 return triangleType;
1300}
1301
1302
1303}//end namespace Dumux
1304#endif
Define some often used mathematical functions.
Simplifies writing multi-file VTK datasets.
Provides methods for transmissibility calculation 2-d.
Defines the properties required for the adaptive sequential 2p2c models.
Finite volume 2p2c pressure model.
void harmonicMeanMatrix(Dune::FieldMatrix< Scalar, m, n > &K, const Dune::FieldMatrix< Scalar, m, n > &Ki, const Dune::FieldMatrix< Scalar, m, n > &Kj)
Calculate the harmonic mean of a fixed-size matrix.
Definition: math.hh:108
Definition: adapt.hh:29
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 permeability() noexcept
I/O name of permeability.
Definition: name.hh:143
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
Provides methods for transmissibility calculation in 2-d.
Definition: 2dtransmissibilitycalculator.hh:44
The finite volume model for the solution of the compositional pressure equation.
Definition: fv2dpressureadaptive.hh:78
const Problem & problem() const
Definition: fv2dpressureadaptive.hh:128
bool enableSecondHalfEdge
Definition: fv2dpressureadaptive.hh:209
int computeTransmissibilities(const IntersectionIterator &, IntersectionIterator &, TransmissivityMatrix &, TransmissivityMatrix &, GlobalPosition &, int &)
Computes the transmissibility coefficients for the MPFA-l method.
Definition: fv2dpressureadaptive.hh:805
bool enableMPFA
Enables mpfa method to calculate the fluxes near hanging nodes.
Definition: fv2dpressureadaptive.hh:208
Problem & problem()
Definition: fv2dpressureadaptive.hh:124
int transmissibilityAdapter_(const IntersectionIterator &, const IntersectionIterator &, IntersectionIterator &, GlobalPosition &, TransmissivityMatrix &)
An adapter to use the traditional 2p implementation for the second interaction region.
Definition: fv2dpressureadaptive.hh:1120
void initializeMatrix()
initializes the matrix to store the system of equations
Definition: fv2dpressureadaptive.hh:223
void adaptPressure()
Adapt primary variables vector after adapting the grid.
Definition: fv2dpressureadaptive.hh:150
FV2dPressure2P2CAdaptive(Problem &problem)
Constructs a FV2dPressure2P2CAdaptive object.
Definition: fv2dpressureadaptive.hh:166
void assemble(bool first)
function which assembles the system of equations to be solved
Definition: fv2dpressureadaptive.hh:364
DimMatrix R_
Matrix for vector rotation used in mpfa.
Definition: fv2dpressureadaptive.hh:206
TransmissibilityCalculator transmissibilityCalculator_
The 2p Mpfa pressure module, that is only used for the calulation of transmissibility of the second i...
Definition: fv2dpressureadaptive.hh:211
bool enableVolumeIntegral
Enables the volume integral of the pressure equation.
Definition: fv2dpressureadaptive.hh:207
void getMpfaFlux(const IntersectionIterator &, const CellData &)
Compute flux over an irregular interface using a mpfa method.
Definition: fv2dpressureadaptive.hh:507
The finite volume model for the solution of the compositional pressure equation.
Definition: fvpressure.hh:73
static constexpr int pressureType
gives kind of pressure used ( , , )
Definition: fvpressure.hh:190
Class including the information of an interaction volume of a MPFA L-method that does not change with...
Definition: linteractionvolume.hh:41
void setSubVolumeElement(const Element &element, int subVolumeIdx)
Store an element of the interaction volume.
Definition: linteractionvolume.hh:137
void setFacePosition(const DimVector &pos, int subVolumeIdx, int subVolumeFaceIdxInInside)
Store position of the center of a flux face.
Definition: linteractionvolume.hh:149
void setIndexOnElement(int indexInInside, int subVolumeIdx, int subVolumeFaceIdxInInside)
Store map from local interaction volume numbering to numbering of the Dune reference element.
Definition: linteractionvolume.hh:234
void setNormal(DimVector &normal, int subVolumeIdx, int subVolumeFaceIdxInInside)
Store a flux face normal.
Definition: linteractionvolume.hh:171
void setCenterPosition(DimVector &centerVertexPos)
Store position of the center vertex of the interaction volume.
Definition: linteractionvolume.hh:127
void setFaceArea(Scalar &faceArea, int subVolumeIdx, int subVolumeFaceIdxInInside)
Store a flux face area.
Definition: linteractionvolume.hh:160
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
PressureSolution & pressure()
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:120
@ rhs
index for the right hand side entry
Definition: sequential/cellcentered/pressure.hh:88
@ matrix
index for the global matrix entry
Definition: sequential/cellcentered/pressure.hh:89