opm-simulators/examples/problems/fingerproblem.hh

597 lines
19 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::FingerProblem
*/
#ifndef EWOMS_FINGER_PROBLEM_HH
#define EWOMS_FINGER_PROBLEM_HH
2019-09-16 03:09:33 -05:00
#include <opm/models/io/structuredgridvanguard.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/ParkerLenhard.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/Air.hpp>
#include <opm/models/immiscible/immiscibleproperties.hh>
#include <opm/models/discretization/common/restrictprolong.hh>
#if HAVE_DUNE_ALUGRID
#include <dune/alugrid/grid.hh>
#endif
#include <dune/common/version.hh>
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <dune/grid/utility/persistentcontainer.hh>
#include <vector>
#include <string>
namespace Opm {
template <class TypeTag>
class FingerProblem;
} // namespace Opm
namespace Opm::Properties {
// Create new type tags
namespace TTag {
struct FingerBaseProblem { using InheritsFrom = std::tuple<StructuredGridVanguard>; };
} // end namespace TTag
#if HAVE_DUNE_ALUGRID
// use dune-alugrid if available
template<class TypeTag>
struct Grid<TypeTag, TTag::FingerBaseProblem>
{ using type = Dune::ALUGrid</*dim=*/2,
/*dimWorld=*/2,
Dune::cube,
Dune::nonconforming>; };
#endif
// declare the properties used by the finger problem
template<class TypeTag, class MyTypeTag>
struct InitialWaterSaturation { using type = UndefinedProperty; };
// Set the problem property
template<class TypeTag>
struct Problem<TypeTag, TTag::FingerBaseProblem> { using type = Opm::FingerProblem<TypeTag>; };
// Set the wetting phase
template<class TypeTag>
struct WettingPhase<TypeTag, TTag::FingerBaseProblem>
{
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::FingerBaseProblem>
{
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::GasPhase<Scalar, Opm::Air<Scalar> >;
};
// Set the material Law
template<class TypeTag>
struct MaterialLaw<TypeTag, TTag::FingerBaseProblem>
{
2020-06-10 06:49:42 -05:00
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
using Traits = Opm::TwoPhaseMaterialTraits<Scalar,
/*wettingPhaseIdx=*/FluidSystem::wettingPhaseIdx,
/*nonWettingPhaseIdx=*/FluidSystem::nonWettingPhaseIdx>;
2013-11-06 07:50:01 -06:00
// use the parker-lenhard hysteresis law
2020-06-10 06:49:42 -05:00
using ParkerLenhard = Opm::ParkerLenhard<Traits>;
using type = ParkerLenhard;
};
// Write the solutions of individual newton iterations?
template<class TypeTag>
struct NewtonWriteConvergence<TypeTag, TTag::FingerBaseProblem> { static constexpr bool value = false; };
// Use forward differences instead of central differences
template<class TypeTag>
struct NumericDifferenceMethod<TypeTag, TTag::FingerBaseProblem> { static constexpr int value = +1; };
// Enable constraints
template<class TypeTag>
struct EnableConstraints<TypeTag, TTag::FingerBaseProblem> { static constexpr int value = true; };
// Enable gravity
template<class TypeTag>
struct EnableGravity<TypeTag, TTag::FingerBaseProblem> { static constexpr bool value = true; };
// define the properties specific for the finger problem
template<class TypeTag>
struct DomainSizeX<TypeTag, TTag::FingerBaseProblem>
{
using type = GetPropType<TypeTag, Scalar>;
static constexpr type value = 0.1;
};
template<class TypeTag>
struct DomainSizeY<TypeTag, TTag::FingerBaseProblem>
{
using type = GetPropType<TypeTag, Scalar>;
static constexpr type value = 0.3;
};
template<class TypeTag>
struct DomainSizeZ<TypeTag, TTag::FingerBaseProblem>
{
using type = GetPropType<TypeTag, Scalar>;
static constexpr type value = 0.1;
};
template<class TypeTag>
struct InitialWaterSaturation<TypeTag, TTag::FingerBaseProblem>
{
using type = GetPropType<TypeTag, Scalar>;
static constexpr type value = 0.01;
};
template<class TypeTag>
struct CellsX<TypeTag, TTag::FingerBaseProblem> { static constexpr int value = 20; };
template<class TypeTag>
struct CellsY<TypeTag, TTag::FingerBaseProblem> { static constexpr int value = 70; };
template<class TypeTag>
struct CellsZ<TypeTag, TTag::FingerBaseProblem> { static constexpr int value = 1; };
// The default for the end time of the simulation
template<class TypeTag>
struct EndTime<TypeTag, TTag::FingerBaseProblem>
{
using type = GetPropType<TypeTag, Scalar>;
static constexpr type value = 215;
};
// The default for the initial time step size of the simulation
template<class TypeTag>
struct InitialTimeStepSize<TypeTag, TTag::FingerBaseProblem>
{
using type = GetPropType<TypeTag, Scalar>;
static constexpr type value = 10;
};
} // namespace Opm::Properties
namespace Opm {
/*!
* \ingroup TestProblems
*
* \brief Two-phase problem featuring some gravity-driven saturation
* fingers.
*
* The domain of this problem is sized 10cm times 1m and is initially
* dry. Water is then injected at three locations on the top of the
* domain which leads to gravity fingering. The boundary conditions
* used are no-flow for the left and right and top of the domain and
* free-flow at the bottom. This problem uses the Parker-Lenhard
* hystersis model which might lead to non-monotonic saturation in the
* fingers if the right material parameters is chosen and the spatial
* discretization is fine enough.
*/
template <class TypeTag>
class FingerProblem : public GetPropType<TypeTag, Properties::BaseProblem>
{
//!\cond SKIP_THIS
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 Constraints = GetPropType<TypeTag, Properties::Constraints>;
using Model = GetPropType<TypeTag, Properties::Model>;
enum {
// number of phases
numPhases = FluidSystem::numPhases,
// phase indices
wettingPhaseIdx = FluidSystem::wettingPhaseIdx,
nonWettingPhaseIdx = FluidSystem::nonWettingPhaseIdx,
// equation indices
contiWettingEqIdx = Indices::conti0EqIdx + wettingPhaseIdx,
// Grid and world dimension
dim = GridView::dimension,
dimWorld = GridView::dimensionworld
};
2020-06-10 06:49:42 -05:00
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
using Stencil = GetPropType<TypeTag, Properties::Stencil> ;
enum { codim = Stencil::Entity::codimension };
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>;
2020-06-10 06:49:42 -05:00
using ParkerLenhard = typename GetProp<TypeTag, Properties::MaterialLaw>::ParkerLenhard;
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>;
using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
2020-06-10 06:49:42 -05:00
using Grid = typename GridView :: Grid;
2020-06-10 06:49:42 -05:00
using MaterialLawParamsContainer = Dune::PersistentContainer< Grid, std::shared_ptr< MaterialLawParams > > ;
//!\endcond
public:
2020-06-10 06:49:42 -05:00
using RestrictProlongOperator = CopyRestrictProlong< Grid, MaterialLawParamsContainer >;
/*!
* \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
FingerProblem(Simulator& simulator)
: ParentType(simulator),
materialParams_( simulator.vanguard().grid(), codim )
{
}
/*!
* \name Auxiliary methods
*/
//! \{
/*!
* \brief \copydoc FvBaseProblem::restrictProlongOperator
*/
RestrictProlongOperator restrictProlongOperator()
{
return RestrictProlongOperator( materialParams_ );
}
/*!
* \copydoc FvBaseProblem::name
*/
std::string name() const
{ return
std::string("finger") +
"_" + Model::name() +
"_" + Model::discretizationName() +
(this->model().enableGridAdaptation()?"_adaptive":"");
}
/*!
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();
2024-04-04 05:49:01 -05:00
Parameters::registerParam<TypeTag, Properties::InitialWaterSaturation>
("The initial saturation in the domain [] of the wetting phase");
}
/*!
* \copydoc FvBaseProblem::finishInit()
*/
void finishInit()
{
ParentType::finishInit();
eps_ = 3e-6;
temperature_ = 273.15 + 20; // -> 20°C
FluidSystem::init();
// parameters for the Van Genuchten law of the main imbibition
// and the main drainage curves.
micParams_.setVgAlpha(0.0037);
micParams_.setVgN(4.7);
2013-11-06 07:50:01 -06:00
micParams_.finalize();
mdcParams_.setVgAlpha(0.0037);
mdcParams_.setVgN(4.7);
2013-11-06 07:50:01 -06:00
mdcParams_.finalize();
// initialize the material parameter objects of the individual
// finite volumes, resize will resize the container to the number of elements
materialParams_.resize();
for (auto it = materialParams_.begin(),
end = materialParams_.end(); it != end; ++it ) {
std::shared_ptr< MaterialLawParams >& materialParams = *it ;
if( ! materialParams )
{
materialParams.reset( new MaterialLawParams() );
materialParams->setMicParams(&micParams_);
materialParams->setMdcParams(&mdcParams_);
materialParams->setSwr(0.0);
materialParams->setSnr(0.1);
materialParams->finalize();
ParkerLenhard::reset(*materialParams);
}
}
K_ = this->toDimMatrix_(4.6e-10);
setupInitialFluidState_();
}
/*!
* \copydoc FvBaseProblem::endTimeStep
*/
void endTimeStep()
{
#ifndef NDEBUG
// checkConservativeness() does not include the effect of constraints, so we
// disable it for this problem...
//this->model().checkConservativeness();
// Calculate storage terms
EqVector storage;
this->model().globalStorage(storage);
// Write mass balance information for rank 0
if (this->gridView().comm().rank() == 0) {
std::cout << "Storage: " << storage << std::endl << std::flush;
}
#endif // NDEBUG
// update the history of the hysteresis law
ElementContext elemCtx(this->simulator());
for (const auto& elem : elements(this->gridView())) {
elemCtx.updateAll(elem);
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
size_t numDofs = elemCtx.numDof(/*timeIdx=*/0);
for (unsigned scvIdx = 0; scvIdx < numDofs; ++scvIdx)
{
MaterialLawParams& materialParam = materialLawParams( elemCtx, scvIdx, /*timeIdx=*/0 );
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 auto& fs = elemCtx.intensiveQuantities(scvIdx, /*timeIdx=*/0).fluidState();
ParkerLenhard::update(materialParam, fs);
}
}
}
//! \}
/*!
* \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::temperature
*/
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
Scalar temperature(const Context& /*context*/, unsigned /*spaceIdx*/, unsigned /*timeIdx*/) const
{ return temperature_; }
/*!
Implement the element centered finite volume spatial discretization This makes eWoms multi-discretization capable. Along the way, this fixes some bugs and does a medium sized reorganization of the source tree. This is a squashed patch of the following commits: -------- 1st commit message: add initial version of the element centered finite volume discretization currently, it is a misnomer as it is just a copy of the vertex centered discretization plus some renames... -------- 2nd commit message: rename [VE]cfvModel -> [VE]cfvDiscretization -------- 3rd commit message: ecfv: prelimary changes required to make it compile but not work yet... -------- 4th commit message: Rename *FvElementGeometry to *Stencil "Stencil" seems to be the standard expression for this concept... (also, it is not specific to finite volume methods and is shorter.) -------- 5th commit message: refactor the stencil class for the element centered finite volume discretization -------- 6th commit message: ECFV: some work on the stencil class -------- 7th commit message: ECFV: make the boundary handling code compile -------- 8th commit message: rename elemContext() to elementContext() -------- 9th commit message: ECFV: make the VTK output modules compile -------- 10th commit message: stencil: introduce the concept of primary DOFs also save an vector of all element pointers in the stencil. -------- 11th commit message: ECFV: try to fix assembly; add missing timeIdx arguments to the num*() methods -------- 12th commit message: ECFV: fix stupid mistake in the assembler -------- 13th commit message: ECFV: remove a few implicit DOF == vertex assumptions the black-oil example now runs without valgrind complaints until it encounters a negative oil mole fraction. -------- 14th commit message: VCFV: make everything compile again all vertex centered FV examples should now work again... -------- 15th commit message: rename [ev]cfvmodel.hh to [ev]cfvdiscretization.hh the classes have already been renamed. -------- 16th commit message: ECFV: make it work to the point where it can write out the initial solution. -------- 17th commit message: ECFV: make it work the local residual/jacobian needed some work in distinguishing primary and secondary DOFs and there was an minor issue with the serialization code. for some reason, it seems still not correct. (-> convergence is too slow.) -------- 18th commit message: VCFV: make it compile for the black oil model again -------- 19th commit message: VCFV: make it compile with the remaining models again -------- 20th commit message: flash model: make it work with ECFV although this breaks its compatibility with VCFV. (-> next commit) -------- 21st commit message: adapt the VCFV to make it compatible with the flash model again -------- 22nd commit message: make all models compile with VCFV again -------- 23rd commit message: VCFV: more cleanups of the stencil VcfvStencil now does not have any public attributes anymore. TODO: do not export attributes in the SubControlVolume and SubControlVolumeFace classes. -------- 24th commit message: VCFV: actually update the element pointer -------- 25th commit message: change the blackoil model back to ECFV -------- 26th commit message: immiscible model: make it compatible with the ECFV discretization -------- 27th commit message: PVS model: make it work with ECFV -------- 28th commit message: NCP model: make it work with ECFV -------- 29th commit message: rename Vcfv*VelocityModule to *VelocityModule -------- 30th commit message: richards model: make it work with ECFV -------- 31st commit message: unify the ECFV and the VCFV VTK output modules and other cleanups -------- 32nd commit message: unify the common code of the VCFV and the ECFV disctretizations -------- 33rd commit message: unify the element contexts between element and vertex centered finite volumes -------- 34th commit message: unify the local jacobian class of the finite volume discretizations -------- 35th commit message: replace [VE]vcf(LocalResidual|ElementContext|BoundaryContext|ConstraintsContext) by generic code -------- 36th commit message: replace the [EV]cfvLocalResidual by generic code -------- 37th commit message: unify the MultiPhaseProblem and Problem classes, introduce NullBorderListCreator -------- 38th commit message: remove the discretization specific boundary context -------- 39th commit message: unify the [EV]cfvDiscretization classes -------- 40th commit message: Unify [EV]cfvMultiPhaseFluxVariables -------- 41st commit message: Unify the [EC]cfvNewton* classes -------- 42nd commit message: Unify [EV]cfvVolumeVariables -------- 43rd commit message: unify [EV]cfvAssembler -------- 44th commit message: unified flux variables: fix stupid mistake when calculating pressure gradients -------- 45th commit message: unify what's to unify for the [EV]CFV properties -------- 46th commit message: make the method to calculate gradients and values at flux approximation points changeable Currently, this is used by the vertex centered finite volume method to be able to use P1-finite element gradients instead of two-point ones... -------- 47th commit message: make the restart code work correctly, use the correct DofMapper for VCFV -------- 48th commit message: actually use the gradient calculator in a model the immiscible model in this case -------- 49th commit message: move some files around to where they belong, use the new gradient calculation code in all models TODO: proper handling of boundary gradients -------- 50th commit message: fix the stokes model currently it only works with the vertex centered finite volume discretization, but the plan is to soon move it to a staggered grid scheme anyway... -------- 51st commit message: move all models back to using the vertex centered finite volume discretization by default -------- 52nd commit message: models: some variable renames and documentation fixes - scv -> dof - vert -> dof - vertex -> dof - replace 'VCFV' - fix some typos -------- 53rd commit message: don't expect UG anymore since it is quite non-free and hard to get. we now use ALUGrid instead! -------- 54th commit message: temporarily disable jacobian recycling -------- 55th commit message: fix writing/reading restart files using the generic code -------- 56th commit message: fix bug where fluxes were only counted once in the stencil this only affected the vertex centered finite volumes discretization... -------- 57th commit message: boundary gradients: use the center of the sub-control volume adjacent to a boundary segment -------- 58th commit message: make it compile on GCC -------- 59th commit message: get rid of most hacks for this, partial reassemble and jacobian recycling was brought back. For the this and the remaining stuff the main trick is the introduction of the GridCommHandleFactory concept which constructs communication handles suited for the respective spatial discretization... -------- 60th commit message: fix a few annoying bugs first, default the convergence criterion for the linear solver did not honor the initial residual which lead to linear solver breakdowns, then some debugging code was left in the discrete fracture model and then there was a bug in the TP gradient approximation class... this has the consequence that we need a new reference solution for the discrete fracture problem... -------- 61st commit message: iterative linear solver: remove the code for the non-default convergence criteria -------- 62nd commit message: provide the FE cache instead of the local FE this fixes a segfault in the stokes model caused by the fact that the local FE was not initialized at this point. -------- 63rd commit message: (Navier-)Stokes: fix bug due to the transition to unit normals now, all tests pass for this branch. The only things which need to be fixed are some annoying performance regressions compared to master and some bug in the splices feature of the property system... -------- 64th commit message: some fix for the local residual of the immiscible model -------- 65th commit message: Navier-Stokes: implement SCV center gradients There seems to be a bug in the previous implementation (the jacobian inverse transposed is evaluated using the local, not the global geometry), so the reference solution for the stokes2c test problem has also been updated... -------- 66th commit message: remove the ALUGrid specialization of the LensGridCreator and the YaspGrid one for the fingerproblem using different grid seems to sometimes cause a different vertex order, which in turn causes the respective test to fail if the reference solution was computed using the other grid... -------- 67th commit message: VCFV: use the correct BorderListCreator this makes MPI parallel computations work again. apart from performance regressions, this branch does not exhibit any known regressions compared to master anymore... -------- 68th commit message: make verything compile with the element centered finite volume discretization except the Navier-Stokes and the two-phase DFM models, of course... -------- 69th commit message: minor fixes - make the navier-stokes model slighly more generic by using the proper (in,ex)teriorIndex() methods on sub-control volumes - make the signature of the calculateValue() template method of the common two-point gradient approximator match the one of the vertex centered finite volume one -------- 70th commit message: fix fallout from the Big Rebase -------- 71st commit message: ECFV: some bugs in the boundary -------- 72nd commit message: make computeFlux() compute area-specific quantities -------- 73rd commit message: fix more bugs in the element centered FV discretization now eWoms should match Dumux pretty closely... -------- 74th commit message: coalesce the common code of the multi phase porous medium models into "MultiPhaseBaseModel" -------- 75th commit message: update reference solutions these were changed because of the screw-up with the area of boundary segments... -------- 76th commit message: rename "ImplicitBase" to "FvBase" because in eWoms, everything is implicit and these are currently the base classes for all finite volume discretizations. -------- 77th commit message: make the spatial discretization selectable using a splice This requires an opm-core with a the patches from https://github.com/OPM/opm-core/pull/446 merged... -------- 78th commit message: rename the properties used for splices to *Splice -------- 79th commit message: move the files in 'tests/models' to 'tests' since 'tests' was empty except for the 'models' subdirectory... -------- 80th commit message: improve and fix the tutorial -------- 81st commit message: remove the -fno-strict-aliasing flag from the provided option files seems like recent versions of Dune have been adapted... -------- 82nd commit message: also compile all CO2 injection simulations using the element centered finite volume discretization -------- 83rd commit message: PVS model: make it work properly with the element-centered finite volume discretiation because DOF != number of vertices
2013-12-12 05:52:44 -06:00
* \copydoc FvBaseMultiPhaseProblem::intrinsicPermeability
*/
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
{ return K_; }
/*!
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>
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
Scalar porosity(const Context& /*context*/, unsigned /*spaceIdx*/, unsigned /*timeIdx*/) const
{ return 0.4; }
/*!
* \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
MaterialLawParams& materialLawParams(const Context& context,
unsigned spaceIdx, unsigned timeIdx)
{
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 auto& entity = context.stencil(timeIdx).entity(spaceIdx);
assert(materialParams_[entity]);
return *materialParams_[entity];
}
/*!
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
{
const auto& entity = context.stencil(timeIdx).entity( spaceIdx );
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
assert(materialParams_[entity]);
return *materialParams_[entity];
}
//! \}
/*!
* \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) || onLowerBoundary_(pos))
values.setNoFlow();
else {
assert(onUpperBoundary_(pos));
values.setFreeFlow(context, spaceIdx, timeIdx, initialFluidState_);
}
// override the value for the liquid phase by forced
// imbibition of water on inlet boundary segments
if (onInlet_(pos)) {
values[contiWettingEqIdx] = -0.001; // [kg/(m^2 s)]
}
}
//! \}
/*!
* \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
{
// assign the primary variables
values.assignNaive(initialFluidState_);
}
/*!
* \copydoc FvBaseProblem::constraints
*/
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 constraints(Constraints& constraints, 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 (onUpperBoundary_(pos) && !onInlet_(pos)) {
constraints.setActive(true);
constraints.assignNaive(initialFluidState_);
}
else if (onLowerBoundary_(pos)) {
constraints.setActive(true);
constraints.assignNaive(initialFluidState_);
}
}
/*!
* \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, 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 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;
if (!onUpperBoundary_(pos))
return false;
Scalar xInject[] = { 0.25, 0.75 };
Scalar injectLen[] = { 0.1, 0.1 };
for (unsigned i = 0; i < sizeof(xInject) / sizeof(Scalar); ++i) {
if (xInject[i] - injectLen[i] / 2 < lambda
&& lambda < xInject[i] + injectLen[i] / 2)
return true;
}
return false;
}
void setupInitialFluidState_()
{
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
auto& fs = initialFluidState_;
fs.setPressure(wettingPhaseIdx, /*pressure=*/1e5);
Scalar Sw = EWOMS_GET_PARAM(TypeTag, Scalar, InitialWaterSaturation);
fs.setSaturation(wettingPhaseIdx, Sw);
fs.setSaturation(nonWettingPhaseIdx, 1 - Sw);
fs.setTemperature(temperature_);
// set the absolute pressures
Scalar pn = 1e5;
fs.setPressure(nonWettingPhaseIdx, pn);
fs.setPressure(wettingPhaseIdx, pn);
typename FluidSystem::template ParameterCache<Scalar> paramCache;
paramCache.updateAll(fs);
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
fs.setDensity(phaseIdx, FluidSystem::density(fs, paramCache, phaseIdx));
fs.setViscosity(phaseIdx, FluidSystem::viscosity(fs, paramCache, phaseIdx));
}
}
DimMatrix K_;
typename MaterialLawParams::VanGenuchtenParams micParams_;
typename MaterialLawParams::VanGenuchtenParams mdcParams_;
MaterialLawParamsContainer materialParams_;
Opm::ImmiscibleFluidState<Scalar, FluidSystem> initialFluidState_;
Scalar temperature_;
Scalar eps_;
};
} // namespace Opm
#endif