CornerPointChopper: add methods opm-parser'y methods

This commit is contained in:
Andreas Lauser
2014-04-11 17:28:49 +02:00
parent a313dd86d6
commit 6a90aa086a

View File

@@ -22,6 +22,14 @@
#include <opm/core/io/eclipse/EclipseGridParser.hpp>
#include <opm/core/utility/parameters/ParameterGroup.hpp>
#include <opm/parser/eclipse/Parser/Parser.hpp>
#include <opm/parser/eclipse/Units/UnitSystem.hpp>
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <opm/parser/eclipse/Deck/DeckKeyword.hpp>
#include <opm/parser/eclipse/Deck/DeckRecord.hpp>
#include <opm/parser/eclipse/Deck/DeckDoubleItem.hpp>
#include <opm/parser/eclipse/Deck/DeckIntItem.hpp>
#include <opm/parser/eclipse/Deck/DeckStringItem.hpp>
#include <iostream>
#include <fstream>
#include <string>
@@ -38,12 +46,18 @@ namespace Opm
CornerPointChopper(const std::string& file)
: parser_(file, false)
{
for (int dd = 0; dd < 3; ++dd) {
dims_[dd] = parser_.getSPECGRID().dimensions[dd];
}
Opm::ParserPtr parser(new Opm::Parser());
deck_ = parser->parseFile(file);
metricUnits_.reset(Opm::UnitSystem::newMETRIC());
Opm::DeckRecordConstPtr specgridRecord = deck_->getKeyword("SPECGRID")->getRecord(0);
dims_[0] = specgridRecord->getItem("NX")->getInt(0);
dims_[1] = specgridRecord->getItem("NY")->getInt(0);
dims_[2] = specgridRecord->getItem("NZ")->getInt(0);
int layersz = 8*dims_[0]*dims_[1];
const std::vector<double>& ZCORN = parser_.getFloatingPointValue("ZCORN");
const std::vector<double>& ZCORN = deck_->getKeyword("ZCORN")->getRawDoubleData();
botmax_ = *std::max_element(ZCORN.begin(), ZCORN.begin() + layersz/2);
topmin_ = *std::min_element(ZCORN.begin() + dims_[2]*layersz - layersz/2,
ZCORN.begin() + dims_[2]*layersz);
@@ -133,7 +147,7 @@ namespace Opm
new_dims_[1] = jmax - jmin;
// Filter the coord field
const std::vector<double>& COORD = parser_.getFloatingPointValue("COORD");
const std::vector<double>& COORD = deck_->getKeyword("COORD")->getRawDoubleData();
int num_coord = COORD.size();
if (num_coord != 6*(dims_[0] + 1)*(dims_[1] + 1)) {
std::cerr << "Error! COORD size (" << COORD.size() << ") not consistent with SPECGRID\n";
@@ -166,7 +180,7 @@ namespace Opm
// coordinate of the bottom surface, while zmax must be less than or
// equal to the lowest coordinate of the top surface.
int layersz = 8*dims_[0]*dims_[1];
const std::vector<double>& ZCORN = parser_.getFloatingPointValue("ZCORN");
const std::vector<double>& ZCORN = deck_->getKeyword("ZCORN")->getRawDoubleData();
int num_zcorn = ZCORN.size();
if (num_zcorn != layersz*dims_[2]) {
std::cerr << "Error! ZCORN size (" << ZCORN.size() << ") not consistent with SPECGRID\n";
@@ -276,8 +290,50 @@ namespace Opm
return sp;
}
/// Return a sub-deck with fields corresponding to the selected subset.
Opm::DeckConstPtr subDeck()
{
Opm::DeckPtr subDeck(new Opm::Deck);
Opm::DeckKeywordPtr specGridKw(new Opm::DeckKeyword("SPECGRID"));
Opm::DeckRecordPtr specGridRecord(new Opm::DeckRecord());
Opm::DeckIntItemPtr nxItem(new Opm::DeckIntItem("NX"));
Opm::DeckIntItemPtr nyItem(new Opm::DeckIntItem("NY"));
Opm::DeckIntItemPtr nzItem(new Opm::DeckIntItem("NZ"));
Opm::DeckIntItemPtr numresItem(new Opm::DeckIntItem("NUMRES"));
Opm::DeckStringItemPtr coordTypeItem(new Opm::DeckStringItem("COORD_TYPE"));
nxItem->push_back(new_dims_[0]);
nyItem->push_back(new_dims_[1]);
nzItem->push_back(new_dims_[2]);
numresItem->push_back(1);
coordTypeItem->push_back("F");
specGridRecord->addItem(nxItem);
specGridRecord->addItem(nyItem);
specGridRecord->addItem(nzItem);
specGridRecord->addItem(numresItem);
specGridRecord->addItem(coordTypeItem);
specGridKw->addRecord(specGridRecord);
subDeck->addKeyword(specGridKw);
addDoubleKeyword_(subDeck, "COORD", /*dimension=*/"Length", new_COORD_);
addDoubleKeyword_(subDeck, "ZCORN", /*dimension=*/"Length", new_ZCORN_);
addIntKeyword_(subDeck, "ACTNUM", new_ACTNUM_);
addDoubleKeyword_(subDeck, "PORO", /*dimension=*/"1", new_PORO_);
addDoubleKeyword_(subDeck, "NTG", /*dimension=*/"1", new_NTG_);
addDoubleKeyword_(subDeck, "SWCR", /*dimension=*/"1", new_SWCR_);
addDoubleKeyword_(subDeck, "SOWCR", /*dimension=*/"1", new_SOWCR_);
addDoubleKeyword_(subDeck, "PERMX", /*dimension=*/"Permeability", new_PERMX_);
addDoubleKeyword_(subDeck, "PERMY", /*dimension=*/"Permeability", new_PERMY_);
addDoubleKeyword_(subDeck, "PERMZ", /*dimension=*/"Permeability", new_PERMZ_);
addIntKeyword_(subDeck, "SATNUM", new_SATNUM_);
return subDeck;
}
void writeGrdecl(const std::string& filename)
{
@@ -307,7 +363,10 @@ namespace Opm
bool hasSOWCR() const {return !new_SOWCR_.empty(); }
private:
EclipseGridParser parser_;
Opm::EclipseGridParser parser_;
Opm::DeckConstPtr deck_;
std::shared_ptr<Opm::UnitSystem> metricUnits_;
double botmax_;
double topmin_;
double abszmin_;
@@ -327,6 +386,49 @@ namespace Opm
int new_dims_[3];
std::vector<int> new_to_old_cell_;
void addDoubleKeyword_(Opm::DeckPtr subDeck,
const std::string& keywordName,
const std::string& dimensionString,
const std::vector<double>& data)
{
if (data.empty())
return;
Opm::DeckKeywordPtr dataKw(new Opm::DeckKeyword(keywordName));
Opm::DeckRecordPtr dataRecord(new Opm::DeckRecord());
Opm::DeckDoubleItemPtr dataItem(new Opm::DeckDoubleItem("DATA"));
for (size_t i = 0; i < data.size(); ++i) {
dataItem->push_back(data[i]);
}
std::shared_ptr<const Dimension> dimension = metricUnits_->parse(dimensionString);
dataItem->push_backDimension(/*active=*/dimension, /*default=*/dimension);
dataRecord->addItem(dataItem);
dataKw->addRecord(dataRecord);
subDeck->addKeyword(dataKw);
}
void addIntKeyword_(Opm::DeckPtr subDeck,
const std::string& keywordName,
const std::vector<int>& data)
{
if (data.empty())
return;
Opm::DeckKeywordPtr dataKw(new Opm::DeckKeyword(keywordName));
Opm::DeckRecordPtr dataRecord(new Opm::DeckRecord());
Opm::DeckIntItemPtr dataItem(new Opm::DeckIntItem("DATA"));
for (size_t i = 0; i < data.size(); ++i) {
dataItem->push_back(data[i]);
}
dataRecord->addItem(dataItem);
dataKw->addRecord(dataRecord);
subDeck->addKeyword(dataKw);
}
template <typename T>
void outputField(std::ostream& os,
@@ -366,16 +468,16 @@ namespace Opm
void filterDoubleField(const std::string& keyword, std::vector<double>& output_field)
{
if (parser_.hasField(keyword)) {
const std::vector<double>& field = parser_.getFloatingPointValue(keyword);
if (deck_->hasKeyword(keyword)) {
const std::vector<double>& field = deck_->getKeyword(keyword)->getRawDoubleData();
filterField(field, output_field);
}
}
void filterIntegerField(const std::string& keyword, std::vector<int>& output_field)
{
if (parser_.hasField(keyword)) {
const std::vector<int>& field = parser_.getIntegerValue(keyword);
if (deck_->hasKeyword(keyword)) {
const std::vector<int>& field = deck_->getKeyword(keyword)->getIntData();
filterField(field, output_field);
}
}