richardslensproblem.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_RICHARDS_LENS_PROBLEM_HH
29#define EWOMS_RICHARDS_LENS_PROBLEM_HH
30
32
33#include <opm/material/components/SimpleH2O.hpp>
34#include <opm/material/fluidsystems/LiquidPhase.hpp>
35#include <opm/material/fluidmatrixinteractions/RegularizedVanGenuchten.hpp>
36#include <opm/material/fluidmatrixinteractions/LinearMaterial.hpp>
37#include <opm/material/fluidmatrixinteractions/EffToAbsLaw.hpp>
38#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
39
40#include <dune/grid/yaspgrid.hh>
41#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
42
43#include <dune/common/version.hh>
44#include <dune/common/fvector.hh>
45#include <dune/common/fmatrix.hh>
46
47namespace Opm {
48template <class TypeTag>
49class RichardsLensProblem;
50
51} // namespace Opm
52
53namespace Opm::Properties {
54
55// Create new type tags
56namespace TTag {
57struct RichardsLensProblem { using InheritsFrom = std::tuple<Richards>; };
58} // end namespace TTag
59
60// Use 2d YaspGrid
61template<class TypeTag>
62struct Grid<TypeTag, TTag::RichardsLensProblem> { using type = Dune::YaspGrid<2>; };
63
64// Set the physical problem to be solved
65template<class TypeTag>
67
68// Set the wetting phase
69template<class TypeTag>
70struct WettingFluid<TypeTag, TTag::RichardsLensProblem>
71{
72private:
74
75public:
76 using type = Opm::LiquidPhase<Scalar, Opm::SimpleH2O<Scalar> >;
77};
78
79// Set the material Law
80template<class TypeTag>
81struct MaterialLaw<TypeTag, TTag::RichardsLensProblem>
82{
83private:
85 enum { wettingPhaseIdx = FluidSystem::wettingPhaseIdx };
86 enum { nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx };
87
89 using Traits = Opm::TwoPhaseMaterialTraits<Scalar,
90 /*wettingPhaseIdx=*/FluidSystem::wettingPhaseIdx,
91 /*nonWettingPhaseIdx=*/FluidSystem::nonWettingPhaseIdx>;
92
93 // define the material law which is parameterized by effective
94 // saturations
95 using EffectiveLaw = Opm::RegularizedVanGenuchten<Traits>;
96
97public:
98 // define the material law parameterized by absolute saturations
99 using type = Opm::EffToAbsLaw<EffectiveLaw>;
100};
101
102} // namespace Opm::Properties
103
104namespace Opm {
105
122template <class TypeTag>
123class RichardsLensProblem : public GetPropType<TypeTag, Properties::BaseProblem>
124{
126
137
139 enum {
140 // copy some indices for convenience
141 pressureWIdx = Indices::pressureWIdx,
142 contiEqIdx = Indices::contiEqIdx,
143 wettingPhaseIdx = FluidSystem::wettingPhaseIdx,
144 nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx,
145 numPhases = FluidSystem::numPhases,
146
147 // Grid and world dimension
148 dimWorld = GridView::dimensionworld
149 };
150
151 // get the material law from the property system
154 using MaterialLawParams = typename MaterialLaw::Params;
155
156 using CoordScalar = typename GridView::ctype;
157 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
158 using PhaseVector = Dune::FieldVector<Scalar, numPhases>;
159 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
160
161public:
165 RichardsLensProblem(Simulator& simulator)
166 : ParentType(simulator)
167 , pnRef_(1e5)
168 {
169 dofIsInLens_.resize(simulator.model().numGridDof());
170 }
171
176 {
177 ParentType::finishInit();
178
179 eps_ = 3e-6;
180 pnRef_ = 1e5;
181
182 lensLowerLeft_[0] = 1.0;
183 lensLowerLeft_[1] = 2.0;
184
185 lensUpperRight_[0] = 4.0;
186 lensUpperRight_[1] = 3.0;
187
188 // parameters for the Van Genuchten law
189 // alpha and n
190 lensMaterialParams_.setVgAlpha(0.00045);
191 lensMaterialParams_.setVgN(7.3);
192 lensMaterialParams_.finalize();
193
194 outerMaterialParams_.setVgAlpha(0.0037);
195 outerMaterialParams_.setVgN(4.7);
196 outerMaterialParams_.finalize();
197
198 // parameters for the linear law
199 // minimum and maximum pressures
200 // lensMaterialParams_.setEntryPC(0);
201 // outerMaterialParams_.setEntryPC(0);
202 // lensMaterialParams_.setMaxPC(0);
203 // outerMaterialParams_.setMaxPC(0);
204
205 lensK_ = this->toDimMatrix_(1e-12);
206 outerK_ = this->toDimMatrix_(5e-12);
207
208 // determine which degrees of freedom are in the lens
209 Stencil stencil(this->gridView(), this->simulator().model().dofMapper() );
210 for (const auto& elem : elements(this->gridView())) {
211 stencil.update(elem);
212 for (unsigned dofIdx = 0; dofIdx < stencil.numPrimaryDof(); ++ dofIdx) {
213 unsigned globalDofIdx = stencil.globalSpaceIndex(dofIdx);
214 const auto& dofPos = stencil.subControlVolume(dofIdx).center();
215 dofIsInLens_[globalDofIdx] = isInLens_(dofPos);
216 }
217 }
218 }
219
223 static void registerParameters()
224 {
225 ParentType::registerParameters();
226
227 Parameters::SetDefault<Parameters::GridFile>("./data/richardslens_24x16.dgf");
228
229 // Use central differences to approximate the Jacobian matrix
231 constexpr bool useFD = std::is_same_v<LLS, Properties::TTag::FiniteDifferenceLocalLinearizer>;
232 if constexpr (useFD) {
233 Parameters::SetDefault<Parameters::NumericDifferenceMethod>(0);
234 }
235
236 Parameters::SetDefault<Parameters::EndTime<Scalar>>(3000.0);
237 Parameters::SetDefault<Parameters::InitialTimeStepSize<Scalar>>(100.0);
238 Parameters::SetDefault<Parameters::NewtonMaxIterations>(28);
239 Parameters::SetDefault<Parameters::NewtonTargetIterations>(18);
240 Parameters::SetDefault<Parameters::EnableGravity>(true);
241 }
242
247
251 std::string name() const
252 {
253 std::ostringstream oss;
254 oss << "lens_richards_"
255 << Model::discretizationName();
256 return oss.str();
257 }
258
263 {
264#ifndef NDEBUG
265 this->model().checkConservativeness();
266
267 // Calculate storage terms
268 EqVector storage;
269 this->model().globalStorage(storage);
270
271 // Write mass balance information for rank 0
272 if (this->gridView().comm().rank() == 0) {
273 std::cout << "Storage: " << storage << std::endl << std::flush;
274 }
275#endif // NDEBUG
276 }
277
281 template <class Context>
282 Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
283 { return temperature(context.globalSpaceIndex(spaceIdx, timeIdx), timeIdx); }
284
285 Scalar temperature(unsigned /*globalSpaceIdx*/, unsigned /*timeIdx*/) const
286 { return 273.15 + 10; } // -> 10°C
287
291 template <class Context>
292 const DimMatrix& intrinsicPermeability(const Context& context,
293 unsigned spaceIdx,
294 unsigned timeIdx) const
295 {
296 const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
297 if (isInLens_(pos))
298 return lensK_;
299 return outerK_;
300 }
301
305 template <class Context>
306 Scalar porosity(const Context& /*context*/,
307 unsigned /*spaceIdx*/,
308 unsigned /*timeIdx*/) const
309 { return 0.4; }
310
314 template <class Context>
315 const MaterialLawParams& materialLawParams(const Context& context,
316 unsigned spaceIdx,
317 unsigned timeIdx) const
318 {
319 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
320 return materialLawParams(globalSpaceIdx, timeIdx);
321 }
322
323 const MaterialLawParams& materialLawParams(unsigned globalSpaceIdx,
324 unsigned /*timeIdx*/) const
325 {
326 if (dofIsInLens_[globalSpaceIdx])
327 return lensMaterialParams_;
328 return outerMaterialParams_;
329 }
330
336 template <class Context>
337 Scalar referencePressure(const Context& context,
338 unsigned spaceIdx,
339 unsigned timeIdx) const
340 { return referencePressure(context.globalSpaceIndex(spaceIdx, timeIdx), timeIdx); }
341
342 // the Richards model does not have an element context available at all places
343 // where the reference pressure is required...
344 Scalar referencePressure(unsigned /*globalSpaceIdx*/,
345 unsigned /*timeIdx*/) const
346 { return pnRef_; }
347
349
354
358 template <class Context>
359 void boundary(BoundaryRateVector& values,
360 const Context& context,
361 unsigned spaceIdx,
362 unsigned timeIdx) const
363 {
364 const auto& pos = context.pos(spaceIdx, timeIdx);
365
366 if (onLeftBoundary_(pos) || onRightBoundary_(pos)) {
367 const auto& materialParams = this->materialLawParams(context, spaceIdx, timeIdx);
368
369 Scalar Sw = 0.0;
370 Opm::ImmiscibleFluidState<Scalar, FluidSystem> fs;
371 fs.setSaturation(wettingPhaseIdx, Sw);
372 fs.setSaturation(nonWettingPhaseIdx, 1.0 - Sw);
373
374 PhaseVector pC;
375 MaterialLaw::capillaryPressures(pC, materialParams, fs);
376 fs.setPressure(wettingPhaseIdx, pnRef_ + pC[wettingPhaseIdx] - pC[nonWettingPhaseIdx]);
377 fs.setPressure(nonWettingPhaseIdx, pnRef_);
378
379 typename FluidSystem::template ParameterCache<Scalar> paramCache;
380 paramCache.updateAll(fs);
381 fs.setDensity(wettingPhaseIdx, FluidSystem::density(fs, paramCache, wettingPhaseIdx));
382 //fs.setDensity(nonWettingPhaseIdx, FluidSystem::density(fs, paramCache, nonWettingPhaseIdx));
383
384 fs.setViscosity(wettingPhaseIdx, FluidSystem::viscosity(fs, paramCache, wettingPhaseIdx));
385 //fs.setViscosity(nonWettingPhaseIdx, FluidSystem::viscosity(fs, paramCache, nonWettingPhaseIdx));
386
387 values.setFreeFlow(context, spaceIdx, timeIdx, fs);
388 }
389 else if (onInlet_(pos)) {
390 RateVector massRate(0.0);
391
392 // inflow of water
393 massRate[contiEqIdx] = -0.04; // kg / (m * s)
394
395 values.setMassRate(massRate);
396 }
397 else
398 values.setNoFlow();
399 }
400
402
407
411 template <class Context>
412 void initial(PrimaryVariables& values,
413 const Context& context,
414 unsigned spaceIdx,
415 unsigned timeIdx) const
416 {
417 const auto& materialParams = this->materialLawParams(context, spaceIdx, timeIdx);
418
419 Scalar Sw = 0.0;
420 Opm::ImmiscibleFluidState<Scalar, FluidSystem> fs;
421 fs.setSaturation(wettingPhaseIdx, Sw);
422 fs.setSaturation(nonWettingPhaseIdx, 1.0 - Sw);
423
424 PhaseVector pC;
425 MaterialLaw::capillaryPressures(pC, materialParams, fs);
426 values[pressureWIdx] = pnRef_ + (pC[wettingPhaseIdx] - pC[nonWettingPhaseIdx]);
427 }
428
435 template <class Context>
436 void source(RateVector& rate,
437 const Context& /*context*/,
438 unsigned /*spaceIdx*/,
439 unsigned /*timeIdx*/) const
440 { rate = Scalar(0.0); }
441
443
444private:
445 bool onLeftBoundary_(const GlobalPosition& pos) const
446 { return pos[0] < this->boundingBoxMin()[0] + eps_; }
447
448 bool onRightBoundary_(const GlobalPosition& pos) const
449 { return pos[0] > this->boundingBoxMax()[0] - eps_; }
450
451 bool onLowerBoundary_(const GlobalPosition& pos) const
452 { return pos[1] < this->boundingBoxMin()[1] + eps_; }
453
454 bool onUpperBoundary_(const GlobalPosition& pos) const
455 { return pos[1] > this->boundingBoxMax()[1] - eps_; }
456
457 bool onInlet_(const GlobalPosition& pos) const
458 {
459 Scalar width = this->boundingBoxMax()[0] - this->boundingBoxMin()[0];
460 Scalar lambda = (this->boundingBoxMax()[0] - pos[0]) / width;
461 return onUpperBoundary_(pos) && 0.5 < lambda && lambda < 2.0 / 3.0;
462 }
463
464 bool isInLens_(const GlobalPosition& pos) const
465 {
466 for (unsigned i = 0; i < dimWorld; ++i) {
467 if (pos[i] < lensLowerLeft_[i] || pos[i] > lensUpperRight_[i])
468 return false;
469 }
470 return true;
471 }
472
473 GlobalPosition lensLowerLeft_;
474 GlobalPosition lensUpperRight_;
475
476 DimMatrix lensK_;
477 DimMatrix outerK_;
478 MaterialLawParams lensMaterialParams_;
479 MaterialLawParams outerMaterialParams_;
480
481 std::vector<bool> dofIsInLens_;
482
483 Scalar eps_;
484 Scalar pnRef_;
485};
486} // namespace Opm
487
488#endif
A water infiltration problem with a low-permeability lens embedded into a high-permeability domain.
Definition: richardslensproblem.hh:124
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition: richardslensproblem.hh:175
static void registerParameters()
Definition: richardslensproblem.hh:223
std::string name() const
The problem name.
Definition: richardslensproblem.hh:251
const MaterialLawParams & materialLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: richardslensproblem.hh:323
const MaterialLawParams & materialLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: richardslensproblem.hh:315
Scalar temperature(unsigned, unsigned) const
Definition: richardslensproblem.hh:285
Scalar referencePressure(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Return the reference pressure [Pa] of the wetting phase.
Definition: richardslensproblem.hh:337
RichardsLensProblem(Simulator &simulator)
Definition: richardslensproblem.hh:165
const DimMatrix & intrinsicPermeability(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: richardslensproblem.hh:292
void source(RateVector &rate, const Context &, unsigned, unsigned) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: richardslensproblem.hh:436
Scalar referencePressure(unsigned, unsigned) const
Definition: richardslensproblem.hh:344
void boundary(BoundaryRateVector &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the boundary conditions for a boundary segment.
Definition: richardslensproblem.hh:359
void endTimeStep()
Called by the simulator after each time integration.
Definition: richardslensproblem.hh:262
Scalar temperature(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: richardslensproblem.hh:282
void initial(PrimaryVariables &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the initial value for a control volume.
Definition: richardslensproblem.hh:412
Scalar porosity(const Context &, unsigned, unsigned) const
Definition: richardslensproblem.hh:306
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
Dune::YaspGrid< 2 > type
Definition: richardslensproblem.hh:62
The type of the DUNE grid.
Definition: basicproperties.hh:100
Opm::EffToAbsLaw< EffectiveLaw > type
Definition: richardslensproblem.hh:99
The material law which ought to be used (extracted from the spatial parameters)
Definition: multiphasebaseproperties.hh:51
The type of the problem.
Definition: fvbaseproperties.hh:81
Definition: richardslensproblem.hh:57
std::tuple< Richards > InheritsFrom
Definition: richardslensproblem.hh:57
Opm::LiquidPhase< Scalar, Opm::SimpleH2O< Scalar > > type
Definition: richardslensproblem.hh:76
Definition: richardsproperties.hh:40