move the code which creates a compressed to cartesian cell index map to a separate free function

This commit is contained in:
Andreas Lauser 2015-09-02 14:39:51 +02:00
parent 9515ec6cfb
commit 73710a01d2
4 changed files with 56 additions and 22 deletions

View File

@ -114,6 +114,7 @@ list (APPEND PUBLIC_HEADER_FILES
opm/autodiff/SolventPropsAdFromDeck.hpp
opm/autodiff/BlackoilPropsAdInterface.hpp
opm/autodiff/CPRPreconditioner.hpp
opm/autodiff/createGlobalCellArray.hpp
opm/autodiff/BlackoilSolventModel.hpp
opm/autodiff/BlackoilSolventModel_impl.hpp
opm/autodiff/BlackoilSolventState.hpp

View File

@ -51,6 +51,7 @@
#include <opm/core/grid/cornerpoint_grid.h>
#include <opm/core/grid/GridManager.hpp>
#include <opm/autodiff/GridHelpers.hpp>
#include <opm/autodiff/createGlobalCellArray.hpp>
#include <opm/core/wells.h>
#include <opm/core/wells/WellsManager.hpp>
@ -251,17 +252,8 @@ try
const PhaseUsage pu = Opm::phaseUsageFromDeck(deck);
Opm::BlackoilOutputWriter outputWriter(grid, param, eclipseState, pu );
int numCells = Opm::UgGridHelpers::numCells(grid);
const auto& globalCell = Opm::UgGridHelpers::globalCell(grid);
std::vector<int> compressedToCartesianIdx(numCells);
for (unsigned cellIdx = 0; cellIdx < numCells; ++cellIdx) {
if (globalCell) {
compressedToCartesianIdx[cellIdx] = globalCell[cellIdx];
}
else {
compressedToCartesianIdx[cellIdx] = cellIdx;
}
}
std::vector<int> compressedToCartesianIdx;
Opm::createGlobalCellArray(grid, compressedToCartesianIdx);
typedef BlackoilPropsAdFromDeck::MaterialLawManager MaterialLawManager;
auto materialLawManager = std::make_shared<MaterialLawManager>();

View File

@ -51,6 +51,7 @@
#include <opm/core/grid/cornerpoint_grid.h>
#include <opm/core/grid/GridManager.hpp>
#include <opm/autodiff/GridHelpers.hpp>
#include <opm/autodiff/createGlobalCellArray.hpp>
#include <opm/core/wells.h>
#include <opm/core/wells/WellsManager.hpp>
@ -247,17 +248,8 @@ try
const PhaseUsage pu = Opm::phaseUsageFromDeck(deck);
Opm::BlackoilOutputWriter outputWriter(grid, param, eclipseState, pu );
int numCells = Opm::UgGridHelpers::numCells(grid);
const auto& globalCell = Opm::UgGridHelpers::globalCell(grid);
std::vector<int> compressedToCartesianIdx(numCells);
for (unsigned cellIdx = 0; cellIdx < numCells; ++cellIdx) {
if (globalCell) {
compressedToCartesianIdx[cellIdx] = globalCell[cellIdx];
}
else {
compressedToCartesianIdx[cellIdx] = cellIdx;
}
}
std::vector<int> compressedToCartesianIdx;
Opm::createGlobalCellArray(grid, compressedToCartesianIdx);
typedef BlackoilPropsAdFromDeck::MaterialLawManager MaterialLawManager;
auto materialLawManager = std::make_shared<MaterialLawManager>();

View File

@ -0,0 +1,49 @@
/*
Copyright 2015 Andreas Lauser
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 OPM_CREATE_GLOBAL_CELL_ARRAY_HPP
#define OPM_CREATE_GLOBAL_CELL_ARRAY_HPP
#include <opm/autodiff/GridHelpers.hpp>
namespace Opm
{
/*!
* \brief Create a mapping from a global cell index of a grid to the logically
* Cartesian index of the ECL deck.
*/
template <class Grid>
void createGlobalCellArray(const Grid &grid, std::vector<int>& dest)
{
int numCells = Opm::UgGridHelpers::numCells(grid);
dest.resize(numCells);
const auto& globalCell = Opm::UgGridHelpers::globalCell(grid);
std::vector<int> compressedToCartesianIdx(numCells);
for (unsigned cellIdx = 0; cellIdx < numCells; ++cellIdx) {
if (globalCell) {
dest[cellIdx] = globalCell[cellIdx];
}
else {
dest[cellIdx] = cellIdx;
}
}
}
}
#endif