Merge pull request #909 from atgeirr/silence-warning-pinchprocessor

Remove unused arguments in PinchProcessor methods.
This commit is contained in:
Atgeirr Flø Rasmussen 2015-10-21 11:53:44 +02:00
commit d94e0fc0f7
2 changed files with 19 additions and 29 deletions

View File

@ -53,10 +53,9 @@ namespace Opm
/// \param[in] htrans half cell transmissibility, size is number of cellfaces.
/// \param[in] multz Z+ transmissibility multiplier for all active cells
/// \param[in] pv pore volume for all the cartesian cells
/// \param[in] dz dz for all the cartesian cells
/// \param[in] nnc non-neighbor connection class
/// Algorithm:
/// 1. Mark all the cells which pv and dz less than minpvValue and thickness.
/// 1. Mark all the cells which pv less than minpvValue.
/// 2. Find out proper pinchouts column and associate top and bottom cells.
/// 3. Compute transmissibility for nncs.
/// 4. Apply multz due to different multz options.
@ -65,7 +64,6 @@ namespace Opm
const std::vector<int>& actnum,
const std::vector<double>& multz,
const std::vector<double>& pv,
const std::vector<double>& dz,
NNC& nnc);
private:
@ -75,10 +73,8 @@ namespace Opm
PinchMode::ModeEnum multzMode_;
/// Mark minpved cells.
std::vector<int> getMinpvCells_(const Grid& grid,
const std::vector<int>& actnum,
const std::vector<double>& pv,
const std::vector<double>& dz);
std::vector<int> getMinpvCells_(const std::vector<int>& actnum,
const std::vector<double>& pv);
/// Get the interface for two cells.
int interface_(const Grid& grid,
@ -94,8 +90,7 @@ namespace Opm
std::vector<std::vector<int> >
getPinchoutsColumn_(const Grid& grid,
const std::vector<int>& actnum,
const std::vector<double>& pv,
const std::vector<double>& dz);
const std::vector<double>& pv);
/// Get global cell index.
int getGlobalIndex_(const int i, const int j, const int k, const int* dims);
@ -108,8 +103,7 @@ namespace Opm
std::vector<double> transCompute_(const Grid& grid,
const std::vector<double>& htrans,
const std::vector<int>& pinCells,
const std::vector<int>& pinFaces,
const std::vector<double>& multz);
const std::vector<int>& pinFaces);
/// Get map between half-trans index and the pair of face index and cell index.
std::vector<int> getHfIdxMap_(const Grid& grid);
@ -124,7 +118,6 @@ namespace Opm
const std::vector<int>& actnum,
const std::vector<double>& multz,
const std::vector<double>& pv,
const std::vector<double>& dz,
NNC& nnc);
/// Item 5 in PINCH keyword.
@ -239,10 +232,8 @@ namespace Opm
template<class Grid>
inline std::vector<int> PinchProcessor<Grid>::getMinpvCells_(const Grid& grid,
const std::vector<int>& actnum,
const std::vector<double>& pv,
const std::vector<double>& dz)
inline std::vector<int> PinchProcessor<Grid>::getMinpvCells_(const std::vector<int>& actnum,
const std::vector<double>& pv)
{
std::vector<int> minpvCells(pv.size(), 0);
for (int idx = 0; idx < static_cast<int>(pv.size()); ++idx) {
@ -297,8 +288,7 @@ namespace Opm
inline std::vector<double> PinchProcessor<Grid>::transCompute_(const Grid& grid,
const std::vector<double>& htrans,
const std::vector<int>& pinCells,
const std::vector<int>& pinFaces,
const std::vector<double>& multz)
const std::vector<int>& pinFaces)
{
const int nc = Opm::UgGridHelpers::numCells(grid);
const int nf = Opm::UgGridHelpers::numFaces(grid);
@ -342,11 +332,10 @@ namespace Opm
template<class Grid>
inline std::vector<std::vector<int>> PinchProcessor<Grid>::getPinchoutsColumn_(const Grid& grid,
const std::vector<int>& actnum,
const std::vector<double>& pv,
const std::vector<double>& dz)
const std::vector<double>& pv)
{
const int* dims = Opm::UgGridHelpers::cartDims(grid);
std::vector<int> minpvCells = getMinpvCells_(grid, actnum, pv, dz);
std::vector<int> minpvCells = getMinpvCells_(actnum, pv);
std::vector<std::vector<int>> segment;
for (int z = 0; z < dims[2]; ++z) {
for (int y = 0; y < dims[1]; ++y) {
@ -382,14 +371,13 @@ namespace Opm
const std::vector<int>& actnum,
const std::vector<double>& multz,
const std::vector<double>& pv,
const std::vector<double>& dz,
NNC& nnc)
{
const int* dims = Opm::UgGridHelpers::cartDims(grid);
std::vector<int> pinFaces;
std::vector<int> pinCells;
std::vector<std::vector<int> > newSeg;
auto minpvSeg = getPinchoutsColumn_(grid, actnum, pv, dz);
auto minpvSeg = getPinchoutsColumn_(grid, actnum, pv);
for (auto& seg : minpvSeg) {
std::array<int, 3> ijk1 = getCartIndex_(seg.front(), dims);
std::array<int, 3> ijk2 = getCartIndex_(seg.back(), dims);
@ -436,7 +424,7 @@ namespace Opm
}
}
auto faceTrans = transCompute_(grid, htrans, pinCells, pinFaces, multz);
auto faceTrans = transCompute_(grid, htrans, pinCells, pinFaces);
auto multzmap = multzOptions_(grid, pinCells, pinFaces, multz, newSeg);
applyMultz_(faceTrans, multzmap);
for (int i = 0; i < static_cast<int>(pinCells.size())/2; ++i) {
@ -499,10 +487,9 @@ namespace Opm
const std::vector<int>& actnum,
const std::vector<double>& multz,
const std::vector<double>& pv,
const std::vector<double>& dz,
NNC& nnc)
{
transTopbot_(grid, htrans, actnum, multz, pv, dz, nnc);
transTopbot_(grid, htrans, actnum, multz, pv, nnc);
}
} // namespace Opm

View File

@ -26,7 +26,11 @@
#define NVERBOSE // Suppress own messages when throw()in
#define BOOST_TEST_MODULE PinchProcessorTest
#include <opm/common/utility/platform_dependent/disable_warnings.h>
#include <boost/test/unit_test.hpp>
#include <opm/common/utility/platform_dependent/reenable_warnings.h>
#include <vector>
#include <opm/core/grid.h>
#include <opm/core/grid/GridManager.hpp>
@ -76,7 +80,7 @@ BOOST_AUTO_TEST_CASE(Processing)
std::vector<int> actnum;
eclgrid->exportACTNUM(actnum);
if (actnum.empty() && (nc == eclgrid->getCartesianSize())) {
if (actnum.empty() && (nc == int(eclgrid->getCartesianSize()))) {
actnum.assign(nc, 1);
}
std::vector<double> htrans(Opm::UgGridHelpers::numCellFaces(grid));
@ -84,12 +88,11 @@ BOOST_AUTO_TEST_CASE(Processing)
tpfa_htrans_compute(ug, rock.permeability(), htrans.data());
auto transMult = eclstate->getTransMult();
std::vector<double> multz(nc, 0.0);
std::vector<double> dz(porv.size(), 0.0);
for (int i = 0; i < nc; ++i) {
multz[i] = transMult->getMultiplier(global_cell[i], Opm::FaceDir::ZPlus);
}
Opm::NNC nnc(deck, eclgrid);
pinch.process(grid, htrans, actnum, multz, porv, dz, nnc);
pinch.process(grid, htrans, actnum, multz, porv, nnc);
auto nnc1 = nnc.nnc1();
auto nnc2 = nnc.nnc2();
auto trans = nnc.trans();