Files
opm-common/opm/parser/eclipse/Utility/WelopenWrapper.hpp
Andreas Lauser 8d2ea3e0cd add wrapper classes for more eclipse keywords
these are convenient to convert everything required to rip out the old
parser of sim_fibo_ad. the concrete keywords are:

 - COMPDAT
 - EQUIL
 - GCONINJE
 - GCONPROD
 - GRUPTREE
 - WCONINJE
 - WCONINJ
 - WCONPROD
 - WELOPEN
 - WELSPECS
 - WGRUPCON

in the medium term it would be nice if these wrapper classes could be
automatically generated from the JSON keyword descriptions.
2014-01-15 15:44:38 +01:00

105 lines
3.4 KiB
C++

/*
Copyright (C) 2014 by 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_PARSER_WELOPEN_WRAPPER_HPP
#define OPM_PARSER_WELOPEN_WRAPPER_HPP
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <vector>
#include <algorithm>
namespace Opm {
class WelopenWrapper {
public:
/*!
* \brief A wrapper class to provide convenient access to the
* data of the 'WELOPEN' keyword.
*/
WelopenWrapper(Opm::DeckKeywordConstPtr keyword)
: m_keyword(keyword)
{
}
/*!
* \brief Return the number of injection wells
*/
int numWells() const
{ return m_keyword->size(); }
/*!
* \brief Return the human-readable name of the well with a given index
*/
std::string wellName(int wellIdx) const
{ return m_keyword->getRecord(wellIdx)->getItem(0)->getString(0); }
/*!
* \brief Return whether a well is open or closed
*
* This is one of:
* - OPEN: Well injects
* - STOP: Well does not reach the reservoir, but it
* injects. (and some of this fluid reaches the reservoir
* via crossflow)
* - SHUT: Well does not influence the reservoir
* - AUTO: Simulation selects one of the above depending in the
* well parameters and reservoir conditions at the well.
*/
std::string wellStatus(int wellIdx) const
{ return m_keyword->getRecord(wellIdx)->getItem(1)->getString(0); }
/*!
* \brief Return the I-coordinate of the connection grid block
*/
int coordinateI(int wellIdx) const
{ return m_keyword->getRecord(wellIdx)->getItem(2)->getInt(0); }
/*!
* \brief Return the J-coordinate of the connection grid block
*/
int coordinateJ(int wellIdx) const
{ return m_keyword->getRecord(wellIdx)->getItem(3)->getInt(0); }
/*!
* \brief Return the K-coordinate of the connection grid block
*/
int coordinateK(int wellIdx) const
{ return m_keyword->getRecord(wellIdx)->getItem(4)->getInt(0); }
/*!
* \brief Return the index of the first well completion for
* which this data applies
*/
int firstCompletionNumber(int wellIdx) const
{ return m_keyword->getRecord(wellIdx)->getItem(5)->getInt(0); }
/*!
* \brief Return the index of the last well completion for
* which this data applies
*/
int lastCompletionNumber(int wellIdx) const
{ return m_keyword->getRecord(wellIdx)->getItem(6)->getInt(0); }
private:
Opm::DeckKeywordConstPtr m_keyword;
};
}
#endif // OPM_PARSER_WELOPEN_WRAPPER_HPP