2015-06-18 06:43:59 -05:00
|
|
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
|
|
|
|
// vi: set et ts=4 sw=4 sts=4:
|
2013-12-02 08:48:04 -06:00
|
|
|
/*
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 2 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
2016-03-14 07:21:47 -05:00
|
|
|
|
|
|
|
Consult the COPYING file in the top-level source directory of this
|
|
|
|
module for the precise wording of the license and the list of
|
|
|
|
copyright holders.
|
2013-12-02 08:48:04 -06:00
|
|
|
*/
|
2013-09-23 11:56:30 -05:00
|
|
|
/*!
|
|
|
|
* \file
|
|
|
|
*
|
2019-09-05 09:21:10 -05:00
|
|
|
* \copydoc Opm::ReservoirProblem
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
#ifndef EWOMS_RESERVOIR_PROBLEM_HH
|
|
|
|
#define EWOMS_RESERVOIR_PROBLEM_HH
|
|
|
|
|
2019-09-16 04:15:31 -05:00
|
|
|
#include <opm/models/blackoil/blackoilproperties.hh>
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2013-11-06 07:50:01 -06:00
|
|
|
#include <opm/material/fluidmatrixinteractions/LinearMaterial.hpp>
|
|
|
|
#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
|
2013-09-23 11:56:30 -05:00
|
|
|
#include <opm/material/fluidstates/CompositionalFluidState.hpp>
|
2015-12-31 06:20:56 -06:00
|
|
|
#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
|
|
|
|
#include <opm/material/constraintsolvers/ComputeFromReferencePhase.hpp>
|
2015-02-03 06:58:18 -06:00
|
|
|
#include <opm/material/fluidsystems/blackoilpvt/DryGasPvt.hpp>
|
|
|
|
#include <opm/material/fluidsystems/blackoilpvt/LiveOilPvt.hpp>
|
|
|
|
#include <opm/material/fluidsystems/blackoilpvt/ConstantCompressibilityWaterPvt.hpp>
|
|
|
|
|
2013-12-01 04:51:05 -06:00
|
|
|
#include <dune/grid/yaspgrid.hh>
|
2014-12-16 05:39:32 -06:00
|
|
|
#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2013-12-01 04:51:05 -06:00
|
|
|
#include <dune/common/version.hh>
|
2013-09-23 11:56:30 -05:00
|
|
|
#include <dune/common/fvector.hh>
|
|
|
|
#include <dune/common/fmatrix.hh>
|
|
|
|
|
|
|
|
#include <vector>
|
|
|
|
#include <string>
|
|
|
|
|
2019-09-05 09:21:10 -05:00
|
|
|
namespace Opm {
|
2013-09-23 11:56:30 -05:00
|
|
|
template <class TypeTag>
|
|
|
|
class ReservoirProblem;
|
|
|
|
|
2019-09-05 09:21:10 -05:00
|
|
|
} // namespace Opm
|
2018-06-14 09:06:05 -05:00
|
|
|
|
2020-06-08 10:11:48 -05:00
|
|
|
namespace Opm::Properties {
|
2018-06-14 09:06:05 -05:00
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2020-06-10 06:07:19 -05:00
|
|
|
namespace TTag {
|
|
|
|
|
|
|
|
struct ReservoirBaseProblem {};
|
|
|
|
|
|
|
|
} // namespace TTag
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// Maximum depth of the reservoir
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct MaxDepth { using type = UndefinedProperty; };
|
2013-09-23 11:56:30 -05:00
|
|
|
// The temperature inside the reservoir
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct Temperature { using type = UndefinedProperty; };
|
2015-12-31 06:20:56 -06:00
|
|
|
// The width of producer/injector wells as a fraction of the width of the spatial domain
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct WellWidth { using type = UndefinedProperty; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// Set the grid type
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct Grid<TypeTag, TTag::ReservoirBaseProblem> { using type = Dune::YaspGrid<2>; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// Set the problem property
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct Problem<TypeTag, TTag::ReservoirBaseProblem> { using type = Opm::ReservoirProblem<TypeTag>; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// Set the material Law
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct MaterialLaw<TypeTag, TTag::ReservoirBaseProblem>
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
|
|
|
private:
|
2020-06-10 06:49:42 -05:00
|
|
|
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
|
|
|
|
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using Traits = Opm::
|
2013-11-29 09:33:46 -06:00
|
|
|
ThreePhaseMaterialTraits<Scalar,
|
2014-04-03 09:59:48 -05:00
|
|
|
/*wettingPhaseIdx=*/FluidSystem::waterPhaseIdx,
|
|
|
|
/*nonWettingPhaseIdx=*/FluidSystem::oilPhaseIdx,
|
2020-06-10 06:49:42 -05:00
|
|
|
/*gasPhaseIdx=*/FluidSystem::gasPhaseIdx>;
|
2013-11-06 07:50:01 -06:00
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
public:
|
2020-06-10 06:49:42 -05:00
|
|
|
using type = Opm::LinearMaterial<Traits>;
|
2013-09-23 11:56:30 -05:00
|
|
|
};
|
|
|
|
|
|
|
|
// Write the Newton convergence behavior to disk?
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct NewtonWriteConvergence<TypeTag, TTag::ReservoirBaseProblem> { static constexpr bool value = false; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// Enable gravity
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct EnableGravity<TypeTag, TTag::ReservoirBaseProblem> { static constexpr bool value = true; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// Enable constraint DOFs?
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct EnableConstraints<TypeTag, TTag::ReservoirBaseProblem> { static constexpr bool value = true; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// set the defaults for some problem specific properties
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct MaxDepth<TypeTag, TTag::ReservoirBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 2500;
|
|
|
|
};
|
|
|
|
template<class TypeTag>
|
|
|
|
struct Temperature<TypeTag, TTag::ReservoirBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 293.15;
|
|
|
|
};
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2015-12-31 06:20:56 -06:00
|
|
|
//! The default for the end time of the simulation [s].
|
|
|
|
//!
|
|
|
|
//! By default this problem spans 1000 days (100 "settle down" days and 900 days of
|
|
|
|
//! production)
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct EndTime<TypeTag, TTag::ReservoirBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 1000.0*24*60*60;
|
|
|
|
};
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// The default for the initial time step size of the simulation [s]
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct InitialTimeStepSize<TypeTag, TTag::ReservoirBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 100e3;
|
|
|
|
};
|
2015-12-31 06:20:56 -06:00
|
|
|
|
|
|
|
// The width of producer/injector wells as a fraction of the width of the spatial domain
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct WellWidth<TypeTag, TTag::ReservoirBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 0.01;
|
|
|
|
};
|
2015-12-31 06:20:56 -06:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Explicitly set the fluid system to the black-oil fluid system
|
|
|
|
*
|
|
|
|
* If the black oil model is used, this is superfluous because that model already sets
|
|
|
|
* the FluidSystem property. Setting it explictly for the problem is a good idea anyway,
|
|
|
|
* though because other models are more generic and thus do not assume a particular fluid
|
|
|
|
* system.
|
|
|
|
*/
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct FluidSystem<TypeTag, TTag::ReservoirBaseProblem>
|
2015-12-31 06:20:56 -06:00
|
|
|
{
|
|
|
|
private:
|
2020-06-10 06:49:42 -05:00
|
|
|
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
|
2015-12-31 06:20:56 -06:00
|
|
|
|
|
|
|
public:
|
2020-06-10 06:49:42 -05:00
|
|
|
using type = Opm::BlackOilFluidSystem<Scalar>;
|
2015-12-31 06:20:56 -06:00
|
|
|
};
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// The default DGF file to load
|
2020-06-09 04:15:16 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct GridFile<TypeTag, TTag::ReservoirBaseProblem> { static constexpr auto value = "data/reservoir.dgf"; };
|
2015-12-31 06:20:56 -06:00
|
|
|
|
|
|
|
// increase the tolerance for this problem to get larger time steps
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct NewtonTolerance<TypeTag, TTag::ReservoirBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 1e-6;
|
|
|
|
};
|
2018-06-14 09:06:05 -05:00
|
|
|
|
2020-06-08 10:11:48 -05:00
|
|
|
} // namespace Opm::Properties
|
2018-06-14 09:06:05 -05:00
|
|
|
|
2019-09-05 09:21:10 -05:00
|
|
|
namespace Opm {
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \ingroup TestProblems
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
|
|
|
* \brief Some simple test problem for the black-oil VCVF discretization
|
|
|
|
* inspired by an oil reservoir.
|
|
|
|
*
|
2015-12-31 06:20:56 -06:00
|
|
|
* The domain is two-dimensional and exhibits a size of 6000m times 60m. Initially, the
|
|
|
|
* reservoir is assumed by oil with a bubble point pressure of 20 MPa, which also the
|
|
|
|
* initial pressure in the domain. No-flow boundaries are used for all boundaries. The
|
|
|
|
* permeability of the lower 10 m is reduced compared to the upper 10 m of the domain
|
|
|
|
* witch capillary pressure always being neglected. Three wells are approximated using
|
|
|
|
* constraints: Two water-injector wells, one at the lower-left boundary one at the
|
|
|
|
* lower-right boundary and one producer well in the upper part of the center of the
|
|
|
|
* domain. The pressure for the producer is assumed to be 2/3 of the reservoir pressure,
|
|
|
|
* the injector wells use a pressure which is 50% above the reservoir pressure.
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
template <class TypeTag>
|
2020-06-08 09:41:02 -05:00
|
|
|
class ReservoirProblem : public GetPropType<TypeTag, Properties::BaseProblem>
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
2020-06-10 06:49:42 -05:00
|
|
|
using ParentType = GetPropType<TypeTag, Properties::BaseProblem>;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using GridView = GetPropType<TypeTag, Properties::GridView>;
|
|
|
|
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
|
|
|
|
using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
|
|
|
|
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2013-11-29 09:33:46 -06:00
|
|
|
// Grid and world dimension
|
|
|
|
enum { dim = GridView::dimension };
|
|
|
|
enum { dimWorld = GridView::dimensionworld };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// copy some indices for convenience
|
2013-11-29 09:33:46 -06:00
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
|
|
|
enum { numComponents = FluidSystem::numComponents };
|
2014-04-03 09:59:48 -05:00
|
|
|
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
|
|
|
|
enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
|
|
|
|
enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
|
|
|
|
enum { gasCompIdx = FluidSystem::gasCompIdx };
|
|
|
|
enum { oilCompIdx = FluidSystem::oilCompIdx };
|
|
|
|
enum { waterCompIdx = FluidSystem::waterCompIdx };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using Model = GetPropType<TypeTag, Properties::Model>;
|
|
|
|
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
|
|
|
|
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
|
|
|
|
using EqVector = GetPropType<TypeTag, Properties::EqVector>;
|
|
|
|
using RateVector = GetPropType<TypeTag, Properties::RateVector>;
|
|
|
|
using BoundaryRateVector = GetPropType<TypeTag, Properties::BoundaryRateVector>;
|
|
|
|
using Constraints = GetPropType<TypeTag, Properties::Constraints>;
|
|
|
|
using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
|
|
|
|
using Simulator = GetPropType<TypeTag, Properties::Simulator>;
|
|
|
|
using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
|
|
|
|
|
|
|
|
using CoordScalar = typename GridView::ctype;
|
|
|
|
using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
|
|
|
|
using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
|
|
|
|
using PhaseVector = Dune::FieldVector<Scalar, numPhases>;
|
|
|
|
|
|
|
|
using InitialFluidState = Opm::CompositionalFluidState<Scalar,
|
|
|
|
FluidSystem,
|
|
|
|
/*enableEnthalpy=*/true>;
|
2015-12-31 06:20:56 -06:00
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
public:
|
|
|
|
/*!
|
|
|
|
* \copydoc Doxygen::defaultProblemConstructor
|
|
|
|
*/
|
2016-11-07 08:14:07 -06:00
|
|
|
ReservoirProblem(Simulator& simulator)
|
2014-04-14 13:32:30 -05:00
|
|
|
: ParentType(simulator)
|
2014-08-06 09:31:48 -05:00
|
|
|
{ }
|
2014-07-09 04:57:48 -05:00
|
|
|
|
2014-08-06 09:31:48 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::finishInit
|
|
|
|
*/
|
2014-07-09 04:57:48 -05:00
|
|
|
void finishInit()
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
2014-07-09 04:57:48 -05:00
|
|
|
ParentType::finishInit();
|
|
|
|
|
2013-09-23 13:25:58 -05:00
|
|
|
temperature_ = EWOMS_GET_PARAM(TypeTag, Scalar, Temperature);
|
|
|
|
maxDepth_ = EWOMS_GET_PARAM(TypeTag, Scalar, MaxDepth);
|
2015-12-31 06:20:56 -06:00
|
|
|
wellWidth_ = EWOMS_GET_PARAM(TypeTag, Scalar, WellWidth);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2015-02-03 06:58:18 -06:00
|
|
|
std::vector<std::pair<Scalar, Scalar> > Bo = {
|
2014-07-09 04:57:48 -05:00
|
|
|
{ 101353, 1.062 },
|
|
|
|
{ 1.82504e+06, 1.15 },
|
|
|
|
{ 3.54873e+06, 1.207 },
|
|
|
|
{ 6.99611e+06, 1.295 },
|
|
|
|
{ 1.38909e+07, 1.435 },
|
|
|
|
{ 1.73382e+07, 1.5 },
|
|
|
|
{ 2.07856e+07, 1.565 },
|
|
|
|
{ 2.76804e+07, 1.695 },
|
|
|
|
{ 3.45751e+07, 1.827 }
|
|
|
|
};
|
|
|
|
std::vector<std::pair<Scalar, Scalar> > muo = {
|
|
|
|
{ 101353, 0.00104 },
|
|
|
|
{ 1.82504e+06, 0.000975 },
|
|
|
|
{ 3.54873e+06, 0.00091 },
|
|
|
|
{ 6.99611e+06, 0.00083 },
|
|
|
|
{ 1.38909e+07, 0.000695 },
|
|
|
|
{ 1.73382e+07, 0.000641 },
|
|
|
|
{ 2.07856e+07, 0.000594 },
|
|
|
|
{ 2.76804e+07, 0.00051 },
|
|
|
|
{ 3.45751e+07, 0.000449 }
|
|
|
|
};
|
|
|
|
std::vector<std::pair<Scalar, Scalar> > Rs = {
|
|
|
|
{ 101353, 0.178108 },
|
|
|
|
{ 1.82504e+06, 16.1187 },
|
|
|
|
{ 3.54873e+06, 32.0594 },
|
|
|
|
{ 6.99611e+06, 66.0779 },
|
|
|
|
{ 1.38909e+07, 113.276 },
|
|
|
|
{ 1.73382e+07, 138.033 },
|
|
|
|
{ 2.07856e+07, 165.64 },
|
|
|
|
{ 2.76804e+07, 226.197 },
|
|
|
|
{ 3.45751e+07, 288.178 }
|
|
|
|
};
|
|
|
|
std::vector<std::pair<Scalar, Scalar> > Bg = {
|
|
|
|
{ 101353, 0.93576 },
|
|
|
|
{ 1.82504e+06, 0.0678972 },
|
|
|
|
{ 3.54873e+06, 0.0352259 },
|
|
|
|
{ 6.99611e+06, 0.0179498 },
|
|
|
|
{ 1.38909e+07, 0.00906194 },
|
|
|
|
{ 1.73382e+07, 0.00726527 },
|
|
|
|
{ 2.07856e+07, 0.00606375 },
|
|
|
|
{ 2.76804e+07, 0.00455343 },
|
|
|
|
{ 3.45751e+07, 0.00364386 },
|
|
|
|
{ 6.21542e+07, 0.00216723 }
|
|
|
|
};
|
|
|
|
std::vector<std::pair<Scalar, Scalar> > mug = {
|
|
|
|
{ 101353, 8e-06 },
|
|
|
|
{ 1.82504e+06, 9.6e-06 },
|
|
|
|
{ 3.54873e+06, 1.12e-05 },
|
|
|
|
{ 6.99611e+06, 1.4e-05 },
|
|
|
|
{ 1.38909e+07, 1.89e-05 },
|
|
|
|
{ 1.73382e+07, 2.08e-05 },
|
|
|
|
{ 2.07856e+07, 2.28e-05 },
|
|
|
|
{ 2.76804e+07, 2.68e-05 },
|
|
|
|
{ 3.45751e+07, 3.09e-05 },
|
|
|
|
{ 6.21542e+07, 4.7e-05 }
|
|
|
|
};
|
|
|
|
|
2015-09-29 04:48:30 -05:00
|
|
|
Scalar rhoRefO = 786.0; // [kg]
|
|
|
|
Scalar rhoRefG = 0.97; // [kg]
|
|
|
|
Scalar rhoRefW = 1037.0; // [kg]
|
2015-12-31 06:20:56 -06:00
|
|
|
FluidSystem::initBegin(/*numPvtRegions=*/1);
|
|
|
|
FluidSystem::setEnableDissolvedGas(true);
|
|
|
|
FluidSystem::setEnableVaporizedOil(false);
|
2015-09-29 04:48:30 -05:00
|
|
|
FluidSystem::setReferenceDensities(rhoRefO, rhoRefW, rhoRefG, /*regionIdx=*/0);
|
|
|
|
|
|
|
|
Opm::GasPvtMultiplexer<Scalar> *gasPvt = new Opm::GasPvtMultiplexer<Scalar>;
|
2021-05-05 14:59:36 -05:00
|
|
|
gasPvt->setApproach(GasPvtApproach::DryGasPvt);
|
|
|
|
auto& dryGasPvt = gasPvt->template getRealPvt<GasPvtApproach::DryGasPvt>();
|
2015-09-29 04:48:30 -05:00
|
|
|
dryGasPvt.setNumRegions(/*numPvtRegion=*/1);
|
|
|
|
dryGasPvt.setReferenceDensities(/*regionIdx=*/0, rhoRefO, rhoRefG, rhoRefW);
|
|
|
|
dryGasPvt.setGasFormationVolumeFactor(/*regionIdx=*/0, Bg);
|
|
|
|
dryGasPvt.setGasViscosity(/*regionIdx=*/0, mug);
|
|
|
|
|
|
|
|
Opm::OilPvtMultiplexer<Scalar> *oilPvt = new Opm::OilPvtMultiplexer<Scalar>;
|
2021-05-05 14:59:36 -05:00
|
|
|
oilPvt->setApproach(OilPvtApproach::LiveOilPvt);
|
|
|
|
auto& liveOilPvt = oilPvt->template getRealPvt<OilPvtApproach::LiveOilPvt>();
|
2015-09-29 04:48:30 -05:00
|
|
|
liveOilPvt.setNumRegions(/*numPvtRegion=*/1);
|
|
|
|
liveOilPvt.setReferenceDensities(/*regionIdx=*/0, rhoRefO, rhoRefG, rhoRefW);
|
|
|
|
liveOilPvt.setSaturatedOilGasDissolutionFactor(/*regionIdx=*/0, Rs);
|
|
|
|
liveOilPvt.setSaturatedOilFormationVolumeFactor(/*regionIdx=*/0, Bo);
|
|
|
|
liveOilPvt.setSaturatedOilViscosity(/*regionIdx=*/0, muo);
|
|
|
|
|
|
|
|
Opm::WaterPvtMultiplexer<Scalar> *waterPvt = new Opm::WaterPvtMultiplexer<Scalar>;
|
2021-05-05 14:59:36 -05:00
|
|
|
waterPvt->setApproach(WaterPvtApproach::ConstantCompressibilityWaterPvt);
|
|
|
|
auto& ccWaterPvt = waterPvt->template getRealPvt<WaterPvtApproach::ConstantCompressibilityWaterPvt>();
|
2015-09-29 04:48:30 -05:00
|
|
|
ccWaterPvt.setNumRegions(/*numPvtRegions=*/1);
|
|
|
|
ccWaterPvt.setReferenceDensities(/*regionIdx=*/0, rhoRefO, rhoRefG, rhoRefW);
|
|
|
|
ccWaterPvt.setViscosity(/*regionIdx=*/0, 9.6e-4);
|
|
|
|
ccWaterPvt.setCompressibility(/*regionIdx=*/0, 1.450377e-10);
|
|
|
|
|
2016-01-04 08:32:55 -06:00
|
|
|
gasPvt->initEnd();
|
|
|
|
oilPvt->initEnd();
|
2015-09-29 04:48:30 -05:00
|
|
|
waterPvt->initEnd();
|
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using GasPvtSharedPtr = std::shared_ptr<Opm::GasPvtMultiplexer<Scalar> >;
|
2015-02-03 06:58:18 -06:00
|
|
|
FluidSystem::setGasPvt(GasPvtSharedPtr(gasPvt));
|
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using OilPvtSharedPtr = std::shared_ptr<Opm::OilPvtMultiplexer<Scalar> >;
|
2015-02-03 06:58:18 -06:00
|
|
|
FluidSystem::setOilPvt(OilPvtSharedPtr(oilPvt));
|
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using WaterPvtSharedPtr = std::shared_ptr<Opm::WaterPvtMultiplexer<Scalar> >;
|
2015-02-03 06:58:18 -06:00
|
|
|
FluidSystem::setWaterPvt(WaterPvtSharedPtr(waterPvt));
|
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
FluidSystem::initEnd();
|
|
|
|
|
2014-07-09 04:57:48 -05:00
|
|
|
pReservoir_ = 330e5;
|
2013-09-23 11:56:30 -05:00
|
|
|
layerBottom_ = 22.0;
|
|
|
|
|
|
|
|
// intrinsic permeabilities
|
|
|
|
fineK_ = this->toDimMatrix_(1e-12);
|
|
|
|
coarseK_ = this->toDimMatrix_(1e-11);
|
|
|
|
|
|
|
|
// porosities
|
|
|
|
finePorosity_ = 0.2;
|
|
|
|
coarsePorosity_ = 0.3;
|
|
|
|
|
2015-11-18 04:54:35 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2013-09-23 11:56:30 -05:00
|
|
|
fineMaterialParams_.setPcMinSat(phaseIdx, 0.0);
|
|
|
|
fineMaterialParams_.setPcMaxSat(phaseIdx, 0.0);
|
|
|
|
|
|
|
|
coarseMaterialParams_.setPcMinSat(phaseIdx, 0.0);
|
|
|
|
coarseMaterialParams_.setPcMaxSat(phaseIdx, 0.0);
|
|
|
|
}
|
|
|
|
|
2013-11-06 07:50:01 -06:00
|
|
|
// wrap up the initialization of the material law's parameters
|
|
|
|
fineMaterialParams_.finalize();
|
|
|
|
coarseMaterialParams_.finalize();
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2016-07-12 12:12:22 -05:00
|
|
|
materialParams_.resize(this->model().numGridDof());
|
|
|
|
ElementContext elemCtx(this->simulator());
|
|
|
|
auto eIt = this->simulator().gridView().template begin<0>();
|
|
|
|
const auto& eEndIt = this->simulator().gridView().template end<0>();
|
|
|
|
for (; eIt != eEndIt; ++eIt) {
|
|
|
|
elemCtx.updateStencil(*eIt);
|
2016-11-07 08:14:07 -06:00
|
|
|
size_t nDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
|
|
|
|
for (unsigned dofIdx = 0; dofIdx < nDof; ++ dofIdx) {
|
|
|
|
unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
|
|
|
const GlobalPosition& pos = elemCtx.pos(dofIdx, /*timeIdx=*/0);
|
2016-07-12 12:12:22 -05:00
|
|
|
|
|
|
|
if (isFineMaterial_(pos))
|
|
|
|
materialParams_[globalDofIdx] = &fineMaterialParams_;
|
|
|
|
else
|
|
|
|
materialParams_[globalDofIdx] = &coarseMaterialParams_;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
initFluidState_();
|
2015-12-31 06:20:56 -06:00
|
|
|
|
|
|
|
// start the first ("settle down") episode for 100 days
|
|
|
|
this->simulator().startNextEpisode(100.0*24*60*60);
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
Implement the element centered finite volume spatial discretization
This makes eWoms multi-discretization capable. Along the way, this
fixes some bugs and does a medium sized reorganization of the source tree.
This is a squashed patch of the following commits:
--------
1st commit message:
add initial version of the element centered finite volume discretization
currently, it is a misnomer as it is just a copy of the vertex
centered discretization plus some renames...
--------
2nd commit message:
rename [VE]cfvModel -> [VE]cfvDiscretization
--------
3rd commit message:
ecfv: prelimary changes required to make it compile
but not work yet...
--------
4th commit message:
Rename *FvElementGeometry to *Stencil
"Stencil" seems to be the standard expression for this concept...
(also, it is not specific to finite volume methods and is shorter.)
--------
5th commit message:
refactor the stencil class for the element centered finite volume discretization
--------
6th commit message:
ECFV: some work on the stencil class
--------
7th commit message:
ECFV: make the boundary handling code compile
--------
8th commit message:
rename elemContext() to elementContext()
--------
9th commit message:
ECFV: make the VTK output modules compile
--------
10th commit message:
stencil: introduce the concept of primary DOFs
also save an vector of all element pointers in the stencil.
--------
11th commit message:
ECFV: try to fix assembly; add missing timeIdx arguments to the num*() methods
--------
12th commit message:
ECFV: fix stupid mistake in the assembler
--------
13th commit message:
ECFV: remove a few implicit DOF == vertex assumptions
the black-oil example now runs without valgrind complaints until it encounters
a negative oil mole fraction.
--------
14th commit message:
VCFV: make everything compile again
all vertex centered FV examples should now work again...
--------
15th commit message:
rename [ev]cfvmodel.hh to [ev]cfvdiscretization.hh
the classes have already been renamed.
--------
16th commit message:
ECFV: make it work to the point where it can write out the initial solution.
--------
17th commit message:
ECFV: make it work
the local residual/jacobian needed some work in distinguishing primary
and secondary DOFs and there was an minor issue with the serialization
code.
for some reason, it seems still not correct. (-> convergence is too slow.)
--------
18th commit message:
VCFV: make it compile for the black oil model again
--------
19th commit message:
VCFV: make it compile with the remaining models again
--------
20th commit message:
flash model: make it work with ECFV
although this breaks its compatibility with VCFV. (-> next commit)
--------
21st commit message:
adapt the VCFV to make it compatible with the flash model again
--------
22nd commit message:
make all models compile with VCFV again
--------
23rd commit message:
VCFV: more cleanups of the stencil
VcfvStencil now does not have any public attributes anymore. TODO: do
not export attributes in the SubControlVolume and SubControlVolumeFace
classes.
--------
24th commit message:
VCFV: actually update the element pointer
--------
25th commit message:
change the blackoil model back to ECFV
--------
26th commit message:
immiscible model: make it compatible with the ECFV discretization
--------
27th commit message:
PVS model: make it work with ECFV
--------
28th commit message:
NCP model: make it work with ECFV
--------
29th commit message:
rename Vcfv*VelocityModule to *VelocityModule
--------
30th commit message:
richards model: make it work with ECFV
--------
31st commit message:
unify the ECFV and the VCFV VTK output modules
and other cleanups
--------
32nd commit message:
unify the common code of the VCFV and the ECFV disctretizations
--------
33rd commit message:
unify the element contexts between element and vertex centered finite volumes
--------
34th commit message:
unify the local jacobian class of the finite volume discretizations
--------
35th commit message:
replace [VE]vcf(LocalResidual|ElementContext|BoundaryContext|ConstraintsContext) by generic code
--------
36th commit message:
replace the [EV]cfvLocalResidual by generic code
--------
37th commit message:
unify the MultiPhaseProblem and Problem classes, introduce NullBorderListCreator
--------
38th commit message:
remove the discretization specific boundary context
--------
39th commit message:
unify the [EV]cfvDiscretization classes
--------
40th commit message:
Unify [EV]cfvMultiPhaseFluxVariables
--------
41st commit message:
Unify the [EC]cfvNewton* classes
--------
42nd commit message:
Unify [EV]cfvVolumeVariables
--------
43rd commit message:
unify [EV]cfvAssembler
--------
44th commit message:
unified flux variables: fix stupid mistake when calculating pressure gradients
--------
45th commit message:
unify what's to unify for the [EV]CFV properties
--------
46th commit message:
make the method to calculate gradients and values at flux approximation points changeable
Currently, this is used by the vertex centered finite volume method to
be able to use P1-finite element gradients instead of two-point
ones...
--------
47th commit message:
make the restart code work correctly, use the correct DofMapper for VCFV
--------
48th commit message:
actually use the gradient calculator in a model
the immiscible model in this case
--------
49th commit message:
move some files around to where they belong, use the new gradient calculation code in all models
TODO: proper handling of boundary gradients
--------
50th commit message:
fix the stokes model
currently it only works with the vertex centered finite volume
discretization, but the plan is to soon move it to a staggered grid
scheme anyway...
--------
51st commit message:
move all models back to using the vertex centered finite volume discretization by default
--------
52nd commit message:
models: some variable renames and documentation fixes
- scv -> dof
- vert -> dof
- vertex -> dof
- replace 'VCFV'
- fix some typos
--------
53rd commit message:
don't expect UG anymore
since it is quite non-free and hard to get. we now use ALUGrid instead!
--------
54th commit message:
temporarily disable jacobian recycling
--------
55th commit message:
fix writing/reading restart files using the generic code
--------
56th commit message:
fix bug where fluxes were only counted once in the stencil
this only affected the vertex centered finite volumes discretization...
--------
57th commit message:
boundary gradients: use the center of the sub-control volume adjacent to a boundary segment
--------
58th commit message:
make it compile on GCC
--------
59th commit message:
get rid of most hacks
for this, partial reassemble and jacobian recycling was brought
back. For the this and the remaining stuff the main trick is the
introduction of the GridCommHandleFactory concept which constructs
communication handles suited for the respective spatial
discretization...
--------
60th commit message:
fix a few annoying bugs
first, default the convergence criterion for the linear solver did not
honor the initial residual which lead to linear solver breakdowns,
then some debugging code was left in the discrete fracture model and
then there was a bug in the TP gradient approximation class...
this has the consequence that we need a new reference solution for the
discrete fracture problem...
--------
61st commit message:
iterative linear solver: remove the code for the non-default convergence criteria
--------
62nd commit message:
provide the FE cache instead of the local FE
this fixes a segfault in the stokes model caused by the fact that the
local FE was not initialized at this point.
--------
63rd commit message:
(Navier-)Stokes: fix bug due to the transition to unit normals
now, all tests pass for this branch. The only things which need to be
fixed are some annoying performance regressions compared to master and
some bug in the splices feature of the property system...
--------
64th commit message:
some fix for the local residual of the immiscible model
--------
65th commit message:
Navier-Stokes: implement SCV center gradients
There seems to be a bug in the previous implementation (the jacobian
inverse transposed is evaluated using the local, not the global
geometry), so the reference solution for the stokes2c test problem has
also been updated...
--------
66th commit message:
remove the ALUGrid specialization of the LensGridCreator and the YaspGrid one for the fingerproblem
using different grid seems to sometimes cause a different vertex
order, which in turn causes the respective test to fail if the
reference solution was computed using the other grid...
--------
67th commit message:
VCFV: use the correct BorderListCreator
this makes MPI parallel computations work again. apart from
performance regressions, this branch does not exhibit any known
regressions compared to master anymore...
--------
68th commit message:
make verything compile with the element centered finite volume discretization
except the Navier-Stokes and the two-phase DFM models, of course...
--------
69th commit message:
minor fixes
- make the navier-stokes model slighly more generic by using the
proper (in,ex)teriorIndex() methods on sub-control volumes
- make the signature of the calculateValue() template method of the
common two-point gradient approximator match the one of the vertex
centered finite volume one
--------
70th commit message:
fix fallout from the Big Rebase
--------
71st commit message:
ECFV: some bugs in the boundary
--------
72nd commit message:
make computeFlux() compute area-specific quantities
--------
73rd commit message:
fix more bugs in the element centered FV discretization
now eWoms should match Dumux pretty closely...
--------
74th commit message:
coalesce the common code of the multi phase porous medium models into "MultiPhaseBaseModel"
--------
75th commit message:
update reference solutions
these were changed because of the screw-up with the area of boundary
segments...
--------
76th commit message:
rename "ImplicitBase" to "FvBase"
because in eWoms, everything is implicit and these are currently the
base classes for all finite volume discretizations.
--------
77th commit message:
make the spatial discretization selectable using a splice
This requires an opm-core with a the patches from
https://github.com/OPM/opm-core/pull/446 merged...
--------
78th commit message:
rename the properties used for splices to *Splice
--------
79th commit message:
move the files in 'tests/models' to 'tests'
since 'tests' was empty except for the 'models' subdirectory...
--------
80th commit message:
improve and fix the tutorial
--------
81st commit message:
remove the -fno-strict-aliasing flag from the provided option files
seems like recent versions of Dune have been adapted...
--------
82nd commit message:
also compile all CO2 injection simulations using the element centered finite volume discretization
--------
83rd commit message:
PVS model: make it work properly with the element-centered finite volume discretiation
because DOF != number of vertices
2013-12-12 05:52:44 -06:00
|
|
|
* \copydoc FvBaseMultiPhaseProblem::registerParameters
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
static void registerParameters()
|
|
|
|
{
|
|
|
|
ParentType::registerParameters();
|
|
|
|
|
2013-11-29 09:33:46 -06:00
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, Scalar, Temperature,
|
|
|
|
"The temperature [K] in the reservoir");
|
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, Scalar, MaxDepth,
|
|
|
|
"The maximum depth [m] of the reservoir");
|
2015-12-31 06:20:56 -06:00
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, Scalar, WellWidth,
|
|
|
|
"The width of producer/injector wells as a fraction of the width"
|
|
|
|
" of the spatial domain");
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
|
|
|
|
2014-07-22 05:41:56 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::name
|
|
|
|
*/
|
|
|
|
std::string name() const
|
2015-12-31 06:20:56 -06:00
|
|
|
{ return std::string("reservoir_") + Model::name() + "_" + Model::discretizationName(); }
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::endEpisode
|
|
|
|
*/
|
|
|
|
void endEpisode()
|
|
|
|
{
|
|
|
|
// in the second episode, the actual work is done (the first is "settle down"
|
|
|
|
// episode). we need to use a pretty short initial time step here as the change
|
|
|
|
// in conditions is quite abrupt.
|
|
|
|
this->simulator().startNextEpisode(1e100);
|
|
|
|
this->simulator().setTimeStepSize(5.0);
|
|
|
|
}
|
2014-07-22 05:41:56 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::endTimeStep
|
|
|
|
*/
|
|
|
|
void endTimeStep()
|
|
|
|
{
|
|
|
|
#ifndef NDEBUG
|
|
|
|
// checkConservativeness() does not include the effect of constraints, so we
|
|
|
|
// disable it for this problem...
|
|
|
|
//this->model().checkConservativeness();
|
|
|
|
|
|
|
|
// Calculate storage terms
|
|
|
|
EqVector storage;
|
|
|
|
this->model().globalStorage(storage);
|
|
|
|
|
|
|
|
// Write mass balance information for rank 0
|
|
|
|
if (this->gridView().comm().rank() == 0) {
|
|
|
|
std::cout << "Storage: " << storage << std::endl << std::flush;
|
|
|
|
}
|
|
|
|
#endif // NDEBUG
|
|
|
|
}
|
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
/*!
|
Implement the element centered finite volume spatial discretization
This makes eWoms multi-discretization capable. Along the way, this
fixes some bugs and does a medium sized reorganization of the source tree.
This is a squashed patch of the following commits:
--------
1st commit message:
add initial version of the element centered finite volume discretization
currently, it is a misnomer as it is just a copy of the vertex
centered discretization plus some renames...
--------
2nd commit message:
rename [VE]cfvModel -> [VE]cfvDiscretization
--------
3rd commit message:
ecfv: prelimary changes required to make it compile
but not work yet...
--------
4th commit message:
Rename *FvElementGeometry to *Stencil
"Stencil" seems to be the standard expression for this concept...
(also, it is not specific to finite volume methods and is shorter.)
--------
5th commit message:
refactor the stencil class for the element centered finite volume discretization
--------
6th commit message:
ECFV: some work on the stencil class
--------
7th commit message:
ECFV: make the boundary handling code compile
--------
8th commit message:
rename elemContext() to elementContext()
--------
9th commit message:
ECFV: make the VTK output modules compile
--------
10th commit message:
stencil: introduce the concept of primary DOFs
also save an vector of all element pointers in the stencil.
--------
11th commit message:
ECFV: try to fix assembly; add missing timeIdx arguments to the num*() methods
--------
12th commit message:
ECFV: fix stupid mistake in the assembler
--------
13th commit message:
ECFV: remove a few implicit DOF == vertex assumptions
the black-oil example now runs without valgrind complaints until it encounters
a negative oil mole fraction.
--------
14th commit message:
VCFV: make everything compile again
all vertex centered FV examples should now work again...
--------
15th commit message:
rename [ev]cfvmodel.hh to [ev]cfvdiscretization.hh
the classes have already been renamed.
--------
16th commit message:
ECFV: make it work to the point where it can write out the initial solution.
--------
17th commit message:
ECFV: make it work
the local residual/jacobian needed some work in distinguishing primary
and secondary DOFs and there was an minor issue with the serialization
code.
for some reason, it seems still not correct. (-> convergence is too slow.)
--------
18th commit message:
VCFV: make it compile for the black oil model again
--------
19th commit message:
VCFV: make it compile with the remaining models again
--------
20th commit message:
flash model: make it work with ECFV
although this breaks its compatibility with VCFV. (-> next commit)
--------
21st commit message:
adapt the VCFV to make it compatible with the flash model again
--------
22nd commit message:
make all models compile with VCFV again
--------
23rd commit message:
VCFV: more cleanups of the stencil
VcfvStencil now does not have any public attributes anymore. TODO: do
not export attributes in the SubControlVolume and SubControlVolumeFace
classes.
--------
24th commit message:
VCFV: actually update the element pointer
--------
25th commit message:
change the blackoil model back to ECFV
--------
26th commit message:
immiscible model: make it compatible with the ECFV discretization
--------
27th commit message:
PVS model: make it work with ECFV
--------
28th commit message:
NCP model: make it work with ECFV
--------
29th commit message:
rename Vcfv*VelocityModule to *VelocityModule
--------
30th commit message:
richards model: make it work with ECFV
--------
31st commit message:
unify the ECFV and the VCFV VTK output modules
and other cleanups
--------
32nd commit message:
unify the common code of the VCFV and the ECFV disctretizations
--------
33rd commit message:
unify the element contexts between element and vertex centered finite volumes
--------
34th commit message:
unify the local jacobian class of the finite volume discretizations
--------
35th commit message:
replace [VE]vcf(LocalResidual|ElementContext|BoundaryContext|ConstraintsContext) by generic code
--------
36th commit message:
replace the [EV]cfvLocalResidual by generic code
--------
37th commit message:
unify the MultiPhaseProblem and Problem classes, introduce NullBorderListCreator
--------
38th commit message:
remove the discretization specific boundary context
--------
39th commit message:
unify the [EV]cfvDiscretization classes
--------
40th commit message:
Unify [EV]cfvMultiPhaseFluxVariables
--------
41st commit message:
Unify the [EC]cfvNewton* classes
--------
42nd commit message:
Unify [EV]cfvVolumeVariables
--------
43rd commit message:
unify [EV]cfvAssembler
--------
44th commit message:
unified flux variables: fix stupid mistake when calculating pressure gradients
--------
45th commit message:
unify what's to unify for the [EV]CFV properties
--------
46th commit message:
make the method to calculate gradients and values at flux approximation points changeable
Currently, this is used by the vertex centered finite volume method to
be able to use P1-finite element gradients instead of two-point
ones...
--------
47th commit message:
make the restart code work correctly, use the correct DofMapper for VCFV
--------
48th commit message:
actually use the gradient calculator in a model
the immiscible model in this case
--------
49th commit message:
move some files around to where they belong, use the new gradient calculation code in all models
TODO: proper handling of boundary gradients
--------
50th commit message:
fix the stokes model
currently it only works with the vertex centered finite volume
discretization, but the plan is to soon move it to a staggered grid
scheme anyway...
--------
51st commit message:
move all models back to using the vertex centered finite volume discretization by default
--------
52nd commit message:
models: some variable renames and documentation fixes
- scv -> dof
- vert -> dof
- vertex -> dof
- replace 'VCFV'
- fix some typos
--------
53rd commit message:
don't expect UG anymore
since it is quite non-free and hard to get. we now use ALUGrid instead!
--------
54th commit message:
temporarily disable jacobian recycling
--------
55th commit message:
fix writing/reading restart files using the generic code
--------
56th commit message:
fix bug where fluxes were only counted once in the stencil
this only affected the vertex centered finite volumes discretization...
--------
57th commit message:
boundary gradients: use the center of the sub-control volume adjacent to a boundary segment
--------
58th commit message:
make it compile on GCC
--------
59th commit message:
get rid of most hacks
for this, partial reassemble and jacobian recycling was brought
back. For the this and the remaining stuff the main trick is the
introduction of the GridCommHandleFactory concept which constructs
communication handles suited for the respective spatial
discretization...
--------
60th commit message:
fix a few annoying bugs
first, default the convergence criterion for the linear solver did not
honor the initial residual which lead to linear solver breakdowns,
then some debugging code was left in the discrete fracture model and
then there was a bug in the TP gradient approximation class...
this has the consequence that we need a new reference solution for the
discrete fracture problem...
--------
61st commit message:
iterative linear solver: remove the code for the non-default convergence criteria
--------
62nd commit message:
provide the FE cache instead of the local FE
this fixes a segfault in the stokes model caused by the fact that the
local FE was not initialized at this point.
--------
63rd commit message:
(Navier-)Stokes: fix bug due to the transition to unit normals
now, all tests pass for this branch. The only things which need to be
fixed are some annoying performance regressions compared to master and
some bug in the splices feature of the property system...
--------
64th commit message:
some fix for the local residual of the immiscible model
--------
65th commit message:
Navier-Stokes: implement SCV center gradients
There seems to be a bug in the previous implementation (the jacobian
inverse transposed is evaluated using the local, not the global
geometry), so the reference solution for the stokes2c test problem has
also been updated...
--------
66th commit message:
remove the ALUGrid specialization of the LensGridCreator and the YaspGrid one for the fingerproblem
using different grid seems to sometimes cause a different vertex
order, which in turn causes the respective test to fail if the
reference solution was computed using the other grid...
--------
67th commit message:
VCFV: use the correct BorderListCreator
this makes MPI parallel computations work again. apart from
performance regressions, this branch does not exhibit any known
regressions compared to master anymore...
--------
68th commit message:
make verything compile with the element centered finite volume discretization
except the Navier-Stokes and the two-phase DFM models, of course...
--------
69th commit message:
minor fixes
- make the navier-stokes model slighly more generic by using the
proper (in,ex)teriorIndex() methods on sub-control volumes
- make the signature of the calculateValue() template method of the
common two-point gradient approximator match the one of the vertex
centered finite volume one
--------
70th commit message:
fix fallout from the Big Rebase
--------
71st commit message:
ECFV: some bugs in the boundary
--------
72nd commit message:
make computeFlux() compute area-specific quantities
--------
73rd commit message:
fix more bugs in the element centered FV discretization
now eWoms should match Dumux pretty closely...
--------
74th commit message:
coalesce the common code of the multi phase porous medium models into "MultiPhaseBaseModel"
--------
75th commit message:
update reference solutions
these were changed because of the screw-up with the area of boundary
segments...
--------
76th commit message:
rename "ImplicitBase" to "FvBase"
because in eWoms, everything is implicit and these are currently the
base classes for all finite volume discretizations.
--------
77th commit message:
make the spatial discretization selectable using a splice
This requires an opm-core with a the patches from
https://github.com/OPM/opm-core/pull/446 merged...
--------
78th commit message:
rename the properties used for splices to *Splice
--------
79th commit message:
move the files in 'tests/models' to 'tests'
since 'tests' was empty except for the 'models' subdirectory...
--------
80th commit message:
improve and fix the tutorial
--------
81st commit message:
remove the -fno-strict-aliasing flag from the provided option files
seems like recent versions of Dune have been adapted...
--------
82nd commit message:
also compile all CO2 injection simulations using the element centered finite volume discretization
--------
83rd commit message:
PVS model: make it work properly with the element-centered finite volume discretiation
because DOF != number of vertices
2013-12-12 05:52:44 -06:00
|
|
|
* \copydoc FvBaseMultiPhaseProblem::intrinsicPermeability
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
|
|
|
* For this problem, a layer with high permability is located
|
|
|
|
* above one with low permeability.
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
const DimMatrix& intrinsicPermeability(const Context& context, unsigned spaceIdx,
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned timeIdx) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
2016-11-07 08:14:07 -06:00
|
|
|
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
|
2013-09-23 11:56:30 -05:00
|
|
|
if (isFineMaterial_(pos))
|
|
|
|
return fineK_;
|
|
|
|
return coarseK_;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
Implement the element centered finite volume spatial discretization
This makes eWoms multi-discretization capable. Along the way, this
fixes some bugs and does a medium sized reorganization of the source tree.
This is a squashed patch of the following commits:
--------
1st commit message:
add initial version of the element centered finite volume discretization
currently, it is a misnomer as it is just a copy of the vertex
centered discretization plus some renames...
--------
2nd commit message:
rename [VE]cfvModel -> [VE]cfvDiscretization
--------
3rd commit message:
ecfv: prelimary changes required to make it compile
but not work yet...
--------
4th commit message:
Rename *FvElementGeometry to *Stencil
"Stencil" seems to be the standard expression for this concept...
(also, it is not specific to finite volume methods and is shorter.)
--------
5th commit message:
refactor the stencil class for the element centered finite volume discretization
--------
6th commit message:
ECFV: some work on the stencil class
--------
7th commit message:
ECFV: make the boundary handling code compile
--------
8th commit message:
rename elemContext() to elementContext()
--------
9th commit message:
ECFV: make the VTK output modules compile
--------
10th commit message:
stencil: introduce the concept of primary DOFs
also save an vector of all element pointers in the stencil.
--------
11th commit message:
ECFV: try to fix assembly; add missing timeIdx arguments to the num*() methods
--------
12th commit message:
ECFV: fix stupid mistake in the assembler
--------
13th commit message:
ECFV: remove a few implicit DOF == vertex assumptions
the black-oil example now runs without valgrind complaints until it encounters
a negative oil mole fraction.
--------
14th commit message:
VCFV: make everything compile again
all vertex centered FV examples should now work again...
--------
15th commit message:
rename [ev]cfvmodel.hh to [ev]cfvdiscretization.hh
the classes have already been renamed.
--------
16th commit message:
ECFV: make it work to the point where it can write out the initial solution.
--------
17th commit message:
ECFV: make it work
the local residual/jacobian needed some work in distinguishing primary
and secondary DOFs and there was an minor issue with the serialization
code.
for some reason, it seems still not correct. (-> convergence is too slow.)
--------
18th commit message:
VCFV: make it compile for the black oil model again
--------
19th commit message:
VCFV: make it compile with the remaining models again
--------
20th commit message:
flash model: make it work with ECFV
although this breaks its compatibility with VCFV. (-> next commit)
--------
21st commit message:
adapt the VCFV to make it compatible with the flash model again
--------
22nd commit message:
make all models compile with VCFV again
--------
23rd commit message:
VCFV: more cleanups of the stencil
VcfvStencil now does not have any public attributes anymore. TODO: do
not export attributes in the SubControlVolume and SubControlVolumeFace
classes.
--------
24th commit message:
VCFV: actually update the element pointer
--------
25th commit message:
change the blackoil model back to ECFV
--------
26th commit message:
immiscible model: make it compatible with the ECFV discretization
--------
27th commit message:
PVS model: make it work with ECFV
--------
28th commit message:
NCP model: make it work with ECFV
--------
29th commit message:
rename Vcfv*VelocityModule to *VelocityModule
--------
30th commit message:
richards model: make it work with ECFV
--------
31st commit message:
unify the ECFV and the VCFV VTK output modules
and other cleanups
--------
32nd commit message:
unify the common code of the VCFV and the ECFV disctretizations
--------
33rd commit message:
unify the element contexts between element and vertex centered finite volumes
--------
34th commit message:
unify the local jacobian class of the finite volume discretizations
--------
35th commit message:
replace [VE]vcf(LocalResidual|ElementContext|BoundaryContext|ConstraintsContext) by generic code
--------
36th commit message:
replace the [EV]cfvLocalResidual by generic code
--------
37th commit message:
unify the MultiPhaseProblem and Problem classes, introduce NullBorderListCreator
--------
38th commit message:
remove the discretization specific boundary context
--------
39th commit message:
unify the [EV]cfvDiscretization classes
--------
40th commit message:
Unify [EV]cfvMultiPhaseFluxVariables
--------
41st commit message:
Unify the [EC]cfvNewton* classes
--------
42nd commit message:
Unify [EV]cfvVolumeVariables
--------
43rd commit message:
unify [EV]cfvAssembler
--------
44th commit message:
unified flux variables: fix stupid mistake when calculating pressure gradients
--------
45th commit message:
unify what's to unify for the [EV]CFV properties
--------
46th commit message:
make the method to calculate gradients and values at flux approximation points changeable
Currently, this is used by the vertex centered finite volume method to
be able to use P1-finite element gradients instead of two-point
ones...
--------
47th commit message:
make the restart code work correctly, use the correct DofMapper for VCFV
--------
48th commit message:
actually use the gradient calculator in a model
the immiscible model in this case
--------
49th commit message:
move some files around to where they belong, use the new gradient calculation code in all models
TODO: proper handling of boundary gradients
--------
50th commit message:
fix the stokes model
currently it only works with the vertex centered finite volume
discretization, but the plan is to soon move it to a staggered grid
scheme anyway...
--------
51st commit message:
move all models back to using the vertex centered finite volume discretization by default
--------
52nd commit message:
models: some variable renames and documentation fixes
- scv -> dof
- vert -> dof
- vertex -> dof
- replace 'VCFV'
- fix some typos
--------
53rd commit message:
don't expect UG anymore
since it is quite non-free and hard to get. we now use ALUGrid instead!
--------
54th commit message:
temporarily disable jacobian recycling
--------
55th commit message:
fix writing/reading restart files using the generic code
--------
56th commit message:
fix bug where fluxes were only counted once in the stencil
this only affected the vertex centered finite volumes discretization...
--------
57th commit message:
boundary gradients: use the center of the sub-control volume adjacent to a boundary segment
--------
58th commit message:
make it compile on GCC
--------
59th commit message:
get rid of most hacks
for this, partial reassemble and jacobian recycling was brought
back. For the this and the remaining stuff the main trick is the
introduction of the GridCommHandleFactory concept which constructs
communication handles suited for the respective spatial
discretization...
--------
60th commit message:
fix a few annoying bugs
first, default the convergence criterion for the linear solver did not
honor the initial residual which lead to linear solver breakdowns,
then some debugging code was left in the discrete fracture model and
then there was a bug in the TP gradient approximation class...
this has the consequence that we need a new reference solution for the
discrete fracture problem...
--------
61st commit message:
iterative linear solver: remove the code for the non-default convergence criteria
--------
62nd commit message:
provide the FE cache instead of the local FE
this fixes a segfault in the stokes model caused by the fact that the
local FE was not initialized at this point.
--------
63rd commit message:
(Navier-)Stokes: fix bug due to the transition to unit normals
now, all tests pass for this branch. The only things which need to be
fixed are some annoying performance regressions compared to master and
some bug in the splices feature of the property system...
--------
64th commit message:
some fix for the local residual of the immiscible model
--------
65th commit message:
Navier-Stokes: implement SCV center gradients
There seems to be a bug in the previous implementation (the jacobian
inverse transposed is evaluated using the local, not the global
geometry), so the reference solution for the stokes2c test problem has
also been updated...
--------
66th commit message:
remove the ALUGrid specialization of the LensGridCreator and the YaspGrid one for the fingerproblem
using different grid seems to sometimes cause a different vertex
order, which in turn causes the respective test to fail if the
reference solution was computed using the other grid...
--------
67th commit message:
VCFV: use the correct BorderListCreator
this makes MPI parallel computations work again. apart from
performance regressions, this branch does not exhibit any known
regressions compared to master anymore...
--------
68th commit message:
make verything compile with the element centered finite volume discretization
except the Navier-Stokes and the two-phase DFM models, of course...
--------
69th commit message:
minor fixes
- make the navier-stokes model slighly more generic by using the
proper (in,ex)teriorIndex() methods on sub-control volumes
- make the signature of the calculateValue() template method of the
common two-point gradient approximator match the one of the vertex
centered finite volume one
--------
70th commit message:
fix fallout from the Big Rebase
--------
71st commit message:
ECFV: some bugs in the boundary
--------
72nd commit message:
make computeFlux() compute area-specific quantities
--------
73rd commit message:
fix more bugs in the element centered FV discretization
now eWoms should match Dumux pretty closely...
--------
74th commit message:
coalesce the common code of the multi phase porous medium models into "MultiPhaseBaseModel"
--------
75th commit message:
update reference solutions
these were changed because of the screw-up with the area of boundary
segments...
--------
76th commit message:
rename "ImplicitBase" to "FvBase"
because in eWoms, everything is implicit and these are currently the
base classes for all finite volume discretizations.
--------
77th commit message:
make the spatial discretization selectable using a splice
This requires an opm-core with a the patches from
https://github.com/OPM/opm-core/pull/446 merged...
--------
78th commit message:
rename the properties used for splices to *Splice
--------
79th commit message:
move the files in 'tests/models' to 'tests'
since 'tests' was empty except for the 'models' subdirectory...
--------
80th commit message:
improve and fix the tutorial
--------
81st commit message:
remove the -fno-strict-aliasing flag from the provided option files
seems like recent versions of Dune have been adapted...
--------
82nd commit message:
also compile all CO2 injection simulations using the element centered finite volume discretization
--------
83rd commit message:
PVS model: make it work properly with the element-centered finite volume discretiation
because DOF != number of vertices
2013-12-12 05:52:44 -06:00
|
|
|
* \copydoc FvBaseMultiPhaseProblem::porosity
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
Scalar porosity(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
2016-11-07 08:14:07 -06:00
|
|
|
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
|
2013-09-23 11:56:30 -05:00
|
|
|
if (isFineMaterial_(pos))
|
|
|
|
return finePorosity_;
|
|
|
|
return coarsePorosity_;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
Implement the element centered finite volume spatial discretization
This makes eWoms multi-discretization capable. Along the way, this
fixes some bugs and does a medium sized reorganization of the source tree.
This is a squashed patch of the following commits:
--------
1st commit message:
add initial version of the element centered finite volume discretization
currently, it is a misnomer as it is just a copy of the vertex
centered discretization plus some renames...
--------
2nd commit message:
rename [VE]cfvModel -> [VE]cfvDiscretization
--------
3rd commit message:
ecfv: prelimary changes required to make it compile
but not work yet...
--------
4th commit message:
Rename *FvElementGeometry to *Stencil
"Stencil" seems to be the standard expression for this concept...
(also, it is not specific to finite volume methods and is shorter.)
--------
5th commit message:
refactor the stencil class for the element centered finite volume discretization
--------
6th commit message:
ECFV: some work on the stencil class
--------
7th commit message:
ECFV: make the boundary handling code compile
--------
8th commit message:
rename elemContext() to elementContext()
--------
9th commit message:
ECFV: make the VTK output modules compile
--------
10th commit message:
stencil: introduce the concept of primary DOFs
also save an vector of all element pointers in the stencil.
--------
11th commit message:
ECFV: try to fix assembly; add missing timeIdx arguments to the num*() methods
--------
12th commit message:
ECFV: fix stupid mistake in the assembler
--------
13th commit message:
ECFV: remove a few implicit DOF == vertex assumptions
the black-oil example now runs without valgrind complaints until it encounters
a negative oil mole fraction.
--------
14th commit message:
VCFV: make everything compile again
all vertex centered FV examples should now work again...
--------
15th commit message:
rename [ev]cfvmodel.hh to [ev]cfvdiscretization.hh
the classes have already been renamed.
--------
16th commit message:
ECFV: make it work to the point where it can write out the initial solution.
--------
17th commit message:
ECFV: make it work
the local residual/jacobian needed some work in distinguishing primary
and secondary DOFs and there was an minor issue with the serialization
code.
for some reason, it seems still not correct. (-> convergence is too slow.)
--------
18th commit message:
VCFV: make it compile for the black oil model again
--------
19th commit message:
VCFV: make it compile with the remaining models again
--------
20th commit message:
flash model: make it work with ECFV
although this breaks its compatibility with VCFV. (-> next commit)
--------
21st commit message:
adapt the VCFV to make it compatible with the flash model again
--------
22nd commit message:
make all models compile with VCFV again
--------
23rd commit message:
VCFV: more cleanups of the stencil
VcfvStencil now does not have any public attributes anymore. TODO: do
not export attributes in the SubControlVolume and SubControlVolumeFace
classes.
--------
24th commit message:
VCFV: actually update the element pointer
--------
25th commit message:
change the blackoil model back to ECFV
--------
26th commit message:
immiscible model: make it compatible with the ECFV discretization
--------
27th commit message:
PVS model: make it work with ECFV
--------
28th commit message:
NCP model: make it work with ECFV
--------
29th commit message:
rename Vcfv*VelocityModule to *VelocityModule
--------
30th commit message:
richards model: make it work with ECFV
--------
31st commit message:
unify the ECFV and the VCFV VTK output modules
and other cleanups
--------
32nd commit message:
unify the common code of the VCFV and the ECFV disctretizations
--------
33rd commit message:
unify the element contexts between element and vertex centered finite volumes
--------
34th commit message:
unify the local jacobian class of the finite volume discretizations
--------
35th commit message:
replace [VE]vcf(LocalResidual|ElementContext|BoundaryContext|ConstraintsContext) by generic code
--------
36th commit message:
replace the [EV]cfvLocalResidual by generic code
--------
37th commit message:
unify the MultiPhaseProblem and Problem classes, introduce NullBorderListCreator
--------
38th commit message:
remove the discretization specific boundary context
--------
39th commit message:
unify the [EV]cfvDiscretization classes
--------
40th commit message:
Unify [EV]cfvMultiPhaseFluxVariables
--------
41st commit message:
Unify the [EC]cfvNewton* classes
--------
42nd commit message:
Unify [EV]cfvVolumeVariables
--------
43rd commit message:
unify [EV]cfvAssembler
--------
44th commit message:
unified flux variables: fix stupid mistake when calculating pressure gradients
--------
45th commit message:
unify what's to unify for the [EV]CFV properties
--------
46th commit message:
make the method to calculate gradients and values at flux approximation points changeable
Currently, this is used by the vertex centered finite volume method to
be able to use P1-finite element gradients instead of two-point
ones...
--------
47th commit message:
make the restart code work correctly, use the correct DofMapper for VCFV
--------
48th commit message:
actually use the gradient calculator in a model
the immiscible model in this case
--------
49th commit message:
move some files around to where they belong, use the new gradient calculation code in all models
TODO: proper handling of boundary gradients
--------
50th commit message:
fix the stokes model
currently it only works with the vertex centered finite volume
discretization, but the plan is to soon move it to a staggered grid
scheme anyway...
--------
51st commit message:
move all models back to using the vertex centered finite volume discretization by default
--------
52nd commit message:
models: some variable renames and documentation fixes
- scv -> dof
- vert -> dof
- vertex -> dof
- replace 'VCFV'
- fix some typos
--------
53rd commit message:
don't expect UG anymore
since it is quite non-free and hard to get. we now use ALUGrid instead!
--------
54th commit message:
temporarily disable jacobian recycling
--------
55th commit message:
fix writing/reading restart files using the generic code
--------
56th commit message:
fix bug where fluxes were only counted once in the stencil
this only affected the vertex centered finite volumes discretization...
--------
57th commit message:
boundary gradients: use the center of the sub-control volume adjacent to a boundary segment
--------
58th commit message:
make it compile on GCC
--------
59th commit message:
get rid of most hacks
for this, partial reassemble and jacobian recycling was brought
back. For the this and the remaining stuff the main trick is the
introduction of the GridCommHandleFactory concept which constructs
communication handles suited for the respective spatial
discretization...
--------
60th commit message:
fix a few annoying bugs
first, default the convergence criterion for the linear solver did not
honor the initial residual which lead to linear solver breakdowns,
then some debugging code was left in the discrete fracture model and
then there was a bug in the TP gradient approximation class...
this has the consequence that we need a new reference solution for the
discrete fracture problem...
--------
61st commit message:
iterative linear solver: remove the code for the non-default convergence criteria
--------
62nd commit message:
provide the FE cache instead of the local FE
this fixes a segfault in the stokes model caused by the fact that the
local FE was not initialized at this point.
--------
63rd commit message:
(Navier-)Stokes: fix bug due to the transition to unit normals
now, all tests pass for this branch. The only things which need to be
fixed are some annoying performance regressions compared to master and
some bug in the splices feature of the property system...
--------
64th commit message:
some fix for the local residual of the immiscible model
--------
65th commit message:
Navier-Stokes: implement SCV center gradients
There seems to be a bug in the previous implementation (the jacobian
inverse transposed is evaluated using the local, not the global
geometry), so the reference solution for the stokes2c test problem has
also been updated...
--------
66th commit message:
remove the ALUGrid specialization of the LensGridCreator and the YaspGrid one for the fingerproblem
using different grid seems to sometimes cause a different vertex
order, which in turn causes the respective test to fail if the
reference solution was computed using the other grid...
--------
67th commit message:
VCFV: use the correct BorderListCreator
this makes MPI parallel computations work again. apart from
performance regressions, this branch does not exhibit any known
regressions compared to master anymore...
--------
68th commit message:
make verything compile with the element centered finite volume discretization
except the Navier-Stokes and the two-phase DFM models, of course...
--------
69th commit message:
minor fixes
- make the navier-stokes model slighly more generic by using the
proper (in,ex)teriorIndex() methods on sub-control volumes
- make the signature of the calculateValue() template method of the
common two-point gradient approximator match the one of the vertex
centered finite volume one
--------
70th commit message:
fix fallout from the Big Rebase
--------
71st commit message:
ECFV: some bugs in the boundary
--------
72nd commit message:
make computeFlux() compute area-specific quantities
--------
73rd commit message:
fix more bugs in the element centered FV discretization
now eWoms should match Dumux pretty closely...
--------
74th commit message:
coalesce the common code of the multi phase porous medium models into "MultiPhaseBaseModel"
--------
75th commit message:
update reference solutions
these were changed because of the screw-up with the area of boundary
segments...
--------
76th commit message:
rename "ImplicitBase" to "FvBase"
because in eWoms, everything is implicit and these are currently the
base classes for all finite volume discretizations.
--------
77th commit message:
make the spatial discretization selectable using a splice
This requires an opm-core with a the patches from
https://github.com/OPM/opm-core/pull/446 merged...
--------
78th commit message:
rename the properties used for splices to *Splice
--------
79th commit message:
move the files in 'tests/models' to 'tests'
since 'tests' was empty except for the 'models' subdirectory...
--------
80th commit message:
improve and fix the tutorial
--------
81st commit message:
remove the -fno-strict-aliasing flag from the provided option files
seems like recent versions of Dune have been adapted...
--------
82nd commit message:
also compile all CO2 injection simulations using the element centered finite volume discretization
--------
83rd commit message:
PVS model: make it work properly with the element-centered finite volume discretiation
because DOF != number of vertices
2013-12-12 05:52:44 -06:00
|
|
|
* \copydoc FvBaseMultiPhaseProblem::materialLawParams
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
const MaterialLawParams& materialLawParams(const Context& context,
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned spaceIdx, unsigned timeIdx) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
2016-07-12 12:12:22 -05:00
|
|
|
unsigned globalIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
|
|
|
return *materialParams_[globalIdx];
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
const MaterialLawParams& materialLawParams(unsigned globalIdx) const
|
2016-07-12 12:12:22 -05:00
|
|
|
{ return *materialParams_[globalIdx]; }
|
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
/*!
|
|
|
|
* \name Problem parameters
|
|
|
|
*/
|
|
|
|
//! \{
|
|
|
|
|
|
|
|
|
|
|
|
/*!
|
Implement the element centered finite volume spatial discretization
This makes eWoms multi-discretization capable. Along the way, this
fixes some bugs and does a medium sized reorganization of the source tree.
This is a squashed patch of the following commits:
--------
1st commit message:
add initial version of the element centered finite volume discretization
currently, it is a misnomer as it is just a copy of the vertex
centered discretization plus some renames...
--------
2nd commit message:
rename [VE]cfvModel -> [VE]cfvDiscretization
--------
3rd commit message:
ecfv: prelimary changes required to make it compile
but not work yet...
--------
4th commit message:
Rename *FvElementGeometry to *Stencil
"Stencil" seems to be the standard expression for this concept...
(also, it is not specific to finite volume methods and is shorter.)
--------
5th commit message:
refactor the stencil class for the element centered finite volume discretization
--------
6th commit message:
ECFV: some work on the stencil class
--------
7th commit message:
ECFV: make the boundary handling code compile
--------
8th commit message:
rename elemContext() to elementContext()
--------
9th commit message:
ECFV: make the VTK output modules compile
--------
10th commit message:
stencil: introduce the concept of primary DOFs
also save an vector of all element pointers in the stencil.
--------
11th commit message:
ECFV: try to fix assembly; add missing timeIdx arguments to the num*() methods
--------
12th commit message:
ECFV: fix stupid mistake in the assembler
--------
13th commit message:
ECFV: remove a few implicit DOF == vertex assumptions
the black-oil example now runs without valgrind complaints until it encounters
a negative oil mole fraction.
--------
14th commit message:
VCFV: make everything compile again
all vertex centered FV examples should now work again...
--------
15th commit message:
rename [ev]cfvmodel.hh to [ev]cfvdiscretization.hh
the classes have already been renamed.
--------
16th commit message:
ECFV: make it work to the point where it can write out the initial solution.
--------
17th commit message:
ECFV: make it work
the local residual/jacobian needed some work in distinguishing primary
and secondary DOFs and there was an minor issue with the serialization
code.
for some reason, it seems still not correct. (-> convergence is too slow.)
--------
18th commit message:
VCFV: make it compile for the black oil model again
--------
19th commit message:
VCFV: make it compile with the remaining models again
--------
20th commit message:
flash model: make it work with ECFV
although this breaks its compatibility with VCFV. (-> next commit)
--------
21st commit message:
adapt the VCFV to make it compatible with the flash model again
--------
22nd commit message:
make all models compile with VCFV again
--------
23rd commit message:
VCFV: more cleanups of the stencil
VcfvStencil now does not have any public attributes anymore. TODO: do
not export attributes in the SubControlVolume and SubControlVolumeFace
classes.
--------
24th commit message:
VCFV: actually update the element pointer
--------
25th commit message:
change the blackoil model back to ECFV
--------
26th commit message:
immiscible model: make it compatible with the ECFV discretization
--------
27th commit message:
PVS model: make it work with ECFV
--------
28th commit message:
NCP model: make it work with ECFV
--------
29th commit message:
rename Vcfv*VelocityModule to *VelocityModule
--------
30th commit message:
richards model: make it work with ECFV
--------
31st commit message:
unify the ECFV and the VCFV VTK output modules
and other cleanups
--------
32nd commit message:
unify the common code of the VCFV and the ECFV disctretizations
--------
33rd commit message:
unify the element contexts between element and vertex centered finite volumes
--------
34th commit message:
unify the local jacobian class of the finite volume discretizations
--------
35th commit message:
replace [VE]vcf(LocalResidual|ElementContext|BoundaryContext|ConstraintsContext) by generic code
--------
36th commit message:
replace the [EV]cfvLocalResidual by generic code
--------
37th commit message:
unify the MultiPhaseProblem and Problem classes, introduce NullBorderListCreator
--------
38th commit message:
remove the discretization specific boundary context
--------
39th commit message:
unify the [EV]cfvDiscretization classes
--------
40th commit message:
Unify [EV]cfvMultiPhaseFluxVariables
--------
41st commit message:
Unify the [EC]cfvNewton* classes
--------
42nd commit message:
Unify [EV]cfvVolumeVariables
--------
43rd commit message:
unify [EV]cfvAssembler
--------
44th commit message:
unified flux variables: fix stupid mistake when calculating pressure gradients
--------
45th commit message:
unify what's to unify for the [EV]CFV properties
--------
46th commit message:
make the method to calculate gradients and values at flux approximation points changeable
Currently, this is used by the vertex centered finite volume method to
be able to use P1-finite element gradients instead of two-point
ones...
--------
47th commit message:
make the restart code work correctly, use the correct DofMapper for VCFV
--------
48th commit message:
actually use the gradient calculator in a model
the immiscible model in this case
--------
49th commit message:
move some files around to where they belong, use the new gradient calculation code in all models
TODO: proper handling of boundary gradients
--------
50th commit message:
fix the stokes model
currently it only works with the vertex centered finite volume
discretization, but the plan is to soon move it to a staggered grid
scheme anyway...
--------
51st commit message:
move all models back to using the vertex centered finite volume discretization by default
--------
52nd commit message:
models: some variable renames and documentation fixes
- scv -> dof
- vert -> dof
- vertex -> dof
- replace 'VCFV'
- fix some typos
--------
53rd commit message:
don't expect UG anymore
since it is quite non-free and hard to get. we now use ALUGrid instead!
--------
54th commit message:
temporarily disable jacobian recycling
--------
55th commit message:
fix writing/reading restart files using the generic code
--------
56th commit message:
fix bug where fluxes were only counted once in the stencil
this only affected the vertex centered finite volumes discretization...
--------
57th commit message:
boundary gradients: use the center of the sub-control volume adjacent to a boundary segment
--------
58th commit message:
make it compile on GCC
--------
59th commit message:
get rid of most hacks
for this, partial reassemble and jacobian recycling was brought
back. For the this and the remaining stuff the main trick is the
introduction of the GridCommHandleFactory concept which constructs
communication handles suited for the respective spatial
discretization...
--------
60th commit message:
fix a few annoying bugs
first, default the convergence criterion for the linear solver did not
honor the initial residual which lead to linear solver breakdowns,
then some debugging code was left in the discrete fracture model and
then there was a bug in the TP gradient approximation class...
this has the consequence that we need a new reference solution for the
discrete fracture problem...
--------
61st commit message:
iterative linear solver: remove the code for the non-default convergence criteria
--------
62nd commit message:
provide the FE cache instead of the local FE
this fixes a segfault in the stokes model caused by the fact that the
local FE was not initialized at this point.
--------
63rd commit message:
(Navier-)Stokes: fix bug due to the transition to unit normals
now, all tests pass for this branch. The only things which need to be
fixed are some annoying performance regressions compared to master and
some bug in the splices feature of the property system...
--------
64th commit message:
some fix for the local residual of the immiscible model
--------
65th commit message:
Navier-Stokes: implement SCV center gradients
There seems to be a bug in the previous implementation (the jacobian
inverse transposed is evaluated using the local, not the global
geometry), so the reference solution for the stokes2c test problem has
also been updated...
--------
66th commit message:
remove the ALUGrid specialization of the LensGridCreator and the YaspGrid one for the fingerproblem
using different grid seems to sometimes cause a different vertex
order, which in turn causes the respective test to fail if the
reference solution was computed using the other grid...
--------
67th commit message:
VCFV: use the correct BorderListCreator
this makes MPI parallel computations work again. apart from
performance regressions, this branch does not exhibit any known
regressions compared to master anymore...
--------
68th commit message:
make verything compile with the element centered finite volume discretization
except the Navier-Stokes and the two-phase DFM models, of course...
--------
69th commit message:
minor fixes
- make the navier-stokes model slighly more generic by using the
proper (in,ex)teriorIndex() methods on sub-control volumes
- make the signature of the calculateValue() template method of the
common two-point gradient approximator match the one of the vertex
centered finite volume one
--------
70th commit message:
fix fallout from the Big Rebase
--------
71st commit message:
ECFV: some bugs in the boundary
--------
72nd commit message:
make computeFlux() compute area-specific quantities
--------
73rd commit message:
fix more bugs in the element centered FV discretization
now eWoms should match Dumux pretty closely...
--------
74th commit message:
coalesce the common code of the multi phase porous medium models into "MultiPhaseBaseModel"
--------
75th commit message:
update reference solutions
these were changed because of the screw-up with the area of boundary
segments...
--------
76th commit message:
rename "ImplicitBase" to "FvBase"
because in eWoms, everything is implicit and these are currently the
base classes for all finite volume discretizations.
--------
77th commit message:
make the spatial discretization selectable using a splice
This requires an opm-core with a the patches from
https://github.com/OPM/opm-core/pull/446 merged...
--------
78th commit message:
rename the properties used for splices to *Splice
--------
79th commit message:
move the files in 'tests/models' to 'tests'
since 'tests' was empty except for the 'models' subdirectory...
--------
80th commit message:
improve and fix the tutorial
--------
81st commit message:
remove the -fno-strict-aliasing flag from the provided option files
seems like recent versions of Dune have been adapted...
--------
82nd commit message:
also compile all CO2 injection simulations using the element centered finite volume discretization
--------
83rd commit message:
PVS model: make it work properly with the element-centered finite volume discretiation
because DOF != number of vertices
2013-12-12 05:52:44 -06:00
|
|
|
* \copydoc FvBaseMultiPhaseProblem::temperature
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
|
|
|
* The black-oil model assumes constant temperature to define its
|
|
|
|
* parameters. Although temperature is thus not really used by the
|
|
|
|
* model, it gets written to the VTK output. Who nows, maybe we
|
|
|
|
* will need it one day?
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2022-08-02 03:40:48 -05:00
|
|
|
Scalar temperature(const Context& /*context*/,
|
|
|
|
unsigned /*spaceIdx*/,
|
|
|
|
unsigned /*timeIdx*/) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{ return temperature_; }
|
|
|
|
|
|
|
|
// \}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \name Boundary conditions
|
|
|
|
*/
|
|
|
|
//! \{
|
|
|
|
|
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \copydoc FvBaseProblem::boundary
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
|
|
|
* The reservoir problem uses constraints to approximate
|
|
|
|
* extraction and production wells, so all boundaries are no-flow.
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
void boundary(BoundaryRateVector& values,
|
2022-08-02 03:40:48 -05:00
|
|
|
const Context& /*context*/,
|
|
|
|
unsigned /*spaceIdx*/,
|
|
|
|
unsigned /*timeIdx*/) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
|
|
|
// no flow on top and bottom
|
|
|
|
values.setNoFlow();
|
|
|
|
}
|
|
|
|
|
|
|
|
//! \}
|
|
|
|
|
|
|
|
/*!
|
2014-06-24 07:54:32 -05:00
|
|
|
* \name Volumetric terms
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
//! \{
|
|
|
|
|
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \copydoc FvBaseProblem::initial
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
|
|
|
* The reservoir problem uses a constant boundary condition for
|
|
|
|
* the whole domain.
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
void initial(PrimaryVariables& values,
|
2022-08-02 03:40:48 -05:00
|
|
|
const Context& /*context*/,
|
|
|
|
unsigned /*spaceIdx*/,
|
|
|
|
unsigned /*timeIdx*/) const
|
2015-09-29 04:48:30 -05:00
|
|
|
{
|
|
|
|
values.assignNaive(initialFluidState_);
|
|
|
|
|
|
|
|
#ifndef NDEBUG
|
|
|
|
for (unsigned pvIdx = 0; pvIdx < values.size(); ++ pvIdx)
|
|
|
|
assert(std::isfinite(values[pvIdx]));
|
|
|
|
#endif
|
|
|
|
}
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \copydoc FvBaseProblem::constraints
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
2015-12-31 06:20:56 -06:00
|
|
|
* The reservoir problem places two water-injection wells on the lower-left and
|
|
|
|
* lower-right of the domain and a production well in the middle. The injection wells
|
|
|
|
* are fully water saturated with a higher pressure, the producer is fully oil
|
|
|
|
* saturated with a lower pressure than the remaining reservoir.
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
template <class Context>
|
2018-07-02 02:16:26 -05:00
|
|
|
void constraints(Constraints& constraintValues,
|
|
|
|
const Context& context,
|
|
|
|
unsigned spaceIdx,
|
|
|
|
unsigned timeIdx) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
2015-12-31 06:20:56 -06:00
|
|
|
if (this->simulator().episodeIndex() == 1)
|
|
|
|
return; // no constraints during the "settle down" episode
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
const auto& pos = context.pos(spaceIdx, timeIdx);
|
2015-12-31 06:20:56 -06:00
|
|
|
if (isInjector_(pos)) {
|
2018-07-02 02:16:26 -05:00
|
|
|
constraintValues.setActive(true);
|
|
|
|
constraintValues.assignNaive(injectorFluidState_);
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
2015-12-31 06:20:56 -06:00
|
|
|
else if (isProducer_(pos)) {
|
2018-07-02 02:16:26 -05:00
|
|
|
constraintValues.setActive(true);
|
|
|
|
constraintValues.assignNaive(producerFluidState_);
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \copydoc FvBaseProblem::source
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
|
|
|
* For this problem, the source term of all components is 0 everywhere.
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
void source(RateVector& rate,
|
2022-08-02 03:40:48 -05:00
|
|
|
const Context& /*context*/,
|
|
|
|
unsigned /*spaceIdx*/,
|
|
|
|
unsigned /*timeIdx*/) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{ rate = Scalar(0.0); }
|
|
|
|
|
|
|
|
//! \}
|
|
|
|
|
|
|
|
private:
|
|
|
|
void initFluidState_()
|
|
|
|
{
|
2016-11-07 08:14:07 -06:00
|
|
|
auto& fs = initialFluidState_;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
//////
|
|
|
|
// set temperatures
|
|
|
|
//////
|
|
|
|
fs.setTemperature(temperature_);
|
|
|
|
|
|
|
|
//////
|
|
|
|
// set saturations
|
|
|
|
//////
|
2014-04-03 09:59:48 -05:00
|
|
|
fs.setSaturation(FluidSystem::oilPhaseIdx, 1.0);
|
|
|
|
fs.setSaturation(FluidSystem::waterPhaseIdx, 0.0);
|
|
|
|
fs.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
//////
|
|
|
|
// set pressures
|
|
|
|
//////
|
|
|
|
Scalar pw = pReservoir_;
|
|
|
|
|
|
|
|
PhaseVector pC;
|
2016-11-07 08:14:07 -06:00
|
|
|
const auto& matParams = fineMaterialParams_;
|
2013-09-23 11:56:30 -05:00
|
|
|
MaterialLaw::capillaryPressures(pC, matParams, fs);
|
|
|
|
|
2014-04-03 09:59:48 -05:00
|
|
|
fs.setPressure(oilPhaseIdx, pw + (pC[oilPhaseIdx] - pC[waterPhaseIdx]));
|
|
|
|
fs.setPressure(waterPhaseIdx, pw + (pC[waterPhaseIdx] - pC[waterPhaseIdx]));
|
|
|
|
fs.setPressure(gasPhaseIdx, pw + (pC[gasPhaseIdx] - pC[waterPhaseIdx]));
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// reset all mole fractions to 0
|
2015-11-18 04:54:35 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
|
|
|
|
for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
|
2013-09-23 11:56:30 -05:00
|
|
|
fs.setMoleFraction(phaseIdx, compIdx, 0.0);
|
|
|
|
|
|
|
|
//////
|
|
|
|
// set composition of the gas and water phases
|
|
|
|
//////
|
2014-04-03 09:59:48 -05:00
|
|
|
fs.setMoleFraction(waterPhaseIdx, waterCompIdx, 1.0);
|
|
|
|
fs.setMoleFraction(gasPhaseIdx, gasCompIdx, 1.0);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
//////
|
|
|
|
// set composition of the oil phase
|
|
|
|
//////
|
2016-01-04 08:32:55 -06:00
|
|
|
Scalar RsSat =
|
|
|
|
FluidSystem::saturatedDissolutionFactor(fs, oilPhaseIdx, /*pvtRegionIdx=*/0);
|
|
|
|
Scalar XoGSat = FluidSystem::convertRsToXoG(RsSat, /*pvtRegionIdx=*/0);
|
|
|
|
Scalar xoGSat = FluidSystem::convertXoGToxoG(XoGSat, /*pvtRegionIdx=*/0);
|
|
|
|
Scalar xoG = 0.95*xoGSat;
|
|
|
|
Scalar xoO = 1.0 - xoG;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// finally set the oil-phase composition
|
2014-04-03 09:59:48 -05:00
|
|
|
fs.setMoleFraction(oilPhaseIdx, gasCompIdx, xoG);
|
|
|
|
fs.setMoleFraction(oilPhaseIdx, oilCompIdx, xoO);
|
2015-12-31 06:20:56 -06:00
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using CFRP = Opm::ComputeFromReferencePhase<Scalar, FluidSystem>;
|
2016-04-13 06:32:15 -05:00
|
|
|
typename FluidSystem::template ParameterCache<Scalar> paramCache;
|
2015-12-31 06:20:56 -06:00
|
|
|
CFRP::solve(fs,
|
|
|
|
paramCache,
|
|
|
|
/*refPhaseIdx=*/oilPhaseIdx,
|
|
|
|
/*setViscosities=*/false,
|
|
|
|
/*setEnthalpies=*/false);
|
|
|
|
|
|
|
|
// set up the fluid state used for the injectors
|
|
|
|
auto& injFs = injectorFluidState_;
|
|
|
|
injFs = initialFluidState_;
|
|
|
|
|
|
|
|
Scalar pInj = pReservoir_ * 1.5;
|
|
|
|
injFs.setPressure(waterPhaseIdx, pInj);
|
|
|
|
injFs.setPressure(oilPhaseIdx, pInj);
|
|
|
|
injFs.setPressure(gasPhaseIdx, pInj);
|
|
|
|
injFs.setSaturation(waterPhaseIdx, 1.0);
|
|
|
|
injFs.setSaturation(oilPhaseIdx, 0.0);
|
|
|
|
injFs.setSaturation(gasPhaseIdx, 0.0);
|
|
|
|
|
|
|
|
// set the composition of the phases to immiscible
|
2016-11-07 08:14:07 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
|
|
|
|
for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
|
2015-12-31 06:20:56 -06:00
|
|
|
injFs.setMoleFraction(phaseIdx, compIdx, 0.0);
|
|
|
|
|
|
|
|
injFs.setMoleFraction(gasPhaseIdx, gasCompIdx, 1.0);
|
|
|
|
injFs.setMoleFraction(oilPhaseIdx, oilCompIdx, 1.0);
|
|
|
|
injFs.setMoleFraction(waterPhaseIdx, waterCompIdx, 1.0);
|
|
|
|
|
|
|
|
CFRP::solve(injFs,
|
|
|
|
paramCache,
|
|
|
|
/*refPhaseIdx=*/waterPhaseIdx,
|
2018-01-22 03:33:55 -06:00
|
|
|
/*setViscosities=*/true,
|
2015-12-31 06:20:56 -06:00
|
|
|
/*setEnthalpies=*/false);
|
|
|
|
|
|
|
|
// set up the fluid state used for the producer
|
|
|
|
auto& prodFs = producerFluidState_;
|
|
|
|
prodFs = initialFluidState_;
|
|
|
|
|
|
|
|
Scalar pProd = pReservoir_ / 1.5;
|
|
|
|
prodFs.setPressure(waterPhaseIdx, pProd);
|
|
|
|
prodFs.setPressure(oilPhaseIdx, pProd);
|
|
|
|
prodFs.setPressure(gasPhaseIdx, pProd);
|
|
|
|
prodFs.setSaturation(waterPhaseIdx, 0.0);
|
|
|
|
prodFs.setSaturation(oilPhaseIdx, 1.0);
|
|
|
|
prodFs.setSaturation(gasPhaseIdx, 0.0);
|
|
|
|
|
|
|
|
CFRP::solve(prodFs,
|
|
|
|
paramCache,
|
|
|
|
/*refPhaseIdx=*/oilPhaseIdx,
|
2018-01-22 03:33:55 -06:00
|
|
|
/*setViscosities=*/true,
|
2015-12-31 06:20:56 -06:00
|
|
|
/*setEnthalpies=*/false);
|
|
|
|
}
|
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
bool isProducer_(const GlobalPosition& pos) const
|
2015-12-31 06:20:56 -06:00
|
|
|
{
|
|
|
|
Scalar x = pos[0] - this->boundingBoxMin()[0];
|
|
|
|
Scalar y = pos[dim - 1] - this->boundingBoxMin()[dim - 1];
|
|
|
|
Scalar width = this->boundingBoxMax()[0] - this->boundingBoxMin()[0];
|
|
|
|
Scalar height = this->boundingBoxMax()[dim - 1] - this->boundingBoxMin()[dim - 1];
|
|
|
|
|
|
|
|
// only the upper half of the center section of the spatial domain is assumed to
|
|
|
|
// be the producer
|
|
|
|
if (y <= height/2.0)
|
|
|
|
return false;
|
|
|
|
|
|
|
|
return width/2.0 - width*1e-5 < x && x < width/2.0 + width*(wellWidth_ + 1e-5);
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
bool isInjector_(const GlobalPosition& pos) const
|
2015-12-31 06:20:56 -06:00
|
|
|
{
|
|
|
|
Scalar x = pos[0] - this->boundingBoxMin()[0];
|
|
|
|
Scalar y = pos[dim - 1] - this->boundingBoxMin()[dim - 1];
|
|
|
|
Scalar width = this->boundingBoxMax()[0] - this->boundingBoxMin()[0];
|
|
|
|
Scalar height = this->boundingBoxMax()[dim - 1] - this->boundingBoxMin()[dim - 1];
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2015-12-31 06:20:56 -06:00
|
|
|
// only the lower half of the leftmost and rightmost part of the spatial domain
|
|
|
|
// are assumed to be the water injectors
|
|
|
|
if (y > height/2.0)
|
|
|
|
return false;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2015-12-31 06:20:56 -06:00
|
|
|
return x < width*wellWidth_ - width*1e-5 || x > width*(1.0 - wellWidth_) + width*1e-5;
|
|
|
|
}
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
bool isFineMaterial_(const GlobalPosition& pos) const
|
2013-11-29 09:33:46 -06:00
|
|
|
{ return pos[dim - 1] > layerBottom_; }
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
DimMatrix fineK_;
|
|
|
|
DimMatrix coarseK_;
|
|
|
|
Scalar layerBottom_;
|
|
|
|
Scalar pReservoir_;
|
|
|
|
|
|
|
|
Scalar finePorosity_;
|
|
|
|
Scalar coarsePorosity_;
|
|
|
|
|
|
|
|
MaterialLawParams fineMaterialParams_;
|
|
|
|
MaterialLawParams coarseMaterialParams_;
|
2016-07-12 12:12:22 -05:00
|
|
|
std::vector<const MaterialLawParams*> materialParams_;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2015-12-31 06:20:56 -06:00
|
|
|
InitialFluidState initialFluidState_;
|
|
|
|
InitialFluidState injectorFluidState_;
|
|
|
|
InitialFluidState producerFluidState_;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
Scalar temperature_;
|
|
|
|
Scalar maxDepth_;
|
2015-12-31 06:20:56 -06:00
|
|
|
Scalar wellWidth_;
|
2013-09-23 11:56:30 -05:00
|
|
|
};
|
2019-09-05 09:21:10 -05:00
|
|
|
} // namespace Opm
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
#endif
|