3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fvpressuremultiphysics.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_FVPRESSURE2P2C_MULTIPHYSICS_HH
25#define DUMUX_FVPRESSURE2P2C_MULTIPHYSICS_HH
26
27#include <dune/common/float_cmp.hh>
28
29// dumux environment
33
34namespace Dumux {
68template<class TypeTag>
70{
76
78 using MaterialLaw = typename SpatialParams::MaterialLaw;
79
82
85
87 enum
88 {
89 dim = GridView::dimension, dimWorld = GridView::dimensionworld
90 };
91 enum
92 {
93 pw = Indices::pressureW
94 };
95 enum
96 {
97 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
98 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
99 contiWEqIdx = Indices::contiWEqIdx, contiNEqIdx = Indices::contiNEqIdx
100 };
101
102 // using declarations to abbreviate several dune classes...
103 using Element = typename GridView::Traits::template Codim<0>::Entity;
104 using Grid = typename GridView::Grid;
105 using Intersection = typename GridView::Intersection;
106
107 // convenience shortcuts for Vectors/Matrices
108 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
109 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
110 using PhaseVector = Dune::FieldVector<Scalar, getPropValue<TypeTag, Properties::NumPhases>()>;
112
114 using EntryType = Dune::FieldVector<Scalar, 2>;
116 Problem& problem()
117 { return this->problem_; }
118 const Problem& problem() const
119 { return this->problem_; }
120public:
121 //function which assembles the system of equations to be solved
122 void assemble(bool first);
123
124 void get1pSource(EntryType& sourceEntry, const Element& elementI, const CellData& cellDataI);
125
126 void get1pStorage(EntryType& storageEntry, const Element& elementI, CellData& cellDataI);
127
128 void get1pFlux(EntryType& entries, const Intersection& intersection,
129 const CellData& cellDataI);
130
131 void get1pFluxOnBoundary(EntryType& entries,
132 const Intersection& intersection,
133 const CellData& cellDataI);
134
135 //initialize multi-physics-specific pressure model constituents
136 void initialize(bool solveTwice = false)
137 {
138 // assign whole domain to most complex subdomain, => 2p
139 int size = this->problem().gridView().size(0);
140 for (int i = 0; i < size; i++)
141 {
142 CellData& cellData = this->problem().variables().cellData(i);
143 cellData.subdomain() = 2;
144 }
145 nextSubdomain.resize(size);
147 }
148
157 void serializeEntity(std::ostream &outstream, const Element &element)
158 {
159 ParentType::serializeEntity(outstream,element);
160 int eIdxGlobal = problem().variables().index(element);
161 CellData& cellData = problem().variables().cellData(eIdxGlobal);
162 outstream <<" "<< cellData.subdomain();
163 }
164
173 void deserializeEntity(std::istream &instream, const Element &element)
174 {
175 ParentType::deserializeEntity(instream,element);
176
177 int eIdxGlobal = problem().variables().index(element);
178 CellData& cellData = problem().variables().cellData(eIdxGlobal);
179 int subdomainIdx;
180 instream >> subdomainIdx;
181 cellData.setSubdomainAndFluidStateType(subdomainIdx);
182 }
183
184
185 //constitutive functions are initialized and stored in the variables object
186 void updateMaterialLaws(bool postTimeStep = false);
187 //updates singlephase secondary variables for one cell and stores in the variables object
188 void update1pMaterialLawsInElement(const Element& elementI, CellData& cellData, bool postTimeStep);
189
194 template<class MultiWriter>
195 void addOutputVtkFields(MultiWriter &writer)
196 {
198
199 if(problem().vtkOutputLevel()>=1)
200 {
201 int size = problem().gridView().size(0);
202 // add multiphysics stuff
203 Dune::BlockVector<Dune::FieldVector<int,1> >* subdomainPtr = writer.template allocateManagedBuffer<int, 1> (size);
204 for (int i = 0; i < size; i++)
205 {
206 CellData& cellData = problem().variables().cellData(i);
207 (*subdomainPtr)[i] = cellData.subdomain();
208 }
209 writer.attachCellData(*subdomainPtr, "subdomain");
210 }
211
212 return;
213 }
214
219 FVPressure2P2CMultiPhysics(Problem& problem) : FVPressure2P2C<TypeTag>(problem),
220 gravity_(problem.gravity()), timer_(false)
221 {}
222
223protected:
224 #if HAVE_MPI
225 using ElementMapper = typename SolutionTypes::ElementMapper;
227 #endif
228
229 // subdomain map
230 Dune::BlockVector<Dune::FieldVector<int,1> > nextSubdomain;
231 const GlobalPosition& gravity_;
233 static constexpr int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
234 Dune::Timer timer_;
235
243 enum
244 {
245 rhs = 1,
246 matrix = 0
247
248 };
249};
250
263template<class TypeTag>
265{
266 if(first)
267 {
268 ParentType::assemble(true);
269 return;
270 }
271 // initialization: set matrix this->A_ to zero
272 this->A_ = 0;
273 this->f_ = 0;
274
275 for (const auto& element : elements(problem().gridView()))
276 {
277 // get the global index of the cell
278 int eIdxGlobalI = problem().variables().index(element);
279
280 // assemble interior element contributions
281 if (element.partitionType() == Dune::InteriorEntity)
282 {
283 // get the cell data
284 CellData& cellDataI = problem().variables().cellData(eIdxGlobalI);
285
286 Dune::FieldVector<Scalar, 2> entries(0.);
287
288 /***** source term ***********/
289 if(cellDataI.subdomain() != 2)
290 problem().pressureModel().get1pSource(entries,element, cellDataI);
291 else
292 problem().pressureModel().getSource(entries,element, cellDataI, first);
293
294 this->f_[eIdxGlobalI] = entries[rhs];
295
296 /***** flux term ***********/
297 // iterate over all faces of the cell
298 for (const auto& intersection : intersections(problem().gridView(), element))
299 {
300 /************* handle interior face *****************/
301 if (intersection.neighbor())
302 {
303 int eIdxGlobalJ = problem().variables().index(intersection.outside());
304
305 if (cellDataI.subdomain() != 2
306 or problem().variables().cellData(eIdxGlobalJ).subdomain() != 2) // cell in the 1p domain
307 get1pFlux(entries, intersection, cellDataI);
308 else
309 problem().pressureModel().getFlux(entries, intersection, cellDataI, first);
310
311 //set right hand side
312 this->f_[eIdxGlobalI] -= entries[rhs];
313 // set diagonal entry
314 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
315 // set off-diagonal entry
316 this->A_[eIdxGlobalI][eIdxGlobalJ] = -entries[matrix];
317 } // end neighbor
318
319
320 /************* boundary face ************************/
321 else
322 {
323 if (cellDataI.subdomain() != 2) //the current cell in the 1p domain
324 problem().pressureModel().get1pFluxOnBoundary(entries, intersection, cellDataI);
325 else
326 problem().pressureModel().getFluxOnBoundary(entries, intersection, cellDataI, first);
327
328 //set right hand side
329 this->f_[eIdxGlobalI] += entries[rhs];
330 // set diagonal entry
331 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
332 }
333 } //end interfaces loop
334 // printmatrix(std::cout, this->A_, "global stiffness matrix", "row", 11, 3);
335
336 /***** storage term ***********/
337 if (cellDataI.subdomain() != 2) //the current cell in the 1p domain
338 problem().pressureModel().get1pStorage(entries, element, cellDataI);
339 else
340 problem().pressureModel().getStorage(entries, element, cellDataI, first);
341
342 this->f_[eIdxGlobalI] += entries[rhs];
343 // set diagonal entry
344 this->A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
345 }
346 // assemble overlap and ghost element contributions
347 else
348 {
349 this->A_[eIdxGlobalI] = 0.0;
350 this->A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
351 this->f_[eIdxGlobalI] = this->pressure()[eIdxGlobalI];
352 }
353 } // end grid traversal
354// printmatrix(std::cout, this->A_, "global stiffness matrix after assempling", "row", 11,3);
355// printvector(std::cout, this->f_, "right hand side", "row", 10);
356 return;
357}
358
370template<class TypeTag>
371void FVPressure2P2CMultiPhysics<TypeTag>::get1pSource(Dune::FieldVector<Scalar, 2>& sourceEntry,
372 const Element& elementI, const CellData& cellDataI)
373{
374 sourceEntry=0.;
375
376 // cell volume & perimeter, assume linear map here
377 Scalar volume = elementI.geometry().volume();
378 int subdomainIdx = cellDataI.subdomain();
379
380 /****************implement source************************/
381 PrimaryVariables source(NAN);
382 problem().source(source, elementI);
383 source[1+subdomainIdx] /= cellDataI.density(subdomainIdx);
384
385 sourceEntry[1] = volume * source[1+subdomainIdx];
386
387 return;
388}
389
403template<class TypeTag>
404void FVPressure2P2CMultiPhysics<TypeTag>::get1pStorage(Dune::FieldVector<Scalar, 2>& storageEntry,
405 const Element& elementI,
406 CellData& cellDataI)
407{
408 storageEntry = 0.;
409 // cell index
410 int eIdxGlobalI = problem().variables().index(elementI);
411 int presentPhaseIdx = cellDataI.subdomain();
412 Scalar volume = elementI.geometry().volume();
413
414 // determine maximum error to scale error-term
415 Scalar timestep_ = problem().timeManager().timeStepSize();
416
417 // compressibility term: 1p domain, so no dv_dp calculated
418 if (true)
419 {
420 Scalar& incp = this->incp_;
421
422 // numerical derivative of fluid volume with respect to pressure
423 PhaseVector p_(incp);
424 p_[nPhaseIdx] += cellDataI.pressure(nPhaseIdx);
425 p_[wPhaseIdx] += cellDataI.pressure(wPhaseIdx);
426
427 Scalar sumC = (cellDataI.massConcentration(wCompIdx) + cellDataI.massConcentration(nCompIdx));
428 Scalar Z0 = cellDataI.massConcentration(wCompIdx) / sumC;
429 // initialize simple fluidstate object
432 flashSolver.concentrationFlash1p2c(pseudoFluidState, Z0, p_, cellDataI.subdomain(),
433 cellDataI.temperature(wPhaseIdx));
434
435 Scalar v_ = 1. / pseudoFluidState.density(presentPhaseIdx);
436 cellDataI.dv_dp() = (sumC * ( v_ - (1. /cellDataI.density(presentPhaseIdx)))) /incp;
437
438 if (cellDataI.dv_dp()>0)
439 {
440 // dV_dp > 0 is unphysical: Try inverse increment for secant
441 Dune::dinfo << "dv_dp larger 0 at Idx " << eIdxGlobalI << " , try and invert secant"<< std::endl;
442
443 p_ -= 2*incp;
444 flashSolver.concentrationFlash1p2c(pseudoFluidState, Z0, p_, cellDataI.subdomain(),
445 cellDataI.temperature(wPhaseIdx));
446 v_ = 1. / pseudoFluidState.density(presentPhaseIdx);
447 cellDataI.dv_dp() = (sumC * ( v_ - (1. /cellDataI.density(presentPhaseIdx)))) /incp;
448 // dV_dp > 0 is unphysical: Try inverse increment for secant
449 if (cellDataI.dv_dp()>0)
450 {
451 Dune::dinfo <<__FILE__<< "dv_dp still larger 0 after inverting secant. regularize"<< std::endl;
452 cellDataI.dv_dp() *= -1;
453 }
454 }
455
456
457 Scalar compress_term = cellDataI.dv_dp() / timestep_;
458
459 storageEntry[matrix] -= compress_term*volume;
460 storageEntry[rhs] -= cellDataI.pressure(pressureType) * compress_term * volume;
461
462 using std::isnan;
463 using std::isinf;
464 if (isnan(compress_term) ||isinf(compress_term))
465 DUNE_THROW(Dune::MathError, "Compressibility term leads to NAN matrix entry at index " << eIdxGlobalI);
466
467 if(!getPropValue<TypeTag, Properties::EnableCompressibility>())
468 DUNE_THROW(Dune::NotImplemented, "Compressibility is switched off???");
469 }
470
471 // Abort error damping if there will be a possibly tiny timestep compared with last one
472 // This might be the case if the episode or simulation comes to an end.
473 if( problem().timeManager().episodeWillBeFinished()
474 || problem().timeManager().willBeFinished())
475 {
476 problem().variables().cellData(eIdxGlobalI).errorCorrection() = 0.;
477 return;
478 }
479
480 // error reduction routine: volumetric error is damped and inserted to right hand side
481 // if damping is not done, the solution method gets unstable!
482 problem().variables().cellData(eIdxGlobalI).volumeError() /= timestep_;
483 Scalar maxError = this->maxError_;
484 Scalar erri = fabs(cellDataI.volumeError());
485 Scalar x_lo = this->ErrorTermLowerBound_;
486 Scalar x_mi = this->ErrorTermUpperBound_;
487 Scalar fac = this->ErrorTermFactor_;
488 if (pressureType == pw)
489 fac = 0.1*this->ErrorTermFactor_;
490 Scalar lofac = 0.;
491 Scalar hifac = 0.;
492
493 if ((erri*timestep_ > 5e-5) && (erri > x_lo * maxError) && (!problem().timeManager().willBeFinished()))
494 {
495 if (erri <= x_mi * maxError)
496 storageEntry[rhs] +=
497 problem().variables().cellData(eIdxGlobalI).errorCorrection() =
498 fac* (1-x_mi*(lofac-1)/(x_lo-x_mi) + (lofac-1)/(x_lo-x_mi)*erri/maxError)
499 * cellDataI.volumeError() * volume;
500 else
501 storageEntry[rhs] +=
502 problem().variables().cellData(eIdxGlobalI).errorCorrection() =
503 fac * (1 + x_mi - hifac*x_mi/(1-x_mi) + (hifac/(1-x_mi)-1)*erri/maxError)
504 * cellDataI.volumeError() * volume;
505 }
506 else
507 problem().variables().cellData(eIdxGlobalI).errorCorrection()=0 ;
508
509 return;
510}
511
526template<class TypeTag>
527void FVPressure2P2CMultiPhysics<TypeTag>::get1pFlux(Dune::FieldVector<Scalar, 2>& entries,
528 const Intersection& intersection, const CellData& cellDataI)
529{
530 entries = 0.;
531 auto elementI = intersection.inside();
532
533 // get global coordinate of cell center
534 const GlobalPosition& globalPos = elementI.geometry().center();
535
536 // get absolute permeability
537 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
538
539 // get normal vector
540 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
541
542 // get face volume
543 Scalar faceArea = intersection.geometry().volume();
544
545 // access neighbor
546 auto neighbor = intersection.outside();
547 int eIdxGlobalJ = problem().variables().index(neighbor);
548 CellData& cellDataJ = problem().variables().cellData(eIdxGlobalJ);
549
550 // gemotry info of neighbor
551 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
552
553 // distance vector between barycenters
554 GlobalPosition distVec = globalPosNeighbor - globalPos;
555
556 // compute distance between cell centers
557 Scalar dist = distVec.two_norm();
558
559 GlobalPosition unitDistVec(distVec);
560 unitDistVec /= dist;
561
562 DimMatrix permeabilityJ
563 = problem().spatialParams().intrinsicPermeability(neighbor);
564
565 // compute vectorized permeabilities
566 DimMatrix meanPermeability(0);
567 harmonicMeanMatrix(meanPermeability, permeabilityI, permeabilityJ);
568
569 Dune::FieldVector<Scalar, dim> permeability(0);
570 meanPermeability.mv(unitDistVec, permeability);
571
572 // due to "safety cell" around subdomain, both cells I and J
573 // have single-phase conditions, although one is in 2p domain.
574 using std::min;
575 int phaseIdx = min(cellDataI.subdomain(), cellDataJ.subdomain());
576
577 Scalar rhoMean = 0.5 * (cellDataI.density(phaseIdx) + cellDataJ.density(phaseIdx));
578
579 // 1p => no pc => only 1 pressure, potential
580 Scalar potential = (cellDataI.pressure(phaseIdx) - cellDataJ.pressure(phaseIdx)) / dist;
581
582 potential += rhoMean * (unitDistVec * gravity_);
583
584 Scalar lambda;
585
586 if (potential > 0.)
587 {
588 lambda = cellDataI.mobility(phaseIdx);
589 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx, false); // store in cellJ since cellI is const
590 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx, false); // store in cellJ since cellI is const
591 }
592 else if (potential < 0.)
593 {
594 lambda = cellDataJ.mobility(phaseIdx);
595 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx, true);
596 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx, true);
597 }
598 else
599 {
600 lambda = harmonicMean(cellDataI.mobility(phaseIdx) , cellDataJ.mobility(phaseIdx));
601 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiWEqIdx, false);
602 cellDataJ.setUpwindCell(intersection.indexInOutside(), contiNEqIdx, false);
603 }
604
605 entries[0] = lambda * faceArea * fabs(permeability * unitOuterNormal) / (dist);
606 entries[1] = rhoMean * lambda;
607 entries[1] *= faceArea * fabs(permeability * unitOuterNormal) * (unitDistVec * gravity_);
608
609 return;
610}
611
629template<class TypeTag>
630void FVPressure2P2CMultiPhysics<TypeTag>::get1pFluxOnBoundary(Dune::FieldVector<Scalar, 2>& entries,
631 const Intersection& intersection, const CellData& cellDataI)
632{
633 entries = 0.;
634 // get global coordinate of cell center
635 auto elementI = intersection.inside();
636 const GlobalPosition& globalPos = elementI.geometry().center();
637// int eIdxGlobalI = problem().variables().index(elementI);
638 int phaseIdx = cellDataI.subdomain();
639
640 // get normal vector
641 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
642 // get face volume
643 Scalar faceArea = intersection.geometry().volume();
644
645 // center of face in global coordinates
646 const GlobalPosition& globalPosFace = intersection.geometry().center();
647
648 // geometrical information
649 GlobalPosition distVec(globalPosFace - globalPos);
650 Scalar dist = distVec.two_norm();
651 GlobalPosition unitDistVec(distVec);
652 unitDistVec /= dist;
653
654 //get boundary condition for boundary face center
655 BoundaryTypes bcType;
656 problem().boundaryTypes(bcType, intersection);
657
658 // prepare pressure boundary condition
659 PhaseVector pressBC(0.);
660
661 /********** Dirichlet Boundary *************/
662 if (bcType.isDirichlet(Indices::pressureEqIdx))
663 {
664 // get absolute permeability
665 DimMatrix permeabilityI(problem().spatialParams().intrinsicPermeability(elementI));
666 if(this->regulateBoundaryPermeability)
667 {
668 int axis = intersection.indexInInside() / 2;
669 if(permeabilityI[axis][axis] < this->minimalBoundaryPermeability)
670 permeabilityI[axis][axis] = this->minimalBoundaryPermeability;
671 }
672 // get mobilities and fractional flow factors
673 Scalar lambdaI = cellDataI.mobility(phaseIdx);
674
675 //permeability vector at boundary
676 Dune::FieldVector<Scalar, dim> permeability(0);
677 permeabilityI.mv(unitDistVec, permeability);
678
679 // create a fluid state for the boundary
680 FluidState BCfluidState;
681
682 //read boundary values
683 PrimaryVariables primaryVariablesOnBoundary(NAN);
684 problem().dirichlet(primaryVariablesOnBoundary, intersection);
685
686 {
687 // read boundary values
688 problem().transportModel().evalBoundary(globalPosFace,
689 intersection,
690 BCfluidState,
691 pressBC);
692
693 // determine fluid properties at the boundary
694 Scalar densityBound =
695 FluidSystem::density(BCfluidState, phaseIdx);
696 Scalar viscosityBound =
697 FluidSystem::viscosity(BCfluidState, phaseIdx);
698
699 // mobility at the boundary
700 Scalar lambdaBound = 0.;
701 switch (getPropValue<TypeTag, Properties::BoundaryMobility>())
702 {
703 case Indices::satDependent:
704 {
705 lambdaBound = BCfluidState.saturation(phaseIdx)
706 / viscosityBound;
707 break;
708 }
709 case Indices::permDependent:
710 {
711 if (phaseIdx == wPhaseIdx)
712 lambdaBound = MaterialLaw::krw(
713 problem().spatialParams().materialLawParams(elementI), BCfluidState.saturation(wPhaseIdx))
714 / viscosityBound;
715 else
716 lambdaBound = MaterialLaw::krn(
717 problem().spatialParams().materialLawParams(elementI), BCfluidState.saturation(wPhaseIdx))
718 / viscosityBound;
719 break;
720 }
721 }
722 Scalar rhoMean = 0.5 * (cellDataI.density(phaseIdx) + densityBound);
723
724 Scalar potential = 0;
725
726 //calculate potential gradient: pc = 0;
727 potential = (cellDataI.pressure(phaseIdx) - pressBC[phaseIdx]) / dist;
728
729 potential += rhoMean * (unitDistVec * gravity_);
730
731 Scalar lambda(0.);
732
733 if (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(potential, 0.0, 1.0e-30))
734 {
735 lambda = 0.5*(lambdaI + lambdaBound);
736 }
737 else if (potential > 0.)
738 {
739 lambda = lambdaI;
740 }
741 else
742 {
743 lambda = lambdaBound;
744 }
745
746 //calculate current matrix entry
747 Scalar entry(0.), rightEntry(0.);
748 entry = lambda * (fabs(permeability * unitOuterNormal) / dist) * faceArea;
749
750 //calculate right hand side
751 rightEntry = lambda * rhoMean * fabs(permeability * unitOuterNormal)
752 * faceArea ;
753
754 // set diagonal entry and right hand side entry
755 entries[0] += entry;
756 entries[1] += entry * primaryVariablesOnBoundary[Indices::pressureEqIdx];
757 entries[1] -= rightEntry * (gravity_ * unitDistVec);
758 } //end of if(first) ... else{...
759 } // end dirichlet
760
761 /**********************************
762 * set neumann boundary condition
763 **********************************/
764 else if(bcType.isNeumann(Indices::pressureEqIdx))
765 {
766 PrimaryVariables J(NAN);
767 problem().neumann(J, intersection);
768 J[1+phaseIdx] /= cellDataI.density(phaseIdx);
769
770 entries[1] -= J[1+phaseIdx] * faceArea;
771 }
772 else
773 DUNE_THROW(Dune::NotImplemented, "Boundary Condition neither Dirichlet nor Neumann!");
774
775
776 return;
777}
778
779
791template<class TypeTag>
793{
794 //get timestep for error term
795 Scalar maxError = 0.;
796
797 // next subdomain map
798 if (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(problem().timeManager().time(), 0.0, 1.0e-30))
799 nextSubdomain = 2; // start with complicated sub in initialization
800 else
801 nextSubdomain = -1; // reduce complexity after first TS
802
803 // Loop A) through leaf grid
804 for (const auto& element : elements(problem().gridView()))
805 {
806 // get global coordinate of cell center
807 int eIdxGlobal = problem().variables().index(element);
808 CellData& cellData = problem().variables().cellData(eIdxGlobal);
809
810 if(cellData.subdomain() == 2) // complex
811 {
812 this->updateMaterialLawsInElement(element, postTimeStep);
813
814 // check subdomain consistency
815 timer_.start();
816 // enshure we are not at source
817 // get sources from problem
818 PrimaryVariables source(NAN);
819 problem().source(source, element);
820
821 if ((cellData.saturation(wPhaseIdx) > 0.0 && cellData.saturation(wPhaseIdx) < 1.0)
822 || Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(source.one_norm(), 0.0, 1.0e-30)) // cell still 2p
823 {
824 // mark this element
825 nextSubdomain[eIdxGlobal] = 2;
826
827 // mark neighbors
828 for (const auto& intersection : intersections(problem().gridView(), element))
829 {
830 if (intersection.neighbor())
831 {
832 int eIdxGlobalJ = problem().variables().index(intersection.outside());
833 // mark neighbor Element
834 nextSubdomain[eIdxGlobalJ] = 2;
835 }
836 }
837 }
838 else if(nextSubdomain[eIdxGlobal] != 2)// update next subdomain if possible
839 {
840 if(Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellData.saturation(wPhaseIdx), 0.0, 1.0e-30))
841 nextSubdomain[eIdxGlobal] = wPhaseIdx;
842 else if (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(cellData.saturation(nPhaseIdx), 0.0, 1.0e-30))
843 nextSubdomain[eIdxGlobal] = nPhaseIdx;
844 }
845 timer_.stop();
846 // end subdomain check
847 }// end complex domain
848 else if (nextSubdomain[eIdxGlobal] != 2) //check if cell remains in simple subdomain
849 nextSubdomain[eIdxGlobal] = cellData.subdomain();
850
851 } //end define complex area of next subdomain
852
853 timer_.start();
854 //communicate next subdomain if parallel
855 #if HAVE_MPI
856 // communicate updated values
857 DataHandle dataHandle(problem().variables().elementMapper(), nextSubdomain);
858 problem().gridView().template communicate<DataHandle>(dataHandle,
859 Dune::InteriorBorder_All_Interface,
860 Dune::ForwardCommunication);
861 #endif
862
863 // Loop B) thorugh leaf grid
864 // investigate cells that were "simple" in current TS
865 for (const auto& element : elements(problem().gridView()))
866 {
867 int eIdxGlobal = problem().variables().index(element);
868 CellData& cellData = problem().variables().cellData(eIdxGlobal);
869
870 // store old subdomain information and assign new info
871 int oldSubdomainI = cellData.subdomain();
872 cellData.subdomain() = nextSubdomain[eIdxGlobal];
873
874 //first check if simple will become complicated
875 if(oldSubdomainI != 2
876 && nextSubdomain[eIdxGlobal] == 2)
877 {
878 // use complex update of the fluidstate
879 timer_.stop();
880 this->updateMaterialLawsInElement(element, postTimeStep);
881 timer_.start();
882 }
883 else if(oldSubdomainI != 2
884 && nextSubdomain[eIdxGlobal] != 2) // will be simple and was simple
885 {
886 // perform simple update
887 this->update1pMaterialLawsInElement(element, cellData, postTimeStep);
888 }
889 //else
890 // a) will remain complex -> everything already done in loop A
891 // b) will be simple and was complex: complex FS available, so next TS
892 // will use comlex FS, next updateMaterialLaw will transform to simple
893
894 using std::max;
895 maxError = max(maxError, fabs(cellData.volumeError()));
896 }// end grid traversal
897 this->maxError_ = maxError/problem().timeManager().timeStepSize();
898
899 timer_.stop();
900
901 if(problem().timeManager().willBeFinished() or problem().timeManager().episodeWillBeFinished())
902 Dune::dinfo << "Subdomain routines took " << timer_.elapsed() << " seconds" << std::endl;
903
904 return;
905}
906
917template<class TypeTag>
918void FVPressure2P2CMultiPhysics<TypeTag>::update1pMaterialLawsInElement(const Element& elementI, CellData& cellData, bool postTimeStep)
919{
920 // get global coordinate of cell center
921 GlobalPosition globalPos = elementI.geometry().center();
922 int eIdxGlobal = problem().variables().index(elementI);
923
924 // determine which phase should be present
925 int presentPhaseIdx = cellData.subdomain(); // this is already =nextSubomainIdx
926
927 // reset to calculate new timeStep if we are at the end of a time step
928 if(postTimeStep)
929 cellData.reset();
930
931 // acess the simple fluid state and prepare for manipulation
932 auto& pseudoFluidState = cellData.manipulateSimpleFluidState();
933
934 // prepare phase pressure for fluid state
935 // both phase pressures are necessary for the case 1p domain is assigned for
936 // the next 2p subdomain
937 PhaseVector pressure(0.);
938 Scalar pc = 0;
939 if(getPropValue<TypeTag, Properties::EnableCapillarity>())
940 pc = MaterialLaw::pc(problem().spatialParams().materialLawParams(elementI),
941 ((presentPhaseIdx == wPhaseIdx) ? 1. : 0.)); // assign sw = 1 if wPhase present, else 0
942 if(pressureType == wPhaseIdx)
943 {
944 pressure[wPhaseIdx] = this->pressure(eIdxGlobal);
945 pressure[nPhaseIdx] = this->pressure(eIdxGlobal)+pc;
946 }
947 else
948 {
949 pressure[wPhaseIdx] = this->pressure(eIdxGlobal)-pc;
950 pressure[nPhaseIdx] = this->pressure(eIdxGlobal);
951 }
952
953 // get the overall mass of first component: Z0 = C^0 / (C^0+C^1) [-]
954 Scalar sumConc = cellData.massConcentration(wCompIdx)
955 + cellData.massConcentration(nCompIdx);
956 Scalar Z0 = cellData.massConcentration(wCompIdx)/ sumConc;
957
959 flashSolver.concentrationFlash1p2c(pseudoFluidState, Z0, pressure, presentPhaseIdx, problem().temperatureAtPos(globalPos));
960
961 // write stuff in fluidstate
962 assert(presentPhaseIdx == pseudoFluidState.presentPhaseIdx());
963
964// cellData.setSimpleFluidState(pseudoFluidState);
965
966 // initialize viscosities
967 cellData.setViscosity(presentPhaseIdx, FluidSystem::viscosity(pseudoFluidState, presentPhaseIdx));
968
969 // initialize mobilities
970 if(presentPhaseIdx == wPhaseIdx)
971 {
972 cellData.setMobility(wPhaseIdx,
973 MaterialLaw::krw(problem().spatialParams().materialLawParams(elementI), pseudoFluidState.saturation(wPhaseIdx))
974 / cellData.viscosity(wPhaseIdx));
975 cellData.setMobility(nPhaseIdx, 0.);
976 }
977 else
978 {
979 cellData.setMobility(nPhaseIdx,
980 MaterialLaw::krn(problem().spatialParams().materialLawParams(elementI), pseudoFluidState.saturation(wPhaseIdx))
981 / cellData.viscosity(nPhaseIdx));
982 cellData.setMobility(wPhaseIdx, 0.);
983 }
984
985 // error term handling
986 Scalar vol(0.);
987 vol = sumConc / pseudoFluidState.density(presentPhaseIdx);
988
989 if (Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(problem().timeManager().timeStepSize(), 0.0, 1.0e-30))
990 cellData.volumeError() = (vol - problem().spatialParams().porosity(elementI));
991 return;
992}
993
994}//end namespace Dumux
995#endif
Contains a class to exchange entries of a vector.
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Finite volume 2p2c pressure model.
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
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 viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
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
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Flash calculation routines for compositional sequential models.
Definition: compositionalflash.hh:51
static void concentrationFlash1p2c(FluidState1p2c &fluidState, const Scalar &Z0, const Dune::FieldVector< Scalar, numPhases > phasePressure, const int presentPhaseIdx, const Scalar &temperature)
The simplest possible update routine for 1p2c "flash" calculations.
Definition: compositionalflash.hh:181
Container for compositional variables in a 1p2c situation.
Definition: pseudo1p2c.hh:44
Scalar density(int phaseIdx) const
Set the density of a phase .
Definition: pseudo1p2c.hh:109
Scalar temperature(int phaseIdx) const
Returns the temperature of the fluids .
Definition: pseudo1p2c.hh:189
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:79
The finite volume model for the solution of the compositional pressure equation.
Definition: fvpressure.hh:73
Problem & problem_
Definition: fvpressure.hh:182
Dune::FieldVector< Scalar, 2 > EntryType
Definition: fvpressure.hh:136
void addOutputVtkFields(MultiWriter &writer)
Write data files.
Definition: fvpressurecompositional.hh:172
The finite volume model for the solution of the compositional pressure equation.
Definition: fvpressuremultiphysics.hh:70
Dune::BlockVector< Dune::FieldVector< int, 1 > > nextSubdomain
vector holding next subdomain
Definition: fvpressuremultiphysics.hh:230
void get1pFluxOnBoundary(EntryType &entries, const Intersection &intersection, const CellData &cellDataI)
The compositional single-phase flux in the multiphysics framework.
Definition: fvpressuremultiphysics.hh:630
const GlobalPosition & gravity_
Definition: fvpressuremultiphysics.hh:231
void update1pMaterialLawsInElement(const Element &elementI, CellData &cellData, bool postTimeStep)
updates secondary variables of one single phase cell
Definition: fvpressuremultiphysics.hh:918
void assemble(bool first)
function which assembles the system of equations to be solved
Definition: fvpressuremultiphysics.hh:264
void addOutputVtkFields(MultiWriter &writer)
Write data files.
Definition: fvpressuremultiphysics.hh:195
FVPressure2P2CMultiPhysics(Problem &problem)
Constructs a FVPressure2P2CPC object.
Definition: fvpressuremultiphysics.hh:219
void get1pStorage(EntryType &storageEntry, const Element &elementI, CellData &cellDataI)
Assembles the storage term for a 1p cell in a multiphysics framework.
Definition: fvpressuremultiphysics.hh:404
void updateMaterialLaws(bool postTimeStep=false)
constitutive functions are updated once if new concentrations are calculated and stored in the variab...
Definition: fvpressuremultiphysics.hh:792
@ rhs
index for the right hand side entry
Definition: fvpressuremultiphysics.hh:245
@ matrix
index for the global matrix entry
Definition: fvpressuremultiphysics.hh:246
void get1pSource(EntryType &sourceEntry, const Element &elementI, const CellData &cellDataI)
Assembles the source term.
Definition: fvpressuremultiphysics.hh:371
Dune::Timer timer_
A timer for the time spent on the multiphysics framework.
Definition: fvpressuremultiphysics.hh:234
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the pressure field.
Definition: fvpressuremultiphysics.hh:157
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the pressure field.
Definition: fvpressuremultiphysics.hh:173
static constexpr int pressureType
gives kind of pressure used ( , , )
Definition: fvpressuremultiphysics.hh:233
void initialize(bool solveTwice=false)
Definition: fvpressuremultiphysics.hh:136
typename SolutionTypes::ElementMapper ElementMapper
Definition: fvpressuremultiphysics.hh:225
void get1pFlux(EntryType &entries, const Intersection &intersection, const CellData &cellDataI)
The compositional single-phase flux in the multiphysics framework.
Definition: fvpressuremultiphysics.hh:527
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
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the pressure field.
Definition: sequential/cellcentered/pressure.hh:264
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the pressure field.
Definition: sequential/cellcentered/pressure.hh:251