mirror of
https://github.com/OPM/opm-simulators.git
synced 2024-11-21 08:54:08 -06:00
isNumericalAquiferCell: put in separate struct for easier reuse
This commit is contained in:
parent
e89c43e54d
commit
4d35ec26de
@ -229,6 +229,7 @@ endif()
|
|||||||
# find tests -name '*.cpp' -a ! -wholename '*/not-unit/*' -printf '\t%p\n' | sort
|
# find tests -name '*.cpp' -a ! -wholename '*/not-unit/*' -printf '\t%p\n' | sort
|
||||||
list (APPEND TEST_SOURCE_FILES
|
list (APPEND TEST_SOURCE_FILES
|
||||||
tests/test_ALQState.cpp
|
tests/test_ALQState.cpp
|
||||||
|
tests/test_aquifergridutils.cpp
|
||||||
tests/test_blackoil_amg.cpp
|
tests/test_blackoil_amg.cpp
|
||||||
tests/test_convergenceoutputconfiguration.cpp
|
tests/test_convergenceoutputconfiguration.cpp
|
||||||
tests/test_convergencereport.cpp
|
tests/test_convergencereport.cpp
|
||||||
@ -379,8 +380,9 @@ list (APPEND PUBLIC_HEADER_FILES
|
|||||||
opm/simulators/wells/WellContainer.hpp
|
opm/simulators/wells/WellContainer.hpp
|
||||||
opm/simulators/aquifers/AquiferAnalytical.hpp
|
opm/simulators/aquifers/AquiferAnalytical.hpp
|
||||||
opm/simulators/aquifers/AquiferCarterTracy.hpp
|
opm/simulators/aquifers/AquiferCarterTracy.hpp
|
||||||
opm/simulators/aquifers/AquiferFetkovich.hpp
|
|
||||||
opm/simulators/aquifers/AquiferConstantFlux.hpp
|
opm/simulators/aquifers/AquiferConstantFlux.hpp
|
||||||
|
opm/simulators/aquifers/AquiferFetkovich.hpp
|
||||||
|
opm/simulators/aquifers/AquiferGridUtils.hpp
|
||||||
opm/simulators/aquifers/AquiferInterface.hpp
|
opm/simulators/aquifers/AquiferInterface.hpp
|
||||||
opm/simulators/aquifers/AquiferNumerical.hpp
|
opm/simulators/aquifers/AquiferNumerical.hpp
|
||||||
opm/simulators/aquifers/BlackoilAquiferModel.hpp
|
opm/simulators/aquifers/BlackoilAquiferModel.hpp
|
||||||
|
68
opm/simulators/aquifers/AquiferGridUtils.hpp
Normal file
68
opm/simulators/aquifers/AquiferGridUtils.hpp
Normal file
@ -0,0 +1,68 @@
|
|||||||
|
/*
|
||||||
|
Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
|
||||||
|
Copyright 2014, 2015 Dr. Blatt - HPC-Simulation-Software & Services
|
||||||
|
Copyright 2014, 2015 Statoil ASA.
|
||||||
|
Copyright 2015 NTNU
|
||||||
|
Copyright 2015, 2016, 2017 IRIS AS
|
||||||
|
|
||||||
|
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 3 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/>.
|
||||||
|
*/
|
||||||
|
|
||||||
|
#ifndef AQUIFER_GRID_UTILS_HEADER_INCLUDED
|
||||||
|
#define AQUIFER_GRID_UTILS_HEADER_INCLUDED
|
||||||
|
|
||||||
|
#include <opm/grid/CpGrid.hpp>
|
||||||
|
|
||||||
|
#include <algorithm>
|
||||||
|
#include <vector>
|
||||||
|
|
||||||
|
namespace Opm {
|
||||||
|
|
||||||
|
template<class Grid>
|
||||||
|
struct IsNumericalAquiferCell {
|
||||||
|
IsNumericalAquiferCell(const Grid&)
|
||||||
|
{}
|
||||||
|
|
||||||
|
template<class T>
|
||||||
|
bool operator()(const T&) const { return false; }
|
||||||
|
};
|
||||||
|
|
||||||
|
template<>
|
||||||
|
struct IsNumericalAquiferCell<Dune::CpGrid> {
|
||||||
|
IsNumericalAquiferCell(const Dune::CpGrid& grid)
|
||||||
|
: grid_(grid)
|
||||||
|
{}
|
||||||
|
|
||||||
|
template<class T>
|
||||||
|
bool operator()(const T& elem) const
|
||||||
|
{
|
||||||
|
const auto& aquiferCells = grid_.sortedNumAquiferCells();
|
||||||
|
if (aquiferCells.empty())
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
auto candidate = std::lower_bound(aquiferCells.begin(),
|
||||||
|
aquiferCells.end(), elem.index());
|
||||||
|
return candidate != aquiferCells.end() && *candidate == elem.index();
|
||||||
|
}
|
||||||
|
|
||||||
|
private:
|
||||||
|
const Dune::CpGrid& grid_;
|
||||||
|
};
|
||||||
|
|
||||||
|
} // namespace Opm
|
||||||
|
|
||||||
|
#endif
|
@ -40,6 +40,7 @@
|
|||||||
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
|
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
|
||||||
#include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
|
#include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
|
||||||
|
|
||||||
|
#include <opm/simulators/aquifers/AquiferGridUtils.hpp>
|
||||||
#include <opm/simulators/aquifers/BlackoilAquiferModel.hpp>
|
#include <opm/simulators/aquifers/BlackoilAquiferModel.hpp>
|
||||||
#include <opm/simulators/flow/countGlobalCells.hpp>
|
#include <opm/simulators/flow/countGlobalCells.hpp>
|
||||||
#include <opm/simulators/flow/partitionCells.hpp>
|
#include <opm/simulators/flow/partitionCells.hpp>
|
||||||
@ -1023,6 +1024,7 @@ namespace Opm {
|
|||||||
|
|
||||||
ElementContext elemCtx(ebosSimulator_);
|
ElementContext elemCtx(ebosSimulator_);
|
||||||
const auto& gridView = ebosSimulator().gridView();
|
const auto& gridView = ebosSimulator().gridView();
|
||||||
|
IsNumericalAquiferCell isNumericalAquiferCell(gridView.grid());
|
||||||
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
||||||
for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
|
for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
|
||||||
elemCtx.updatePrimaryStencil(elem);
|
elemCtx.updatePrimaryStencil(elem);
|
||||||
@ -1036,7 +1038,7 @@ namespace Opm {
|
|||||||
ebosModel.dofTotalVolume(cell_idx);
|
ebosModel.dofTotalVolume(cell_idx);
|
||||||
pvSumLocal += pvValue;
|
pvSumLocal += pvValue;
|
||||||
|
|
||||||
if (isNumericalAquiferCell(gridView.grid(), elem))
|
if (isNumericalAquiferCell(elem))
|
||||||
{
|
{
|
||||||
numAquiferPvSumLocal += pvValue;
|
numAquiferPvSumLocal += pvValue;
|
||||||
}
|
}
|
||||||
@ -1075,6 +1077,7 @@ namespace Opm {
|
|||||||
ElementContext elemCtx(ebosSimulator_);
|
ElementContext elemCtx(ebosSimulator_);
|
||||||
const auto& gridView = domain.view;
|
const auto& gridView = domain.view;
|
||||||
const auto& elemEndIt = gridView.template end</*codim=*/0>();
|
const auto& elemEndIt = gridView.template end</*codim=*/0>();
|
||||||
|
IsNumericalAquiferCell isNumericalAquiferCell(gridView.grid());
|
||||||
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
||||||
for (auto elemIt = gridView.template begin</*codim=*/0>();
|
for (auto elemIt = gridView.template begin</*codim=*/0>();
|
||||||
elemIt != elemEndIt;
|
elemIt != elemEndIt;
|
||||||
@ -1095,7 +1098,7 @@ namespace Opm {
|
|||||||
ebosModel.dofTotalVolume(cell_idx);
|
ebosModel.dofTotalVolume(cell_idx);
|
||||||
pvSumLocal += pvValue;
|
pvSumLocal += pvValue;
|
||||||
|
|
||||||
if (isNumericalAquiferCell(gridView.grid(), elem))
|
if (isNumericalAquiferCell(elem))
|
||||||
{
|
{
|
||||||
numAquiferPvSumLocal += pvValue;
|
numAquiferPvSumLocal += pvValue;
|
||||||
}
|
}
|
||||||
@ -1157,13 +1160,14 @@ namespace Opm {
|
|||||||
const auto& ebosResid = ebosSimulator_.model().linearizer().residual();
|
const auto& ebosResid = ebosSimulator_.model().linearizer().residual();
|
||||||
const auto& gridView = ebosSimulator().gridView();
|
const auto& gridView = ebosSimulator().gridView();
|
||||||
ElementContext elemCtx(ebosSimulator_);
|
ElementContext elemCtx(ebosSimulator_);
|
||||||
|
IsNumericalAquiferCell isNumericalAquiferCell(gridView.grid());
|
||||||
|
|
||||||
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
||||||
|
|
||||||
for (const auto& elem : elements(gridView, Dune::Partitions::interiorBorder))
|
for (const auto& elem : elements(gridView, Dune::Partitions::interiorBorder))
|
||||||
{
|
{
|
||||||
// Skip cells of numerical Aquifer
|
// Skip cells of numerical Aquifer
|
||||||
if (isNumericalAquiferCell(gridView.grid(), elem))
|
if (isNumericalAquiferCell(elem))
|
||||||
{
|
{
|
||||||
continue;
|
continue;
|
||||||
}
|
}
|
||||||
@ -1557,26 +1561,6 @@ namespace Opm {
|
|||||||
}
|
}
|
||||||
|
|
||||||
private:
|
private:
|
||||||
template<class T>
|
|
||||||
bool isNumericalAquiferCell(const Dune::CpGrid& grid, const T& elem)
|
|
||||||
{
|
|
||||||
const auto& aquiferCells = grid.sortedNumAquiferCells();
|
|
||||||
if (aquiferCells.empty())
|
|
||||||
{
|
|
||||||
return false;
|
|
||||||
}
|
|
||||||
auto candidate = std::lower_bound(aquiferCells.begin(), aquiferCells.end(),
|
|
||||||
elem.index());
|
|
||||||
return candidate != aquiferCells.end() && *candidate == elem.index();
|
|
||||||
}
|
|
||||||
|
|
||||||
template<class G, class T>
|
|
||||||
typename std::enable_if<!std::is_same<G,Dune::CpGrid>::value, bool>::type
|
|
||||||
isNumericalAquiferCell(const G&, const T&)
|
|
||||||
{
|
|
||||||
return false;
|
|
||||||
}
|
|
||||||
|
|
||||||
template<class FluidState, class Residual>
|
template<class FluidState, class Residual>
|
||||||
void getMaxCoeff(const unsigned cell_idx,
|
void getMaxCoeff(const unsigned cell_idx,
|
||||||
const IntensiveQuantities& intQuants,
|
const IntensiveQuantities& intQuants,
|
||||||
|
141
tests/test_aquifergridutils.cpp
Normal file
141
tests/test_aquifergridutils.cpp
Normal file
@ -0,0 +1,141 @@
|
|||||||
|
/*
|
||||||
|
Copyright 2018 SINTEF Digital, Mathematics and Cybernetics.
|
||||||
|
Copyright 2018 Equinor.
|
||||||
|
|
||||||
|
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 3 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/>.
|
||||||
|
*/
|
||||||
|
|
||||||
|
#include <config.h>
|
||||||
|
#define BOOST_TEST_MODULE AquiferGridUtils
|
||||||
|
#define BOOST_TEST_NO_MAIN
|
||||||
|
#include <boost/test/unit_test.hpp>
|
||||||
|
|
||||||
|
#include <dune/grid/common/rangegenerators.hh>
|
||||||
|
|
||||||
|
#include <opm/grid/polyhedralgrid.hh>
|
||||||
|
|
||||||
|
#include <opm/input/eclipse/Deck/Deck.hpp>
|
||||||
|
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
|
||||||
|
#include <opm/input/eclipse/Parser/Parser.hpp>
|
||||||
|
|
||||||
|
#include <opm/simulators/aquifers/AquiferGridUtils.hpp>
|
||||||
|
|
||||||
|
namespace {
|
||||||
|
|
||||||
|
Opm::Deck createNumericalAquiferDeck()
|
||||||
|
{
|
||||||
|
const char *deckData = R"(
|
||||||
|
DIMENS
|
||||||
|
8 15 3 /
|
||||||
|
AQUDIMS
|
||||||
|
3 2 1* 1* 1* 50 1* 1* /
|
||||||
|
GRID
|
||||||
|
|
||||||
|
DX
|
||||||
|
360*10./
|
||||||
|
DY
|
||||||
|
360*10./
|
||||||
|
DZ
|
||||||
|
360*1./
|
||||||
|
TOPS
|
||||||
|
360*100./
|
||||||
|
|
||||||
|
PORO
|
||||||
|
0. 0.25 0. 357*0.25/
|
||||||
|
PERMX
|
||||||
|
360*1000./
|
||||||
|
PERMY
|
||||||
|
360*1000./
|
||||||
|
PERMZ
|
||||||
|
360*10./
|
||||||
|
|
||||||
|
BOX
|
||||||
|
1 8 15 15 3 3 /
|
||||||
|
|
||||||
|
MULTY
|
||||||
|
1e9 1e-9 1.0 2.0 3.0 4.0 5.0 6.0/
|
||||||
|
|
||||||
|
ENDBOX
|
||||||
|
|
||||||
|
-- setting the three cells for numerical aquifer to be inactive
|
||||||
|
ACTNUM
|
||||||
|
0 1 0 0 356*1 /
|
||||||
|
|
||||||
|
AQUNUM
|
||||||
|
--AQnr. I J K Area Length PHI K Depth Initial.Pr PVTNUM SATNUM
|
||||||
|
1 1 1 1 1000000.0 10000 0.25 400 2585.00 285.00 2 2 /
|
||||||
|
1 3 1 1 1500000.0 20000 0.24 600 2585.00 285.00 3 * /
|
||||||
|
1 4 1 1 2000000.0 30000 * 700 2585.00 285.00 * 3 /
|
||||||
|
/
|
||||||
|
AQUCON
|
||||||
|
-- Connect numerical aquifer to the reservoir
|
||||||
|
-- Id.nr I1 I2 J1 J2 K1 K2 Face Trans.mult. Trans.opt.
|
||||||
|
1 1 8 15 15 3 3 'J+' 1.00 1 /
|
||||||
|
/
|
||||||
|
)";
|
||||||
|
|
||||||
|
Opm::Parser parser;
|
||||||
|
return parser.parseString(deckData);
|
||||||
|
}
|
||||||
|
|
||||||
|
} // Anonymous namespace
|
||||||
|
|
||||||
|
struct Fixture {
|
||||||
|
Fixture()
|
||||||
|
: numaquifer_deck(createNumericalAquiferDeck())
|
||||||
|
, ecl_state(numaquifer_deck)
|
||||||
|
, ecl_grid(ecl_state.getInputGrid())
|
||||||
|
{
|
||||||
|
}
|
||||||
|
|
||||||
|
Opm::Deck numaquifer_deck;
|
||||||
|
Opm::EclipseState ecl_state;
|
||||||
|
const Opm::EclipseGrid& ecl_grid;
|
||||||
|
};
|
||||||
|
|
||||||
|
BOOST_FIXTURE_TEST_CASE(NumericalAquiferCellUnsupported, Fixture)
|
||||||
|
{
|
||||||
|
Dune::PolyhedralGrid<3,3,double> grid(ecl_grid);
|
||||||
|
|
||||||
|
Opm::IsNumericalAquiferCell isNumericalAquiferCell(grid);
|
||||||
|
|
||||||
|
BOOST_CHECK_EQUAL(isNumericalAquiferCell(1), false);
|
||||||
|
}
|
||||||
|
|
||||||
|
BOOST_FIXTURE_TEST_CASE(NumericalAquiferCpGrid, Fixture)
|
||||||
|
{
|
||||||
|
Dune::CpGrid grid;
|
||||||
|
grid.processEclipseFormat(&ecl_grid, &ecl_state, false, false, false);
|
||||||
|
|
||||||
|
Opm::IsNumericalAquiferCell isNumericalAquiferCell(grid);
|
||||||
|
|
||||||
|
for (const auto& elem : elements(grid.leafGridView())) {
|
||||||
|
BOOST_CHECK_EQUAL(isNumericalAquiferCell(elem), elem.index() == 0 ||
|
||||||
|
elem.index() == 2 ||
|
||||||
|
elem.index() == 3);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
bool init_unit_test_func()
|
||||||
|
{
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
int main(int argc, char** argv)
|
||||||
|
{
|
||||||
|
Dune::MPIHelper::instance(argc, argv);
|
||||||
|
return boost::unit_test::unit_test_main(&init_unit_test_func, argc, argv);
|
||||||
|
}
|
Loading…
Reference in New Issue
Block a user