3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fv3dpressureadaptive.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_FV3DPRESSURE2P2C_ADAPTIVE_HH
25#define DUMUX_FV3DPRESSURE2P2C_ADAPTIVE_HH
26
27// dune environent:
28#include <dune/istl/bvector.hh>
29#include <dune/istl/operators.hh>
30#include <dune/istl/solvers.hh>
31#include <dune/istl/preconditioners.hh>
32
33// dumux environment
35#include <dumux/common/math.hh>
38
39// include pressure model from Markus
43
44namespace Dumux {
45namespace Properties {
46template<class TypeTag>
47struct MPFAInteractionVolume<TypeTag, TTag::SequentialTwoPTwoCAdaptive> { using type = FvMpfaL3dInteractionVolumeAdaptive<TypeTag>; };
48template<class TypeTag>
49struct MPFAInteractionVolumeContainer<TypeTag, TTag::SequentialTwoPTwoCAdaptive> { using type = FvMpfaL3d2P2CInteractionVolumeContainerAdaptive<TypeTag>; };
50} // end namespace Properties
51
78template<class TypeTag> class FV3dPressure2P2CAdaptive
79: public FVPressure2P2CMultiPhysics<TypeTag>
80{
81 //the model implementation
85
90
92 using MaterialLaw = typename SpatialParams::MaterialLaw;
93
96
99
101 enum
102 {
103 dim = GridView::dimension, dimWorld = GridView::dimensionworld,
104 NumPhases = getPropValue<TypeTag, Properties::NumPhases>(), NumComponents = getPropValue<TypeTag, Properties::NumComponents>()
105 };
106 enum
107 {
108 pw = Indices::pressureW,
109 pn = Indices::pressureN,
110 pGlobal = Indices::pressureGlobal,
111 Sw = Indices::saturationW,
112 Sn = Indices::saturationN
113 };
114 enum
115 {
116 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
117 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
118 contiWEqIdx = Indices::contiWEqIdx, contiNEqIdx = Indices::contiNEqIdx
119 };
120 enum
121 {
122 rhs = BaseType::rhs, matrix = BaseType::matrix,
123 };
124
125 // using declarations to abbreviate several dune classes...
126 using Vertex = typename GridView::Traits::template Codim<dim>::Entity;
127 using Element = typename GridView::Traits::template Codim<0>::Entity;
128 using ReferenceElementContainer = Dune::ReferenceElements<Scalar, dim>;
129
130 using Grid = typename GridView::Grid;
131 using Intersection = typename GridView::Intersection;
132 using IntersectionIterator = typename GridView::IntersectionIterator;
133
134 // convenience shortcuts for Vectors/Matrices
135 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
136 using TransmissivityMatrix = Dune::FieldVector<Scalar,dim+1>;
137 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
138 using PhaseVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumPhases>()>;
139 using ComponentVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumComponents>()>;
141
142 // the typenames used for the stiffness matrix and solution vector
145
146 // Dumux MPFA types
148 using InteractionVolume = typename InteractionVolumeContainer::InteractionVolume;
149
150protected:
152 Problem& problem()
153 {
154 return this->problem_;
155 }
156 const Problem& problem() const
157 {
158 return this->problem_;
159 }
161
162public:
163 void update()
164 {
168 {
171 new InteractionVolumeContainer(problem());
172
174 }
175 asImp_().initializeMatrix();
177 }
178 //initializes the matrix to store the system of equations
179 void initializeMatrix();
182
183 void initialize(bool solveTwice = false){
185 {
188 new InteractionVolumeContainer(problem());
189
191 }
192 ParentType::initialize(solveTwice);
193 }
194
195 //function which assembles the system of equations to be solved
196 void assemble(bool first);
197
198 void getMpfaFlux(const IntersectionIterator&, const CellData&);
199
200 void get1pMpfaFlux(const IntersectionIterator&, const CellData&);
201
202 //constitutive functions are initialized and stored in the variables object
203 void updateMaterialLaws(bool fromPostTimestep = false);
204
205 // mpfa transmissibilities
206 int computeTransmissibilities(const IntersectionIterator&,
207 TransmissivityMatrix&,
208 GlobalPosition&,
209 int&,
210 GlobalPosition&,
211 int&);
212
215 {
216 int gridSize = problem().gridView().size(0);
217 this->pressure().resize(gridSize);
218
219 for (const auto& element : elements(problem().gridView()))
220 {
221 // get the global index of the cell
222 int eIdxGlobalI = problem().variables().index(element);
223
224 // assemble interior element contributions
225 if (element.partitionType() == Dune::InteriorEntity)
226 {
227 this->pressure()[eIdxGlobalI]
228 = problem().variables().cellData(eIdxGlobalI).pressure(this->pressureType);
229 }
230 }
231#if HAVE_MPI
232 // communicate updated values
234 using ElementMapper = typename SolutionTypes::ElementMapper;
237
238 DataHandle dataHandle(problem().variables().elementMapper(), this->pressure());
239 problem().gridView().template communicate<DataHandle>(dataHandle,
240 Dune::InteriorBorder_All_Interface,
241 Dune::ForwardCommunication);
242#endif
243 }
244
249 FV3dPressure2P2CAdaptive(Problem& problem) : FVPressure2P2CMultiPhysics<TypeTag>(problem),
250 enableMPFA(false),
252 {
254 enableMPFA = getParam<bool>("GridAdapt.EnableMultiPointFluxApproximation");
255
256 maxInteractionVolumes = getParam<int>("GridAdapt.MaxInteractionVolumes");
257 }
258
259private:
261 Implementation &asImp_()
262 { return *static_cast<Implementation *>(this);}
263
265 const Implementation &asImp_() const
266 { return *static_cast<const Implementation *>(this);}
267
268 int searchCommonVertex_(const Intersection& is, Vertex& vertex)
269 {
270 /******* get corner of interest ************/
271 // search through corners of large cell with isIt
272 int localIdxLarge = 0;
273 for(localIdxLarge = 0; localIdxLarge < is.inside().subEntities(dim); ++localIdxLarge)
274 {
275 auto vLarge = is.inside().template subEntity<dim>(localIdxLarge);
276
277 // search through corners of small cell with isIt
278 for(int verticeSmall = 0; verticeSmall<is.outside().subEntities(dim); ++verticeSmall)
279 {
280 auto vSmall = is.outside().template subEntity<dim>(verticeSmall);
281
282 if(problem().variables().index(vSmall) == problem().variables().index(vLarge) )
283 {
284 vertex = vSmall;
285 return localIdxLarge;
286 }
287 }
288 }
289 return -1;
290 }
291
292protected:
293 int transmissibilityAdapter_(const IntersectionIterator& isIt,
294 InteractionVolume& interactionVolume,
295 const int& subVolumeFaceIdx,
296 bool properFluxDirection,
297 Element& additional2,
298 Element& additional3,
299 TransmissivityMatrix& additionalT);
300
301 std::map<int, std::vector<int> > irregularCellMap_;
305
307 InteractionVolumeContainer* interactionVolumesContainer_;
310};
311
313template<class TypeTag>
315{
316 int gridSize_ = problem().gridView().size(0);
317 // update RHS vector, matrix
318 this->A_.setSize (gridSize_,gridSize_); //
319 this->f_.resize(gridSize_);
320 irregularCellMap_.clear();
321
322 this->initializeMatrixRowSize();
323 this->A_.endrowsizes();
324 this->initializeMatrixIndices();
325 this->A_.endindices();
326}
327
329template<class TypeTag>
331{
332 // determine matrix row sizes
333 for (const auto& element : elements(problem().gridView()))
334 {
335 // cell index
336 int eIdxGlobalI = problem().variables().index(element);
337 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
338
339 // initialize row size
340 int rowSize = 1;
341
342 // set perimeter to zero
343 cellDataI.perimeter() = 0;
344
345 // prepare storage for all found 3rd cells
346 std::vector<int> foundAdditionals;
347
348 int numberOfIntersections = 0;
349 // run through all intersections with neighbors
350 for (const auto& intersection : intersections(problem().gridView(), element))
351 {
352 cellDataI.perimeter() += intersection.geometry().volume();
353 numberOfIntersections++;
354 if (intersection.neighbor())
355 {
356 rowSize++;
357
358 // special treatment for hanging nodes in the mpfa case
359 if (enableMPFA && (element.level() < intersection.outside().level()))
360 {
361 // each cell might have 4 hanging nodes with 4 irregular neighbors each
362 // get global ID of Interface from larger cell
363 int intersectionID = problem().grid().localIdSet().subId(element,
364 intersection.indexInInside(), 1);
365 //index outside
366 int eIdxGlobalJ = problem().variables().index(intersection.outside());
367
368 // add Entry of current neighbor cell to the IS seen from large cell
369 irregularCellMap_[intersectionID].push_back(eIdxGlobalJ);
370 }
371 }
372 }
373 cellDataI.fluxData().resize(numberOfIntersections);
374 this->A_.incrementrowsize(eIdxGlobalI, rowSize);
375 } //end first loop (that already reserved enough space for the MPFA connections on hanging nodes)
376
377 // second loop to determine which matrix entries will be occupied
378 if(enableMPFA && maxInteractionVolumes>1)
379 {
380 //prepare map for additional cell-connection through mpfa
381 std::multimap<int, int> addionalRelations;
382 using IntPair = std::pair<int,int>;
383 std::pair<std::multimap<int,int>::iterator,std::multimap<int,int>::iterator> range;
384 std::multimap<int,int>::iterator rangeIt;
385
386 // second loop for further sub-faces
387 for (const auto& element : elements(problem().gridView()))
388 {
389 // cell index
390 int eIdxGlobalI = problem().variables().index(element);
391 // run through all intersections with neighbors
392 const auto isEndIt = problem().gridView().iend(element);
393 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
394 {
395 const auto& intersection = *isIt;
396
397 if (intersection.neighbor())
398 {
399 //index outside
400 int eIdxGlobalJ = problem().variables().index(intersection.outside());
401
402 // if mpfa is used, more entries might be needed if all interactionRegions are regarded
403 if (intersection.outside().level() > element.level()) //look from larger cell
404 {
405 // Prepare MPFA
406 /* get geometric Info, transmissibility matrix */
407 GlobalPosition globalPos3(0.);
408 int eIdxGlobal3=-1;
409 GlobalPosition globalPos4(0.);
410 int eIdxGlobal4=-1;
411 TransmissivityMatrix T(0.);
412
413 int interactionRegions
414 = problem().variables().getMpfaData3D(intersection, T,
415 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
416 if (interactionRegions == 0)
417 interactionRegions = problem().pressureModel().computeTransmissibilities(isIt,T,
418 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
419 if(!interactionRegions)
420 Dune::dgrave << "something went wrong getting mpfa data on cell " << eIdxGlobalI << std::endl;
421 if (interactionRegions == 1) // no second subface
422 continue;
423
424 //loop over all found interaction regions
425 for (int cocumber=1; cocumber<interactionRegions; cocumber++ )
426 {
427 problem().variables().getMpfaData3D(intersection, T,
428 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4, cocumber);
429 // indices
430 int additionalIdx2 = eIdxGlobal3;
431 int additionalIdx3 = eIdxGlobal4;
432
433 bool addIndex = true;
434
435 // check if additionals are "normal" neighbors of I
436 bool additional2isNeighbor(false), additional3isNeighbor(false);
437 // run through all intersections with neighbors if eIt
438 for (const auto& checkIntersection
439 : intersections(problem().gridView(), element))
440 {
441 if (checkIntersection.neighbor())
442 {
443 if(additionalIdx2==problem().variables().index(checkIntersection.outside()))
444 additional2isNeighbor = true;
445 if(additionalIdx3 == problem().variables().index(checkIntersection.outside()))
446 additional3isNeighbor = true;
447 }
448
449 }
450 // if not "normal" neighbor, increase row size
451 if(!additional2isNeighbor)
452 {
453 // check if relation not already added
454 IntPair intPair(eIdxGlobalI,additionalIdx2);
455 if(eIdxGlobalI > additionalIdx2)
456 {
457 using std::swap;
458 swap(intPair.first, intPair.second);
459 }
460 range = addionalRelations.equal_range(intPair.first);
461 for (rangeIt=range.first; range.first!=range.second
462 && rangeIt!=range.second; ++rangeIt)
463 if((*rangeIt).second == intPair.second)
464 addIndex = false;
465 if(addIndex)
466 {
467 this->A_.incrementrowsize(eIdxGlobalI);
468 // add space for additional itsself
469 this->A_.incrementrowsize(additionalIdx2);
470 // mark relation as added
471 addionalRelations.insert(intPair);
472 }
473 }
474
475 // if not "normal" neighbor, increase row size
476 if(!additional3isNeighbor)
477 {
478 addIndex = true;
479 // check if relation not already added
480 IntPair intPair(eIdxGlobalI,additionalIdx3);
481 if(eIdxGlobalI > additionalIdx3)
482 {
483 using std::swap;
484 swap(intPair.first, intPair.second);
485 }
486 range = addionalRelations.equal_range(intPair.first);
487 for (rangeIt=range.first; range.first!=range.second
488 && rangeIt!=range.second; ++rangeIt)
489 if((*rangeIt).second == intPair.second)
490 addIndex = false;
491 if(addIndex)
492 {
493 this->A_.incrementrowsize(eIdxGlobalI);
494 // add space for additional itsself
495 this->A_.incrementrowsize(additionalIdx3);
496 // mark relation as added
497 addionalRelations.insert(intPair);
498 }
499 }
500
501 //reset bools for J
502 additional2isNeighbor = additional3isNeighbor = false;
503 // run through all intersections with neighbors of J
504 for (const auto& checkIntersection
505 : intersections(problem().gridView(), intersection.outside()))
506 {
507 if (checkIntersection.neighbor())
508 {
509 if(additionalIdx2 == problem().variables().index(checkIntersection.outside()))
510 additional2isNeighbor = true;
511 if(additionalIdx3 == problem().variables().index(checkIntersection.outside()))
512 additional3isNeighbor = true;
513 }
514 }
515
516 // if not "normal" neighbor, increase row size
517 if(!additional2isNeighbor)
518 {
519 addIndex = true;
520 // check if relation not already added
521 IntPair intPair(eIdxGlobalJ,additionalIdx2);
522 if(eIdxGlobalJ > additionalIdx2)
523 {
524 using std::swap;
525 swap(intPair.first, intPair.second);
526 }
527 range = addionalRelations.equal_range(intPair.first);
528 for (rangeIt=range.first; range.first!=range.second
529 && rangeIt!=range.second; ++rangeIt)
530 if((*rangeIt).second == intPair.second)
531 addIndex = false;
532 if(addIndex)
533 {
534 this->A_.incrementrowsize(eIdxGlobalJ);
535 // add space for additional itsself
536 this->A_.incrementrowsize(additionalIdx2);
537 // mark relation as added
538 addionalRelations.insert(intPair);
539 }
540 }
541
542 // if not "normal" neighbor, increase row size
543 if(!additional3isNeighbor)
544 {
545 addIndex = true;
546 // check if relation not already added
547 IntPair intPair(eIdxGlobalJ,additionalIdx3);
548 if(eIdxGlobalJ > additionalIdx3)
549 {
550 using std::swap;
551 swap(intPair.first, intPair.second);
552 }
553 range = addionalRelations.equal_range(intPair.first);
554 for (rangeIt=range.first; range.first!=range.second
555 && rangeIt!=range.second; ++rangeIt)
556 if((*rangeIt).second == intPair.second)
557 addIndex = false;
558 if(addIndex)
559 {
560 this->A_.incrementrowsize(eIdxGlobalJ);
561 // add space for additional itsself
562 this->A_.incrementrowsize(additionalIdx3);
563 // mark relation as added
564 addionalRelations.insert(intPair);
565 }
566 }
567 }
568 }
569 }
570 }
571 } //end second loop
572 }
573}
574
576/* Method adds TPFA and MPFA matrix entries
577 */
578template<class TypeTag>
580{
581 // determine position of matrix entries
582 for (const auto& element : elements(problem().gridView()))
583 {
584 // cell index
585 int eIdxGlobalI = problem().variables().index(element);
586
587 // add diagonal index
588 this->A_.addindex(eIdxGlobalI, eIdxGlobalI);
589
590 // run through all intersections with neighbors
591 const auto isEndIt = problem().gridView().iend(element);
592 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
593 {
594 const auto& intersection = *isIt;
595
596 if (intersection.neighbor())
597 {
598 // access neighbor
599 int eIdxGlobalJ = problem().variables().index(intersection.outside());
600
601 // add off diagonal index
602 this->A_.addindex(eIdxGlobalI, eIdxGlobalJ);
603
604 // special treatment for hanging nodes in the mpfa case
605 if (enableMPFA && (element.level() < intersection.outside().level()))
606 {
607 // prepare stuff to enter transmissibility calculation
608 GlobalPosition globalPos3(0.);
609 int eIdxGlobal3=-1;
610 GlobalPosition globalPos4(0.);
611 int eIdxGlobal4=-1;
612 TransmissivityMatrix T(0.);
613 TransmissivityMatrix additionalT(0.);
614
615 int interactionRegions
616 = problem().variables().getMpfaData3D(intersection, T, globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
617 if (interactionRegions == 0)
618 interactionRegions = problem().pressureModel().computeTransmissibilities(isIt,T,
619 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
620
621 for (int cocumber=1; cocumber<interactionRegions; cocumber++ )
622 {
623 problem().variables().getMpfaData3D(intersection, T,
624 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4, cocumber);
625
626 // add off diagonal index in both directions!!
627 this->A_.addindex(eIdxGlobalI, eIdxGlobal3);
628 this->A_.addindex(eIdxGlobal3, eIdxGlobalI);
629 this->A_.addindex(eIdxGlobalI, eIdxGlobal4);
630 this->A_.addindex(eIdxGlobal4, eIdxGlobalI);
631 this->A_.addindex(eIdxGlobalJ, eIdxGlobal3);
632 this->A_.addindex(eIdxGlobal3, eIdxGlobalJ);
633 this->A_.addindex(eIdxGlobalJ, eIdxGlobal4);
634 this->A_.addindex(eIdxGlobal4, eIdxGlobalJ);
635 }
636 }
637 }
638 }
639 }
640}
641
643template<class TypeTag>
645{
646 if(first)
647 {
648 BaseType::assemble(true);
649 return;
650 }
651
652 // initialization: set matrix A_ to zero
653 this->A_ = 0;
654 this->f_ = 0;
655
656 for (const auto& element : elements(problem().gridView()))
657 {
658 // get the global index of the cell
659 int eIdxGlobalI = problem().variables().index(element);
660
661 // assemble interior element contributions
662 if (element.partitionType() == Dune::InteriorEntity)
663 {
664 // get the cell data
665 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
666
667 Dune::FieldVector<Scalar, 2> entries(0.);
668
669 /***** source term ***********/
670#ifndef noMultiphysics
671 if(cellDataI.subdomain() != 2)
672 problem().pressureModel().get1pSource(entries,element, cellDataI);
673 else
674#endif
675 problem().pressureModel().getSource(entries,element, cellDataI, first);
676
677 this->f_[eIdxGlobalI] += entries[rhs];
678
679 /***** flux term ***********/
680 // iterate over all faces of the cell
681 const auto isEndIt = problem().gridView().iend(element);
682 for (auto isIt = problem().gridView().ibegin(element); isIt != isEndIt; ++isIt)
683 {
684 const auto& intersection = *isIt;
685
686 /************* handle interior face *****************/
687 if (intersection.neighbor())
688 {
689 auto neighbor = intersection.outside();
690 int eIdxGlobalJ = problem().variables().index(neighbor);
691 //check for hanging nodes
692 //take a hanging node never from the element with smaller level!
693 bool haveSameLevel = (element.level() == neighbor.level());
694 // calculate only from one side, but add matrix entries for both sides
695 // the last condition is needed to properly assemble in the presence
696 // of ghost elements
697 if (getPropValue<TypeTag, Properties::VisitFacesOnlyOnce>()
698 && (eIdxGlobalI > eIdxGlobalJ) && haveSameLevel
699 && neighbor.partitionType() == Dune::InteriorEntity)
700 continue;
701
702 entries = 0;
703 //check for hanging nodes
704 if(!haveSameLevel && enableMPFA)
705 {
706 if (cellDataI.subdomain() != 2
707 || problem().variables().cellData(eIdxGlobalJ).subdomain() != 2) // cell in the 1p domain
708 {
709 asImp_().get1pMpfaFlux(isIt, cellDataI);
710 }
711 else
712 {
713 asImp_().getMpfaFlux(isIt, cellDataI);
714 }
715 }
716 else
717 {
718 CellData cellDataJ = problem().variables().cellData(eIdxGlobalJ);
719 if (cellDataI.subdomain() != 2
720 || problem().variables().cellData(eIdxGlobalJ).subdomain() != 2) // cell in the 1p domain
721 {
722 asImp_().get1pFlux(entries, intersection, cellDataI);
723 }
724 else
725 {
726 asImp_().getFlux(entries, intersection, cellDataI, first);
727 }
728
729 //set right hand side
730 this->f_[eIdxGlobalI] -= entries[rhs];
731
732 // set diagonal entry
733 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
734
735 // set off-diagonal entry
736 this->A_[eIdxGlobalI][eIdxGlobalJ] -= entries[matrix];
737
738 // The second condition is needed to not spoil the ghost element entries
739 if (getPropValue<TypeTag, Properties::VisitFacesOnlyOnce>()
740 && neighbor.partitionType() == Dune::InteriorEntity)
741 {
742 this->f_[eIdxGlobalJ] += entries[rhs];
743 this->A_[eIdxGlobalJ][eIdxGlobalJ] += entries[matrix];
744 this->A_[eIdxGlobalJ][eIdxGlobalI] -= entries[matrix];
745 }
746 }
747 } // end neighbor
748
749 /************* boundary face ************************/
750 else
751 {
752 entries = 0;
753 if (cellDataI.subdomain() != 2) //the current cell in the 1p domain
754 asImp_().get1pFluxOnBoundary(entries, intersection, cellDataI);
755 else
756 asImp_().getFluxOnBoundary(entries, intersection, cellDataI, first);
757
758 //set right hand side
759 this->f_[eIdxGlobalI] += entries[rhs];
760 // set diagonal entry
761 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
762 }
763 } //end interfaces loop
764
765 /***** storage term ***********/
766 if (cellDataI.subdomain() != 2) //the current cell in the 1p domain
767 asImp_().get1pStorage(entries, element, cellDataI);
768 else
769 asImp_().getStorage(entries, element, cellDataI, first);
770
771 this->f_[eIdxGlobalI] += entries[rhs];
772 // set diagonal entry
773 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
774 }
775 // assemble overlap and ghost element contributions
776 else
777 {
778 this->A_[eIdxGlobalI] = 0.0;
779 this->A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
780 this->f_[eIdxGlobalI] = this->pressure()[eIdxGlobalI];
781 }
782 } // end grid traversal
783}
784
786
804template<class TypeTag>
805void FV3dPressure2P2CAdaptive<TypeTag>::getMpfaFlux(const IntersectionIterator& isIt,
806 const CellData& cellDataI)
807{
808 const auto& intersection = *isIt;
809
810 // acess Cell I
811 auto elementI = intersection.inside();
812 int eIdxGlobalI = problem().variables().index(elementI);
813
814 // get global coordinate of cell center
815 const GlobalPosition& globalPos = elementI.geometry().center();
816
817 // cell volume & perimeter, assume linear map here
818 Scalar volume = elementI.geometry().volume();
819 Scalar perimeter = cellDataI.perimeter();
820
821 // get absolute permeability
822 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
823
824 // access neighbor
825 auto neighbor = intersection.outside();
826 int eIdxGlobalJ = problem().variables().index(neighbor);
827 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
828
829 // gemotry info of neighbor
830 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
831
832 // distance vector between barycenters
833 GlobalPosition distVec = globalPosNeighbor - globalPos;
834
835 // compute distance between cell centers
836 Scalar dist = distVec.two_norm();
837
838 GlobalPosition unitDistVec(distVec);
839 unitDistVec /= dist;
840
841 DimMatrix permeabilityJ
842 = problem().spatialParams().intrinsicPermeability(neighbor);
843
844 // compute vectorized permeabilities
845 DimMatrix meanPermeability(0);
846 harmonicMeanMatrix(meanPermeability, permeabilityI, permeabilityJ);
847
848 Dune::FieldVector<Scalar, dim> permeability(0);
849 meanPermeability.mv(unitDistVec, permeability);
850
851 // get average density for gravity flux
852 PhaseVector rhoMean(0.);
853 for (int phaseIdx=0; phaseIdx<NumPhases; phaseIdx++)
854 rhoMean[phaseIdx] =0.5 * (cellDataI.density(phaseIdx)+cellDataJ.density(phaseIdx));
855
856 // reset potential gradients
857 Dune::FieldVector<Scalar,2> potential(0.);
858
859 // determine volume derivatives in neighbor
860 if (!cellDataJ.hasVolumeDerivatives())
861 asImp_().volumeDerivatives(globalPosNeighbor, neighbor);
862
863 ComponentVector dv_dC(0.), graddv_dC(0.);
864 for (int compIdx = 0; compIdx < NumComponents; ++compIdx)
865 {
866 dv_dC[compIdx]= (cellDataJ.dv(compIdx) // dV/dm1= dv/dC^1
867 + cellDataI.dv(compIdx)) * 0.5;
868 graddv_dC[compIdx] = (cellDataJ.dv(compIdx)
869 - cellDataI.dv(compIdx)) / dist;
870 }
871// potentialW = problem().variables().potentialWetting(eIdxGlobalI, isIndex);
872// potentialNW = problem().variables().potentialNonwetting(eIdxGlobalI, isIndex);
873//
874// densityW = (potentialW > 0.) ? densityWI : densityWJ;
875// densityNW = (potentialNW > 0.) ? densityNWI : densityNWJ;
876//
877// densityW = (potentialW == 0.) ? rhoMeanW : densityW;
878// densityNW = (potentialNW == 0.) ? rhoMeanNW : densityNW;
879 //jochen: central weighting for gravity term
880
881 // Prepare MPFA
882 /* get geometrical Info, transmissibility matrix */
883 GlobalPosition globalPos3(0.);
884 int eIdxGlobal3=-1;
885 GlobalPosition globalPos4(0.);
886 int eIdxGlobal4=-1;
887 TransmissivityMatrix T(0.);
888
889 // prepare second interaction region
890 GlobalPosition globalPosAdditional3(0.);
891 int eIdxGlobalAdditional3=-1;
892 GlobalPosition globalPosAdditional4(0.);
893 int eIdxGlobalAdditional4=-1;
894
895 TransmissivityMatrix additionalT(0.);
896
897 int interactionRegions
898 = problem().variables().getMpfaData3D(intersection, T,
899 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
900 if (interactionRegions == 0)
901 interactionRegions = problem().pressureModel().computeTransmissibilities(isIt,T,
902 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
903 if(!interactionRegions)
904 Dune::dgrave << "something went wrong getting mpfa data on cell " << eIdxGlobalI << std::endl;
905
906 // shortcurts mpfa case
907 CellData& cellData3 = problem().variables().cellData(eIdxGlobal3);
908 CellData& cellData4 = problem().variables().cellData(eIdxGlobal4);
909 Scalar temp1 = globalPos * this->gravity_;
910 Scalar temp2 = globalPosNeighbor * this->gravity_;
911 Scalar temp3 = globalPos3 * this->gravity_;
912 Scalar temp4 = globalPos4 * this->gravity_;
913
914 for (int phaseIdx = 0; phaseIdx < NumPhases; ++phaseIdx)
915 {
916 potential[phaseIdx] = (cellDataI.pressure(phaseIdx)-temp1*rhoMean[phaseIdx]) * T[0]
917 + (cellDataJ.pressure(phaseIdx)-temp2*rhoMean[phaseIdx]) * T[1]
918 + (cellData3.pressure(phaseIdx)-temp3*rhoMean[phaseIdx]) * T[2]
919 + (cellData4.pressure(phaseIdx)-temp4*rhoMean[phaseIdx]) * T[3];
920 }
921
922 // regard more interaction regions, if there are more
923 if(interactionRegions != 1)
924 {
925 for(int banana = 1; banana < interactionRegions; banana ++)
926 {
927 // get data for second interaction region
928 problem().variables().getMpfaData3D(intersection, additionalT,
929 globalPosAdditional3, eIdxGlobalAdditional3,
930 globalPosAdditional4, eIdxGlobalAdditional4 ,
931 banana); // offset for second interaction region
932
933 Scalar gravityContributionAdditonal
934 = temp1 * additionalT[0] + temp2 * additionalT[1]
935 + globalPosAdditional3*this->gravity_ * additionalT[2]
936 + globalPosAdditional4*this->gravity_ * additionalT[3];
937 CellData& cellDataA3 = problem().variables().cellData(eIdxGlobalAdditional3);
938 CellData& cellDataA4 = problem().variables().cellData(eIdxGlobalAdditional4);
939
940 for (int phaseIdx = 0; phaseIdx < NumPhases; ++phaseIdx)
941 {
942 potential[phaseIdx] += cellDataI.pressure(phaseIdx) * additionalT[0]
943 + cellDataJ.pressure(phaseIdx) * additionalT[1]
944 + cellDataA3.pressure(phaseIdx) * additionalT[2]
945 + cellDataA4.pressure(phaseIdx) * additionalT[3];
946 potential[phaseIdx] -= gravityContributionAdditonal * rhoMean[phaseIdx];
947 }
948 }
949 }
950
951 // initialize convenience shortcuts
952 PhaseVector lambda(0.), dV(0.), gV(0.);
953
954 //do the upwinding of the mobility depending on the phase potentials
955 std::vector<const CellData*> upwindCellData(NumPhases);
956
957 for (int phaseIdx = 0; phaseIdx < NumPhases; ++phaseIdx)
958 {
959 int eqIdx = phaseIdx + 1;
960 if (potential[phaseIdx] > 0.)
961 upwindCellData[phaseIdx] = &cellDataI;
962 else if (potential[phaseIdx] < 0.)
963 upwindCellData[phaseIdx] = &cellDataJ;
964 else
965 {
966 if(cellDataI.isUpwindCell(intersection.indexInInside(), eqIdx))
967 upwindCellData[phaseIdx] = &cellDataI;
968 else if(cellDataJ.isUpwindCell(intersection.indexInOutside(), eqIdx))
969 upwindCellData[phaseIdx] = &cellDataJ;
970 //else
971 // upwinding is not done!
972 }
973
974 //perform upwinding if desired
975 if(!upwindCellData[phaseIdx])
976 {
977 Scalar averagedMassFraction[NumComponents];
978 for (int compIdx = 0; compIdx < NumComponents; ++compIdx)
979 averagedMassFraction[compIdx]
980 = harmonicMean(cellDataI.massFraction(phaseIdx, compIdx), cellDataJ.massFraction(phaseIdx, compIdx));
981 Scalar averageDensity = harmonicMean(cellDataI.density(phaseIdx), cellDataJ.density(phaseIdx));
982
983 for (int compIdx = 0; compIdx < NumComponents; ++compIdx)
984 {
985 dV[phaseIdx] += dv_dC[compIdx] * averagedMassFraction[compIdx];
986 gV[phaseIdx] += graddv_dC[compIdx] * averagedMassFraction[compIdx];
987 }
988 dV[phaseIdx] *= averageDensity;
989 gV[phaseIdx] *= averageDensity;
990 lambda[phaseIdx] = harmonicMean(cellDataI.mobility(phaseIdx), cellDataJ.mobility(phaseIdx));
991 }
992 else //perform upwinding
993 {
994 for (int compIdx = 0; compIdx < NumComponents; ++compIdx)
995 {
996 dV[phaseIdx] += dv_dC[compIdx] * upwindCellData[phaseIdx]->massFraction(phaseIdx, compIdx);
997 gV[phaseIdx] += graddv_dC[compIdx] * upwindCellData[phaseIdx]->massFraction(phaseIdx, compIdx);
998 }
999 lambda[phaseIdx] = upwindCellData[phaseIdx]->mobility(phaseIdx);
1000 dV[phaseIdx] *= upwindCellData[phaseIdx]->density(phaseIdx);
1001 gV[phaseIdx] *= upwindCellData[phaseIdx]->density(phaseIdx);
1002 }
1003 }
1004
1005 /* compute matrix entry: advective fluxes */
1006 /* extend T with other matrix entries and assemble to A_ */
1007 this->A_[eIdxGlobalI][eIdxGlobalI] += (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * T[0];
1008 this->A_[eIdxGlobalI][eIdxGlobalJ] += (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * T[1];
1009 this->A_[eIdxGlobalI][eIdxGlobal3] += (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * T[2];
1010 this->A_[eIdxGlobalI][eIdxGlobal4] += (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * T[3];
1011
1012
1013 // add gravity to RHS vector
1014 Scalar gravityContribution = temp1 * T[0] + temp2 * T[1] + temp3 * T[2] + temp4 * T[3];
1015 this->f_[eIdxGlobalI] += (rhoMean[wPhaseIdx] * lambda[wPhaseIdx] * dV[wPhaseIdx]
1016 + rhoMean[nPhaseIdx] * lambda[nPhaseIdx] * dV[nPhaseIdx]) * gravityContribution;
1017
1018 // weithing accounts for the fraction of the subcontrol volume
1019 Scalar weightingFactor = volume / perimeter; // transforms flux through area A -> V * A/perimeter
1020 if(enableVolumeIntegral_) // switch off volume integral for mpfa case
1021 {
1022 // correct for area integral
1023 this->A_[eIdxGlobalI][eIdxGlobalI] -=
1024 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * T[0];
1025 this->A_[eIdxGlobalI][eIdxGlobalJ] -=
1026 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * T[1];
1027 this->A_[eIdxGlobalI][eIdxGlobal3] -=
1028 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * T[2];
1029 this->A_[eIdxGlobalI][eIdxGlobal4] -=
1030 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * T[3];
1031
1032 // add gravity to RHS vector
1033 this->f_[eIdxGlobalI] -= weightingFactor * gravityContribution *
1034 (rhoMean[wPhaseIdx] * lambda[wPhaseIdx] * gV[wPhaseIdx] + rhoMean[nPhaseIdx] * lambda[nPhaseIdx] * gV[nPhaseIdx]);
1035 }
1036
1037 // capillary pressure flux
1038 Scalar pcGradient = cellDataI.capillaryPressure() * T[0]
1039 + cellDataJ.capillaryPressure() * T[1]
1040 + cellData3.capillaryPressure() * T[2]
1041 + cellData4.capillaryPressure() * T[3];
1042
1043 if (this->pressureType == pw)
1044 pcGradient *= + lambda[nPhaseIdx] * dV[nPhaseIdx]
1045 - enableVolumeIntegral_ * weightingFactor * lambda[nPhaseIdx] * gV[nPhaseIdx];
1046 else if (this->pressureType == pn)
1047 pcGradient *= - lambda[wPhaseIdx] * dV[wPhaseIdx]
1048 + enableVolumeIntegral_ * weightingFactor * lambda[wPhaseIdx] * gV[wPhaseIdx];
1049
1050 this->f_[eIdxGlobalI] += pcGradient;
1051
1052 // regard more interaction regions, if there are more
1053 if(interactionRegions != 1)
1054 {
1055 for(int banana = 1; banana < interactionRegions; banana ++)
1056 {
1057 // get data for second interaction region
1058 problem().variables().getMpfaData3D(intersection, additionalT,
1059 globalPosAdditional3, eIdxGlobalAdditional3,
1060 globalPosAdditional4, eIdxGlobalAdditional4 ,
1061 banana); // offset for second interaction region
1062
1063 Scalar gravityContributionAdditonal
1064 = temp1 * additionalT[0] + temp2 * additionalT[1]
1065 + globalPosAdditional3*this->gravity_ * additionalT[2]
1066 + globalPosAdditional4*this->gravity_ * additionalT[3];
1067 CellData& cellDataA3 = problem().variables().cellData(eIdxGlobalAdditional3);
1068 CellData& cellDataA4 = problem().variables().cellData(eIdxGlobalAdditional4);
1069
1070 /* compute matrix entry: advective fluxes */
1071 /* extend T with other matrix entries and assemble to A_ */
1072 this->A_[eIdxGlobalI][eIdxGlobalI] +=
1073 (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * additionalT[0];
1074 this->A_[eIdxGlobalI][eIdxGlobalJ] +=
1075 (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * additionalT[1];
1076 this->A_[eIdxGlobalI][eIdxGlobalAdditional3] +=
1077 (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * additionalT[2];
1078 this->A_[eIdxGlobalI][eIdxGlobalAdditional4] +=
1079 (lambda[wPhaseIdx] * dV[wPhaseIdx] + lambda[nPhaseIdx] * dV[nPhaseIdx]) * additionalT[3];
1080
1081
1082 // add gravity to RHS vector
1083 this->f_[eIdxGlobalI] += (rhoMean[wPhaseIdx] * lambda[wPhaseIdx] * dV[wPhaseIdx]
1084 + rhoMean[nPhaseIdx] * lambda[nPhaseIdx] * dV[nPhaseIdx]) * gravityContributionAdditonal;
1085
1086 // weithing accounts for the fraction of the subcontrol volume
1087 if(enableVolumeIntegral_) // switch off volume integral for mpfa case
1088 {
1089 // correct for area integral
1090 this->A_[eIdxGlobalI][eIdxGlobalI] -=
1091 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * additionalT[0];
1092 this->A_[eIdxGlobalI][eIdxGlobalJ] -=
1093 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * additionalT[1];
1094 this->A_[eIdxGlobalI][eIdxGlobalAdditional3] -=
1095 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * additionalT[2];
1096 this->A_[eIdxGlobalI][eIdxGlobalAdditional4] -=
1097 weightingFactor * (lambda[wPhaseIdx] * gV[wPhaseIdx] + lambda[nPhaseIdx] * gV[nPhaseIdx]) * additionalT[3];
1098
1099 // add gravity to RHS vector
1100 this->f_[eIdxGlobalI] -= weightingFactor * gravityContribution *
1101 (rhoMean[wPhaseIdx] * lambda[wPhaseIdx] * gV[wPhaseIdx] + rhoMean[nPhaseIdx] * lambda[nPhaseIdx] * gV[nPhaseIdx]);
1102 }
1103
1104 // capillary pressure flux
1105 Scalar addPcGradient = cellDataI.capillaryPressure() * additionalT[0]
1106 + cellDataJ.capillaryPressure() * additionalT[1]
1107 + cellDataA3.capillaryPressure() * additionalT[2]
1108 + cellDataA4.capillaryPressure() * additionalT[3];
1109
1110 if (this->pressureType == pw)
1111 addPcGradient *= + lambda[nPhaseIdx] * dV[nPhaseIdx]
1112 - enableVolumeIntegral_ * weightingFactor * lambda[nPhaseIdx] * gV[nPhaseIdx];
1113 else if (this->pressureType == pn)
1114 addPcGradient *= - lambda[wPhaseIdx] * dV[wPhaseIdx]
1115 + enableVolumeIntegral_ * weightingFactor * lambda[wPhaseIdx] * gV[wPhaseIdx];
1116
1117 this->f_[eIdxGlobalI] += addPcGradient;
1118 }
1119 }
1120}
1121
1122
1141template<class TypeTag>
1142void FV3dPressure2P2CAdaptive<TypeTag>::get1pMpfaFlux(const IntersectionIterator& isIt,
1143 const CellData& cellDataI)
1144{
1145 const auto& intersection = *isIt;
1146
1147 // acess Cell I
1148 auto elementI = intersection.inside();
1149 int eIdxGlobalI = problem().variables().index(elementI);
1150
1151 // get global coordinate of cell center
1152 const GlobalPosition& globalPos = elementI.geometry().center();
1153
1154 // access neighbor
1155 auto neighbor = intersection.outside();
1156 int eIdxGlobalJ = problem().variables().index(neighbor);
1157 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
1158
1159 // gemotry info of neighbor
1160 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
1161
1162 // due to "safety cell" around subdomain, both cells I and J
1163 // have single-phase conditions, although one is in 2p domain.
1164 using std::min;
1165 int phaseIdx = min(cellDataI.subdomain(), cellDataJ.subdomain());
1166
1167 // get average density for gravity flux
1168 Scalar rhoMean = 0.5 * (cellDataI.density(phaseIdx) + cellDataJ.density(phaseIdx));
1169
1170 // 1p => no pC => only 1 pressure, potential
1171 Scalar potential = 0.;
1172
1173 // Prepare MPFA
1175 GlobalPosition globalPos3(0.);
1176 int eIdxGlobal3=-1;
1177 GlobalPosition globalPos4(0.);
1178 int eIdxGlobal4=-1;
1179 TransmissivityMatrix T(0.);
1180
1181 // prepare second interaction region
1182 GlobalPosition globalPosAdditional3(0.);
1183 int eIdxGlobalAdditional3=-1;
1184 GlobalPosition globalPosAdditional4(0.);
1185 int eIdxGlobalAdditional4=-1;
1186
1187 TransmissivityMatrix additionalT(0.);
1188
1189 int interactionRegions
1190 = problem().variables().getMpfaData3D(intersection, T,
1191 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
1192 if (interactionRegions == 0)
1193 interactionRegions = problem().pressureModel().computeTransmissibilities(isIt,T,
1194 globalPos3, eIdxGlobal3, globalPos4, eIdxGlobal4 );
1195 if(!interactionRegions)
1196 Dune::dgrave << "something went wrong getting mpfa data on cell " << eIdxGlobalI << std::endl;
1197
1198 // shortcurts mpfa case
1199 CellData& cellData3 = problem().variables().cellData(eIdxGlobal3);
1200 CellData& cellData4 = problem().variables().cellData(eIdxGlobal4);
1201 Scalar temp1 = globalPos * this->gravity_;
1202 Scalar temp2 = globalPosNeighbor * this->gravity_;
1203 Scalar temp3 = globalPos3 * this->gravity_;
1204 Scalar temp4 = globalPos4 * this->gravity_;
1205
1206 potential = (cellDataI.pressure(phaseIdx)-temp1*rhoMean) * T[0]
1207 + (cellDataJ.pressure(phaseIdx)-temp2*rhoMean) * T[1]
1208 + (cellData3.pressure(phaseIdx)-temp3*rhoMean) * T[2]
1209 + (cellData4.pressure(phaseIdx)-temp4*rhoMean) * T[3];
1210
1211 // regard more interaction regions, if there are more
1212 if(interactionRegions != 1)
1213 {
1214 for(int banana = 1; banana < interactionRegions; banana ++)
1215 {
1216 // get data for second interaction region
1217 problem().variables().getMpfaData3D(intersection, additionalT,
1218 globalPosAdditional3, eIdxGlobalAdditional3,
1219 globalPosAdditional4, eIdxGlobalAdditional4 ,
1220 banana); // offset for second interaction region
1221
1222 Scalar gravityContributionAdditonal
1223 = temp1 * additionalT[0] + temp2 * additionalT[1]
1224 + globalPosAdditional3*this->gravity_ * additionalT[2]
1225 + globalPosAdditional4*this->gravity_ * additionalT[3];
1226 CellData& cellDataA3 = problem().variables().cellData(eIdxGlobalAdditional3);
1227 CellData& cellDataA4 = problem().variables().cellData(eIdxGlobalAdditional4);
1228
1229 potential += cellDataI.pressure(phaseIdx) * additionalT[0]
1230 + cellDataJ.pressure(phaseIdx) * additionalT[1]
1231 + cellDataA3.pressure(phaseIdx) * additionalT[2]
1232 + cellDataA4.pressure(phaseIdx) * additionalT[3];
1233 potential -= gravityContributionAdditonal * rhoMean;
1234 }
1235 }
1236
1237 // initialize convenience shortcuts
1238 Scalar lambda(0.);
1239
1240 //do the upwinding of the mobility depending on the phase potentials
1241 if (potential >= 0.)
1242 lambda = cellDataI.mobility(phaseIdx);
1243 else
1244 lambda = cellDataJ.mobility(phaseIdx);
1245
1247 /* extend T with other matrix entries and assemble to A_ */
1248 this->A_[eIdxGlobalI][eIdxGlobalI] += lambda * T[0];
1249 this->A_[eIdxGlobalI][eIdxGlobalJ] += lambda * T[1];
1250 this->A_[eIdxGlobalI][eIdxGlobal3] += lambda * T[2];
1251 this->A_[eIdxGlobalI][eIdxGlobal4] += lambda * T[3];
1252
1253 // add gravity to RHS vector
1254 Scalar gravityContribution = temp1 * T[0] + temp2 * T[1] + temp3 * T[2] + temp4 * T[3];
1255 this->f_[eIdxGlobalI] += lambda * rhoMean * gravityContribution;
1256
1257 // regard more interaction regions, if there are more
1258 if(interactionRegions != 1)
1259 {
1260 for(int banana = 1; banana < interactionRegions; banana ++)
1261 {
1262 // get data for second interaction region
1263 problem().variables().getMpfaData3D(intersection, additionalT,
1264 globalPosAdditional3, eIdxGlobalAdditional3,
1265 globalPosAdditional4, eIdxGlobalAdditional4 ,
1266 banana); // offset for second interaction region
1267
1268 Scalar gravityContributionAdditonal
1269 = temp1 * additionalT[0] + temp2 * additionalT[1]
1270 + globalPosAdditional3*this->gravity_ * additionalT[2]
1271 + globalPosAdditional4*this->gravity_ * additionalT[3];
1272
1274 /* extend T with other matrix entries and assemble to A_ */
1275 this->A_[eIdxGlobalI][eIdxGlobalI] += lambda * additionalT[0];
1276 this->A_[eIdxGlobalI][eIdxGlobalJ] += lambda * additionalT[1];
1277 this->A_[eIdxGlobalI][eIdxGlobalAdditional3] += lambda * additionalT[2];
1278 this->A_[eIdxGlobalI][eIdxGlobalAdditional4] += lambda * additionalT[3];
1279
1280 // add gravity to RHS vector
1281 this->f_[eIdxGlobalI] += lambda * rhoMean* gravityContributionAdditonal;
1282 }
1283 }
1284}
1285
1288
1294template<class TypeTag>
1296{
1297//#define noMultiphysics
1298#ifdef noMultiphysics
1300 return;
1301#endif
1302
1303 if(!fromPostTimestep) // this happens after grid update
1304 {
1305 Scalar maxError = 0.;
1306 // iterate through leaf grid an evaluate c0 at cell center
1307 for (const auto& element : elements(problem().gridView()))
1308 {
1309 int eIdxGlobal = problem().variables().index(element);
1310 CellData& cellData = problem().variables().cellData(eIdxGlobal);
1311
1312 if(cellData.fluidStateType() == 0) // i.e. it is complex
1313 problem().pressureModel().updateMaterialLawsInElement(element, fromPostTimestep);
1314 else
1315 problem().pressureModel().update1pMaterialLawsInElement(element, cellData, fromPostTimestep);
1316 using std::max;
1317 maxError = max(maxError, fabs(cellData.volumeError()));
1318 }
1319 this->maxError_ = maxError/problem().timeManager().timeStepSize();
1320 }
1321 else
1322 {
1323 // resize multiphysics vectors
1324 FVPressure2P2CMultiPhysics<TypeTag>::nextSubdomain.resize(problem().gridView().size(0));
1325
1327 }
1328}
1329
1348template <class TypeTag>
1350 TransmissivityMatrix& T,
1351 GlobalPosition& globalPos4,
1352 int& eIdxGlobal4,
1353 GlobalPosition& globalPos6,
1354 int& eIdxGlobal6)
1355{
1356 const auto& intersection = *isIt;
1357
1358 // get geometry information of cellI = cell1, cellJ = cell2
1359 auto element = intersection.inside();
1360 auto neighbor = intersection.outside();
1361 GlobalPosition globalPos1 = element.geometry().center();
1362 GlobalPosition globalPos2 = neighbor.geometry().center();
1363 DimMatrix K1(problem().spatialParams().intrinsicPermeability(element));
1364 DimMatrix K2(problem().spatialParams().intrinsicPermeability(neighbor));
1365
1366 // determine ID of intersection seen from larger cell
1367 int intersectionID = 0;
1368 if(intersection.inside().level() < intersection.outside().level())
1369 intersectionID = problem().grid().localIdSet().subId(element,
1370 intersection.indexInInside(), 1);
1371 else
1372 DUNE_THROW(Dune::NotImplemented, " ABORT, transmiss calculated from wrong side!!");
1373
1374 std::vector<int> localIrregularCells = irregularCellMap_[intersectionID];
1375
1377 // geometry and Data of face IJ in nomenclature of mpfa
1378 GlobalPosition globalPosFace12 = intersection.geometry().center(); // = 'x'1
1379 GlobalPosition outerNormaln12 = intersection.centerUnitOuterNormal();
1380
1381 /**** a) search for intersection 24 and 26 ************/
1382 auto face24=isIt; // as long as face24 = isIt, it is still not found!
1383 auto face26=isIt; // as long as face26 = isIt, it is still not found!
1384
1385 const auto isEndIt = problem().gridView().iend(neighbor);
1386 for (auto isIt2 = problem().gridView().ibegin(neighbor); isIt2 != isEndIt; ++isIt2)
1387 {
1388 const auto& intersection2 = *isIt2;
1389
1390 // continue if no neighbor or arrived at intersection
1391 if(!(intersection2.neighbor()) || intersection2.outside() == element)
1392 continue;
1393
1394 int currentNeighbor = problem().variables().index(intersection2.outside());
1395
1396 // have we found face24?
1397 if (find(localIrregularCells.begin(), localIrregularCells.end(),
1398 currentNeighbor) != localIrregularCells.end() && face24==isIt)
1399 face24 = isIt2;
1400 else if (find(localIrregularCells.begin(), localIrregularCells.end(),
1401 currentNeighbor) != localIrregularCells.end() && face26==isIt)
1402 {
1403 // we now found both intersections, but we have to investigate orientation:
1404 // Calculate the vector product of the normals in current orientation
1405 GlobalPosition vectorProduct = crossProduct(face24->centerUnitOuterNormal(), intersection2.centerUnitOuterNormal());
1406 // the vector product of the two isNormals for the desired configuration points
1407 // in the direction of outerNormaln12 => ScalarProduct > 0.
1408 if((vectorProduct * outerNormaln12) > 0.)
1409 face26 = isIt2;
1410 else
1411 {
1412 // revert orientation: switch face24 and face26
1413 face26 = face24;
1414 face24 = isIt2;
1415 }
1416 }
1417 }
1418
1419 // get information of cell4
1420 globalPos4 = face24->outside().geometry().center();
1421 eIdxGlobal4 = problem().variables().index(face24->outside());
1422 GlobalPosition outerNormaln24 = face24->centerUnitOuterNormal();
1423 // get absolute permeability of neighbor cell 3
1424 DimMatrix K4(problem().spatialParams().intrinsicPermeability(face24->outside()));
1425
1426 // get information of cell6
1427 globalPos6 = face26->outside().geometry().center();
1428 eIdxGlobal6 = problem().variables().index(face26->outside());
1429 GlobalPosition outerNormaln26 = face26->centerUnitOuterNormal();
1430 // get absolute permeability of neighbor cell 3
1431 DimMatrix K6(problem().spatialParams().intrinsicPermeability(face26->outside()));
1432
1433 /**** b) Get Points on the edges (in plane of isIt), 'x'4 and 'x'5 ***********/
1434 int localFace12 = intersection.indexInOutside(); // isIt is face12
1435 int localFace24 = face24->indexInInside();
1436 int localFace26 = face26->indexInInside();
1437
1438 const auto referenceElement = ReferenceElementContainer::general(neighbor.type());
1439
1440 //find 'x'5 = edgeCoord1226
1441 int edge1226;
1442 // search through edges of face 12
1443 for(int nectarine=0; nectarine < referenceElement.size(localFace12, 1, dim-1); nectarine++)
1444 {
1445 // get local Idx of edge on face 12
1446 int localEdgeOn12 = referenceElement.subEntity(localFace12, 1, nectarine, dim-1);
1447 // search through edges of face 26
1448 for(int plum = 0; plum < referenceElement.size(localFace26, 1,dim-1); plum++)
1449 {
1450// int localEdge26 = referenceElement.subEntity(localFace26, 1, plum, dim-1);
1451 if(referenceElement.subEntity(localFace12, 1, nectarine, dim-1)
1452 == referenceElement.subEntity(localFace26, 1, plum, dim-1))
1453 {
1454 edge1226 = localEdgeOn12;
1455 break;
1456 }
1457 }
1458 }
1459 GlobalPosition edgeCoord1226 = // 'x'5
1460 neighbor.geometry().global(referenceElement.position(edge1226, dim-1));
1461
1462 //find 'x'4 = edgeCoord1224
1463 int edge1224;
1464 // search through edges of face 12
1465 for(int nectarine=0; nectarine < referenceElement.size(localFace12, 1, dim-1); nectarine++)
1466 {
1467 // get local Idx of edge on face 12
1468 int localEdgeOn12 = referenceElement.subEntity(localFace12, 1, nectarine, dim-1);
1469 // search through edges of face 24
1470 for(int plum = 0; plum < referenceElement.size(localFace24, 1, dim-1); plum++)
1471 {
1472 int localEdge24 = referenceElement.subEntity(localFace24, 1, plum, dim-1);
1473 if(localEdgeOn12 == localEdge24)
1474 {
1475 edge1224 = localEdgeOn12;
1476 break;
1477 }
1478 }
1479 }
1480 GlobalPosition edgeCoord1224 = // 'x'4
1481 neighbor.geometry().global(referenceElement.position(edge1224, dim-1));
1482
1483 //find 'x'6 = edgeCoord2426
1484 int edge2426;
1485 // search through edges of face 24
1486 for(int nectarine=0; nectarine < referenceElement.size(localFace24, 1, dim-1); nectarine++)
1487 {
1488 // get local Idx of edge on face 24
1489 int localEdgeOn24 = referenceElement.subEntity(localFace24, 1, nectarine, dim-1);
1490 // search through edges of face 26
1491 for(int plum = 0; plum < referenceElement.size(localFace26, 1, dim-1); plum++)
1492 {
1493 int localEdge26 = referenceElement.subEntity(localFace26, 1, plum, dim-1);
1494 if(localEdgeOn24 == localEdge26)
1495 {
1496 edge2426 = localEdgeOn24;
1497 break;
1498 }
1499 }
1500 }
1501 GlobalPosition edgeCoord2426 = // 'x'6
1502 neighbor.geometry().global(referenceElement.position(edge2426, dim-1));
1503
1505 // center of face in global coordinates, i.e., the midpoint of face 'isIt24'
1506 GlobalPosition globalPosFace24 = face24->geometry().center();
1507 GlobalPosition globalPosFace26 = face26->geometry().center();
1508
1509 // get face volumes
1510 //Area12 = | 'x'1->'x'5 x 'x'1->'x'4 |
1511 Scalar subFaceArea12 = crossProduct(edgeCoord1226-globalPosFace12, edgeCoord1224-globalPosFace12).two_norm();
1512
1513 //Area24 = | 'x'2->'x'4 x 'x'2->'x'6 |
1514 Scalar subFaceArea24 = crossProduct(edgeCoord1224-globalPosFace24, edgeCoord2426-globalPosFace24).two_norm();
1515
1516 //Area26 = | 'x'3->'x'5 x 'x'3->'x'6 |
1517 Scalar subFaceArea26 = crossProduct(edgeCoord1226-globalPosFace26, edgeCoord2426-globalPosFace26).two_norm();
1518
1519 // compute normal vectors for case 2 (idx1, idx2, idx4, idx6)
1520 GlobalPosition nu11C2 = crossProduct(edgeCoord1226-globalPos1, globalPosFace12-globalPos1);
1521 GlobalPosition nu12C2 = crossProduct(edgeCoord1224-globalPos1, edgeCoord1226-globalPos1);
1522 GlobalPosition nu13C2 = crossProduct(globalPosFace12-globalPos1, edgeCoord1224-globalPos1);
1523
1524 GlobalPosition nu21C2 = crossProduct(globalPosFace26-globalPos2, globalPosFace24-globalPos2);
1525 GlobalPosition nu22C2 = crossProduct(globalPosFace12-globalPos2, globalPosFace26-globalPos2);
1526 GlobalPosition nu23C2 = crossProduct(globalPosFace24-globalPos2, globalPosFace12-globalPos2);
1527
1528 GlobalPosition nu41C2 = crossProduct(edgeCoord2426-globalPos4, edgeCoord1224-globalPos4);
1529 GlobalPosition nu42C2 = crossProduct(globalPosFace24-globalPos4, edgeCoord2426-globalPos4);
1530 GlobalPosition nu43C2 = crossProduct(edgeCoord1224-globalPos4, globalPosFace24-globalPos4);
1531
1532 GlobalPosition nu61C2 = crossProduct(globalPosFace26-globalPos6, edgeCoord1226-globalPos6);
1533 GlobalPosition nu62C2 = crossProduct(edgeCoord2426-globalPos6, globalPosFace26-globalPos6);
1534 GlobalPosition nu63C2 = crossProduct(edgeCoord1226-globalPos6, edgeCoord2426-globalPos6);
1535
1536 // compute T, i.e., the volume of the subtetrahedron
1537 Scalar T1C2 = (globalPosFace12-globalPos1) * crossProduct(edgeCoord1224-globalPos1, edgeCoord1226-globalPos1);
1538 Scalar T2C2 = (globalPosFace12-globalPos2) * crossProduct(globalPosFace26-globalPos2, globalPosFace24-globalPos2);
1539 Scalar T4C2 = (globalPosFace24-globalPos4) * crossProduct(edgeCoord2426-globalPos4, edgeCoord1224-globalPos4);
1540 Scalar T6C2 = (globalPosFace26-globalPos6) * crossProduct(edgeCoord1226-globalPos6, edgeCoord2426-globalPos6);
1541
1542 // compute components needed for flux calculation, denoted as 'omega' and 'r'
1543 GlobalPosition K1nu11C2(0);
1544 K1.mv(nu11C2, K1nu11C2);
1545 GlobalPosition K1nu12C2(0);
1546 K1.mv(nu12C2, K1nu12C2);
1547 GlobalPosition K1nu13C2(0);
1548 K1.mv(nu13C2, K1nu13C2);
1549
1550 GlobalPosition K2nu21C2(0);
1551 K2.mv(nu21C2, K2nu21C2);
1552 GlobalPosition K2nu22C2(0);
1553 K2.mv(nu22C2, K2nu22C2);
1554 GlobalPosition K2nu23C2(0);
1555 K2.mv(nu23C2, K2nu23C2);
1556
1557 GlobalPosition K4nu41C2(0);
1558 K4.mv(nu41C2, K4nu41C2);
1559 GlobalPosition K4nu42C2(0);
1560 K4.mv(nu42C2, K4nu42C2);
1561 GlobalPosition K4nu43C2(0);
1562 K4.mv(nu43C2, K4nu43C2);
1563
1564 GlobalPosition K6nu61C2(0);
1565 K6.mv(nu61C2, K6nu61C2);
1566 GlobalPosition K6nu62C2(0);
1567 K6.mv(nu62C2, K6nu62C2);
1568 GlobalPosition K6nu63C2(0);
1569 K6.mv(nu63C2, K6nu63C2);
1570
1571 // upwinding of lambda is not included here to "recycle" it-> only possible if lambda does
1572 // not vary harshly. If it would, refinement indicator would not allow different grid levels here.
1573 Scalar omega111C2 = /*lambda12C2 */ (outerNormaln12 * K1nu11C2) * subFaceArea12/T1C2;
1574 Scalar omega112C2 = /*lambda12C2 */ (outerNormaln12 * K1nu12C2) * subFaceArea12/T1C2;
1575 Scalar omega113C2 = /*lambda12C2 */ (outerNormaln12 * K1nu13C2) * subFaceArea12/T1C2;
1576
1577 Scalar omega121C2 = /*lambda21C2 */ (outerNormaln12 * K2nu21C2) * subFaceArea12/T2C2;
1578 Scalar omega122C2 = /*lambda21C2 */ (outerNormaln12 * K2nu22C2) * subFaceArea12/T2C2;
1579 Scalar omega123C2 = /*lambda21C2 */ (outerNormaln12 * K2nu23C2) * subFaceArea12/T2C2;
1580
1581 Scalar omega221C2 = /*lambda24C2 */ (outerNormaln24 * K2nu21C2) * subFaceArea24/T2C2;
1582 Scalar omega222C2 = /*lambda24C2 */ (outerNormaln24 * K2nu22C2) * subFaceArea24/T2C2;
1583 Scalar omega223C2 = /*lambda24C2 */ (outerNormaln24 * K2nu23C2) * subFaceArea24/T2C2;
1584
1585 Scalar omega241C2 = /*lambda42C2 */ (outerNormaln24 * K4nu41C2) * subFaceArea24/T4C2;
1586 Scalar omega242C2 = /*lambda42C2 */ (outerNormaln24 * K4nu42C2) * subFaceArea24/T4C2;
1587 Scalar omega243C2 = /*lambda42C2 */ (outerNormaln24 * K4nu43C2) * subFaceArea24/T4C2;
1588
1589 Scalar omega321C2 = /*lambda26C2 */ (outerNormaln26 * K2nu21C2) * subFaceArea26/T2C2;
1590 Scalar omega322C2 = /*lambda26C2 */ (outerNormaln26 * K2nu22C2) * subFaceArea26/T2C2;
1591 Scalar omega323C2 = /*lambda26C2 */ (outerNormaln26 * K2nu23C2) * subFaceArea26/T2C2;
1592
1593 Scalar omega361C2 = /*lambda62C2 */ (outerNormaln26 * K6nu61C2) * subFaceArea26/T6C2;
1594 Scalar omega362C2 = /*lambda62C2 */ (outerNormaln26 * K6nu62C2) * subFaceArea26/T6C2;
1595 Scalar omega363C2 = /*lambda62C2 */ (outerNormaln26 * K6nu63C2) * subFaceArea26/T6C2;
1596
1597 Scalar r211C2 = (nu21C2 * (edgeCoord1224-globalPos2))/T2C2;
1598 Scalar r212C2 = (nu21C2 * (edgeCoord1226-globalPos2))/T2C2;
1599 Scalar r213C2 = (nu21C2 * (edgeCoord2426-globalPos2))/T2C2;
1600
1601 Scalar r221C2 = (nu22C2 * (edgeCoord1224-globalPos2))/T2C2;
1602 Scalar r222C2 = (nu22C2 * (edgeCoord1226-globalPos2))/T2C2;
1603 Scalar r223C2 = (nu22C2 * (edgeCoord2426-globalPos2))/T2C2;
1604
1605 Scalar r231C2 = (nu23C2 * (edgeCoord1224-globalPos2))/T2C2;
1606 Scalar r232C2 = (nu23C2 * (edgeCoord1226-globalPos2))/T2C2;
1607 Scalar r233C2 = (nu23C2 * (edgeCoord2426-globalPos2))/T2C2;
1608
1609
1611 // compute transmissibility matrix T = CA^{-1}B+D
1612 DimMatrix C(0), A(0);
1613 Dune::FieldMatrix<Scalar,dim,dim+1> D(0), B(0);
1614
1615 // evaluate matrix C, D, A, B
1616 C[0][0] = -omega121C2;
1617 C[0][1] = -omega122C2;
1618 C[0][2] = -omega123C2;
1619 C[1][0] = -omega221C2;
1620 C[1][1] = -omega222C2;
1621 C[1][2] = -omega223C2;
1622 C[2][0] = -omega321C2;
1623 C[2][1] = -omega322C2;
1624 C[2][2] = -omega323C2;
1625
1626 D[0][1] = omega121C2 + omega122C2 + omega123C2;
1627 D[1][1] = omega221C2 + omega222C2 + omega223C2;
1628 D[2][1] = omega321C2 + omega322C2 + omega323C2;
1629
1630 A[0][0] = omega121C2 - omega112C2 - omega111C2*r211C2 - omega113C2*r212C2;
1631 A[0][1] = omega122C2 - omega111C2*r221C2 - omega113C2*r222C2;
1632 A[0][2] = omega123C2 - omega111C2*r231C2 - omega113C2*r232C2;
1633 A[1][0] = omega221C2 - omega242C2*r211C2 - omega243C2*r213C2;
1634 A[1][1] = omega222C2 - omega241C2 - omega242C2*r221C2 - omega243C2*r223C2;
1635 A[1][2] = omega223C2 - omega242C2*r231C2 - omega243C2*r233C2;
1636 A[2][0] = omega321C2 - omega361C2*r213C2 - omega362C2*r212C2;
1637 A[2][1] = omega322C2 - omega361C2*r223C2 - omega362C2*r222C2;
1638 A[2][2] = omega323C2 - omega363C2 - omega361C2*r233C2 - omega362C2*r232C2;
1639
1640 B[0][0] = -omega111C2 - omega112C2 - omega113C2;
1641 B[0][1] = omega121C2 + omega122C2 + omega123C2 + omega111C2*(1.0 - r211C2 - r221C2 -r231C2)
1642 + omega113C2*(1.0 - r212C2 - r222C2 - r232C2);
1643 B[1][1] = omega221C2 + omega222C2 + omega223C2 + omega242C2*(1.0 - r211C2 - r221C2 - r231C2)
1644 + omega243C2*(1.0 - r213C2 - r223C2 - r233C2);
1645 B[1][2] = -omega241C2 - omega242C2 - omega243C2;
1646 B[2][1] = omega321C2 + omega322C2 + omega323C2 + omega361C2*(1.0 - r213C2 - r223C2 - r233C2)
1647 + omega362C2*(1.0 - r212C2 -r222C2 -r232C2);
1648 B[2][3] = -omega361C2 - omega362C2 -omega363C2;
1649
1650// printmatrix(std::cout, A, "matrix A ", "row", 11, 3);
1651// printmatrix(std::cout, B, "matrix B ", "row", 11, 3);
1652// printmatrix(std::cout, C, "matrix C ", "row", 11, 3);
1653// printmatrix(std::cout, D, "matrix D ", "row", 11, 3);
1654
1655 // compute T
1656 A.invert();
1657 D += B.leftmultiply(C.rightmultiply(A));
1658 T = D[0];
1659
1661 if(maxInteractionVolumes ==1)
1662 {
1663 T *= intersection.geometry().volume()/subFaceArea12 ;
1664
1665 // set your map entry
1666 problem().variables().storeMpfaData3D(intersection, T, globalPos4, eIdxGlobal4, globalPos6, eIdxGlobal6);
1667 return 1; // indicates that only 1 interaction region was regarded
1668 }
1669 else
1670 {
1671 int countInteractionRegions = 1;
1672 //initialize additional transmissitivity matrix
1673 TransmissivityMatrix additionalT(0.);
1674
1675 auto outerCorner = intersection.inside().template subEntity<dim>(0); //initialize with rubbish
1676 // prepare additonal pointer to cells
1677 auto additional2 = intersection.inside(); //initialize with something wrong!
1678 auto additional3 = intersection.inside();
1679 int caseL = -2;
1680
1681 /**** 2nd interaction region: get corner of interest ************/
1682 // search through corners of large cell with isIt
1683 int localIdxLarge = searchCommonVertex_(intersection, outerCorner);
1684
1685 //in the parallel case, skip all border entities
1686 #if HAVE_MPI
1687 if (problem().gridView().comm().size() > 1)
1688 if(outerCorner.partitionType() != Dune::InteriorEntity)
1689 caseL = -1; // abort this specific interaction volume
1690 #endif
1691
1692 // get Interaction Volume object
1693 int vIdxGlobal = problem().variables().index(outerCorner);
1694 InteractionVolume& interactionVolume
1695 = interactionVolumesContainer_->interactionVolume(vIdxGlobal);
1696
1697 // abort if we are on boundary
1698 if(interactionVolume.isBoundaryInteractionVolume())
1699 caseL = -1;
1700
1701 // use Markus mpfa-l implementation
1702 int subVolumeFaceIdx = -1;
1703 bool properFluxDirection = true;
1704 if(caseL == -2)
1705 {
1706 //get Hanging Node type: no hanging node would mean -1
1707 int hangingNodeType = interactionVolume.getHangingNodeType();
1708
1709 if(hangingNodeType == InteractionVolume::noHangingNode)
1710 subVolumeFaceIdx = interactionVolumesContainer_->getMpfaCase8cells(isIt, localIdxLarge, interactionVolume, properFluxDirection);
1711 else if(hangingNodeType == InteractionVolume::sixSmallCells)
1712 subVolumeFaceIdx = interactionVolumesContainer_->getMpfaCase6cells(isIt, interactionVolume, properFluxDirection);
1713 else
1714 Dune::dgrave << "HangingType " << hangingNodeType << " not implemented " << std::endl;
1715
1716 caseL = this->transmissibilityAdapter_(isIt, interactionVolume, subVolumeFaceIdx,
1717 properFluxDirection, additional2, additional3, additionalT);
1718 }
1719
1720 //store everything as second interaction region
1721 if(caseL != -1) //check if we regard 2 interaction regions
1722 {
1723 problem().variables().storeMpfaData3D(intersection, additionalT,
1724 additional2.geometry().center(), problem().variables().index(additional2),
1725 additional3.geometry().center(), problem().variables().index(additional3),
1726 1); // offset for second interaction region
1727 countInteractionRegions++;
1728 }
1729
1730 /***** 3rd and 4th interaction region *******/
1731 if(maxInteractionVolumes>2)
1732 {
1733 // loop through remaining 2 points
1734 std::vector<int> diagonal;
1735 for(int verticeSmall = 0; verticeSmall < intersection.outside().subEntities(dim); ++verticeSmall)
1736 {
1737 auto vSmall = intersection.outside().template subEntity<dim>(verticeSmall);
1738
1739 //in the parallel case, skip all border entities
1740 #if HAVE_MPI
1741 if (problem().gridView().comm().size() > 1)
1742 if(vSmall.partitionType() != Dune::InteriorEntity)
1743 continue;
1744 #endif
1745
1746 // get postion as seen from element
1747 GlobalPosition vertexOnElement
1748 = referenceElement.position(verticeSmall, dim);
1749
1750 for (int indexOnFace = 0; indexOnFace < 4; indexOnFace++)
1751 {
1752 // get position as seen from interface
1753 GlobalPosition vertexOnInterface
1754 = intersection.geometryInOutside().corner(indexOnFace);
1755
1756 if(vSmall != outerCorner
1757 && ((vertexOnInterface - vertexOnElement).two_norm()<1e-5))
1758 {
1759 int vIdxGlobal2 = problem().variables().index(vSmall);
1760 // acess interactionVolume
1761 InteractionVolume& interactionVolume2
1762 = interactionVolumesContainer_->interactionVolume(vIdxGlobal2);
1763 if(interactionVolume2.isBoundaryInteractionVolume())
1764 continue;
1765
1766 int hangingNodeType = interactionVolume2.getHangingNodeType();
1767 // reset flux direction indicator
1768 properFluxDirection = true;
1769
1770 if(hangingNodeType != InteractionVolume::fourSmallCellsFace)
1771 {
1772 diagonal.push_back(problem().variables().index(vSmall));
1773 // a) take interaction volume and determine fIdx
1774 if(hangingNodeType == InteractionVolume::noHangingNode)
1775 {
1776 //TODO determine current localIdxLarge!!!!
1777 Dune::dgrave << " args, noHangingNode on additional interaction region";
1778 // subVolumeFaceIdx = getMpfaCase8cells_(isIt, localIdxLarge, interactionVolume2, properFluxDirection);
1779 }
1780 else if(hangingNodeType == InteractionVolume::sixSmallCells)
1781 subVolumeFaceIdx = interactionVolumesContainer_->getMpfaCase6cells(isIt,
1782 interactionVolume2, properFluxDirection);
1783 else
1784 subVolumeFaceIdx = interactionVolumesContainer_->getMpfaCase2or4cells(isIt,
1785 interactionVolume2, properFluxDirection);
1786
1787 // b) calculate T, eIdxGlobal3+4
1788 caseL = this->transmissibilityAdapter_(isIt, interactionVolume2, subVolumeFaceIdx,
1789 properFluxDirection, additional2, additional3, additionalT);
1790
1791 // c) store it
1792 //store everything as third/4th interaction region
1793 if(caseL != -1) //check if we regard this interaction region
1794 {
1795 problem().variables().storeMpfaData3D(intersection, additionalT,
1796 additional2.geometry().center(), problem().variables().index(additional2),
1797 additional3.geometry().center(), problem().variables().index(additional3),
1798 countInteractionRegions); // offset for this interaction region
1799 countInteractionRegions++;
1800 }
1801 }
1802 //else we would have type 1, which is handled explicitly in the beginning of this method
1803 }
1804 }
1805 }
1806 }
1807
1808 // set your map entry
1809 problem().variables().storeMpfaData3D(intersection, T, globalPos4, eIdxGlobal4, globalPos6, eIdxGlobal6);
1810
1811 // determine weights
1812 Scalar weight = intersection.geometry().volume()/subFaceArea12; // =4 for 1 interaction region
1813 weight /= static_cast<Scalar>(countInteractionRegions);
1814
1815 // perform weighting of transmissibilies
1816 problem().variables().performTransmissitivityWeighting(intersection, weight);
1817
1818 return countInteractionRegions;
1819 }
1820
1821 return 0;
1822}
1823
1841template<class TypeTag>
1843 InteractionVolume& interactionVolume,
1844 const int& subVolumeFaceIdx,
1845 bool properFluxDirection,
1846 Element& additional2,
1847 Element& additional3,
1848 TransmissivityMatrix& additionalT)
1849{
1850 // abort if we are on boundary
1851 if(interactionVolume.isBoundaryInteractionVolume())
1852 return -1;
1853
1854 //get Hanging Node type: no hanging node would mean -1
1855 int hangingNodeType = interactionVolume.getHangingNodeType();
1856
1857 // prepare necessary variables etc
1858 Dune::FieldMatrix<Scalar, dim, 2 * dim - dim + 1> T(0);
1859 int caseL = -1;
1860 Dune::FieldVector<Scalar, dim> unity(1.);
1861 std::vector<Dune::FieldVector<Scalar, dim> > lambda
1862 = {unity, unity, unity, unity, unity, unity, unity, unity};
1863 Dune::FieldVector<bool, 4> useCases(false);
1864
1865 /*********** calculate transmissibility for that case ******/
1866 switch(subVolumeFaceIdx)
1867 {
1868 case 0:{
1869 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1870 lambda, 0, 1, 2, 3, 4, 5);
1871
1872 additionalT = T[0];
1873
1874 // determine cells for T-matrix entry 2 and 3
1875 switch (caseL)
1876 {
1877 case 1:{
1878 additional2 = interactionVolume.getSubVolumeElement(2);
1879 additional3 = interactionVolume.getSubVolumeElement(4);
1880 break; }
1881 case 2:{
1882 additional2 = interactionVolume.getSubVolumeElement(3);
1883 additional3 = interactionVolume.getSubVolumeElement(5);
1884 break; }
1885 case 3:{
1886 additional2 = interactionVolume.getSubVolumeElement(3);
1887 additional3 = interactionVolume.getSubVolumeElement(4);
1888 break; }
1889 case 4:{
1890 additional2 = interactionVolume.getSubVolumeElement(2);
1891 additional3 = interactionVolume.getSubVolumeElement(5);
1892 break; }
1893 }
1894 break;
1895 }
1896
1897 case 1:{
1898 if (hangingNodeType == InteractionVolume::twoSmallCells
1899 || hangingNodeType == InteractionVolume::fourSmallCellsDiag )
1900 {
1901 useCases[0] = true;
1902 useCases[1] = false;
1903 useCases[2] = false;
1904 if(hangingNodeType != InteractionVolume::twoSmallCells)
1905 useCases[3] = true; //differs from 2p Implementation
1906 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1907 lambda, 1, 3, 0, 2, 5, 7, useCases);
1908 }
1909 else
1910 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1911 lambda, 1, 3, 0, 2, 5, 7);
1912
1913 additionalT = T[0];
1914
1915 // determine cells for T-matrix entry 2 and 3
1916 switch (caseL)
1917 {
1918 case 1:
1919 {
1920 additional2 = interactionVolume.getSubVolumeElement(0);
1921 additional3 = interactionVolume.getSubVolumeElement(5);
1922 break;}
1923 case 2:
1924 {
1925 additional2 = interactionVolume.getSubVolumeElement(2);
1926 additional3 = interactionVolume.getSubVolumeElement(7);
1927 break;}
1928 case 3:
1929 {
1930 additional2 = interactionVolume.getSubVolumeElement(2);
1931 additional3 = interactionVolume.getSubVolumeElement(5);
1932 break;}
1933 case 4:
1934 {
1935 additional2 = interactionVolume.getSubVolumeElement(0);
1936 additional3 = interactionVolume.getSubVolumeElement(7);
1937 break;}
1938 }
1939 break;
1940 }
1941
1942 case 2:
1943 {
1944 assert (hangingNodeType != InteractionVolume::twoSmallCells
1945 && hangingNodeType != InteractionVolume::fourSmallCellsDiag);
1946
1947 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1948 lambda, 3, 2, 1, 0, 7, 6);
1949
1950 additionalT = T[0];
1951
1952 // determine cells for T-matrix entry 2 and 3
1953 switch (caseL)
1954 {
1955 case 1:
1956 {
1957 additional2 = interactionVolume.getSubVolumeElement(1);
1958 additional3 = interactionVolume.getSubVolumeElement(7);
1959 break;}
1960 case 2:
1961 {
1962 additional2 = interactionVolume.getSubVolumeElement(0);
1963 additional3 = interactionVolume.getSubVolumeElement(6);
1964 break;}
1965 case 3:
1966 {
1967 additional2 = interactionVolume.getSubVolumeElement(0);
1968 additional3 = interactionVolume.getSubVolumeElement(7);
1969 break;}
1970 case 4:
1971 {
1972 additional2 = interactionVolume.getSubVolumeElement(1);
1973 additional3 = interactionVolume.getSubVolumeElement(6);
1974 break;}
1975 }
1976 break;
1977 }
1978
1979 case 3:{
1980 if (hangingNodeType == InteractionVolume::twoSmallCells
1981 || hangingNodeType == InteractionVolume::fourSmallCellsDiag)
1982 {
1983 useCases[0] = false;
1984 useCases[1] = true;
1985 if(hangingNodeType != InteractionVolume::twoSmallCells)
1986 useCases[2] = true; //TODO: warning, hacked from markus!!
1987 useCases[3] = false;
1988 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T,
1989 interactionVolume, lambda, 2, 0, 3, 1, 6, 4, useCases);
1990 }
1991 else
1992 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
1993 lambda, 2, 0, 3, 1, 6, 4);
1994
1995 additionalT = T[0];
1996
1997 // determine cells for T-matrix entry 2 and 3
1998 switch (caseL)
1999 {
2000 case 1:
2001 {
2002 additional2 = interactionVolume.getSubVolumeElement(3);
2003 additional3 = interactionVolume.getSubVolumeElement(6);
2004 break;}
2005 case 2:
2006 {
2007 additional2 = interactionVolume.getSubVolumeElement(1);
2008 additional3 = interactionVolume.getSubVolumeElement(4);
2009 break;}
2010 case 3:
2011 {
2012 additional2 = interactionVolume.getSubVolumeElement(1);
2013 additional3 = interactionVolume.getSubVolumeElement(6);
2014 break;}
2015 case 4:
2016 {
2017 additional2 = interactionVolume.getSubVolumeElement(3);
2018 additional3 = interactionVolume.getSubVolumeElement(4);
2019 break;}
2020 }
2021 break;
2022 }
2023
2024 case 4:{
2025 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge) // should never happen, cause TPFA should be applied
2026 {
2027 assert(problem().variables().index(interactionVolume.getSubVolumeElement(4))
2028 != problem().variables().index(interactionVolume.getSubVolumeElement(5))); // else there would not be a subVolFaceIdx 4
2029 Dune::dgrave << " SubVolumeFace4 in hanging node type 3 should be modelled by"
2030 << " Tpfa!!" <<std::endl; // TODO: or use 1,5,3,4 and case 4
2031 return -1;
2032 }
2033 else
2034 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2035 lambda, 5, 4, 7, 6, 1, 0);
2036
2037 additionalT = T[0];
2038
2039 // determine cells for T-matrix entry 2 and 3
2040 switch (caseL)
2041 {
2042 case 1:{
2043 additional2 = interactionVolume.getSubVolumeElement(7);
2044 additional3 = interactionVolume.getSubVolumeElement(1);
2045 break; }
2046 case 2:{
2047 additional2 = interactionVolume.getSubVolumeElement(6);
2048 additional3 = interactionVolume.getSubVolumeElement(0);
2049 break; }
2050 case 3:{
2051 additional2 = interactionVolume.getSubVolumeElement(6);
2052 additional3 = interactionVolume.getSubVolumeElement(1);
2053 break; }
2054 case 4:{
2055 additional2 = interactionVolume.getSubVolumeElement(7);
2056 additional3 = interactionVolume.getSubVolumeElement(0);
2057 break; }
2058 }
2059 break;
2060 }
2061
2062 case 5:{
2063 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2064 {
2065 assert (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2066 != problem().variables().index(interactionVolume.getSubVolumeElement(6))); // else there would not be a subVolFaceIdx 5
2067 Dune::dgrave << " SubVolumeFace5 in hanging node type 3 should be modelled by"
2068 << " Tpfa!!" <<std::endl; // TODO: or use 1,5,7,0 and case 3
2069 return -1;
2070 }
2071 else if (hangingNodeType == InteractionVolume::sixSmallCells)
2072 {
2073 useCases[0] = false;
2074 useCases[1] = true;
2075 useCases[2] = true;
2076 useCases[3] = false;
2077
2078 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2079 lambda, 7, 5, 6, 4, 3, 1, useCases);
2080 }
2081 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2082 {
2083 useCases[0] = true;
2084 useCases[1] = false;
2085 useCases[2] = false;
2086 useCases[3] = true;
2087 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2088 lambda, 7, 5, 6, 4, 3, 1, useCases);
2089 }
2090 else
2091 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2092 lambda, 7, 5, 6, 4, 3, 1);
2093
2094 additionalT = T[0];
2095
2096 // determine cells for T-matrix entry 2 and 3
2097 switch (caseL)
2098 {
2099 case 1:
2100 {
2101 additional2 = interactionVolume.getSubVolumeElement(6);
2102 additional3 = interactionVolume.getSubVolumeElement(3);
2103 break;}
2104 case 2:
2105 {
2106 additional2 = interactionVolume.getSubVolumeElement(4);
2107 additional3 = interactionVolume.getSubVolumeElement(1);
2108 break;}
2109 case 3:
2110 {
2111 additional2 = interactionVolume.getSubVolumeElement(4);
2112 additional3 = interactionVolume.getSubVolumeElement(3);
2113 break;}
2114 case 4:
2115 {
2116 additional2 = interactionVolume.getSubVolumeElement(6);
2117 additional3 = interactionVolume.getSubVolumeElement(1);
2118 break;}
2119 }
2120 break;
2121 }
2122
2123 case 6:{
2124 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2125 {
2126 assert (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2127 != problem().variables().index(interactionVolume.getSubVolumeElement(5))); // else there would not be a subVolFaceIdx 4
2128 Dune::dgrave << " SubVolumeFace6 in hanging node type 3 should be modelled by"
2129 << " Tpfa!!" <<std::endl; // TODO: or use 2,6,0,7 and case 4
2130 return -1;
2131 }
2132 else
2133 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2134 lambda, 6, 7, 4, 5, 2, 3);
2135
2136 additionalT = T[0];
2137
2138 // determine cells for T-matrix entry 2 and 3
2139 switch (caseL)
2140 {
2141 case 1:
2142 {
2143 additional2 = interactionVolume.getSubVolumeElement(4);
2144 additional3 = interactionVolume.getSubVolumeElement(2);
2145 break;}
2146 case 2:
2147 {
2148 additional2 = interactionVolume.getSubVolumeElement(5);
2149 additional3 = interactionVolume.getSubVolumeElement(3);
2150 break;}
2151 case 3:
2152 {
2153 additional2 = interactionVolume.getSubVolumeElement(5);
2154 additional3 = interactionVolume.getSubVolumeElement(2);
2155 break;}
2156 case 4:
2157 {
2158 additional2 = interactionVolume.getSubVolumeElement(4);
2159 additional3 = interactionVolume.getSubVolumeElement(3);
2160 break;}
2161 }
2162 break;
2163 }
2164
2165 case 7:
2166 {
2167 if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2168 {
2169 assert (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2170 != problem().variables().index(interactionVolume.getSubVolumeElement(6))); // else there would not be a subVolFaceIdx 5
2171 Dune::dgrave << " SubVolumeFace5 in hanging node type 3 should be modelled by"
2172 << " Tpfa!!" <<std::endl; // TODO: or use 4,0,6,1 and case 4
2173 return -1;
2174 }
2175 else if (hangingNodeType == InteractionVolume::sixSmallCells)
2176 {
2177 useCases[0] = true;
2178 useCases[1] = false;
2179 useCases[2] = false;
2180 useCases[3] = true;
2181
2182 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2183 lambda, 4, 6, 5, 7, 0, 2, useCases);
2184 }
2185 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2186 {
2187 useCases[0] = false;
2188 useCases[1] = true;
2189 useCases[2] = true;
2190 useCases[3] = false;
2191 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2192 lambda, 4, 6, 5, 7, 0, 2, useCases);
2193 }
2194 else
2195 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2196 lambda, 4, 6, 5, 7, 0, 2);
2197
2198 additionalT = T[0];
2199
2200 // determine cells for T-matrix entry 2 and 3
2201 switch (caseL)
2202 {
2203 case 1:
2204 {
2205 additional2 = interactionVolume.getSubVolumeElement(5);
2206 additional3 = interactionVolume.getSubVolumeElement(0);
2207 break;}
2208 case 2:
2209 {
2210 additional2 = interactionVolume.getSubVolumeElement(7);
2211 additional3 = interactionVolume.getSubVolumeElement(2);
2212 break;}
2213 case 3:
2214 {
2215 additional2 = interactionVolume.getSubVolumeElement(7);
2216 additional3 = interactionVolume.getSubVolumeElement(0);
2217 break;}
2218 case 4:
2219 {
2220 additional2 = interactionVolume.getSubVolumeElement(5);
2221 additional3 = interactionVolume.getSubVolumeElement(2);
2222 break;}
2223 }
2224 break;
2225 }
2226
2227 case 8:
2228 {
2229 if(hangingNodeType == InteractionVolume::noHangingNode
2230 || hangingNodeType == InteractionVolume::sixSmallCells)
2231 {
2232 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2233 lambda, 4, 0, 6, 2, 5, 1);
2234 }
2235 else if (hangingNodeType == InteractionVolume::twoSmallCells
2236 || hangingNodeType == InteractionVolume::fourSmallCellsFace)
2237 {
2238 caseL = mpfal3DTransmissibilityCalculator_.transmissibilityCaseTwo(T, interactionVolume,
2239 lambda, 4, 0, 2, 1);
2240 }
2241 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
2242 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge
2243 && (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2244 != problem().variables().index(interactionVolume.getSubVolumeElement(6)))) )
2245 {
2246 useCases[0] = false;
2247 useCases[1] = true;
2248 useCases[2] = false;
2249 useCases[3] = true;
2250
2251 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2252 lambda, 4, 0, 6, 2, 5, 1, useCases);
2253 }
2254 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2255 {
2256 useCases[0] = false;
2257 useCases[1] = true;
2258 useCases[2] = true;
2259 useCases[3] = false;
2260
2261 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2262 lambda, 4, 0, 6, 2, 5, 1, useCases);
2263 }
2264
2265 additionalT = T[0];
2266
2267 // determine cells for T-matrix entry 2 and 3
2268 switch (caseL)
2269 {
2270 case 1:
2271 {
2272 additional2 = interactionVolume.getSubVolumeElement(6);
2273 additional3 = interactionVolume.getSubVolumeElement(5);
2274 break;}
2275 case 2:
2276 {
2277 additional2 = interactionVolume.getSubVolumeElement(2);
2278 additional3 = interactionVolume.getSubVolumeElement(1);
2279 break;}
2280 case 3:
2281 {
2282 additional2 = interactionVolume.getSubVolumeElement(2);
2283 additional3 = interactionVolume.getSubVolumeElement(5);
2284 break;}
2285 case 4:
2286 {
2287 additional2 = interactionVolume.getSubVolumeElement(6);
2288 additional3 = interactionVolume.getSubVolumeElement(1);
2289 break;}
2290 }
2291 break;
2292 }
2293
2294 case 9:
2295 {
2296 if(hangingNodeType == InteractionVolume::noHangingNode
2297 || hangingNodeType == InteractionVolume::sixSmallCells)
2298 {
2299 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2300 lambda, 1, 5, 3, 7, 0, 4);
2301 }
2302 else if (hangingNodeType == InteractionVolume::twoSmallCells || hangingNodeType == InteractionVolume::fourSmallCellsFace)
2303 {
2304 caseL = mpfal3DTransmissibilityCalculator_.transmissibilityCaseOne(T, interactionVolume,
2305 lambda, 1, 5, 3, 0);
2306 }
2307 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag
2308 || (hangingNodeType == InteractionVolume::fourSmallCellsEdge
2309 &&(problem().variables().index(interactionVolume.getSubVolumeElement(4))
2310 != problem().variables().index(interactionVolume.getSubVolumeElement(6))) ))
2311 {
2312 useCases[0] = true;
2313 useCases[1] = false;
2314 useCases[2] = true;
2315 useCases[3] = false;
2316
2317 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2318 lambda, 1, 5, 3, 7, 0, 4, useCases);
2319 }
2320 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge
2321 &&(problem().variables().index(interactionVolume.getSubVolumeElement(4))
2322 != problem().variables().index(interactionVolume.getSubVolumeElement(5))) )
2323 {
2324 useCases[0] = true;
2325 useCases[1] = false;
2326 useCases[2] = false;
2327 useCases[3] = true;
2328
2329 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2330 lambda, 1, 5, 3, 7, 0, 4, useCases);
2331 }
2332 else
2333 Dune::dgrave << " Missing case for subVolFaceIdx 9!!" <<std::endl;
2334
2335 additionalT = T[0];
2336
2337 // determine cells for T-matrix entry 2 and 3
2338 switch (caseL)
2339 {
2340 case 1:
2341 {
2342 additional2 = interactionVolume.getSubVolumeElement(3);
2343 additional3 = interactionVolume.getSubVolumeElement(0);
2344 break;}
2345 case 2:
2346 {
2347 additional2 = interactionVolume.getSubVolumeElement(7);
2348 additional3 = interactionVolume.getSubVolumeElement(4);
2349 break;}
2350 case 3:
2351 {
2352 additional2 = interactionVolume.getSubVolumeElement(7);
2353 additional3 = interactionVolume.getSubVolumeElement(0);
2354 break;}
2355 case 4:
2356 {
2357 additional2 = interactionVolume.getSubVolumeElement(3);
2358 additional3 = interactionVolume.getSubVolumeElement(4);
2359 break;}
2360 }
2361 break;
2362 }
2363
2364 case 10:
2365 {
2366 if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
2367 caseL = mpfal3DTransmissibilityCalculator_.transmissibilityCaseTwo(T, interactionVolume,
2368 lambda, 7, 3, 1, 2);
2369 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge
2370 &&(problem().variables().index(interactionVolume.getSubVolumeElement(4))
2371 != problem().variables().index(interactionVolume.getSubVolumeElement(6))) )
2372 || hangingNodeType == InteractionVolume::sixSmallCells)
2373 {
2374 useCases[0] = false;
2375 useCases[1] = true;
2376 useCases[2] = false;
2377 useCases[3] = true;
2378
2379 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2380 lambda, 7, 3, 5, 1, 6, 2, useCases);
2381 }
2382 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge)
2383 {
2384 useCases[0] = false;
2385 useCases[1] = true;
2386 useCases[2] = true;
2387 useCases[3] = false;
2388
2389 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2390 lambda, 7, 3, 5, 1, 6, 2, useCases);
2391 }
2392 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2393 {
2394 useCases[0] = true;
2395 useCases[1] = false;
2396 useCases[2] = true;
2397 useCases[3] = false;
2398
2399 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2400 lambda, 7, 3, 5, 1, 6, 2, useCases);
2401 }
2402 else
2403 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2404 lambda, 7, 3, 5, 1, 6, 2);
2405
2406 additionalT = T[0];
2407
2408 // determine cells for T-matrix entry 2 and 3
2409 switch (caseL)
2410 {
2411 case 1:
2412 {
2413 additional2 = interactionVolume.getSubVolumeElement(5);
2414 additional3 = interactionVolume.getSubVolumeElement(6);
2415 break;}
2416 case 2:
2417 {
2418 additional2 = interactionVolume.getSubVolumeElement(1);
2419 additional3 = interactionVolume.getSubVolumeElement(2);
2420 break;}
2421 case 3:
2422 {
2423 additional2 = interactionVolume.getSubVolumeElement(1);
2424 additional3 = interactionVolume.getSubVolumeElement(6);
2425 break;}
2426 case 4:
2427 {
2428 additional2 = interactionVolume.getSubVolumeElement(5);
2429 additional3 = interactionVolume.getSubVolumeElement(2);
2430 break;}
2431 }
2432 break;
2433 }
2434
2435 case 11:
2436 {
2437 if(!interactionVolume.isHangingNodeVolume())
2438 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2439 lambda, 2, 6, 0, 4, 3, 7);
2440 else if (hangingNodeType == InteractionVolume::fourSmallCellsFace)
2441 caseL = mpfal3DTransmissibilityCalculator_.transmissibilityCaseOne(T, interactionVolume,
2442 lambda, 2, 6, 0, 3);
2443 else if ((hangingNodeType == InteractionVolume::fourSmallCellsEdge
2444 && (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2445 != problem().variables().index(interactionVolume.getSubVolumeElement(6))))
2446 || hangingNodeType == InteractionVolume::sixSmallCells)
2447 {
2448 useCases[0] = true;
2449 useCases[1] = false;
2450 useCases[2] = true;
2451 useCases[3] = false;
2452
2453 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2454 lambda, 2, 6, 0, 4, 3, 7, useCases);
2455 }
2456 else if (hangingNodeType == InteractionVolume::fourSmallCellsEdge
2457 && (problem().variables().index(interactionVolume.getSubVolumeElement(4))
2458 != problem().variables().index(interactionVolume.getSubVolumeElement(5))) )
2459 {
2460 useCases[0] = true;
2461 useCases[1] = false;
2462 useCases[2] = false;
2463 useCases[3] = true;
2464
2465 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2466 lambda, 2, 6, 0, 4, 3, 7, useCases);
2467 }
2468 else if (hangingNodeType == InteractionVolume::fourSmallCellsDiag)
2469 {
2470 useCases[0] = false;
2471 useCases[1] = true;
2472 useCases[2] = false;
2473 useCases[3] = true;
2474
2475 caseL = mpfal3DTransmissibilityCalculator_.transmissibility(T, interactionVolume,
2476 lambda, 2, 6, 0, 4, 3, 7, useCases);
2477 }
2478
2479 additionalT = T[0];
2480
2481 // determine cells for T-matrix entry 2 and 3
2482 switch (caseL)
2483 {
2484 case 1:
2485 {
2486 additional2 = interactionVolume.getSubVolumeElement(0);
2487 additional3 = interactionVolume.getSubVolumeElement(3);
2488 break;}
2489 case 2:
2490 {
2491 additional2 = interactionVolume.getSubVolumeElement(4);
2492 additional3 = interactionVolume.getSubVolumeElement(7);
2493 break;}
2494 case 3:
2495 {
2496 additional2 = interactionVolume.getSubVolumeElement(4);
2497 additional3 = interactionVolume.getSubVolumeElement(3);
2498 break;}
2499 case 4:
2500 {
2501 additional2 = interactionVolume.getSubVolumeElement(0);
2502 additional3 = interactionVolume.getSubVolumeElement(7);
2503 break;}
2504 }
2505 break;}
2506 }
2507
2508 //if necessary, correct for reversely calculated flux direction
2509 if(!properFluxDirection)
2510 {
2511 // mpfal3DTransmissibilityCalculator_.setTransChoiceThreshold(transChoiceThreshold_);
2512 // a) reverse direction
2513 additionalT *= -1;
2514 // b) swap cell I and J
2515 using std::swap;
2516 swap(additionalT[0], additionalT[1]);
2517 }
2518
2519 return caseL;
2520}
2521
2522
2523}//end namespace Dumux
2524#endif
Define some often used mathematical functions.
Simplifies writing multi-file VTK datasets.
Provides methods for transmissibility calculation in 3-d.
Defines the properties required for the adaptive sequential 2p2c models.
Interaction volume container for compositional adaptive 3-d (using MPFA L-method).
Finite volume 2p2c pressure model with multi-physics.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
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
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:631
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 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
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:79
Provides methods for transmissibility calculation in 3-d.
Definition: 3dtransmissibilitycalculator.hh:51
The finite volume model for the solution of the compositional pressure equation.
Definition: fv3dpressureadaptive.hh:80
void update()
Definition: fv3dpressureadaptive.hh:163
void adaptPressure()
Adapt primary variables vector after adapting the grid.
Definition: fv3dpressureadaptive.hh:214
bool enableMPFA
Enables mpfa on hanging nodes (on by default)
Definition: fv3dpressureadaptive.hh:303
int transmissibilityAdapter_(const IntersectionIterator &isIt, InteractionVolume &interactionVolume, const int &subVolumeFaceIdx, bool properFluxDirection, Element &additional2, Element &additional3, TransmissivityMatrix &additionalT)
Adapter to use the general implementation of the mpfa-l for the compositional models.
Definition: fv3dpressureadaptive.hh:1842
void initializeMatrix()
initializes the matrix to store the system of equations
Definition: fv3dpressureadaptive.hh:314
void initializeMatrixRowSize()
Initialize the row sizes of the sparse global matrix.
Definition: fv3dpressureadaptive.hh:330
FvMpfaL3dTransmissibilityCalculator< TypeTag > mpfal3DTransmissibilityCalculator_
The common implementation to calculate the Transmissibility with the mpfa-L-method.
Definition: fv3dpressureadaptive.hh:309
bool enableVolumeIntegral_
Calculates the volume integral (on by default)
Definition: fv3dpressureadaptive.hh:302
void initializeMatrixIndices()
Determine position of matrix entries.
Definition: fv3dpressureadaptive.hh:579
int computeTransmissibilities(const IntersectionIterator &, TransmissivityMatrix &, GlobalPosition &, int &, GlobalPosition &, int &)
Computes the transmissibility coefficients for the MPFA-l method in 3D.
Definition: fv3dpressureadaptive.hh:1349
void getMpfaFlux(const IntersectionIterator &, const CellData &)
Compute flux through an irregular interface using a mpfa method.
Definition: fv3dpressureadaptive.hh:805
void assemble(bool first)
function which assembles the system of equations to be solved
Definition: fv3dpressureadaptive.hh:644
FV3dPressure2P2CAdaptive(Problem &problem)
Constructs a FVPressure2P2C object.
Definition: fv3dpressureadaptive.hh:249
std::map< int, std::vector< int > > irregularCellMap_
Container to store all cell's Indice with a hanging node.
Definition: fv3dpressureadaptive.hh:301
void initialize(bool solveTwice=false)
Definition: fv3dpressureadaptive.hh:183
void get1pMpfaFlux(const IntersectionIterator &, const CellData &)
Compute single-phase flux through an irregular interface using a mpfa method.
Definition: fv3dpressureadaptive.hh:1142
InteractionVolumeContainer * interactionVolumesContainer_
A pointer to the adaptive interaction volumes container.
Definition: fv3dpressureadaptive.hh:307
void updateMaterialLaws(bool fromPostTimestep=false)
Definition: fv3dpressureadaptive.hh:1295
int maxInteractionVolumes
Maximum number of interaction volumes considered (4 by default)
Definition: fv3dpressureadaptive.hh:304
Interaction volume container for compositional adaptive 3-d (using MPFA L-method) Model.
Definition: fvmpfal3dinteractionvolumecontaineradaptive.hh:47
bool enableVolumeIntegral
Enables the volume integral of the pressure equation.
Definition: fvpressure.hh:183
Problem & problem_
Definition: fvpressure.hh:182
void updateMaterialLaws(bool postTimeStep=false)
Updates secondary variables.
Definition: fvpressurecompositional.hh:711
void update()
Compositional pressure solution routine: update estimate for secants, assemble, solve.
Definition: fvpressurecompositional.hh:131
The finite volume model for the solution of the compositional pressure equation.
Definition: fvpressuremultiphysics.hh:70
void updateMaterialLaws(bool postTimeStep=false)
constitutive functions are updated once if new concentrations are calculated and stored in the variab...
Definition: fvpressuremultiphysics.hh:792
static constexpr int pressureType
gives kind of pressure used ( , , )
Definition: fvpressuremultiphysics.hh:233
typename SolutionTypes::ElementMapper ElementMapper
Definition: fvpressuremultiphysics.hh:225
Class including the information of a 3d interaction volume of an adaptive MPFA L-method that does not...
Definition: linteractionvolume3dadaptive.hh:192
Type of the data container for one interaction volume.
Definition: porousmediumflow/sequential/cellcentered/mpfa/properties.hh:102
Type of the data container for one interaction volume.
Definition: porousmediumflow/sequential/cellcentered/mpfa/properties.hh:104
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:213
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
Properties for a MPFA method.