powerinjectionproblem.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_POWER_INJECTION_PROBLEM_HH
29#define EWOMS_POWER_INJECTION_PROBLEM_HH
30
33
34#include <opm/material/fluidmatrixinteractions/RegularizedVanGenuchten.hpp>
35#include <opm/material/fluidmatrixinteractions/LinearMaterial.hpp>
36#include <opm/material/fluidmatrixinteractions/EffToAbsLaw.hpp>
37#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
38#include <opm/material/fluidsystems/TwoPhaseImmiscibleFluidSystem.hpp>
39#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
40#include <opm/material/components/SimpleH2O.hpp>
41#include <opm/material/components/Air.hpp>
42
43#include <dune/grid/yaspgrid.hh>
44
45#include <dune/common/version.hh>
46#include <dune/common/fvector.hh>
47#include <dune/common/fmatrix.hh>
48
49#include <sstream>
50#include <string>
51#include <type_traits>
52#include <iostream>
53
54namespace Opm {
55template <class TypeTag>
56class PowerInjectionProblem;
57}
58
59namespace Opm::Properties {
60
61namespace TTag {
63}
64
65// Set the grid implementation to be used
66template<class TypeTag>
67struct Grid<TypeTag, TTag::PowerInjectionBaseProblem> { using type = Dune::YaspGrid</*dim=*/1>; };
68
69// set the Vanguard property
70template<class TypeTag>
71struct Vanguard<TypeTag, TTag::PowerInjectionBaseProblem> { using type = Opm::CubeGridVanguard<TypeTag>; };
72
73// Set the problem property
74template<class TypeTag>
75struct Problem<TypeTag, TTag::PowerInjectionBaseProblem> { using type = Opm::PowerInjectionProblem<TypeTag>; };
76
77// Set the wetting phase
78template<class TypeTag>
79struct WettingPhase<TypeTag, TTag::PowerInjectionBaseProblem>
80{
81private:
83
84public:
85 using type = Opm::LiquidPhase<Scalar, Opm::SimpleH2O<Scalar> >;
86};
87
88// Set the non-wetting phase
89template<class TypeTag>
90struct NonwettingPhase<TypeTag, TTag::PowerInjectionBaseProblem>
91{
92private:
94
95public:
96 using type = Opm::GasPhase<Scalar, Opm::Air<Scalar> >;
97};
98
99// Set the material Law
100template<class TypeTag>
101struct MaterialLaw<TypeTag, TTag::PowerInjectionBaseProblem>
102{
103private:
105 enum { wettingPhaseIdx = FluidSystem::wettingPhaseIdx };
106 enum { nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx };
107
109 using Traits = Opm::TwoPhaseMaterialTraits<Scalar,
110 /*wettingPhaseIdx=*/FluidSystem::wettingPhaseIdx,
111 /*nonWettingPhaseIdx=*/FluidSystem::nonWettingPhaseIdx>;
112
113 // define the material law which is parameterized by effective
114 // saturations
115 using EffectiveLaw = Opm::RegularizedVanGenuchten<Traits>;
116
117public:
118 // define the material law parameterized by absolute saturations
119 using type = Opm::EffToAbsLaw<EffectiveLaw>;
120};
121
122} // namespace Opm::Properties
123
124namespace Opm {
137template <class TypeTag>
138class PowerInjectionProblem : public GetPropType<TypeTag, Properties::BaseProblem>
139{
141
153
154 enum {
155 // number of phases
156 numPhases = FluidSystem::numPhases,
157
158 // phase indices
159 wettingPhaseIdx = FluidSystem::wettingPhaseIdx,
160 nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx,
161
162 // equation indices
163 contiNEqIdx = Indices::conti0EqIdx + nonWettingPhaseIdx,
164
165 // Grid and world dimension
166 dim = GridView::dimension,
167 dimWorld = GridView::dimensionworld
168 };
169
172
173 using CoordScalar = typename GridView::ctype;
174 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
175
176 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
177
178public:
182 PowerInjectionProblem(Simulator& simulator)
183 : ParentType(simulator)
184 { }
185
190 {
191 ParentType::finishInit();
192
193 eps_ = 3e-6;
194 FluidSystem::init();
195
196 temperature_ = 273.15 + 26.6;
197
198 // parameters for the Van Genuchten law
199 // alpha and n
200 materialParams_.setVgAlpha(0.00045);
201 materialParams_.setVgN(7.3);
202 materialParams_.finalize();
203
204 K_ = this->toDimMatrix_(5.73e-08); // [m^2]
205
206 setupInitialFluidState_();
207 }
208
212 static void registerParameters()
213 {
214 ParentType::registerParameters();
215
216 Parameters::SetDefault<Parameters::CellsX>(250);
217 Parameters::SetDefault<Parameters::DomainSizeX<Scalar>>(100.0);
218
219 if constexpr (dim > 1) {
220 Parameters::SetDefault<Parameters::CellsY>(1);
221 Parameters::SetDefault<Parameters::DomainSizeY<Scalar>>(1.0);
222 }
223 if constexpr (dim == 3) {
224 Parameters::SetDefault<Parameters::CellsZ>(1);
225 Parameters::SetDefault<Parameters::DomainSizeZ<Scalar>>(1.0);
226 }
227
228 Parameters::SetDefault<Parameters::EndTime<Scalar>>(100.0);
229 Parameters::SetDefault<Parameters::InitialTimeStepSize<Scalar>>(1e-3);
230 Parameters::SetDefault<Parameters::VtkWriteFilterVelocities>(true);
231 }
232
237
241 std::string name() const
242 {
243 std::ostringstream oss;
244 oss << "powerinjection_";
247 oss << "darcy";
248 else
249 oss << "forchheimer";
250
253 oss << "_" << "ad";
254 else
255 oss << "_" << "fd";
256
257 return oss.str();
258 }
259
264 {
265#ifndef NDEBUG
266 this->model().checkConservativeness();
267
268 // Calculate storage terms
269 EqVector storage;
270 this->model().globalStorage(storage);
271
272 // Write mass balance information for rank 0
273 if (this->gridView().comm().rank() == 0) {
274 std::cout << "Storage: " << storage << std::endl << std::flush;
275 }
276#endif // NDEBUG
277 }
279
284
288 template <class Context>
289 const DimMatrix& intrinsicPermeability(const Context& /*context*/,
290 unsigned /*spaceIdx*/,
291 unsigned /*timeIdx*/) const
292 { return K_; }
293
297 template <class Context>
298 Scalar ergunCoefficient(const Context& /*context*/,
299 unsigned /*spaceIdx*/,
300 unsigned /*timeIdx*/) const
301 { return 0.3866; }
302
306 template <class Context>
307 Scalar porosity(const Context& /*context*/,
308 unsigned /*spaceIdx*/,
309 unsigned /*timeIdx*/) const
310 { return 0.558; }
311
315 template <class Context>
316 const MaterialLawParams&
317 materialLawParams(const Context& /*context*/,
318 unsigned /*spaceIdx*/,
319 unsigned /*timeIdx*/) const
320 { return materialParams_; }
321
325 template <class Context>
326 Scalar temperature(const Context& /*context*/,
327 unsigned /*spaceIdx*/,
328 unsigned /*timeIdx*/) const
329 { return temperature_; }
330
332
337
344 template <class Context>
345 void boundary(BoundaryRateVector& values,
346 const Context& context,
347 unsigned spaceIdx,
348 unsigned timeIdx) const
349 {
350 const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
351
352 if (onLeftBoundary_(pos)) {
353 RateVector massRate(0.0);
354 massRate = 0.0;
355 massRate[contiNEqIdx] = -1.00; // kg / (m^2 * s)
356
357 // impose a forced flow boundary
358 values.setMassRate(massRate);
359 }
360 else if (onRightBoundary_(pos))
361 // free flow boundary with initial condition on the right
362 values.setFreeFlow(context, spaceIdx, timeIdx, initialFluidState_);
363 else
364 values.setNoFlow();
365 }
366
368
373
377 template <class Context>
378 void initial(PrimaryVariables& values,
379 const Context& /*context*/,
380 unsigned /*spaceIdx*/,
381 unsigned /*timeIdx*/) const
382 {
383 // assign the primary variables
384 values.assignNaive(initialFluidState_);
385 }
386
393 template <class Context>
394 void source(RateVector& rate,
395 const Context& /*context*/,
396 unsigned /*spaceIdx*/,
397 unsigned /*timeIdx*/) const
398 { rate = Scalar(0.0); }
399
401
402private:
403 bool onLeftBoundary_(const GlobalPosition& pos) const
404 { return pos[0] < this->boundingBoxMin()[0] + eps_; }
405
406 bool onRightBoundary_(const GlobalPosition& pos) const
407 { return pos[0] > this->boundingBoxMax()[0] - eps_; }
408
409 void setupInitialFluidState_()
410 {
411 initialFluidState_.setTemperature(temperature_);
412
413 Scalar Sw = 1.0;
414 initialFluidState_.setSaturation(wettingPhaseIdx, Sw);
415 initialFluidState_.setSaturation(nonWettingPhaseIdx, 1 - Sw);
416
417 Scalar p = 1e5;
418 initialFluidState_.setPressure(wettingPhaseIdx, p);
419 initialFluidState_.setPressure(nonWettingPhaseIdx, p);
420
421 typename FluidSystem::template ParameterCache<Scalar> paramCache;
422 paramCache.updateAll(initialFluidState_);
423 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
424 initialFluidState_.setDensity(phaseIdx,
425 FluidSystem::density(initialFluidState_, paramCache, phaseIdx));
426 initialFluidState_.setViscosity(phaseIdx,
427 FluidSystem::viscosity(initialFluidState_, paramCache, phaseIdx));
428 }
429 }
430
431 DimMatrix K_;
432 MaterialLawParams materialParams_;
433
434 Opm::ImmiscibleFluidState<Scalar, FluidSystem> initialFluidState_;
435 Scalar temperature_;
436 Scalar eps_;
437};
438
439} // namespace Opm
440
441#endif
Provides a simulator vanguad which creates a regular grid made of quadrilaterals.
Definition: cubegridvanguard.hh:53
1D Problem with very fast injection of gas on the left.
Definition: powerinjectionproblem.hh:139
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition: powerinjectionproblem.hh:189
Scalar porosity(const Context &, unsigned, unsigned) const
Definition: powerinjectionproblem.hh:307
static void registerParameters()
Definition: powerinjectionproblem.hh:212
Scalar temperature(const Context &, unsigned, unsigned) const
Definition: powerinjectionproblem.hh:326
void endTimeStep()
Called by the simulator after each time integration.
Definition: powerinjectionproblem.hh:263
void initial(PrimaryVariables &values, const Context &, unsigned, unsigned) const
Evaluate the initial value for a control volume.
Definition: powerinjectionproblem.hh:378
Scalar ergunCoefficient(const Context &, unsigned, unsigned) const
Returns the Ergun coefficient.
Definition: powerinjectionproblem.hh:298
const MaterialLawParams & materialLawParams(const Context &, unsigned, unsigned) const
Definition: powerinjectionproblem.hh:317
std::string name() const
The problem name.
Definition: powerinjectionproblem.hh:241
void source(RateVector &rate, const Context &, unsigned, unsigned) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: powerinjectionproblem.hh:394
PowerInjectionProblem(Simulator &simulator)
Definition: powerinjectionproblem.hh:182
void boundary(BoundaryRateVector &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the boundary conditions for a boundary segment.
Definition: powerinjectionproblem.hh:345
const DimMatrix & intrinsicPermeability(const Context &, unsigned, unsigned) const
Definition: powerinjectionproblem.hh:289
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
Specifies a flux module which uses the Darcy relation.
Definition: darcyfluxmodule.hh:61
Dune::YaspGrid< 1 > type
Definition: powerinjectionproblem.hh:67
The type of the DUNE grid.
Definition: basicproperties.hh:100
Opm::EffToAbsLaw< EffectiveLaw > type
Definition: powerinjectionproblem.hh:119
The material law which ought to be used (extracted from the spatial parameters)
Definition: multiphasebaseproperties.hh:51
Opm::GasPhase< Scalar, Opm::Air< Scalar > > type
Definition: powerinjectionproblem.hh:96
The non-wetting phase for two-phase models.
Definition: immiscibleproperties.hh:44
The type of the problem.
Definition: fvbaseproperties.hh:81
Definition: fvbaseadlocallinearizer.hh:53
Definition: powerinjectionproblem.hh:62
Property which provides a Vanguard (manages grids)
Definition: basicproperties.hh:96
Opm::LiquidPhase< Scalar, Opm::SimpleH2O< Scalar > > type
Definition: powerinjectionproblem.hh:85
The wetting phase for two-phase models.
Definition: immiscibleproperties.hh:41