opm-simulators/examples/problems/lensproblem.hh

689 lines
23 KiB
C++
Raw Normal View History

// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
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/>.
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.
*/
/*!
* \file
*
* \copydoc Opm::LensProblem
*/
#ifndef EWOMS_LENS_PROBLEM_HH
#define EWOMS_LENS_PROBLEM_HH
2019-09-16 03:09:33 -05:00
#include <opm/models/io/structuredgridvanguard.hh>
#include <opm/models/immiscible/immiscibleproperties.hh>
#include <opm/models/discretization/common/fvbaseadlocallinearizer.hh>
#include <opm/models/discretization/ecfv/ecfvdiscretization.hh>
#include <opm/models/common/transfluxmodule.hh>
2013-11-06 07:50:01 -06:00
#include <opm/material/fluidmatrixinteractions/RegularizedVanGenuchten.hpp>
#include <opm/material/fluidmatrixinteractions/LinearMaterial.hpp>
#include <opm/material/fluidmatrixinteractions/EffToAbsLaw.hpp>
#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
#include <opm/material/fluidsystems/TwoPhaseImmiscibleFluidSystem.hpp>
#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
#include <opm/material/components/SimpleH2O.hpp>
#include <opm/material/components/Dnapl.hpp>
#include <dune/common/version.hh>
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <sstream>
#include <string>
#include <iostream>
namespace Opm {
template <class TypeTag>
class LensProblem;
}
namespace Opm::Properties {
// Create new type tags
namespace TTag {
struct LensBaseProblem { using InheritsFrom = std::tuple<StructuredGridVanguard>; };
} // end namespace TTag
// Set the problem property
template<class TypeTag>
struct Problem<TypeTag, TTag::LensBaseProblem> { using type = Opm::LensProblem<TypeTag>; };
// Use Dune-grid's YaspGrid
template<class TypeTag>
struct Grid<TypeTag, TTag::LensBaseProblem> { using type = Dune::YaspGrid<2>; };
// Set the wetting phase
template<class TypeTag>
struct WettingPhase<TypeTag, TTag::LensBaseProblem>
{
private:
2020-06-10 06:49:42 -05:00
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
public:
2020-06-10 06:49:42 -05:00
using type = Opm::LiquidPhase<Scalar, Opm::SimpleH2O<Scalar> >;
};
// Set the non-wetting phase
template<class TypeTag>
struct NonwettingPhase<TypeTag, TTag::LensBaseProblem>
{
private:
2020-06-10 06:49:42 -05:00
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
public:
2020-06-10 06:49:42 -05:00
using type = Opm::LiquidPhase<Scalar, Opm::DNAPL<Scalar> >;
};
// Set the material Law
template<class TypeTag>
struct MaterialLaw<TypeTag, TTag::LensBaseProblem>
{
private:
2020-06-10 06:49:42 -05:00
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
enum { wettingPhaseIdx = FluidSystem::wettingPhaseIdx };
enum { nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx };
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::wettingPhaseIdx,
/*nonWettingPhaseIdx=*/FluidSystem::nonWettingPhaseIdx>;
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 EffectiveLaw = Opm::RegularizedVanGenuchten<Traits>;
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<EffectiveLaw>;
};
} // namespace Opm::Properties
namespace Opm::Parameters {
// declare the properties specific for the lens problem
template<class TypeTag, class MyTypeTag>
struct LensLowerLeftX { using type = Properties::UndefinedProperty; };
template<class TypeTag, class MyTypeTag>
struct LensLowerLeftY { using type = Properties::UndefinedProperty; };
template<class TypeTag, class MyTypeTag>
struct LensLowerLeftZ { using type = Properties::UndefinedProperty; };
template<class TypeTag, class MyTypeTag>
struct LensUpperRightX { using type = Properties::UndefinedProperty; };
template<class TypeTag, class MyTypeTag>
struct LensUpperRightY { using type = Properties::UndefinedProperty; };
template<class TypeTag, class MyTypeTag>
struct LensUpperRightZ { using type = Properties::UndefinedProperty; };
// Enable gravity
template<class TypeTag>
struct EnableGravity<TypeTag, Properties::TTag::LensBaseProblem>
{ static constexpr bool value = true; };
// define the properties specific for the lens problem
template<class TypeTag>
struct LensLowerLeftX<TypeTag, Properties::TTag::LensBaseProblem>
{
using type = GetPropType<TypeTag, Properties::Scalar>;
static constexpr type value = 1.0;
};
template<class TypeTag>
struct LensLowerLeftY<TypeTag, Properties::TTag::LensBaseProblem>
{
using type = GetPropType<TypeTag, Properties::Scalar>;
static constexpr type value = 2.0;
};
template<class TypeTag>
struct LensLowerLeftZ<TypeTag, Properties::TTag::LensBaseProblem>
{
using type = GetPropType<TypeTag, Properties::Scalar>;
static constexpr type value = 0.0;
};
template<class TypeTag>
struct LensUpperRightX<TypeTag, Properties::TTag::LensBaseProblem>
{
using type = GetPropType<TypeTag, Properties::Scalar>;
static constexpr type value = 4.0;
};
template<class TypeTag>
struct LensUpperRightY<TypeTag, Properties::TTag::LensBaseProblem>
{
using type = GetPropType<TypeTag, Properties::Scalar>;
static constexpr type value = 3.0;
};
template<class TypeTag>
struct LensUpperRightZ<TypeTag, Properties::TTag::LensBaseProblem>
{
using type = GetPropType<TypeTag, Properties::Scalar>;
static constexpr type value = 1.0;
};
// By default, include the intrinsic permeability tensor to the VTK output files
template<class TypeTag>
struct VtkWriteIntrinsicPermeabilities<TypeTag, Properties::TTag::LensBaseProblem>
{ static constexpr bool value = true; };
} // namespace Opm::Parameters
namespace Opm {
/*!
* \ingroup TestProblems
*
* \brief Soil contamination problem where DNAPL infiltrates a fully
* water saturated medium.
*
* The domain is sized 6m times 4m and features a rectangular lens
* with low permeablility which spans from (1m, 2m) to (4m, 3m)
* and is surrounded by a medium with higher permability. Note that
* this problem is discretized using only two dimensions, so from the
* point of view of the model, the depth of the domain is implicitly
* assumed to be 1 m everywhere.
*
* On the top and the bottom of the domain no-flow boundary conditions
* are used, while free-flow conditions apply on the left and right
* boundaries; DNAPL is injected at the top boundary from 3m to 4m at
* a rate of 0.04 kg/(s m^2).
*
* At the boundary on the left, a free-flow condition using the
* hydrostatic pressure scaled by a factor of 1.125 is imposed, while
* on the right, it is just the hydrostatic pressure. The DNAPL
* saturation on both sides is zero.
*/
template <class TypeTag>
class LensProblem : public GetPropType<TypeTag, Properties::BaseProblem>
{
2020-06-10 06:49:42 -05:00
using ParentType = GetPropType<TypeTag, Properties::BaseProblem>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using GridView = GetPropType<TypeTag, Properties::GridView>;
using Indices = GetPropType<TypeTag, Properties::Indices>;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
using WettingPhase = GetPropType<TypeTag, Properties::WettingPhase>;
using NonwettingPhase = GetPropType<TypeTag, Properties::NonwettingPhase>;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using Simulator = GetPropType<TypeTag, Properties::Simulator>;
using Model = GetPropType<TypeTag, Properties::Model>;
enum {
// number of phases
numPhases = FluidSystem::numPhases,
// phase indices
wettingPhaseIdx = FluidSystem::wettingPhaseIdx,
nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx,
// equation indices
contiNEqIdx = Indices::conti0EqIdx + nonWettingPhaseIdx,
// Grid and world dimension
dim = GridView::dimension,
dimWorld = GridView::dimensionworld
};
2020-06-10 06:49:42 -05:00
using EqVector = GetPropType<TypeTag, Properties::EqVector>;
using RateVector = GetPropType<TypeTag, Properties::RateVector>;
using BoundaryRateVector = GetPropType<TypeTag, Properties::BoundaryRateVector>;
using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
2020-06-10 06:49:42 -05:00
using CoordScalar = typename GridView::ctype;
using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
2020-06-10 06:49:42 -05:00
using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
public:
/*!
* \copydoc Doxygen::defaultProblemConstructor
*/
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
LensProblem(Simulator& simulator)
: ParentType(simulator)
{ }
/*!
* \copydoc FvBaseProblem::finishInit
*/
void finishInit()
{
ParentType::finishInit();
eps_ = 3e-6;
FluidSystem::init();
temperature_ = 273.15 + 20; // -> 20°C
lensLowerLeft_[0] = Parameters::get<TypeTag, Parameters::LensLowerLeftX>();
lensLowerLeft_[1] = Parameters::get<TypeTag, Parameters::LensLowerLeftY>();
lensUpperRight_[0] = Parameters::get<TypeTag, Parameters::LensUpperRightX>();
lensUpperRight_[1] = Parameters::get<TypeTag, Parameters::LensUpperRightY>();
if constexpr (dim == 3) {
lensLowerLeft_[2] = Parameters::get<TypeTag, Parameters::LensLowerLeftZ>();
lensUpperRight_[2] = Parameters::get<TypeTag, Parameters::LensUpperRightZ>();
}
// residual saturations
lensMaterialParams_.setResidualSaturation(wettingPhaseIdx, 0.18);
lensMaterialParams_.setResidualSaturation(nonWettingPhaseIdx, 0.0);
outerMaterialParams_.setResidualSaturation(wettingPhaseIdx, 0.05);
outerMaterialParams_.setResidualSaturation(nonWettingPhaseIdx, 0.0);
// parameters for the Van Genuchten law: alpha and n
lensMaterialParams_.setVgAlpha(0.00045);
lensMaterialParams_.setVgN(7.3);
outerMaterialParams_.setVgAlpha(0.0037);
outerMaterialParams_.setVgN(4.7);
2013-11-06 07:50:01 -06:00
lensMaterialParams_.finalize();
outerMaterialParams_.finalize();
lensK_ = this->toDimMatrix_(9.05e-12);
outerK_ = this->toDimMatrix_(4.6e-10);
if (dimWorld == 3) {
this->gravity_ = 0;
this->gravity_[1] = -9.81;
}
}
/*!
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
*/
static void registerParameters()
{
ParentType::registerParameters();
Parameters::registerParam<TypeTag, Parameters::LensLowerLeftX>
2024-04-04 05:49:01 -05:00
("The x-coordinate of the lens' lower-left corner [m].");
Parameters::registerParam<TypeTag, Parameters::LensLowerLeftY>
2024-04-04 05:49:01 -05:00
("The y-coordinate of the lens' lower-left corner [m].");
Parameters::registerParam<TypeTag, Parameters::LensUpperRightX>
2024-04-04 05:49:01 -05:00
("The x-coordinate of the lens' upper-right corner [m].");
Parameters::registerParam<TypeTag, Parameters::LensUpperRightY>
2024-04-04 05:49:01 -05:00
("The y-coordinate of the lens' upper-right corner [m].");
if constexpr (dim == 3) {
Parameters::registerParam<TypeTag, Parameters::LensLowerLeftZ>
2024-04-04 05:49:01 -05:00
("The z-coordinate of the lens' lower-left corner [m].");
Parameters::registerParam<TypeTag, Parameters::LensUpperRightZ>
2024-04-04 05:49:01 -05:00
("The z-coordinate of the lens' upper-right corner [m].");
}
Parameters::SetDefault<Parameters::CellsX>(48);
Parameters::SetDefault<Parameters::CellsY>(32);
Parameters::SetDefault<Parameters::DomainSizeX<Scalar>>(6.0);
Parameters::SetDefault<Parameters::DomainSizeY<Scalar>>(4.0);
if constexpr (dim == 3) {
Parameters::SetDefault<Parameters::CellsZ>(16);
Parameters::SetDefault<Parameters::DomainSizeZ<Scalar>>(1.0);
}
// Use forward differences
using LLS = GetPropType<TypeTag, Properties::LocalLinearizerSplice>;
constexpr bool useFD = std::is_same_v<LLS, Properties::TTag::FiniteDifferenceLocalLinearizer>;
if constexpr (useFD) {
Parameters::SetDefault<Parameters::NumericDifferenceMethod>(+1);
}
Parameters::SetDefault<Parameters::EndTime<Scalar>>(30e3);
Parameters::SetDefault<Parameters::EnableIntensiveQuantityCache>(true);
Parameters::SetDefault<Parameters::EnableStorageCache>(true);
Parameters::SetDefault<Parameters::InitialTimeStepSize<Scalar>>(250.0);
}
/*!
* \copydoc FvBaseProblem::briefDescription
*/
static std::string briefDescription()
{
std::string thermal = "isothermal";
constexpr bool enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>();
if constexpr (enableEnergy)
thermal = "non-isothermal";
std::string deriv = "finite difference";
2020-06-10 06:49:42 -05:00
using LLS = GetPropType<TypeTag, Properties::LocalLinearizerSplice>;
constexpr bool useAutoDiff = std::is_same_v<LLS, Properties::TTag::AutoDiffLocalLinearizer>;
if constexpr (useAutoDiff) {
deriv = "automatic differentiation";
}
std::string disc = "vertex centered finite volume";
2020-06-10 06:49:42 -05:00
using D = GetPropType<TypeTag, Properties::Discretization>;
constexpr bool useEcfv = std::is_same<D, Opm::EcfvDiscretization<TypeTag>>::value;
if constexpr (useEcfv)
disc = "element centered finite volume";
return std::string("")+
"Ground remediation problem where a dense oil infiltrates "+
"an aquifer with an embedded low-permability lens. " +
"This is the binary for the "+thermal+" variant using "+deriv+
"and the "+disc+" discretization";
}
/*!
* \name Soil parameters
*/
//! \{
/*!
Implement the element centered finite volume spatial discretization This makes eWoms multi-discretization capable. Along the way, this fixes some bugs and does a medium sized reorganization of the source tree. This is a squashed patch of the following commits: -------- 1st commit message: add initial version of the element centered finite volume discretization currently, it is a misnomer as it is just a copy of the vertex centered discretization plus some renames... -------- 2nd commit message: rename [VE]cfvModel -> [VE]cfvDiscretization -------- 3rd commit message: ecfv: prelimary changes required to make it compile but not work yet... -------- 4th commit message: Rename *FvElementGeometry to *Stencil "Stencil" seems to be the standard expression for this concept... (also, it is not specific to finite volume methods and is shorter.) -------- 5th commit message: refactor the stencil class for the element centered finite volume discretization -------- 6th commit message: ECFV: some work on the stencil class -------- 7th commit message: ECFV: make the boundary handling code compile -------- 8th commit message: rename elemContext() to elementContext() -------- 9th commit message: ECFV: make the VTK output modules compile -------- 10th commit message: stencil: introduce the concept of primary DOFs also save an vector of all element pointers in the stencil. -------- 11th commit message: ECFV: try to fix assembly; add missing timeIdx arguments to the num*() methods -------- 12th commit message: ECFV: fix stupid mistake in the assembler -------- 13th commit message: ECFV: remove a few implicit DOF == vertex assumptions the black-oil example now runs without valgrind complaints until it encounters a negative oil mole fraction. -------- 14th commit message: VCFV: make everything compile again all vertex centered FV examples should now work again... -------- 15th commit message: rename [ev]cfvmodel.hh to [ev]cfvdiscretization.hh the classes have already been renamed. -------- 16th commit message: ECFV: make it work to the point where it can write out the initial solution. -------- 17th commit message: ECFV: make it work the local residual/jacobian needed some work in distinguishing primary and secondary DOFs and there was an minor issue with the serialization code. for some reason, it seems still not correct. (-> convergence is too slow.) -------- 18th commit message: VCFV: make it compile for the black oil model again -------- 19th commit message: VCFV: make it compile with the remaining models again -------- 20th commit message: flash model: make it work with ECFV although this breaks its compatibility with VCFV. (-> next commit) -------- 21st commit message: adapt the VCFV to make it compatible with the flash model again -------- 22nd commit message: make all models compile with VCFV again -------- 23rd commit message: VCFV: more cleanups of the stencil VcfvStencil now does not have any public attributes anymore. TODO: do not export attributes in the SubControlVolume and SubControlVolumeFace classes. -------- 24th commit message: VCFV: actually update the element pointer -------- 25th commit message: change the blackoil model back to ECFV -------- 26th commit message: immiscible model: make it compatible with the ECFV discretization -------- 27th commit message: PVS model: make it work with ECFV -------- 28th commit message: NCP model: make it work with ECFV -------- 29th commit message: rename Vcfv*VelocityModule to *VelocityModule -------- 30th commit message: richards model: make it work with ECFV -------- 31st commit message: unify the ECFV and the VCFV VTK output modules and other cleanups -------- 32nd commit message: unify the common code of the VCFV and the ECFV disctretizations -------- 33rd commit message: unify the element contexts between element and vertex centered finite volumes -------- 34th commit message: unify the local jacobian class of the finite volume discretizations -------- 35th commit message: replace [VE]vcf(LocalResidual|ElementContext|BoundaryContext|ConstraintsContext) by generic code -------- 36th commit message: replace the [EV]cfvLocalResidual by generic code -------- 37th commit message: unify the MultiPhaseProblem and Problem classes, introduce NullBorderListCreator -------- 38th commit message: remove the discretization specific boundary context -------- 39th commit message: unify the [EV]cfvDiscretization classes -------- 40th commit message: Unify [EV]cfvMultiPhaseFluxVariables -------- 41st commit message: Unify the [EC]cfvNewton* classes -------- 42nd commit message: Unify [EV]cfvVolumeVariables -------- 43rd commit message: unify [EV]cfvAssembler -------- 44th commit message: unified flux variables: fix stupid mistake when calculating pressure gradients -------- 45th commit message: unify what's to unify for the [EV]CFV properties -------- 46th commit message: make the method to calculate gradients and values at flux approximation points changeable Currently, this is used by the vertex centered finite volume method to be able to use P1-finite element gradients instead of two-point ones... -------- 47th commit message: make the restart code work correctly, use the correct DofMapper for VCFV -------- 48th commit message: actually use the gradient calculator in a model the immiscible model in this case -------- 49th commit message: move some files around to where they belong, use the new gradient calculation code in all models TODO: proper handling of boundary gradients -------- 50th commit message: fix the stokes model currently it only works with the vertex centered finite volume discretization, but the plan is to soon move it to a staggered grid scheme anyway... -------- 51st commit message: move all models back to using the vertex centered finite volume discretization by default -------- 52nd commit message: models: some variable renames and documentation fixes - scv -> dof - vert -> dof - vertex -> dof - replace 'VCFV' - fix some typos -------- 53rd commit message: don't expect UG anymore since it is quite non-free and hard to get. we now use ALUGrid instead! -------- 54th commit message: temporarily disable jacobian recycling -------- 55th commit message: fix writing/reading restart files using the generic code -------- 56th commit message: fix bug where fluxes were only counted once in the stencil this only affected the vertex centered finite volumes discretization... -------- 57th commit message: boundary gradients: use the center of the sub-control volume adjacent to a boundary segment -------- 58th commit message: make it compile on GCC -------- 59th commit message: get rid of most hacks for this, partial reassemble and jacobian recycling was brought back. For the this and the remaining stuff the main trick is the introduction of the GridCommHandleFactory concept which constructs communication handles suited for the respective spatial discretization... -------- 60th commit message: fix a few annoying bugs first, default the convergence criterion for the linear solver did not honor the initial residual which lead to linear solver breakdowns, then some debugging code was left in the discrete fracture model and then there was a bug in the TP gradient approximation class... this has the consequence that we need a new reference solution for the discrete fracture problem... -------- 61st commit message: iterative linear solver: remove the code for the non-default convergence criteria -------- 62nd commit message: provide the FE cache instead of the local FE this fixes a segfault in the stokes model caused by the fact that the local FE was not initialized at this point. -------- 63rd commit message: (Navier-)Stokes: fix bug due to the transition to unit normals now, all tests pass for this branch. The only things which need to be fixed are some annoying performance regressions compared to master and some bug in the splices feature of the property system... -------- 64th commit message: some fix for the local residual of the immiscible model -------- 65th commit message: Navier-Stokes: implement SCV center gradients There seems to be a bug in the previous implementation (the jacobian inverse transposed is evaluated using the local, not the global geometry), so the reference solution for the stokes2c test problem has also been updated... -------- 66th commit message: remove the ALUGrid specialization of the LensGridCreator and the YaspGrid one for the fingerproblem using different grid seems to sometimes cause a different vertex order, which in turn causes the respective test to fail if the reference solution was computed using the other grid... -------- 67th commit message: VCFV: use the correct BorderListCreator this makes MPI parallel computations work again. apart from performance regressions, this branch does not exhibit any known regressions compared to master anymore... -------- 68th commit message: make verything compile with the element centered finite volume discretization except the Navier-Stokes and the two-phase DFM models, of course... -------- 69th commit message: minor fixes - make the navier-stokes model slighly more generic by using the proper (in,ex)teriorIndex() methods on sub-control volumes - make the signature of the calculateValue() template method of the common two-point gradient approximator match the one of the vertex centered finite volume one -------- 70th commit message: fix fallout from the Big Rebase -------- 71st commit message: ECFV: some bugs in the boundary -------- 72nd commit message: make computeFlux() compute area-specific quantities -------- 73rd commit message: fix more bugs in the element centered FV discretization now eWoms should match Dumux pretty closely... -------- 74th commit message: coalesce the common code of the multi phase porous medium models into "MultiPhaseBaseModel" -------- 75th commit message: update reference solutions these were changed because of the screw-up with the area of boundary segments... -------- 76th commit message: rename "ImplicitBase" to "FvBase" because in eWoms, everything is implicit and these are currently the base classes for all finite volume discretizations. -------- 77th commit message: make the spatial discretization selectable using a splice This requires an opm-core with a the patches from https://github.com/OPM/opm-core/pull/446 merged... -------- 78th commit message: rename the properties used for splices to *Splice -------- 79th commit message: move the files in 'tests/models' to 'tests' since 'tests' was empty except for the 'models' subdirectory... -------- 80th commit message: improve and fix the tutorial -------- 81st commit message: remove the -fno-strict-aliasing flag from the provided option files seems like recent versions of Dune have been adapted... -------- 82nd commit message: also compile all CO2 injection simulations using the element centered finite volume discretization -------- 83rd commit message: PVS model: make it work properly with the element-centered finite volume discretiation because DOF != number of vertices
2013-12-12 05:52:44 -06:00
* \copydoc FvBaseMultiPhaseProblem::intrinsicPermeability
*/
template <class Context>
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
const DimMatrix& intrinsicPermeability(const Context& context, unsigned spaceIdx,
unsigned timeIdx) const
{
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
const GlobalPosition& globalPos = context.pos(spaceIdx, timeIdx);
if (isInLens_(globalPos))
return lensK_;
return outerK_;
}
/*!
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
*/
template <class Context>
2022-08-02 03:40:48 -05:00
Scalar porosity(const Context& /*context*/,
unsigned /*spaceIdx*/,
unsigned /*timeIdx*/) const
{ return 0.4; }
/*!
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
*/
template <class Context>
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
const MaterialLawParams& materialLawParams(const Context& context,
unsigned spaceIdx, unsigned timeIdx) const
{
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
const GlobalPosition& globalPos = context.pos(spaceIdx, timeIdx);
if (isInLens_(globalPos))
return lensMaterialParams_;
return outerMaterialParams_;
}
/*!
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
*/
template <class Context>
2022-08-02 03:40:48 -05:00
Scalar temperature(const Context& /*context*/,
unsigned /*spaceIdx*/,
unsigned /*timeIdx*/) const
{ return temperature_; }
//! \}
/*!
* \name Auxiliary methods
*/
//! \{
/*!
* \copydoc FvBaseProblem::name
*/
std::string name() const
{
2020-06-10 06:49:42 -05:00
using LLS = GetPropType<TypeTag, Properties::LocalLinearizerSplice>;
constexpr bool useAutoDiff = std::is_same_v<LLS, Properties::TTag::AutoDiffLocalLinearizer>;
using FM = GetPropType<TypeTag, Properties::FluxModule>;
constexpr bool useTrans = std::is_same_v<FM, Opm::TransFluxModule<TypeTag>>;
std::ostringstream oss;
2014-07-11 09:44:45 -05:00
oss << "lens_" << Model::name()
<< "_" << Model::discretizationName()
<< "_" << (useAutoDiff?"ad":"fd");
if (useTrans)
oss << "_trans";
return oss.str();
}
/*!
* \copydoc FvBaseProblem::beginTimeStep
*/
void beginTimeStep()
{ }
/*!
* \copydoc FvBaseProblem::beginIteration
*/
void beginIteration()
{ }
/*!
* \copydoc FvBaseProblem::endTimeStep
*/
void endTimeStep()
{
#ifndef NDEBUG
//this->model().checkConservativeness();
// Calculate storage terms
EqVector storage;
this->model().globalStorage(storage);
// Write mass balance information for rank 0
if (this->gridView().comm().rank() == 0) {
std::cout << "Storage: " << storage << std::endl << std::flush;
}
#endif // NDEBUG
}
//! \}
/*!
* \name Boundary conditions
*/
//! \{
/*!
* \copydoc FvBaseProblem::boundary
*/
template <class Context>
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
void boundary(BoundaryRateVector& values,
const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
if (onLeftBoundary_(pos) || onRightBoundary_(pos)) {
// free flow boundary. we assume incompressible fluids
Scalar densityW = WettingPhase::density(temperature_, /*pressure=*/Scalar(1e5));
Scalar densityN = NonwettingPhase::density(temperature_, /*pressure=*/Scalar(1e5));
Scalar T = temperature(context, spaceIdx, timeIdx);
Scalar pw, Sw;
// set wetting phase pressure and saturation
if (onLeftBoundary_(pos)) {
2013-12-27 11:25:54 -06:00
Scalar height = this->boundingBoxMax()[1] - this->boundingBoxMin()[1];
Scalar depth = this->boundingBoxMax()[1] - pos[1];
Scalar alpha = (1 + 1.5 / height);
// hydrostatic pressure scaled by alpha
pw = 1e5 - alpha * densityW * this->gravity()[1] * depth;
Sw = 1.0;
}
else {
2013-12-27 11:25:54 -06:00
Scalar depth = this->boundingBoxMax()[1] - pos[1];
// hydrostatic pressure
pw = 1e5 - densityW * this->gravity()[1] * depth;
Sw = 1.0;
}
// specify a full fluid state using pw and Sw
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
const MaterialLawParams& matParams = this->materialLawParams(context, spaceIdx, timeIdx);
Opm::ImmiscibleFluidState<Scalar, FluidSystem,
/*storeEnthalpy=*/false> fs;
fs.setSaturation(wettingPhaseIdx, Sw);
fs.setSaturation(nonWettingPhaseIdx, 1 - Sw);
fs.setTemperature(T);
Scalar pC[numPhases];
MaterialLaw::capillaryPressures(pC, matParams, fs);
fs.setPressure(wettingPhaseIdx, pw);
fs.setPressure(nonWettingPhaseIdx, pw + pC[nonWettingPhaseIdx] - pC[wettingPhaseIdx]);
fs.setDensity(wettingPhaseIdx, densityW);
fs.setDensity(nonWettingPhaseIdx, densityN);
fs.setViscosity(wettingPhaseIdx, WettingPhase::viscosity(temperature_, fs.pressure(wettingPhaseIdx)));
fs.setViscosity(nonWettingPhaseIdx, NonwettingPhase::viscosity(temperature_, fs.pressure(nonWettingPhaseIdx)));
// impose an freeflow boundary condition
values.setFreeFlow(context, spaceIdx, timeIdx, fs);
}
else if (onInlet_(pos)) {
RateVector massRate(0.0);
massRate = 0.0;
massRate[contiNEqIdx] = -0.04; // kg / (m^2 * s)
// impose a forced flow boundary
values.setMassRate(massRate);
}
else {
// no flow boundary
values.setNoFlow();
}
}
//! \}
/*!
* \name Volumetric terms
*/
//! \{
/*!
* \copydoc FvBaseProblem::initial
*/
template <class Context>
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
void initial(PrimaryVariables& values, const Context& context, unsigned spaceIdx, unsigned timeIdx) const
{
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
2013-12-27 11:25:54 -06:00
Scalar depth = this->boundingBoxMax()[1] - pos[1];
Opm::ImmiscibleFluidState<Scalar, FluidSystem> fs;
fs.setPressure(wettingPhaseIdx, /*pressure=*/1e5);
Scalar Sw = 1.0;
fs.setSaturation(wettingPhaseIdx, Sw);
fs.setSaturation(nonWettingPhaseIdx, 1 - Sw);
fs.setTemperature(temperature_);
typename FluidSystem::template ParameterCache<Scalar> paramCache;
paramCache.updatePhase(fs, wettingPhaseIdx);
Scalar densityW = FluidSystem::density(fs, paramCache, wettingPhaseIdx);
// hydrostatic pressure (assuming incompressibility)
Scalar pw = 1e5 - densityW * this->gravity()[1] * depth;
// calculate the capillary pressure
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
const MaterialLawParams& matParams = this->materialLawParams(context, spaceIdx, timeIdx);
Scalar pC[numPhases];
MaterialLaw::capillaryPressures(pC, matParams, fs);
// make a full fluid state
fs.setPressure(wettingPhaseIdx, pw);
fs.setPressure(nonWettingPhaseIdx, pw + (pC[wettingPhaseIdx] - pC[nonWettingPhaseIdx]));
// assign the primary variables
values.assignNaive(fs);
}
/*!
* \copydoc FvBaseProblem::source
*
* For this problem, the source term of all components is 0
* everywhere.
*/
template <class Context>
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
void source(RateVector& rate,
2022-08-02 03:40:48 -05:00
const Context& /*context*/,
unsigned /*spaceIdx*/,
unsigned /*timeIdx*/) const
{ rate = Scalar(0.0); }
//! \}
private:
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
bool isInLens_(const GlobalPosition& pos) const
{
for (unsigned i = 0; i < dim; ++i) {
if (pos[i] < lensLowerLeft_[i] - eps_ || pos[i] > lensUpperRight_[i]
+ eps_)
return false;
}
return true;
}
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
bool onLeftBoundary_(const GlobalPosition& pos) const
2013-12-27 11:25:54 -06:00
{ return pos[0] < this->boundingBoxMin()[0] + eps_; }
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
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_; }
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
bool onLowerBoundary_(const GlobalPosition& pos) const
2013-12-27 11:25:54 -06:00
{ return pos[1] < this->boundingBoxMin()[1] + eps_; }
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
bool onUpperBoundary_(const GlobalPosition& pos) const
2013-12-27 11:25:54 -06:00
{ return pos[1] > this->boundingBoxMax()[1] - eps_; }
fix most pedantic compiler warnings in the basic infrastructure i.e., using clang 3.8 to compile the test suite with the following flags: ``` -Weverything -Wno-documentation -Wno-documentation-unknown-command -Wno-c++98-compat -Wno-c++98-compat-pedantic -Wno-undef -Wno-padded -Wno-global-constructors -Wno-exit-time-destructors -Wno-weak-vtables -Wno-float-equal ``` should not produce any warnings anymore. In my opinion the only flag which would produce beneficial warnings is -Wdocumentation. This has not been fixed in this patch because writing documentation is left for another day (or, more likely, year). note that this patch consists of a heavy dose of the OPM_UNUSED macro and plenty of static_casts (to fix signedness issues). Fixing the singedness issues were quite a nightmare and the fact that the Dune API is quite inconsistent in that regard was not exactly helpful. :/ Finally this patch includes quite a few formatting changes (e.g., all occurences of 'T &t' should be changed to `T& t`) and some fixes for minor issues which I've found during the excercise. I've made sure that all unit tests the test suite still pass successfully and I've made sure that flow_ebos still works for Norne and that it did not regress w.r.t. performance. (Note that this patch does not fix compiler warnings triggered `ebos` and `flow_ebos` but only those caused by the basic infrastructure or the unit tests.) v2: fix the warnings that occur if the dune-localfunctions module is not available. thanks to [at]atgeirr for testing. v3: fix dune 2.3 build issue
2016-11-07 08:14:07 -06:00
bool onInlet_(const GlobalPosition& pos) const
{
2013-12-27 11:25:54 -06:00
Scalar width = this->boundingBoxMax()[0] - this->boundingBoxMin()[0];
Scalar lambda = (this->boundingBoxMax()[0] - pos[0]) / width;
return onUpperBoundary_(pos) && 0.5 < lambda && lambda < 2.0 / 3.0;
}
GlobalPosition lensLowerLeft_;
GlobalPosition lensUpperRight_;
DimMatrix lensK_;
DimMatrix outerK_;
MaterialLawParams lensMaterialParams_;
MaterialLawParams outerMaterialParams_;
Scalar temperature_;
Scalar eps_;
};
} // namespace Opm
#endif