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