3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fvpressurecompositional.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_FVPRESSURECOMPOSITIONAL_HH
25#define DUMUX_FVPRESSURECOMPOSITIONAL_HH
26
27#include <cmath>
28#include <dune/common/float_cmp.hh>
29
30// dumux environment
31#include <dumux/common/math.hh>
36
38
39namespace Dumux {
67template<class TypeTag> class FVPressureCompositional
68: public FVPressure<TypeTag>
69{
70 //the model implementation
73
78
80
83
85 enum
86 {
87 dim = GridView::dimension, dimWorld = GridView::dimensionworld
88 };
89 enum
90 {
91 pw = Indices::pressureW,
92 pn = Indices::pressureN
93 };
94 enum
95 {
96 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
97 wCompIdx = Indices::wCompIdx, nCompIdx = Indices::nCompIdx,
98 contiWEqIdx = Indices::contiWEqIdx, contiNEqIdx = Indices::contiNEqIdx
99 };
100 enum
101 {
102 numPhases = getPropValue<TypeTag, Properties::NumPhases>(),
103 numComponents = getPropValue<TypeTag, Properties::NumComponents>()
104 };
105
106 // using declarations to abbreviate a dune class...
107 using Element = typename GridView::Traits::template Codim<0>::Entity;
108
109 // convenience shortcuts for Vectors/Matrices
110 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
111 using PhaseVector = Dune::FieldVector<Scalar, numPhases>;
112 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
113
114public:
115 //the variables object is initialized, non-compositional before and compositional after first pressure calculation
116 void initialMaterialLaws(bool compositional);
117
118 //initialization routine to prepare first timestep
119 void initialize(bool solveTwice = false);
120
130 void update()
131 {
132 //pre-transport to estimate update vector
133 Scalar dt_estimate = 0.;
134 Dune::dinfo << "secant guess"<< std::endl;
135 problem_.transportModel().update(problem_.timeManager().time(), dt_estimate, updateEstimate_, false);
136 //last argument false in update() makes sure that this is estimate and no "real" transport step
137
138 // if we just started a new episode, the TS size of the update Estimate is a better
139 // estimate then the size of the last time step
140 if(Dune::FloatCmp::eq<Scalar>(problem_.timeManager().time(), problem_.timeManager().episodeStartTime())
141 && problem_.timeManager().episodeIndex() > 1)
142 problem_.timeManager().setTimeStepSize(dt_estimate*getParam<Scalar>("Impet.CFLFactor"));
143
144 updateEstimate_ *= problem_.timeManager().timeStepSize();
145
146 problem_.pressureModel().assemble(false); Dune::dinfo << "pressure calculation"<< std::endl;
147 problem_.pressureModel().solve();
148
149 return;
150 }
151
152 //constitutive functions are initialized and stored in the variables object
153 void updateMaterialLaws(bool postTimeStep = false);
154
155 //numerical volume derivatives wrt changes in mass, pressure
156 void volumeDerivatives(const GlobalPosition&,const Element& ep);
157
170 template<class MultiWriter>
171 void addOutputVtkFields(MultiWriter &writer)
172 {
173 using ScalarSolutionType = typename GetProp<TypeTag, Properties::SolutionTypes>::ScalarSolution;
174 int size = problem_.gridView().size(0);
175 ScalarSolutionType *pressureW = writer.allocateManagedBuffer(size);
176 ScalarSolutionType *pressureN = writer.allocateManagedBuffer(size);
177 ScalarSolutionType *totalConcentration1 = writer.allocateManagedBuffer (size);
178 ScalarSolutionType *totalConcentration2 = writer.allocateManagedBuffer (size);
179 ScalarSolutionType *pc = writer.allocateManagedBuffer(size);
180 ScalarSolutionType *saturationW = writer.allocateManagedBuffer(size);
181 ScalarSolutionType *densityWetting = writer.allocateManagedBuffer(size);
182 ScalarSolutionType *densityNonwetting = writer.allocateManagedBuffer(size);
183 ScalarSolutionType *mobilityW = writer.allocateManagedBuffer (size);
184 ScalarSolutionType *mobilityNW = writer.allocateManagedBuffer (size);
185 ScalarSolutionType *massfraction1W = writer.allocateManagedBuffer (size);
186 ScalarSolutionType *massfraction1NW = writer.allocateManagedBuffer (size);
187 // numerical stuff
188 ScalarSolutionType *volErr = writer.allocateManagedBuffer (size);
189
190
191 for (int i = 0; i < size; i++)
192 {
193 // basic level 0 output
194 CellData& cellData = problem_.variables().cellData(i);
195 (*pressureW)[i] = cellData.pressure(wPhaseIdx);
196 (*pressureN)[i] = cellData.pressure(nPhaseIdx);
197 (*totalConcentration1)[i] = cellData.massConcentration(wCompIdx);
198 (*totalConcentration2)[i] = cellData.massConcentration(nCompIdx);
199 (*saturationW)[i] = cellData.saturation(wPhaseIdx);
200 // output standard secondary variables
201 if(problem_.vtkOutputLevel()>=1)
202 {
203 (*pc)[i] = cellData.capillaryPressure();
204 (*densityWetting)[i] = cellData.density(wPhaseIdx);
205 (*densityNonwetting)[i] = cellData.density(nPhaseIdx);
206 (*mobilityW)[i] = cellData.mobility(wPhaseIdx);
207 (*mobilityNW)[i] = cellData.mobility(nPhaseIdx);
208 (*massfraction1W)[i] = cellData.massFraction(wPhaseIdx,wCompIdx);
209 (*massfraction1NW)[i] = cellData.massFraction(nPhaseIdx,wCompIdx);
210 (*volErr)[i] = cellData.volumeError();
211 }
212 }
213 writer.attachCellData(*pressureW, "wetting pressure");
214 writer.attachCellData(*pressureN, "nonwetting pressure");
215 writer.attachCellData(*saturationW, "wetting saturation");
216 writer.attachCellData(*totalConcentration1, "C^w from cellData");
217 writer.attachCellData(*totalConcentration2, "C^n from cellData");
218 if(problem_.vtkOutputLevel()>=1)
219 {
220 writer.attachCellData(*pc, "capillary pressure");
221 writer.attachCellData(*densityWetting, "wetting density");
222 writer.attachCellData(*densityNonwetting, "nonwetting density");
223 writer.attachCellData(*mobilityW, "mobility w_phase");
224 writer.attachCellData(*mobilityNW, "mobility nw_phase");
225 std::ostringstream oss1, oss2;
226 oss1 << "mass fraction " << FluidSystem::componentName(0) << " in " << FluidSystem::phaseName(0) << "-phase";
227 writer.attachCellData(*massfraction1W, oss1.str());
228 oss2 << "mass fraction " << FluidSystem::componentName(0) << " in " << FluidSystem::phaseName(1) << "-phase";
229 writer.attachCellData(*massfraction1NW, oss2.str());
230 writer.attachCellData(*volErr, "volume Error");
231 }
232 // verbose output with numerical details
233 if(problem_.vtkOutputLevel()>=2)
234 {
235 // add debug stuff
236 ScalarSolutionType *errorCorrPtr = writer.allocateManagedBuffer (size);
237 ScalarSolutionType *dv_dpPtr = writer.allocateManagedBuffer (size);
238 ScalarSolutionType *dV_dC1Ptr = writer.allocateManagedBuffer (size);
239 ScalarSolutionType *dV_dC2Ptr = writer.allocateManagedBuffer (size);
240 ScalarSolutionType *updEstimate1 = writer.allocateManagedBuffer (size);
241 ScalarSolutionType *updEstimate2 = writer.allocateManagedBuffer (size);
242 for (int i = 0; i < size; i++)
243 {
244 CellData& cellData = problem_.variables().cellData(i);
245 (*errorCorrPtr)[i] = cellData.errorCorrection();
246 (*dv_dpPtr)[i] = cellData.dv_dp();
247 (*dV_dC1Ptr)[i] = cellData.dv(wCompIdx);
248 (*dV_dC2Ptr)[i] = cellData.dv(nCompIdx);
249 (*updEstimate1)[i] = updateEstimate_[0][i];
250 (*updEstimate2)[i] = updateEstimate_[1][i];
251 }
252 writer.attachCellData(*errorCorrPtr, "Error Correction");
253 writer.attachCellData(*dv_dpPtr, "dv_dp");
254 writer.attachCellData(*dV_dC1Ptr, "dV_dC1");
255 writer.attachCellData(*dV_dC2Ptr, "dV_dC2");
256 writer.attachCellData(*updEstimate1, "updEstimate comp 1");
257 writer.attachCellData(*updEstimate2, "updEstimate comp 2");
258 }
259
260 // very verbose output
261 if(problem_.vtkOutputLevel()>=3)
262 {
263 ScalarSolutionType *pressurePV = writer.allocateManagedBuffer(size);
264 ScalarSolutionType *viscosityWetting = writer.allocateManagedBuffer(size);
265 ScalarSolutionType *viscosityNonwetting = writer.allocateManagedBuffer(size);
266 // ScalarSolutionType *nun = writer.allocateManagedBuffer(size);
267 // ScalarSolutionType *nuw = writer.allocateManagedBuffer(size);
268 ScalarSolutionType *faceUpwindW = writer.allocateManagedBuffer(size);
269 ScalarSolutionType *faceUpwindN = writer.allocateManagedBuffer(size);
270 for (int i = 0; i < size; i++)
271 {
272 CellData& cellData = problem_.variables().cellData(i);
273 (*viscosityWetting)[i] = cellData.viscosity(wPhaseIdx);
274 (*viscosityNonwetting)[i] = cellData.viscosity(nPhaseIdx);
275 // (*nun)[i] = cellData.phaseMassFraction(nPhaseIdx);
276 // (*nuw)[i] = cellData.phaseMassFraction(wPhaseIdx);
277 (*faceUpwindW)[i] = 0;
278 (*faceUpwindN)[i] = 0;
279 // run thorugh all local face idx and collect upwind information
280 for(int fIdx = 0; fIdx<cellData.fluxData().size(); fIdx++)
281 {
282 if(cellData.isUpwindCell(fIdx, contiWEqIdx))
283 (*faceUpwindW)[i] += pow(10,static_cast<double>(3-fIdx));
284 if(cellData.isUpwindCell(fIdx, contiNEqIdx))
285 (*faceUpwindN)[i] += pow(10,static_cast<double>(3-fIdx));
286 }
287 }
288 // writer.attachCellData(*nun, "phase mass fraction n-phase");
289 // writer.attachCellData(*nuw, "phase mass fraction w-phase");
290 *pressurePV = this->pressure();
291 writer.attachCellData(*faceUpwindW, "isUpwind w-phase");
292 writer.attachCellData(*faceUpwindN, "isUpwind n-phase");
293 writer.attachCellData(*pressurePV, "pressure (Primary Variable");
294 writer.attachCellData(*viscosityWetting, "wetting viscosity");
295 writer.attachCellData(*viscosityNonwetting, "nonwetting viscosity");
296 }
297 return;
298 }
299
308 void initializationOutput(double pseudoTS = 0.)
309 {
310 std::cout << "Writing debug for current time step\n";
311 initializationOutputWriter_.beginWrite(problem_.timeManager().time() + pseudoTS);
312 asImp_().addOutputVtkFields(initializationOutputWriter_);
313
314#if DUNE_MINIMAL_DEBUG_LEVEL <= 2
315 int size_ = problem_.gridView().size(0);
316 // output porosity, permeability
317 Dune::BlockVector<Dune::FieldVector<double,1> > *poroPtr = initializationOutputWriter_.allocateManagedBuffer (size_);
318 Dune::BlockVector<Dune::FieldVector<double,1> > *permPtr = initializationOutputWriter_.allocateManagedBuffer (size_);
319
320 Dune::BlockVector<Dune::FieldVector<double,1> > poro_(0.), perm_(0.);
321 poro_.resize(size_); perm_.resize(size_);
322 // iterate over all elements of domain
323 for (const auto& element : elements(problem_.gridView()))
324 {
325 // get index
326 int eIdxGlobal = problem_.variables().index(element);
327 poro_[eIdxGlobal] = problem_.spatialParams().porosity(element);
328 perm_[eIdxGlobal] = problem_.spatialParams().intrinsicPermeability(element)[0][0];
329 }
330 *poroPtr = poro_;
331 *permPtr = perm_;
332 initializationOutputWriter_.attachCellData(*poroPtr, "porosity");
333 initializationOutputWriter_.attachCellData(*permPtr, "permeability");
334
335 if(problem_.vtkOutputLevel()>=3)
336 {
337 Dune::BlockVector<Dune::FieldVector<double,1> > *permPtrY = initializationOutputWriter_.allocateManagedBuffer (size_);
338 Dune::BlockVector<Dune::FieldVector<double,1> > *permPtrZ = initializationOutputWriter_.allocateManagedBuffer (size_);
339 Dune::BlockVector<Dune::FieldVector<double,1> > permY_(0.), permZ_(0.);
340 permY_.resize(size_); permZ_.resize(size_);
341 // iterate over all elements of domain
342 for (const auto& element : elements(problem_.gridView()))
343 {
344 // get index
345 int eIdxGlobal = problem_.variables().index(element);
346 if(dim >=2)
347 permY_[eIdxGlobal] = problem_.spatialParams().intrinsicPermeability(element)[1][1];
348 if(dim >=3)
349 permZ_[eIdxGlobal] = problem_.spatialParams().intrinsicPermeability(element)[2][2];
350 }
351 if(dim >=2)
352 {
353 *permPtrY = permY_;
354 initializationOutputWriter_.attachCellData(*permPtrY, "permeability Y");
355 }
356 if(dim >=3)
357 {
358 *permPtrZ = permZ_;
359 initializationOutputWriter_.attachCellData(*permPtrZ, "permeability Z");
360 }
361 }
362#endif
363
365 return;
366 }
368
373 FVPressureCompositional(Problem& problem) : FVPressure<TypeTag>(problem),
374 problem_(problem), initializationOutputWriter_(problem.gridView(),"initOutput2p2c"),
375 maxError_(0.0), incp_(1.0e1)
376 {
377 updateEstimate_.resize(getPropValue<TypeTag, Properties::NumPhases>());
378 for (int i=0; i<getPropValue<TypeTag, Properties::NumPhases>(); i++)
379 updateEstimate_[i].resize(problem.gridView().size(0));
380
381 ErrorTermFactor_ = getParam<Scalar>("Impet.ErrorTermFactor");
382 ErrorTermLowerBound_ = getParam<Scalar>("Impet.ErrorTermLowerBound", 0.2);
383 ErrorTermUpperBound_ = getParam<Scalar>("Impet.ErrorTermUpperBound");
384
385 if (pressureType == Indices::pressureGlobal)
386 {
387 DUNE_THROW(Dune::NotImplemented, "Global Pressure type not supported!");
388 }
389 }
390protected:
391 TransportSolutionType updateEstimate_;
392 Problem& problem_;
393
396
397 Scalar maxError_;
398 Scalar incp_;
403 static constexpr int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>();
404
405private:
406 Problem& problem()
407 {
408 return problem_;
409 }
410 const Problem& problem() const
411 {
412 return problem_;
413 }
414
416 Implementation &asImp_()
417 { return *static_cast<Implementation *>(this);}
418
420 const Implementation &asImp_() const
421 { return *static_cast<const Implementation *>(this);}
422};
423
424
433template<class TypeTag>
435{
436 // inizialize matrix, RHS, and condition Matrix
437 ParentType::initialize();
438
439 // initialguess: set saturations, determine visco and mobility for initial pressure equation
440 // at this moment, the pressure is unknown. Hence, dont regard compositional effects.
441 Dune::dinfo << "first saturation guess"<<std::endl; //=J: initialGuess()
442 problem_.pressureModel().initialMaterialLaws(false);
443 #if DUNE_MINIMAL_DEBUG_LEVEL <= 3
444 // this update produces fluxes on init pressures
445 Scalar dummy;
446 problem_.transportModel().update(0.,dummy, updateEstimate_, false);
447 initializationOutput();
448 #endif
449 Dune::dinfo << "first pressure guess"<<std::endl;
450 problem_.pressureModel().assemble(true);
451 problem_.pressureModel().solve();
452 #if DUNE_MINIMAL_DEBUG_LEVEL <= 3
453 initializationOutput(1e-6);
454 #endif
455 // update the compositional variables (hence true)
456 Dune::dinfo << "first guess for mass fractions"<<std::endl;//= J: transportInitial()
457 problem_.pressureModel().initialMaterialLaws(true);
458
459 // perform concentration update to determine secants
460 Dune::dinfo << "secant guess"<< std::endl;
461 Scalar dt_estimate = 0.;
462 problem_.transportModel().update(0., dt_estimate, updateEstimate_, false);
463 using std::min;
464 dt_estimate = min ( problem_.timeManager().timeStepSize(), dt_estimate);
465 //make sure the right time-step is used by all processes in the parallel case
466 if (problem_.gridView().comm().size() > 1)
467 dt_estimate = problem_.gridView().comm().min(dt_estimate);
468
469 updateEstimate_ *= dt_estimate;
470 //communicate in parallel case
471// problem_.variables().communicateUpdateEstimate();
472 #if DUNE_MINIMAL_DEBUG_LEVEL <= 3
473 initializationOutput(2e-6);
474 #endif
475 // pressure calculation
476 Dune::dinfo << "second pressure guess"<<std::endl;
477 problem_.pressureModel().assemble(false);
478 problem_.pressureModel().solve();
479 #if DUNE_MINIMAL_DEBUG_LEVEL <= 3
480 initializationOutput(3e-6);
481 #endif
482
483 // update the compositional variables
484 problem_.pressureModel().initialMaterialLaws(true);
485
486
487 if (solveTwice)
488 {
489 Dune::BlockVector<Dune::FieldVector<Scalar, 1> > pressureOld(this->pressure());
490 Dune::BlockVector<Dune::FieldVector<Scalar, 1> > pressureDiff;
491 Scalar pressureNorm = 1.; //dummy initialization to perform at least 1 iteration
492 int numIter = 1;
493
494 while (pressureNorm > 1e-5 && numIter < 10)
495 {
496 Scalar dt_dummy=0.; // without this dummy, it will never converge!
497 // update for secants
498 problem_.transportModel().update(0., dt_dummy, updateEstimate_, false); Dune::dinfo << "secant guess"<< std::endl;
499 updateEstimate_ *= dt_estimate;
500
501 //communicate in parallel case
502// problem_.variables().communicateUpdateEstimate();
503
504 // pressure calculation
505 problem_.pressureModel().assemble(false); Dune::dinfo << "pressure guess number "<< numIter <<std::endl;
506 problem_.pressureModel().solve();
507 // update the compositional variables
508 problem_.pressureModel().initialMaterialLaws(true);
509
510 pressureDiff = pressureOld;
511 pressureDiff -= this->pressure();
512 pressureNorm = pressureDiff.infinity_norm();
513 pressureOld = this->pressure();
514 pressureNorm /= pressureOld.infinity_norm();
515
516 numIter++;
517 }
518 }
519 return;
520}
521
530template<class TypeTag>
532{
533 // iterate through leaf grid an evaluate c0 at cell center
534 for (const auto& element : elements(problem_.gridView()))
535 {
536 // get global coordinate of cell center
537 GlobalPosition globalPos = element.geometry().center();
538
539 // assign an Index for convenience
540 int eIdxGlobal = problem_.variables().index(element);
541
542 // get the temperature
543 Scalar temperature_ = problem_.temperatureAtPos(globalPos);
544 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
545 // acess the fluid state and prepare for manipulation
546 FluidState& fluidState = cellData.manipulateFluidState();
548
549 // old material law interface is deprecated: Replace this by
550 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
551 // after the release of 3.3, when the deprecated interface is no longer supported
552 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
553
554 // initial conditions
555 PhaseVector pressure(0.);
556 Scalar sat_0=0.;
557
558 typename Indices::BoundaryFormulation icFormulation;
559 problem_.initialFormulation(icFormulation, element); // get type of initial condition
560
561 if(!compositional) //means that we do the first approximate guess without compositions
562 {
563 // phase pressures are unknown, so start with an exemplary
564 Scalar exemplaryPressure = problem_.referencePressure(element);
565 pressure[wPhaseIdx] = pressure[nPhaseIdx] = this->pressure()[eIdxGlobal] = exemplaryPressure;
566 if (icFormulation == Indices::saturation) // saturation initial condition
567 {
568 sat_0 = problem_.initSat(element);
569 flashSolver.saturationFlash2p2c(fluidState, sat_0, pressure, temperature_);
570 }
571 else if (icFormulation == Indices::concentration) // concentration initial condition
572 {
573 Scalar Z0 = problem_.initConcentration(element);
574 flashSolver.concentrationFlash2p2c(fluidState, Z0, pressure, temperature_);
575 }
576 }
577 else if(compositional) //means we regard compositional effects since we know an estimate pressure field
578 {
579 if (icFormulation == Indices::saturation) // saturation initial condition
580 {
581 //get saturation, determine pc
582 sat_0 = problem_.initSat(element);
583 Scalar pc=0.;
584 if(getPropValue<TypeTag, Properties::EnableCapillarity>())
585 {
586 pc = fluidMatrixInteraction.pc(sat_0);
587 }
588 else
589 pc = 0.;
590
591 //determine phase pressures from primary pressure variable
592 switch (pressureType)
593 {
594 case pw:
595 {
596 pressure[wPhaseIdx] = this->pressure()[eIdxGlobal];
597 pressure[nPhaseIdx] = this->pressure()[eIdxGlobal] +pc;
598 break;
599 }
600 case pn:
601 {
602 pressure[wPhaseIdx] = this->pressure()[eIdxGlobal]-pc;
603 pressure[nPhaseIdx] = this->pressure()[eIdxGlobal];
604 break;
605 }
606 }
607 flashSolver.saturationFlash2p2c(fluidState, sat_0, pressure, temperature_);
608 }
609 else if (icFormulation == Indices::concentration) // concentration initial condition
610 {
611 Scalar Z0 = problem_.initConcentration(element);
612 // If total concentrations are given at the boundary, saturation is unknown.
613 // This may affect pc and hence p_alpha and hence again saturation -> iteration.
614
615 // iterations in case of enabled capillary pressure
616 if(getPropValue<TypeTag, Properties::EnableCapillarity>())
617 {
618 //start with pc from last TS
619 Scalar pc(cellData.capillaryPressure());
620
621 int maxiter = 3;
622 //start iteration loop
623 for(int iter=0; iter < maxiter; iter++)
624 {
625 //determine phase pressures from primary pressure variable
626 switch (pressureType)
627 {
628 case pw:
629 {
630 pressure[wPhaseIdx] = this->pressure()[eIdxGlobal];
631 pressure[nPhaseIdx] = this->pressure()[eIdxGlobal] + pc;
632 break;
633 }
634 case pn:
635 {
636 pressure[wPhaseIdx] = this->pressure()[eIdxGlobal] - pc;
637 pressure[nPhaseIdx] = this->pressure()[eIdxGlobal];
638 break;
639 }
640 }
641
642 //store old pc
643 Scalar oldPc = pc;
644 //update with better pressures
645 flashSolver.concentrationFlash2p2c(fluidState, Z0, pressure, problem_.temperatureAtPos(globalPos));
646 pc = fluidMatrixInteraction.pc(fluidState.saturation(wPhaseIdx));
647 // TODO: get right criterion, do output for evaluation
648 //converge criterion
649 using std::abs;
650 if (abs(oldPc - pc) < 10.0)
651 iter = maxiter;
652
653 pc = fluidMatrixInteraction.pc(fluidState.saturation(wPhaseIdx));
654 }
655 }
656 else // capillary pressure neglected
657 {
658 pressure[wPhaseIdx] = pressure[nPhaseIdx]
659 = this->pressure()[eIdxGlobal];
660 flashSolver.concentrationFlash2p2c(fluidState, Z0, pressure, temperature_);
661 }
662 } //end conc initial condition
663 } //end compositional
664
665 cellData.calculateMassConcentration(problem_.spatialParams().porosity(element));
666
667 problem_.transportModel().totalConcentration(wCompIdx,eIdxGlobal) = cellData.massConcentration(wCompIdx);
668 problem_.transportModel().totalConcentration(nCompIdx,eIdxGlobal) = cellData.massConcentration(nCompIdx);
669
670 // initialize mobilities
671 cellData.setMobility(wPhaseIdx, fluidMatrixInteraction.krw(fluidState.saturation(wPhaseIdx))
672 / cellData.viscosity(wPhaseIdx));
673 cellData.setMobility(nPhaseIdx, fluidMatrixInteraction.krn(fluidState.saturation(wPhaseIdx))
674 / cellData.viscosity(nPhaseIdx));
675
676 // calculate perimeter used as weighting factor
677 if(!compositional)
678 {
679 // run through all intersections with neighbors
680 for (const auto& intersection : intersections(problem_.gridView(), element))
681 {
682 cellData.perimeter()
683 += intersection.geometry().volume();
684 }
685 cellData.globalIdx() = eIdxGlobal;
686
687 // set dv to zero to prevent output errors
688 cellData.dv_dp() = 0.;
689 cellData.dv(wPhaseIdx) = 0.;
690 cellData.dv(nPhaseIdx) = 0.;
691 }
692 // all
693 cellData.reset();
694 }
695 return;
696}
697
705template<class TypeTag>
707{
708 Scalar maxError = 0.;
709 // iterate through leaf grid an evaluate c0 at cell center
710 for (const auto& element : elements(problem().gridView()))
711 {
712 int eIdxGlobal = problem().variables().index(element);
713
714 CellData& cellData = problem().variables().cellData(eIdxGlobal);
715
716 asImp_().updateMaterialLawsInElement(element, postTimeStep);
717
718 using std::max;
719 maxError = max(maxError, fabs(cellData.volumeError()));
720 }
721 if (problem_.gridView().comm().size() > 1)
722 maxError_ = problem_.gridView().comm().max(maxError_);
723
724 maxError_ = maxError/problem().timeManager().timeStepSize();
725 return;
726}
727
740template<class TypeTag>
741void FVPressureCompositional<TypeTag>::volumeDerivatives(const GlobalPosition& globalPos, const Element& element)
742{
743 // cell index
744 int eIdxGlobal = problem_.variables().index(element);
745
746 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
747
748 // get cell temperature
749 Scalar temperature_ = cellData.temperature(wPhaseIdx);
750
751 // initialize an Fluid state and a flash solver
752 FluidState updFluidState;
754
755 /**********************************
756 * a) get necessary variables
757 **********************************/
758 //determine phase pressures for flash calculation
759 PhaseVector pressure(0.);
760 for(int phaseIdx = 0; phaseIdx< numPhases; phaseIdx++)
761 pressure[phaseIdx] = cellData.pressure(phaseIdx);
762
763 // mass of components inside the cell
764 ComponentVector mass(0.);
765 for(int compIdx = 0; compIdx< numComponents; compIdx++)
766 mass[compIdx] = cellData.massConcentration(compIdx);
767
768 // actual fluid volume
769 // see Fritz 2011 (Dissertation) eq.3.76
770 Scalar specificVolume(0.); // = \sum_{\alpha} \nu_{\alpha} / \rho_{\alpha}
771 for(int phaseIdx = 0; phaseIdx< numPhases; phaseIdx++)
772 specificVolume += cellData.phaseMassFraction(phaseIdx) / cellData.density(phaseIdx);
773 Scalar volalt = mass.one_norm() * specificVolume;
774// volalt = cellData.volumeError()+problem_.spatialParams().porosity(element);
775 // = \sum_{\kappa} C^{\kappa} + \sum_{\alpha} \nu_{\alpha} / \rho_{\alpha}
776
777 /**********************************
778 * b) define increments
779 **********************************/
780 // increments for numerical derivatives
781 ComponentVector massIncrement(0.);
782 for(int compIdx = 0; compIdx< numComponents; compIdx++)
783 {
784 massIncrement[compIdx] = updateEstimate_[compIdx][eIdxGlobal];
785 if(fabs(massIncrement[compIdx]) < 1e-8 * cellData.density(compIdx))
786 massIncrement[compIdx] = 1e-8* cellData.density(compIdx); // as phaseIdx = compIdx
787 }
788 Scalar& incp = incp_;
789
790 /**********************************
791 * c) Secant method for derivatives
792 **********************************/
793
794 // numerical derivative of fluid volume with respect to pressure
795 PhaseVector p_(incp);
796 p_ += pressure;
797 Scalar Z0 = mass[0] / mass.one_norm();
798 flashSolver.concentrationFlash2p2c(updFluidState, Z0, p_, temperature_);
799
800 specificVolume=0.; // = \sum_{\alpha} \nu_{\alpha} / \rho_{\alpha}
801 for(int phaseIdx = 0; phaseIdx< numPhases; phaseIdx++)
802 specificVolume += updFluidState.phaseMassFraction(phaseIdx) / updFluidState.density(phaseIdx);
803 Scalar dv_dp = ((mass.one_norm() * specificVolume) - volalt) /incp;
804
805 if (dv_dp>0)
806 {
807 // dV_dp > 0 is unphysical: Try inverse increment for secant
808 Dune::dinfo << "dv_dp larger 0 at Idx " << eIdxGlobal << " , try and invert secant"<< std::endl;
809
810 p_ -= 2*incp;
811 flashSolver.concentrationFlash2p2c(updFluidState, Z0, p_, temperature_);
812
813 specificVolume=0.; // = \sum_{\alpha} \nu_{\alpha} / \rho_{\alpha}
814 for(int phaseIdx = 0; phaseIdx< numPhases; phaseIdx++)
815 specificVolume += updFluidState.phaseMassFraction(phaseIdx) / updFluidState.density(phaseIdx);
816 dv_dp = ((mass.one_norm() * specificVolume) - volalt) /incp;
817
818 // dV_dp > 0 is unphysical: Try inverse increment for secant
819 if (dv_dp>0)
820 {
821 Dune::dwarn << "dv_dp still larger 0 after inverting secant at idx"<< eIdxGlobal<< std::endl;
822 p_ += 2*incp;
823 flashSolver.concentrationFlash2p2c(updFluidState, Z0, p_, temperature_);
824 // neglect effects of phase split, only regard changes in phase densities
825 specificVolume=0.; // = \sum_{\alpha} \nu_{\alpha} / \rho_{\alpha}
826 for(int phaseIdx = 0; phaseIdx< numPhases; phaseIdx++)
827 specificVolume += cellData.phaseMassFraction(phaseIdx) / updFluidState.density(phaseIdx);
828 dv_dp = ((mass.one_norm() * specificVolume) - volalt) /incp;
829 if (dv_dp>0)
830 {
831 std::cout << "dv_dp still larger 0 after both inverts at idx " << eIdxGlobal << std::endl;
832 dv_dp = cellData.dv_dp();
833 }
834 }
835 }
836 cellData.dv_dp()=dv_dp;
837
838 // numerical derivative of fluid volume with respect to mass of components
839 for (int compIdx = 0; compIdx<numComponents; compIdx++)
840 {
841 mass[compIdx] += massIncrement[compIdx];
842 Z0 = mass[0] / mass.one_norm();
843
844 flashSolver.concentrationFlash2p2c(updFluidState, Z0, pressure, temperature_);
845
846 specificVolume=0.; // = \sum_{\alpha} \nu_{\alpha} / \rho_{\alpha}
847 for(int phaseIdx = 0; phaseIdx< numPhases; phaseIdx++)
848 specificVolume += updFluidState.phaseMassFraction(phaseIdx) / updFluidState.density(phaseIdx);
849
850 cellData.dv(compIdx) = ((mass.one_norm() * specificVolume) - volalt) / massIncrement[compIdx];
851 mass[compIdx] -= massIncrement[compIdx];
852
853 //check routines if derivatives are meaningful
854 using std::isnan;
855 using std::isinf;
856 if (isnan(cellData.dv(compIdx)) || isinf(cellData.dv(compIdx)) )
857 {
858 DUNE_THROW(Dune::MathError, "NAN/inf of dV_dm. If that happens in first timestep, try smaller firstDt!");
859 }
860 }
861 cellData.volumeDerivativesAvailable(true);
862}
863
864}//end namespace Dumux
865#endif
Helpers for deprecation.
Define some often used mathematical functions.
Simplifies writing multi-file VTK datasets.
Determines the pressures and saturations of all fluid phases given the total mass of all components.
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 saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:146
Flash calculation routines for compositional sequential models.
Definition: compositionalflash.hh:51
static void saturationFlash2p2c(FluidState &fluidState, const Scalar &saturation, const PhaseVector &phasePressure, const Scalar &porosity, const Scalar &temperature)
Definition: compositionalflash.hh:219
static void concentrationFlash2p2c(FluidState &fluidState, const Scalar &Z0, const PhaseVector &phasePressure, const Scalar &porosity, const Scalar &temperature)
Definition: compositionalflash.hh:71
The finite volume model for the solution of the compositional pressure equation.
Definition: fvpressurecompositional.hh:69
Scalar incp_
Increment for the volume derivative w.r.t pressure.
Definition: fvpressurecompositional.hh:398
void volumeDerivatives(const GlobalPosition &, const Element &ep)
Partial derivatives of the volumes w.r.t. changes in total concentration and pressure.
Definition: fvpressurecompositional.hh:741
Scalar ErrorTermLowerBound_
Handling of error term: lower bound for error dampening.
Definition: fvpressurecompositional.hh:400
void initializationOutput(double pseudoTS=0.)
Write additional debug info in a special writer.
Definition: fvpressurecompositional.hh:308
Scalar maxError_
Maximum volume error of all cells.
Definition: fvpressurecompositional.hh:397
void addOutputVtkFields(MultiWriter &writer)
Write data files.
Definition: fvpressurecompositional.hh:171
void updateMaterialLaws(bool postTimeStep=false)
Updates secondary variables.
Definition: fvpressurecompositional.hh:706
Problem & problem_
Definition: fvpressurecompositional.hh:392
FVPressureCompositional(Problem &problem)
Constructs a FVPressureCompositional object.
Definition: fvpressurecompositional.hh:373
Scalar ErrorTermFactor_
Handling of error term: relaxation factor.
Definition: fvpressurecompositional.hh:399
static constexpr int pressureType
gives kind of pressure used ( , , )
Definition: fvpressurecompositional.hh:403
Scalar ErrorTermUpperBound_
Definition: fvpressurecompositional.hh:401
void initialMaterialLaws(bool compositional)
initializes the fluid distribution and hereby the variables container
Definition: fvpressurecompositional.hh:531
VtkMultiWriter< GridView > initializationOutputWriter_
output for the initialization procedure
Definition: fvpressurecompositional.hh:395
void update()
Compositional pressure solution routine: update estimate for secants, assemble, solve.
Definition: fvpressurecompositional.hh:130
TransportSolutionType updateEstimate_
Update estimate for changes in volume for the pressure equation.
Definition: fvpressurecompositional.hh:391
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
Defines the properties required for the sequential 2p2c models.
Finite Volume Diffusion Model.