lensproblem.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_LENS_PROBLEM_HH
29#define EWOMS_LENS_PROBLEM_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/version.hh>
34
35#include <opm/material/components/Dnapl.hpp>
36#include <opm/material/components/SimpleH2O.hpp>
37#include <opm/material/fluidmatrixinteractions/RegularizedVanGenuchten.hpp>
38#include <opm/material/fluidmatrixinteractions/LinearMaterial.hpp>
39#include <opm/material/fluidmatrixinteractions/EffToAbsLaw.hpp>
40#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
41#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
42#include <opm/material/fluidsystems/TwoPhaseImmiscibleFluidSystem.hpp>
43
46
49
51
54
55#include <iostream>
56#include <sstream>
57#include <string>
58
59namespace Opm {
60template <class TypeTag>
61class LensProblem;
62}
63
64namespace Opm::Properties {
65
66// Create new type tags
67namespace TTag {
68struct LensBaseProblem { using InheritsFrom = std::tuple<StructuredGridVanguard>; };
69} // end namespace TTag
70
71// Set the problem property
72template<class TypeTag>
73struct Problem<TypeTag, TTag::LensBaseProblem> { using type = Opm::LensProblem<TypeTag>; };
74
75// Use Dune-grid's YaspGrid
76template<class TypeTag>
77struct Grid<TypeTag, TTag::LensBaseProblem> { using type = Dune::YaspGrid<2>; };
78
79// Set the wetting phase
80template<class TypeTag>
81struct WettingPhase<TypeTag, TTag::LensBaseProblem>
82{
83private:
85
86public:
87 using type = Opm::LiquidPhase<Scalar, Opm::SimpleH2O<Scalar> >;
88};
89
90// Set the non-wetting phase
91template<class TypeTag>
92struct NonwettingPhase<TypeTag, TTag::LensBaseProblem>
93{
94private:
96
97public:
98 using type = Opm::LiquidPhase<Scalar, Opm::DNAPL<Scalar> >;
99};
100
101// Set the material Law
102template<class TypeTag>
103struct MaterialLaw<TypeTag, TTag::LensBaseProblem>
104{
105private:
107 enum { wettingPhaseIdx = FluidSystem::wettingPhaseIdx };
108 enum { nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx };
109
111 using Traits = Opm::TwoPhaseMaterialTraits<Scalar,
112 /*wettingPhaseIdx=*/FluidSystem::wettingPhaseIdx,
113 /*nonWettingPhaseIdx=*/FluidSystem::nonWettingPhaseIdx>;
114
115 // define the material law which is parameterized by effective
116 // saturations
117 using EffectiveLaw = Opm::RegularizedVanGenuchten<Traits>;
118
119public:
120 // define the material law parameterized by absolute saturations
121 using type = Opm::EffToAbsLaw<EffectiveLaw>;
122};
123
124} // namespace Opm::Properties
125
126namespace Opm::Parameters {
127
128// define the properties specific for the lens problem
129template<class Scalar>
130struct LensLowerLeftX { static constexpr Scalar value = 1.0; };
131
132template<class Scalar>
133struct LensLowerLeftY { static constexpr Scalar value = 2.0; };
134
135template<class Scalar>
136struct LensLowerLeftZ { static constexpr Scalar value = 0.0; };
137
138template<class Scalar>
139struct LensUpperRightX { static constexpr Scalar value = 4.0; };
140
141template<class Scalar>
142struct LensUpperRightY { static constexpr Scalar value = 3.0; };
143
144template<class Scalar>
145struct LensUpperRightZ { static constexpr Scalar value = 1.0; };
146
147} // namespace Opm::Parameters
148
149namespace Opm {
150
174template <class TypeTag>
175class LensProblem : public GetPropType<TypeTag, Properties::BaseProblem>
176{
178
188
189 enum {
190 // number of phases
191 numPhases = FluidSystem::numPhases,
192
193 // phase indices
194 wettingPhaseIdx = FluidSystem::wettingPhaseIdx,
195 nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx,
196
197 // equation indices
198 contiNEqIdx = Indices::conti0EqIdx + nonWettingPhaseIdx,
199
200 // Grid and world dimension
201 dim = GridView::dimension,
202 dimWorld = GridView::dimensionworld
203 };
204
210
211 using CoordScalar = typename GridView::ctype;
212 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
213
214 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
215
216public:
220 LensProblem(Simulator& simulator)
221 : ParentType(simulator)
222 { }
223
228 {
229 ParentType::finishInit();
230
231 eps_ = 3e-6;
232 FluidSystem::init();
233
234 temperature_ = 273.15 + 20; // -> 20°C
235 lensLowerLeft_[0] = Parameters::Get<Parameters::LensLowerLeftX<Scalar>>();
236 lensLowerLeft_[1] = Parameters::Get<Parameters::LensLowerLeftY<Scalar>>();
237 lensUpperRight_[0] = Parameters::Get<Parameters::LensUpperRightX<Scalar>>();
238 lensUpperRight_[1] = Parameters::Get<Parameters::LensUpperRightY<Scalar>>();
239
240 if constexpr (dim == 3) {
241 lensLowerLeft_[2] = Parameters::Get<Parameters::LensLowerLeftZ<Scalar>>();
242 lensUpperRight_[2] = Parameters::Get<Parameters::LensUpperRightZ<Scalar>>();
243 }
244
245 // residual saturations
246 lensMaterialParams_.setResidualSaturation(wettingPhaseIdx, 0.18);
247 lensMaterialParams_.setResidualSaturation(nonWettingPhaseIdx, 0.0);
248 outerMaterialParams_.setResidualSaturation(wettingPhaseIdx, 0.05);
249 outerMaterialParams_.setResidualSaturation(nonWettingPhaseIdx, 0.0);
250
251 // parameters for the Van Genuchten law: alpha and n
252 lensMaterialParams_.setVgAlpha(0.00045);
253 lensMaterialParams_.setVgN(7.3);
254 outerMaterialParams_.setVgAlpha(0.0037);
255 outerMaterialParams_.setVgN(4.7);
256
257 lensMaterialParams_.finalize();
258 outerMaterialParams_.finalize();
259
260 lensK_ = this->toDimMatrix_(9.05e-12);
261 outerK_ = this->toDimMatrix_(4.6e-10);
262
263 if (dimWorld == 3) {
264 this->gravity_ = 0;
265 this->gravity_[1] = -9.81;
266 }
267 }
268
272 static void registerParameters()
273 {
274 ParentType::registerParameters();
275
276 Parameters::Register<Parameters::LensLowerLeftX<Scalar>>
277 ("The x-coordinate of the lens' lower-left corner [m].");
278 Parameters::Register<Parameters::LensLowerLeftY<Scalar>>
279 ("The y-coordinate of the lens' lower-left corner [m].");
280 Parameters::Register<Parameters::LensUpperRightX<Scalar>>
281 ("The x-coordinate of the lens' upper-right corner [m].");
282 Parameters::Register<Parameters::LensUpperRightY<Scalar>>
283 ("The y-coordinate of the lens' upper-right corner [m].");
284
285 if constexpr (dim == 3) {
286 Parameters::Register<Parameters::LensLowerLeftZ<Scalar>>
287 ("The z-coordinate of the lens' lower-left corner [m].");
288 Parameters::Register<Parameters::LensUpperRightZ<Scalar>>
289 ("The z-coordinate of the lens' upper-right corner [m].");
290 }
291
292 Parameters::SetDefault<Parameters::CellsX>(48);
293 Parameters::SetDefault<Parameters::CellsY>(32);
294 Parameters::SetDefault<Parameters::DomainSizeX<Scalar>>(6.0);
295 Parameters::SetDefault<Parameters::DomainSizeY<Scalar>>(4.0);
296
297 if constexpr (dim == 3) {
298 Parameters::SetDefault<Parameters::CellsZ>(16);
299 Parameters::SetDefault<Parameters::DomainSizeZ<Scalar>>(1.0);
300 }
301
302 // Use forward differences
304 constexpr bool useFD = std::is_same_v<LLS, Properties::TTag::FiniteDifferenceLocalLinearizer>;
305 if constexpr (useFD) {
306 Parameters::SetDefault<Parameters::NumericDifferenceMethod>(+1);
307 }
308
309 Parameters::SetDefault<Parameters::EndTime<Scalar>>(30e3);
310 Parameters::SetDefault<Parameters::EnableIntensiveQuantityCache>(true);
311 Parameters::SetDefault<Parameters::EnableStorageCache>(true);
312 Parameters::SetDefault<Parameters::InitialTimeStepSize<Scalar>>(250.0);
313 Parameters::SetDefault<Parameters::VtkWriteIntrinsicPermeabilities>(true);
314 Parameters::SetDefault<Parameters::EnableGravity>(true);
315 }
316
320 static std::string briefDescription()
321 {
322 std::string thermal = "isothermal";
323 constexpr bool enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>();
324 if constexpr (enableEnergy)
325 thermal = "non-isothermal";
326
327 std::string deriv = "finite difference";
329 constexpr bool useAutoDiff = std::is_same_v<LLS, Properties::TTag::AutoDiffLocalLinearizer>;
330 if constexpr (useAutoDiff) {
331 deriv = "automatic differentiation";
332 }
333
334 std::string disc = "vertex centered finite volume";
336 constexpr bool useEcfv = std::is_same<D, Opm::EcfvDiscretization<TypeTag>>::value;
337 if constexpr (useEcfv)
338 disc = "element centered finite volume";
339
340 return std::string("")+
341 "Ground remediation problem where a dense oil infiltrates "+
342 "an aquifer with an embedded low-permability lens. " +
343 "This is the binary for the "+thermal+" variant using "+deriv+
344 "and the "+disc+" discretization";
345 }
346
351
355 template <class Context>
356 const DimMatrix& intrinsicPermeability(const Context& context, unsigned spaceIdx,
357 unsigned timeIdx) const
358 {
359 const GlobalPosition& globalPos = context.pos(spaceIdx, timeIdx);
360
361 if (isInLens_(globalPos))
362 return lensK_;
363 return outerK_;
364 }
365
369 template <class Context>
370 Scalar porosity(const Context& /*context*/,
371 unsigned /*spaceIdx*/,
372 unsigned /*timeIdx*/) const
373 { return 0.4; }
374
378 template <class Context>
379 const MaterialLawParams& materialLawParams(const Context& context,
380 unsigned spaceIdx, unsigned timeIdx) const
381 {
382 const GlobalPosition& globalPos = context.pos(spaceIdx, timeIdx);
383
384 if (isInLens_(globalPos))
385 return lensMaterialParams_;
386 return outerMaterialParams_;
387 }
388
392 template <class Context>
393 Scalar temperature(const Context& /*context*/,
394 unsigned /*spaceIdx*/,
395 unsigned /*timeIdx*/) const
396 { return temperature_; }
397
399
404
408 std::string name() const
409 {
411
412 constexpr bool useAutoDiff = std::is_same_v<LLS, Properties::TTag::AutoDiffLocalLinearizer>;
413
415 constexpr bool useTrans = std::is_same_v<FM, Opm::TransFluxModule<TypeTag>>;
416
417 std::ostringstream oss;
418 oss << "lens_" << Model::name()
419 << "_" << Model::discretizationName()
420 << "_" << (useAutoDiff?"ad":"fd");
421 if (useTrans)
422 oss << "_trans";
423
424 return oss.str();
425 }
426
431 { }
432
437 { }
438
443 {
444#ifndef NDEBUG
445 //this->model().checkConservativeness();
446
447 // Calculate storage terms
448 EqVector storage;
449 this->model().globalStorage(storage);
450
451 // Write mass balance information for rank 0
452 if (this->gridView().comm().rank() == 0) {
453 std::cout << "Storage: " << storage << std::endl << std::flush;
454 }
455#endif // NDEBUG
456 }
457
459
464
468 template <class Context>
469 void boundary(BoundaryRateVector& values,
470 const Context& context,
471 unsigned spaceIdx,
472 unsigned timeIdx) const
473 {
474 const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
475
476 if (onLeftBoundary_(pos) || onRightBoundary_(pos)) {
477 // free flow boundary. we assume incompressible fluids
478 Scalar densityW = WettingPhase::density(temperature_, /*pressure=*/Scalar(1e5));
479 Scalar densityN = NonwettingPhase::density(temperature_, /*pressure=*/Scalar(1e5));
480
481 Scalar T = temperature(context, spaceIdx, timeIdx);
482 Scalar pw, Sw;
483
484 // set wetting phase pressure and saturation
485 if (onLeftBoundary_(pos)) {
486 Scalar height = this->boundingBoxMax()[1] - this->boundingBoxMin()[1];
487 Scalar depth = this->boundingBoxMax()[1] - pos[1];
488 Scalar alpha = (1 + 1.5 / height);
489
490 // hydrostatic pressure scaled by alpha
491 pw = 1e5 - alpha * densityW * this->gravity()[1] * depth;
492 Sw = 1.0;
493 }
494 else {
495 Scalar depth = this->boundingBoxMax()[1] - pos[1];
496
497 // hydrostatic pressure
498 pw = 1e5 - densityW * this->gravity()[1] * depth;
499 Sw = 1.0;
500 }
501
502 // specify a full fluid state using pw and Sw
503 const MaterialLawParams& matParams = this->materialLawParams(context, spaceIdx, timeIdx);
504
505 Opm::ImmiscibleFluidState<Scalar, FluidSystem,
506 /*storeEnthalpy=*/false> fs;
507 fs.setSaturation(wettingPhaseIdx, Sw);
508 fs.setSaturation(nonWettingPhaseIdx, 1 - Sw);
509 fs.setTemperature(T);
510
511 Scalar pC[numPhases];
512 MaterialLaw::capillaryPressures(pC, matParams, fs);
513 fs.setPressure(wettingPhaseIdx, pw);
514 fs.setPressure(nonWettingPhaseIdx, pw + pC[nonWettingPhaseIdx] - pC[wettingPhaseIdx]);
515
516 fs.setDensity(wettingPhaseIdx, densityW);
517 fs.setDensity(nonWettingPhaseIdx, densityN);
518
519 fs.setViscosity(wettingPhaseIdx, WettingPhase::viscosity(temperature_, fs.pressure(wettingPhaseIdx)));
520 fs.setViscosity(nonWettingPhaseIdx, NonwettingPhase::viscosity(temperature_, fs.pressure(nonWettingPhaseIdx)));
521
522 // impose an freeflow boundary condition
523 values.setFreeFlow(context, spaceIdx, timeIdx, fs);
524 }
525 else if (onInlet_(pos)) {
526 RateVector massRate(0.0);
527 massRate = 0.0;
528 massRate[contiNEqIdx] = -0.04; // kg / (m^2 * s)
529
530 // impose a forced flow boundary
531 values.setMassRate(massRate);
532 }
533 else {
534 // no flow boundary
535 values.setNoFlow();
536 }
537 }
538
540
545
549 template <class Context>
550 void initial(PrimaryVariables& values, const Context& context, unsigned spaceIdx, unsigned timeIdx) const
551 {
552 const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
553 Scalar depth = this->boundingBoxMax()[1] - pos[1];
554
555 Opm::ImmiscibleFluidState<Scalar, FluidSystem> fs;
556 fs.setPressure(wettingPhaseIdx, /*pressure=*/1e5);
557
558 Scalar Sw = 1.0;
559 fs.setSaturation(wettingPhaseIdx, Sw);
560 fs.setSaturation(nonWettingPhaseIdx, 1 - Sw);
561
562 fs.setTemperature(temperature_);
563
564 typename FluidSystem::template ParameterCache<Scalar> paramCache;
565 paramCache.updatePhase(fs, wettingPhaseIdx);
566 Scalar densityW = FluidSystem::density(fs, paramCache, wettingPhaseIdx);
567
568 // hydrostatic pressure (assuming incompressibility)
569 Scalar pw = 1e5 - densityW * this->gravity()[1] * depth;
570
571 // calculate the capillary pressure
572 const MaterialLawParams& matParams = this->materialLawParams(context, spaceIdx, timeIdx);
573 Scalar pC[numPhases];
574 MaterialLaw::capillaryPressures(pC, matParams, fs);
575
576 // make a full fluid state
577 fs.setPressure(wettingPhaseIdx, pw);
578 fs.setPressure(nonWettingPhaseIdx, pw + (pC[wettingPhaseIdx] - pC[nonWettingPhaseIdx]));
579
580 // assign the primary variables
581 values.assignNaive(fs);
582 }
583
590 template <class Context>
591 void source(RateVector& rate,
592 const Context& /*context*/,
593 unsigned /*spaceIdx*/,
594 unsigned /*timeIdx*/) const
595 { rate = Scalar(0.0); }
596
598
599private:
600 bool isInLens_(const GlobalPosition& pos) const
601 {
602 for (unsigned i = 0; i < dim; ++i) {
603 if (pos[i] < lensLowerLeft_[i] - eps_ || pos[i] > lensUpperRight_[i]
604 + eps_)
605 return false;
606 }
607 return true;
608 }
609
610 bool onLeftBoundary_(const GlobalPosition& pos) const
611 { return pos[0] < this->boundingBoxMin()[0] + eps_; }
612
613 bool onRightBoundary_(const GlobalPosition& pos) const
614 { return pos[0] > this->boundingBoxMax()[0] - eps_; }
615
616 bool onLowerBoundary_(const GlobalPosition& pos) const
617 { return pos[1] < this->boundingBoxMin()[1] + eps_; }
618
619 bool onUpperBoundary_(const GlobalPosition& pos) const
620 { return pos[1] > this->boundingBoxMax()[1] - eps_; }
621
622 bool onInlet_(const GlobalPosition& pos) const
623 {
624 Scalar width = this->boundingBoxMax()[0] - this->boundingBoxMin()[0];
625 Scalar lambda = (this->boundingBoxMax()[0] - pos[0]) / width;
626 return onUpperBoundary_(pos) && 0.5 < lambda && lambda < 2.0 / 3.0;
627 }
628
629 GlobalPosition lensLowerLeft_;
630 GlobalPosition lensUpperRight_;
631
632 DimMatrix lensK_;
633 DimMatrix outerK_;
634 MaterialLawParams lensMaterialParams_;
635 MaterialLawParams outerMaterialParams_;
636
637 Scalar temperature_;
638 Scalar eps_;
639};
640
641} // namespace Opm
642
643#endif
Soil contamination problem where DNAPL infiltrates a fully water saturated medium.
Definition: lensproblem.hh:176
Scalar temperature(const Context &, unsigned, unsigned) const
Definition: lensproblem.hh:393
static void registerParameters()
Definition: lensproblem.hh:272
void beginTimeStep()
Called by the simulator before each time integration.
Definition: lensproblem.hh:430
static std::string briefDescription()
Returns a human readable description of the problem for the help message.
Definition: lensproblem.hh:320
const MaterialLawParams & materialLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: lensproblem.hh:379
void initial(PrimaryVariables &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the initial value for a control volume.
Definition: lensproblem.hh:550
Scalar porosity(const Context &, unsigned, unsigned) const
Definition: lensproblem.hh:370
LensProblem(Simulator &simulator)
Definition: lensproblem.hh:220
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition: lensproblem.hh:227
void boundary(BoundaryRateVector &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the boundary conditions for a boundary segment.
Definition: lensproblem.hh:469
const DimMatrix & intrinsicPermeability(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: lensproblem.hh:356
std::string name() const
The problem name.
Definition: lensproblem.hh:408
void endTimeStep()
Called by the simulator after each time integration.
Definition: lensproblem.hh:442
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition: lensproblem.hh:436
void source(RateVector &rate, const Context &, unsigned, unsigned) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: lensproblem.hh:591
Defines the properties required for the immiscible multi-phase model.
Defines the common parameters for the porous medium multi-phase models.
Definition: blackoilnewtonmethodparams.hpp:31
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
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:235
static constexpr Scalar value
Definition: groundwaterproblem.hh:95
static constexpr Scalar value
Definition: groundwaterproblem.hh:98
static constexpr Scalar value
Definition: groundwaterproblem.hh:101
static constexpr Scalar value
Definition: groundwaterproblem.hh:104
static constexpr Scalar value
Definition: groundwaterproblem.hh:107
static constexpr Scalar value
Definition: groundwaterproblem.hh:110
Dune::YaspGrid< 2 > type
Definition: lensproblem.hh:77
The type of the DUNE grid.
Definition: basicproperties.hh:100
Opm::EffToAbsLaw< EffectiveLaw > type
Definition: lensproblem.hh:121
The material law which ought to be used (extracted from the spatial parameters)
Definition: multiphasebaseproperties.hh:51
Opm::LiquidPhase< Scalar, Opm::DNAPL< Scalar > > type
Definition: lensproblem.hh:98
The non-wetting phase for two-phase models.
Definition: immiscibleproperties.hh:44
The type of the problem.
Definition: fvbaseproperties.hh:81
Definition: lensproblem.hh:68
std::tuple< StructuredGridVanguard > InheritsFrom
Definition: lensproblem.hh:68
Opm::LiquidPhase< Scalar, Opm::SimpleH2O< Scalar > > type
Definition: lensproblem.hh:87
The wetting phase for two-phase models.
Definition: immiscibleproperties.hh:41
This file contains the flux module that uses transmissibilities.