Renamed getNNC() to getInputNNC() avoids confusion

* The NNC in the eclipse state is and will always be the nnc from deck
* updated tests
This commit is contained in:
Pål Grønås Drange
2016-05-13 14:09:29 +02:00
parent 64fe528f6d
commit 607e4e2b9a
5 changed files with 34 additions and 36 deletions

View File

@@ -75,7 +75,7 @@ namespace Opm {
m_schedule( std::make_shared<Schedule>( m_parseContext, m_inputGrid, deck, m_ioConfig ) ),
m_simulationConfig( std::make_shared<SimulationConfig>( m_parseContext, deck, m_eclipseProperties ) ),
m_summaryConfig( deck, *m_schedule, m_eclipseProperties, m_inputGrid->getNXYZ() ),
m_nnc( deck, m_inputGrid ),
m_inputNnc( deck, m_inputGrid ),
m_deckUnitSystem( deck.getActiveUnitSystem() )
{
@@ -168,12 +168,12 @@ namespace Opm {
return m_transMult;
}
const NNC& EclipseState::getNNC() const {
return m_nnc;
const NNC& EclipseState::getInputNNC() const {
return m_inputNnc;
}
bool EclipseState::hasNNC() const {
return m_nnc.hasNNC();
bool EclipseState::hasInputNNC() const {
return m_inputNnc.hasNNC();
}
std::string EclipseState::getTitle() const {

View File

@@ -90,8 +90,11 @@ namespace Opm {
const FaultCollection& getFaults() const;
std::shared_ptr<const TransMult> getTransMult() const;
const NNC& getNNC() const;
bool hasNNC() const;
/// non-neighboring connections
/// the non-standard adjacencies as specified in input deck
const NNC& getInputNNC() const;
bool hasInputNNC() const;
const Eclipse3DProperties& get3DProperties() const;
@@ -130,7 +133,7 @@ namespace Opm {
std::string m_title;
std::shared_ptr<TransMult> m_transMult;
FaultCollection m_faults;
NNC m_nnc;
NNC m_inputNnc;
MessageContainer m_messageContainer;
const UnitSystem& m_deckUnitSystem;

View File

@@ -29,6 +29,7 @@
namespace Opm
{
/// [[deprecated]]
NNC::NNC(std::shared_ptr<const Deck> deck, EclipseGridConstPtr eclipseGrid) :
NNC(*deck, eclipseGrid)
{}

View File

@@ -17,8 +17,11 @@
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef NNC_HPP
#define NNC_HPP
#ifndef OPM_PARSER_NNC_HPP
#define OPM_PARSER_NNC_HPP
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/EclipseGrid.hpp>
#include <cstddef>
#include <memory>
@@ -26,8 +29,6 @@
namespace Opm
{
class Deck;
class EclipseGrid;
struct NNCdata {
size_t cell1;
@@ -35,12 +36,17 @@ struct NNCdata {
double trans;
};
/// Represents non-neighboring connections (non-standard adjacencies).
/// This class is essentially a directed weighted graph.
class NNC
{
public:
/// Construct from input deck.
NNC();
/// [[deprecated]]
NNC(std::shared_ptr<const Deck> deck_ptr, std::shared_ptr< const EclipseGrid > eclipseGrid);
/// Construct from input deck.
NNC(const Deck& deck, std::shared_ptr< const EclipseGrid > eclipseGrid);
void addNNC(const size_t cell1, const size_t cell2, const double trans);
const std::vector<NNCdata>& nncdata() const { return m_nnc; }
@@ -56,4 +62,4 @@ private:
} // namespace Opm
#endif // NNC_HPP
#endif // OPM_PARSER_NNC_HPP

View File

@@ -40,25 +40,16 @@ using namespace Opm;
BOOST_AUTO_TEST_CASE(noNNC)
{
Opm::ParseContext parseContext;
const std::string filename = "testdata/integration_tests/NNC/noNNC.DATA";
Opm::ParserPtr parser(new Opm::Parser());
Opm::DeckConstPtr deck(parser->parseFile(filename, parseContext));
Opm::EclipseStateConstPtr eclipseState(new EclipseState(deck , parseContext));
auto eclGrid = eclipseState->getInputGrid();
Opm::NNC nnc(deck, eclGrid);
auto eclipseState = Parser::parse("testdata/integration_tests/NNC/noNNC.DATA");
const auto& nnc = eclipseState.getInputNNC();
BOOST_CHECK(!eclipseState.hasInputNNC());
BOOST_CHECK(!nnc.hasNNC());
}
BOOST_AUTO_TEST_CASE(readDeck)
{
Opm::ParseContext parseContext;
const std::string filename = "testdata/integration_tests/NNC/NNC.DATA";
Opm::ParserPtr parser(new Opm::Parser());
Opm::DeckConstPtr deck(parser->parseFile(filename, parseContext));
Opm::EclipseStateConstPtr eclipseState(new EclipseState(deck , parseContext));
auto eclGrid = eclipseState->getInputGrid();
Opm::NNC nnc(deck, eclGrid);
auto eclipseState = Parser::parse("testdata/integration_tests/NNC/NNC.DATA");
const auto& nnc = eclipseState.getInputNNC();
BOOST_CHECK(nnc.hasNNC());
const std::vector<NNCdata>& nncdata = nnc.nncdata();
@@ -75,17 +66,14 @@ BOOST_AUTO_TEST_CASE(readDeck)
BOOST_AUTO_TEST_CASE(addNNCfromDeck)
{
Opm::ParseContext parseContext;
const std::string filename = "testdata/integration_tests/NNC/NNC.DATA";
Opm::ParserPtr parser(new Opm::Parser());
Opm::DeckConstPtr deck(parser->parseFile(filename, parseContext));
Opm::EclipseStateConstPtr eclipseState(new EclipseState(deck , parseContext));
auto eclGrid = eclipseState->getInputGrid();
Opm::NNC nnc(deck, eclGrid);
auto eclipseState = Parser::parse("testdata/integration_tests/NNC/NNC.DATA");
auto nnc = eclipseState.getInputNNC();
BOOST_CHECK(nnc.hasNNC());
const std::vector<NNCdata>& nncdata = nnc.nncdata();
BOOST_CHECK_EQUAL(nnc.numNNC(), 4);
// test add NNC
nnc.addNNC(2,2,2.0);
nnc.addNNC(2, 2, 2.0);
BOOST_CHECK_EQUAL(nnc.numNNC(), 5);
BOOST_CHECK_EQUAL(nncdata[4].cell1, 2);
BOOST_CHECK_EQUAL(nncdata[4].cell2, 2);