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::Co2InjectionProblem
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
#ifndef EWOMS_CO2_INJECTION_PROBLEM_HH
|
|
|
|
#define EWOMS_CO2_INJECTION_PROBLEM_HH
|
|
|
|
|
2019-09-16 05:07:20 -05:00
|
|
|
#include <opm/models/immiscible/immisciblemodel.hh>
|
2019-09-09 02:20:08 -05:00
|
|
|
#include <opm/simulators/linalg/parallelamgbackend.hh>
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
#include <opm/material/fluidsystems/H2ON2FluidSystem.hpp>
|
|
|
|
#include <opm/material/fluidsystems/BrineCO2FluidSystem.hpp>
|
|
|
|
#include <opm/material/fluidstates/CompositionalFluidState.hpp>
|
|
|
|
#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
|
|
|
|
#include <opm/material/constraintsolvers/ComputeFromReferencePhase.hpp>
|
2013-11-06 07:50:01 -06:00
|
|
|
#include <opm/material/fluidmatrixinteractions/LinearMaterial.hpp>
|
|
|
|
#include <opm/material/fluidmatrixinteractions/RegularizedBrooksCorey.hpp>
|
|
|
|
#include <opm/material/fluidmatrixinteractions/EffToAbsLaw.hpp>
|
|
|
|
#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
|
2018-01-04 08:26:07 -06:00
|
|
|
#include <opm/material/thermal/SomertonThermalConductionLaw.hpp>
|
2017-12-11 05:58:30 -06:00
|
|
|
#include <opm/material/thermal/ConstantSolidHeatCapLaw.hpp>
|
2013-09-23 11:56:30 -05:00
|
|
|
#include <opm/material/binarycoefficients/Brine_CO2.hpp>
|
2015-04-28 04:58:21 -05:00
|
|
|
#include <opm/material/common/UniformTabulated2DFunction.hpp>
|
2018-02-01 07:40:01 -06:00
|
|
|
#include <opm/material/common/Unused.hpp>
|
2013-09-23 11:56:30 -05:00
|
|
|
|
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 <sstream>
|
|
|
|
#include <iostream>
|
|
|
|
#include <string>
|
|
|
|
|
2019-09-05 09:21:10 -05:00
|
|
|
namespace Opm {
|
2015-06-18 07:18:42 -05:00
|
|
|
//! \cond SKIP_THIS
|
2013-09-23 11:56:30 -05:00
|
|
|
template <class TypeTag>
|
|
|
|
class Co2InjectionProblem;
|
|
|
|
|
|
|
|
namespace Co2Injection {
|
|
|
|
#include <opm/material/components/co2tables.inc>
|
2013-11-29 09:33:46 -06:00
|
|
|
}
|
2015-06-18 07:18:42 -05:00
|
|
|
//! \endcond
|
2018-06-14 09:06:05 -05:00
|
|
|
}
|
|
|
|
|
2020-06-08 10:11:48 -05:00
|
|
|
namespace Opm::Properties {
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2020-06-10 06:07:19 -05:00
|
|
|
namespace TTag {
|
|
|
|
struct Co2InjectionBaseProblem {};
|
|
|
|
}
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// declare the CO2 injection problem specific property tags
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct FluidSystemPressureLow { using type = UndefinedProperty; };
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct FluidSystemPressureHigh { using type = UndefinedProperty; };
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct FluidSystemNumPressure { using type = UndefinedProperty; };
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct FluidSystemTemperatureLow { using type = UndefinedProperty; };
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct FluidSystemTemperatureHigh { using type = UndefinedProperty; };
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct FluidSystemNumTemperature { using type = UndefinedProperty; };
|
|
|
|
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct MaxDepth { using type = UndefinedProperty; };
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct Temperature { using type = UndefinedProperty; };
|
|
|
|
template<class TypeTag, class MyTypeTag>
|
|
|
|
struct SimulationName { 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::Co2InjectionBaseProblem> { using type = Dune::YaspGrid<2>; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// Set the problem property
|
2020-06-10 06:07:19 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct Problem<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{ using type = Opm::Co2InjectionProblem<TypeTag>; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// Set fluid configuration
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct FluidSystem<TypeTag, TTag::Co2InjectionBaseProblem>
|
2013-11-29 09:33:46 -06:00
|
|
|
{
|
|
|
|
private:
|
2020-06-10 06:49:42 -05:00
|
|
|
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
|
|
|
|
using CO2Tables = Opm::Co2Injection::CO2Tables;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
public:
|
2020-06-10 06:49:42 -05:00
|
|
|
using type = Opm::BrineCO2FluidSystem<Scalar, CO2Tables>;
|
|
|
|
//using type = Opm::H2ON2FluidSystem<Scalar, /*useComplexRelations=*/false>;
|
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::Co2InjectionBaseProblem>
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
|
|
|
private:
|
2020-06-10 06:49:42 -05:00
|
|
|
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
|
2014-04-03 09:59:48 -05:00
|
|
|
enum { liquidPhaseIdx = FluidSystem::liquidPhaseIdx };
|
|
|
|
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
|
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
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
|
|
|
|
using Traits = Opm::TwoPhaseMaterialTraits<Scalar,
|
|
|
|
/*wettingPhaseIdx=*/FluidSystem::liquidPhaseIdx,
|
|
|
|
/*nonWettingPhaseIdx=*/FluidSystem::gasPhaseIdx>;
|
2013-11-06 07:50:01 -06:00
|
|
|
|
|
|
|
// define the material law which is parameterized by effective
|
|
|
|
// saturations
|
2020-06-10 06:49:42 -05:00
|
|
|
using EffMaterialLaw = Opm::RegularizedBrooksCorey<Traits>;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
public:
|
2013-11-06 07:50:01 -06:00
|
|
|
// define the material law parameterized by absolute saturations
|
2020-06-10 06:49:42 -05:00
|
|
|
using type = Opm::EffToAbsLaw<EffMaterialLaw>;
|
2013-09-23 11:56:30 -05:00
|
|
|
};
|
|
|
|
|
2018-01-04 08:26:07 -06:00
|
|
|
// Set the thermal conduction law
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct ThermalConductionLaw<TypeTag, TTag::Co2InjectionBaseProblem>
|
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
|
|
|
|
|
|
|
public:
|
|
|
|
// define the material law parameterized by absolute saturations
|
2020-06-10 06:49:42 -05:00
|
|
|
using type = Opm::SomertonThermalConductionLaw<FluidSystem, Scalar>;
|
2013-09-23 11:56:30 -05:00
|
|
|
};
|
|
|
|
|
2018-01-04 08:26:07 -06:00
|
|
|
// set the energy storage law for the solid phase
|
2020-06-10 06:07:19 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct SolidEnergyLaw<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{ using type = Opm::ConstantSolidHeatCapLaw<GetPropType<TypeTag, Properties::Scalar>>; };
|
2017-12-11 05:58:30 -06:00
|
|
|
|
2014-12-21 11:36:59 -06:00
|
|
|
// Use the algebraic multi-grid linear solver for this problem
|
2020-06-09 03:55:25 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct LinearSolverSplice<TypeTag, TTag::Co2InjectionBaseProblem> { using type = TTag::ParallelAmgLinearSolver; };
|
2014-12-21 11:36:59 -06:00
|
|
|
|
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::Co2InjectionBaseProblem> { 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::Co2InjectionBaseProblem> { static constexpr bool value = true; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// set the defaults for the problem specific properties
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct FluidSystemPressureLow<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 3e7;
|
|
|
|
};
|
|
|
|
template<class TypeTag>
|
|
|
|
struct FluidSystemPressureHigh<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 4e7;
|
|
|
|
};
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct FluidSystemNumPressure<TypeTag, TTag::Co2InjectionBaseProblem> { static constexpr int value = 100; };
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct FluidSystemTemperatureLow<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 290;
|
|
|
|
};
|
|
|
|
template<class TypeTag>
|
|
|
|
struct FluidSystemTemperatureHigh<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 500;
|
|
|
|
};
|
2020-06-08 09:41:02 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct FluidSystemNumTemperature<TypeTag, TTag::Co2InjectionBaseProblem> { static constexpr int value = 100; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct MaxDepth<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 2500;
|
|
|
|
};
|
|
|
|
template<class TypeTag>
|
|
|
|
struct Temperature<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 293.15;
|
|
|
|
};
|
2020-06-09 04:15:16 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct SimulationName<TypeTag, TTag::Co2InjectionBaseProblem> { static constexpr auto value = "co2injection"; };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// The default for the end time of the simulation
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct EndTime<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 1e4;
|
|
|
|
};
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// The default for the initial time step size of the simulation
|
2020-06-09 03:43:28 -05:00
|
|
|
template<class TypeTag>
|
|
|
|
struct InitialTimeStepSize<TypeTag, TTag::Co2InjectionBaseProblem>
|
|
|
|
{
|
|
|
|
using type = GetPropType<TypeTag, Scalar>;
|
|
|
|
static constexpr type value = 250;
|
|
|
|
};
|
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::Co2InjectionBaseProblem> { static constexpr auto value = "data/co2injection.dgf"; };
|
2018-06-14 09:06:05 -05:00
|
|
|
|
2020-06-08 10:11:48 -05:00
|
|
|
} // namespace Opm::Properties
|
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
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \ingroup TestProblems
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
|
|
|
* \brief Problem where \f$CO_2\f$ is injected under a low permeable
|
|
|
|
* layer at a depth of 2700m.
|
|
|
|
*
|
|
|
|
* The domain is sized 60m times 40m and consists of two layers, one
|
|
|
|
* which is moderately permeable (\f$K = 10^{-12}\;m^2\f$) for \f$ y >
|
|
|
|
* 22\; m\f$ and one with a lower intrinsic permeablility (\f$
|
|
|
|
* K=10^{-13}\;m^2\f$) in the rest of the domain.
|
|
|
|
*
|
|
|
|
* \f$CO_2\f$ gets injected by means of a forced-flow boundary
|
|
|
|
* condition into water-filled aquifer, which is situated 2700m below
|
|
|
|
* sea level, at the lower-right boundary (\f$5m<y<15m\f$) and
|
|
|
|
* migrates upwards due to buoyancy. It accumulates and eventually
|
|
|
|
* enters the lower permeable aquitard.
|
|
|
|
*
|
|
|
|
* The boundary conditions applied by this problem are no-flow
|
|
|
|
* conditions on the top bottom and right boundaries and a free-flow
|
|
|
|
* boundary condition on the left. For the free-flow condition,
|
|
|
|
* hydrostatic pressure is assumed.
|
|
|
|
*/
|
|
|
|
template <class TypeTag>
|
2020-06-08 09:41:02 -05:00
|
|
|
class Co2InjectionProblem : 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 Scalar = GetPropType<TypeTag, Properties::Scalar>;
|
|
|
|
using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
|
|
|
|
using GridView = GetPropType<TypeTag, Properties::GridView>;
|
|
|
|
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2013-11-29 09:33:46 -06:00
|
|
|
enum { dim = GridView::dimension };
|
|
|
|
enum { dimWorld = GridView::dimensionworld };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// copy some indices for convenience
|
2020-06-10 06:49:42 -05:00
|
|
|
using Indices = GetPropType<TypeTag, Properties::Indices>;
|
2013-11-29 09:33:46 -06:00
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
2014-04-03 09:59:48 -05:00
|
|
|
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
|
|
|
|
enum { liquidPhaseIdx = FluidSystem::liquidPhaseIdx };
|
2013-11-29 09:33:46 -06:00
|
|
|
enum { CO2Idx = FluidSystem::CO2Idx };
|
|
|
|
enum { BrineIdx = FluidSystem::BrineIdx };
|
|
|
|
enum { conti0EqIdx = Indices::conti0EqIdx };
|
|
|
|
enum { contiCO2EqIdx = conti0EqIdx + CO2Idx };
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
|
|
|
|
using RateVector = GetPropType<TypeTag, Properties::RateVector>;
|
|
|
|
using BoundaryRateVector = GetPropType<TypeTag, Properties::BoundaryRateVector>;
|
|
|
|
using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
|
|
|
|
using Simulator = GetPropType<TypeTag, Properties::Simulator>;
|
|
|
|
using Model = GetPropType<TypeTag, Properties::Model>;
|
|
|
|
using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
|
|
|
|
using ThermalConductionLaw = GetPropType<TypeTag, Properties::ThermalConductionLaw>;
|
|
|
|
using SolidEnergyLawParams = GetPropType<TypeTag, Properties::SolidEnergyLawParams>;
|
|
|
|
using ThermalConductionLawParams = typename ThermalConductionLaw::Params;
|
|
|
|
|
|
|
|
using Toolbox = Opm::MathToolbox<Evaluation>;
|
|
|
|
using CoordScalar = typename GridView::ctype;
|
|
|
|
using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
|
|
|
|
using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
public:
|
|
|
|
/*!
|
|
|
|
* \copydoc Doxygen::defaultProblemConstructor
|
|
|
|
*/
|
2016-11-07 08:14:07 -06:00
|
|
|
Co2InjectionProblem(Simulator& simulator)
|
2014-04-14 13:32:30 -05:00
|
|
|
: ParentType(simulator)
|
2014-08-06 09:31:48 -05:00
|
|
|
{ }
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::finishInit
|
|
|
|
*/
|
|
|
|
void finishInit()
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
2014-08-06 09:31:48 -05:00
|
|
|
ParentType::finishInit();
|
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
eps_ = 1e-6;
|
|
|
|
|
2014-02-23 11:25:09 -06:00
|
|
|
temperatureLow_ = EWOMS_GET_PARAM(TypeTag, Scalar, FluidSystemTemperatureLow);
|
|
|
|
temperatureHigh_ = EWOMS_GET_PARAM(TypeTag, Scalar, FluidSystemTemperatureHigh);
|
2016-11-07 08:14:07 -06:00
|
|
|
nTemperature_ = EWOMS_GET_PARAM(TypeTag, unsigned, FluidSystemNumTemperature);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2013-09-23 13:25:58 -05:00
|
|
|
pressureLow_ = EWOMS_GET_PARAM(TypeTag, Scalar, FluidSystemPressureLow);
|
2014-02-23 11:25:09 -06:00
|
|
|
pressureHigh_ = EWOMS_GET_PARAM(TypeTag, Scalar, FluidSystemPressureHigh);
|
2016-11-07 08:14:07 -06:00
|
|
|
nPressure_ = EWOMS_GET_PARAM(TypeTag, unsigned, FluidSystemNumPressure);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2013-09-23 13:25:58 -05:00
|
|
|
maxDepth_ = EWOMS_GET_PARAM(TypeTag, Scalar, MaxDepth);
|
|
|
|
temperature_ = EWOMS_GET_PARAM(TypeTag, Scalar, Temperature);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// initialize the tables of the fluid system
|
2013-11-29 09:33:46 -06:00
|
|
|
// FluidSystem::init();
|
2013-09-23 11:56:30 -05:00
|
|
|
FluidSystem::init(/*Tmin=*/temperatureLow_,
|
|
|
|
/*Tmax=*/temperatureHigh_,
|
|
|
|
/*nT=*/nTemperature_,
|
|
|
|
/*pmin=*/pressureLow_,
|
|
|
|
/*pmax=*/pressureHigh_,
|
|
|
|
/*np=*/nPressure_);
|
|
|
|
|
|
|
|
fineLayerBottom_ = 22.0;
|
|
|
|
|
|
|
|
// intrinsic permeabilities
|
|
|
|
fineK_ = this->toDimMatrix_(1e-13);
|
|
|
|
coarseK_ = this->toDimMatrix_(1e-12);
|
|
|
|
|
|
|
|
// porosities
|
|
|
|
finePorosity_ = 0.3;
|
|
|
|
coarsePorosity_ = 0.3;
|
|
|
|
|
|
|
|
// residual saturations
|
2014-04-03 09:59:48 -05:00
|
|
|
fineMaterialParams_.setResidualSaturation(liquidPhaseIdx, 0.2);
|
|
|
|
fineMaterialParams_.setResidualSaturation(gasPhaseIdx, 0.0);
|
|
|
|
coarseMaterialParams_.setResidualSaturation(liquidPhaseIdx, 0.2);
|
|
|
|
coarseMaterialParams_.setResidualSaturation(gasPhaseIdx, 0.0);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
// parameters for the Brooks-Corey law
|
2013-11-06 07:50:01 -06:00
|
|
|
fineMaterialParams_.setEntryPressure(1e4);
|
|
|
|
coarseMaterialParams_.setEntryPressure(5e3);
|
2013-09-23 11:56:30 -05:00
|
|
|
fineMaterialParams_.setLambda(2.0);
|
|
|
|
coarseMaterialParams_.setLambda(2.0);
|
|
|
|
|
2013-11-06 07:50:01 -06:00
|
|
|
fineMaterialParams_.finalize();
|
|
|
|
coarseMaterialParams_.finalize();
|
|
|
|
|
2018-01-04 08:26:07 -06:00
|
|
|
// parameters for the somerton law thermal conduction
|
|
|
|
computeThermalCondParams_(fineThermalCondParams_, finePorosity_);
|
|
|
|
computeThermalCondParams_(coarseThermalCondParams_, coarsePorosity_);
|
2017-12-11 05:58:30 -06:00
|
|
|
|
2018-01-04 08:26:07 -06:00
|
|
|
// assume constant heat capacity and granite
|
|
|
|
solidEnergyLawParams_.setSolidHeatCapacity(790.0 // specific heat capacity of granite [J / (kg K)]
|
|
|
|
* 2700.0); // density of granite [kg/m^3]
|
|
|
|
solidEnergyLawParams_.finalize();
|
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, FluidSystemTemperatureLow,
|
|
|
|
"The lower temperature [K] for tabulation of the "
|
|
|
|
"fluid system");
|
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, Scalar, FluidSystemTemperatureHigh,
|
|
|
|
"The upper temperature [K] for tabulation of the "
|
|
|
|
"fluid system");
|
2016-11-10 13:17:39 -06:00
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, unsigned, FluidSystemNumTemperature,
|
2013-11-29 09:33:46 -06:00
|
|
|
"The number of intervals between the lower and "
|
|
|
|
"upper temperature");
|
|
|
|
|
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, Scalar, FluidSystemPressureLow,
|
|
|
|
"The lower pressure [Pa] for tabulation of the "
|
|
|
|
"fluid system");
|
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, Scalar, FluidSystemPressureHigh,
|
|
|
|
"The upper pressure [Pa] for tabulation of the "
|
|
|
|
"fluid system");
|
2016-11-10 13:17:39 -06:00
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, unsigned, FluidSystemNumPressure,
|
2013-11-29 09:33:46 -06:00
|
|
|
"The number of intervals between the lower and "
|
|
|
|
"upper pressure");
|
|
|
|
|
|
|
|
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");
|
|
|
|
EWOMS_REGISTER_PARAM(TypeTag, std::string, SimulationName,
|
|
|
|
"The name of the simulation used for the output "
|
|
|
|
"files");
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \name Problem parameters
|
|
|
|
*/
|
|
|
|
//! \{
|
|
|
|
|
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \copydoc FvBaseProblem::name
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
2014-04-25 10:22:28 -05:00
|
|
|
std::string name() const
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
|
|
|
std::ostringstream oss;
|
2014-03-07 04:34:55 -06:00
|
|
|
oss << EWOMS_GET_PARAM(TypeTag, std::string, SimulationName)
|
|
|
|
<< "_" << Model::name();
|
2020-06-08 09:41:02 -05:00
|
|
|
if (getPropValue<TypeTag, Properties::EnableEnergy>())
|
2013-09-23 11:56:30 -05:00
|
|
|
oss << "_ni";
|
2014-03-07 04:34:55 -06:00
|
|
|
oss << "_" << Model::discretizationName();
|
2013-09-23 11:56:30 -05:00
|
|
|
return oss.str();
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \copydoc FvBaseProblem::endTimeStep
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
2014-07-17 09:54:05 -05:00
|
|
|
void endTimeStep()
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
2014-07-22 05:41:56 -05:00
|
|
|
#ifndef NDEBUG
|
2016-11-02 09:45:19 -05:00
|
|
|
Scalar tol = this->model().newtonMethod().tolerance()*1e5;
|
2014-07-22 05:41:56 -05:00
|
|
|
this->model().checkConservativeness(tol);
|
|
|
|
|
2013-09-23 11:56:30 -05:00
|
|
|
// Calculate storage terms
|
|
|
|
PrimaryVariables storageL, storageG;
|
|
|
|
this->model().globalPhaseStorage(storageL, /*phaseIdx=*/0);
|
|
|
|
this->model().globalPhaseStorage(storageG, /*phaseIdx=*/1);
|
|
|
|
|
|
|
|
// Write mass balance information for rank 0
|
|
|
|
if (this->gridView().comm().rank() == 0) {
|
2013-11-29 09:33:46 -06:00
|
|
|
std::cout << "Storage: liquid=[" << storageL << "]"
|
2014-03-06 12:32:04 -06:00
|
|
|
<< " gas=[" << storageG << "]\n" << std::flush;
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
2014-07-22 05:41:56 -05:00
|
|
|
#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::temperature
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
Scalar temperature(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 auto& pos = context.pos(spaceIdx, timeIdx);
|
2013-09-23 11:56:30 -05:00
|
|
|
if (inHighTemperatureRegion_(pos))
|
|
|
|
return temperature_ + 100;
|
|
|
|
return temperature_;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
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
|
|
|
*/
|
|
|
|
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-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 fineMaterialParams_;
|
|
|
|
return coarseMaterialParams_;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
2017-12-11 05:58:30 -06:00
|
|
|
* \brief Return the parameters for the heat storage law of the rock
|
2013-09-23 11:56:30 -05:00
|
|
|
*
|
|
|
|
* In this case, we assume the rock-matrix to be granite.
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2017-12-11 05:58:30 -06:00
|
|
|
const SolidEnergyLawParams&
|
2018-01-04 08:26:07 -06:00
|
|
|
solidEnergyLawParams(const Context& context OPM_UNUSED,
|
|
|
|
unsigned spaceIdx OPM_UNUSED,
|
|
|
|
unsigned timeIdx OPM_UNUSED) const
|
|
|
|
{ return solidEnergyLawParams_; }
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
/*!
|
2018-01-04 08:26:07 -06:00
|
|
|
* \copydoc FvBaseMultiPhaseProblem::thermalConductionParams
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
template <class Context>
|
2018-01-04 08:26:07 -06:00
|
|
|
const ThermalConductionLawParams &
|
|
|
|
thermalConductionLawParams(const Context& context,
|
2017-12-11 05:58:30 -06:00
|
|
|
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))
|
2018-01-04 08:26:07 -06:00
|
|
|
return fineThermalCondParams_;
|
|
|
|
return coarseThermalCondParams_;
|
2013-09-23 11:56:30 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
//! \}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \name Boundary conditions
|
|
|
|
*/
|
|
|
|
//! \{
|
|
|
|
|
|
|
|
/*!
|
2014-07-17 10:24:48 -05:00
|
|
|
* \copydoc FvBaseProblem::boundary
|
2013-09-23 11:56:30 -05:00
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
void boundary(BoundaryRateVector& values, 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-11-07 08:14:07 -06:00
|
|
|
const auto& pos = context.pos(spaceIdx, timeIdx);
|
2013-09-23 11:56:30 -05:00
|
|
|
if (onLeftBoundary_(pos)) {
|
|
|
|
Opm::CompositionalFluidState<Scalar, FluidSystem> fs;
|
|
|
|
initialFluidState_(fs, context, spaceIdx, timeIdx);
|
|
|
|
fs.checkDefined();
|
|
|
|
|
|
|
|
// impose an freeflow boundary condition
|
|
|
|
values.setFreeFlow(context, spaceIdx, timeIdx, fs);
|
|
|
|
}
|
|
|
|
else if (onInlet_(pos)) {
|
|
|
|
RateVector massRate(0.0);
|
|
|
|
massRate[contiCO2EqIdx] = -1e-3; // [kg/(m^3 s)]
|
|
|
|
|
2020-06-10 06:49:42 -05:00
|
|
|
using FluidState = Opm::ImmiscibleFluidState<Scalar, FluidSystem>;
|
2016-04-13 06:32:15 -05:00
|
|
|
FluidState fs;
|
2014-04-03 09:59:48 -05:00
|
|
|
fs.setSaturation(gasPhaseIdx, 1.0);
|
2015-05-21 09:18:40 -05:00
|
|
|
const auto& pg =
|
|
|
|
context.intensiveQuantities(spaceIdx, timeIdx).fluidState().pressure(gasPhaseIdx);
|
|
|
|
fs.setPressure(gasPhaseIdx, Toolbox::value(pg));
|
2013-09-23 11:56:30 -05:00
|
|
|
fs.setTemperature(temperature(context, spaceIdx, timeIdx));
|
2016-04-13 06:32:15 -05:00
|
|
|
|
|
|
|
typename FluidSystem::template ParameterCache<Scalar> paramCache;
|
2014-04-03 09:59:48 -05:00
|
|
|
paramCache.updatePhase(fs, gasPhaseIdx);
|
2016-04-13 06:32:15 -05:00
|
|
|
Scalar h = FluidSystem::template enthalpy<FluidState, Scalar>(fs, paramCache, gasPhaseIdx);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2015-05-21 09:18:40 -05:00
|
|
|
// impose an forced inflow boundary condition for pure CO2
|
2013-09-23 11:56:30 -05:00
|
|
|
values.setMassRate(massRate);
|
|
|
|
values.setEnthalpyRate(massRate[contiCO2EqIdx] * h);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
// 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
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
void initial(PrimaryVariables& values, const Context& context, unsigned spaceIdx,
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned timeIdx) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
|
|
|
Opm::CompositionalFluidState<Scalar, FluidSystem> fs;
|
|
|
|
initialFluidState_(fs, context, spaceIdx, timeIdx);
|
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
// const auto& matParams = this->materialLawParams(context, spaceIdx,
|
2013-11-29 09:33:46 -06:00
|
|
|
// timeIdx);
|
|
|
|
// values.assignMassConservative(fs, matParams, /*inEquilibrium=*/true);
|
2013-09-23 11:56:30 -05:00
|
|
|
values.assignNaive(fs);
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
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,
|
2017-01-17 06:21:16 -06:00
|
|
|
const Context& context OPM_UNUSED,
|
|
|
|
unsigned spaceIdx OPM_UNUSED,
|
|
|
|
unsigned timeIdx OPM_UNUSED) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{ rate = Scalar(0.0); }
|
|
|
|
|
|
|
|
//! \}
|
|
|
|
|
|
|
|
private:
|
|
|
|
template <class Context, class FluidState>
|
2016-11-07 08:14:07 -06:00
|
|
|
void initialFluidState_(FluidState& fs,
|
|
|
|
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
|
|
|
|
|
|
|
//////
|
|
|
|
// set temperature
|
|
|
|
//////
|
|
|
|
fs.setTemperature(temperature(context, spaceIdx, timeIdx));
|
|
|
|
|
|
|
|
//////
|
|
|
|
// set saturations
|
|
|
|
//////
|
2014-04-03 09:59:48 -05:00
|
|
|
fs.setSaturation(FluidSystem::liquidPhaseIdx, 1.0);
|
|
|
|
fs.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
//////
|
|
|
|
// set pressures
|
|
|
|
//////
|
2015-05-21 09:18:40 -05:00
|
|
|
Scalar densityL = FluidSystem::Brine::liquidDensity(temperature_, Scalar(1e5));
|
2013-11-29 09:33:46 -06:00
|
|
|
Scalar depth = maxDepth_ - pos[dim - 1];
|
|
|
|
Scalar pl = 1e5 - densityL * this->gravity()[dim - 1] * depth;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
Scalar pC[numPhases];
|
2016-11-07 08:14:07 -06:00
|
|
|
const auto& matParams = this->materialLawParams(context, spaceIdx, timeIdx);
|
2013-09-23 11:56:30 -05:00
|
|
|
MaterialLaw::capillaryPressures(pC, matParams, fs);
|
|
|
|
|
2014-04-03 09:59:48 -05:00
|
|
|
fs.setPressure(liquidPhaseIdx, pl + (pC[liquidPhaseIdx] - pC[liquidPhaseIdx]));
|
|
|
|
fs.setPressure(gasPhaseIdx, pl + (pC[gasPhaseIdx] - pC[liquidPhaseIdx]));
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
//////
|
|
|
|
// set composition of the liquid phase
|
|
|
|
//////
|
2014-04-03 09:59:48 -05:00
|
|
|
fs.setMoleFraction(liquidPhaseIdx, CO2Idx, 0.005);
|
|
|
|
fs.setMoleFraction(liquidPhaseIdx, BrineIdx,
|
|
|
|
1.0 - fs.moleFraction(liquidPhaseIdx, CO2Idx));
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2016-04-13 06:32:15 -05:00
|
|
|
typename FluidSystem::template ParameterCache<Scalar> paramCache;
|
2020-06-10 06:49:42 -05:00
|
|
|
using CFRP = Opm::ComputeFromReferencePhase<Scalar, FluidSystem>;
|
2013-11-29 09:33:46 -06:00
|
|
|
CFRP::solve(fs, paramCache,
|
2014-04-03 09:59:48 -05:00
|
|
|
/*refPhaseIdx=*/liquidPhaseIdx,
|
2013-09-23 11:56:30 -05:00
|
|
|
/*setViscosity=*/true,
|
|
|
|
/*setEnthalpy=*/true);
|
|
|
|
}
|
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
bool onLeftBoundary_(const GlobalPosition& pos) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{ return pos[0] < eps_; }
|
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
bool onRightBoundary_(const GlobalPosition& pos) const
|
2013-12-27 11:25:54 -06:00
|
|
|
{ return pos[0] > this->boundingBoxMax()[0] - eps_; }
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
bool onInlet_(const GlobalPosition& pos) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{ return onRightBoundary_(pos) && (5 < pos[1]) && (pos[1] < 15); }
|
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
bool inHighTemperatureRegion_(const GlobalPosition& pos) const
|
2013-09-23 11:56:30 -05:00
|
|
|
{ return (pos[0] > 20) && (pos[0] < 30) && (pos[1] > 5) && (pos[1] < 35); }
|
|
|
|
|
2018-01-04 08:26:07 -06:00
|
|
|
void computeThermalCondParams_(ThermalConductionLawParams& params, Scalar poro)
|
2013-09-23 11:56:30 -05:00
|
|
|
{
|
|
|
|
Scalar lambdaWater = 0.6;
|
|
|
|
Scalar lambdaGranite = 2.8;
|
|
|
|
|
2013-11-29 09:33:46 -06:00
|
|
|
Scalar lambdaWet = std::pow(lambdaGranite, (1 - poro))
|
|
|
|
* std::pow(lambdaWater, poro);
|
|
|
|
Scalar lambdaDry = std::pow(lambdaGranite, (1 - poro));
|
2013-09-23 11:56:30 -05:00
|
|
|
|
2014-04-03 09:59:48 -05:00
|
|
|
params.setFullySaturatedLambda(gasPhaseIdx, lambdaDry);
|
|
|
|
params.setFullySaturatedLambda(liquidPhaseIdx, lambdaWet);
|
2013-09-23 11:56:30 -05:00
|
|
|
params.setVacuumLambda(lambdaDry);
|
|
|
|
}
|
|
|
|
|
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] > fineLayerBottom_; }
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
DimMatrix fineK_;
|
|
|
|
DimMatrix coarseK_;
|
|
|
|
Scalar fineLayerBottom_;
|
|
|
|
|
|
|
|
Scalar finePorosity_;
|
|
|
|
Scalar coarsePorosity_;
|
|
|
|
|
|
|
|
MaterialLawParams fineMaterialParams_;
|
|
|
|
MaterialLawParams coarseMaterialParams_;
|
|
|
|
|
2018-01-04 08:26:07 -06:00
|
|
|
ThermalConductionLawParams fineThermalCondParams_;
|
|
|
|
ThermalConductionLawParams coarseThermalCondParams_;
|
|
|
|
SolidEnergyLawParams solidEnergyLawParams_;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
Scalar temperature_;
|
|
|
|
Scalar maxDepth_;
|
|
|
|
Scalar eps_;
|
|
|
|
|
2016-11-07 08:14:07 -06:00
|
|
|
unsigned nTemperature_;
|
|
|
|
unsigned nPressure_;
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
Scalar pressureLow_, pressureHigh_;
|
|
|
|
Scalar temperatureLow_, temperatureHigh_;
|
|
|
|
};
|
2019-09-05 09:21:10 -05:00
|
|
|
} // namespace Opm
|
2013-09-23 11:56:30 -05:00
|
|
|
|
|
|
|
#endif
|