20 #ifndef OPM_ROCKFROMDECK_HEADER_INCLUDED
21 #define OPM_ROCKFROMDECK_HEADER_INCLUDED
23 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
27 struct UnstructuredGrid;
36 friend class BlackoilPropsDataHandle;
50 void init(
const Opm::EclipseState& eclState,
51 int number_of_cells,
const int* global_cell,
52 const int* cart_dims);
63 return porosity_.size();
77 return &permeability_[0];
91 const int number_of_cells,
92 const int* global_cell,
94 const double perm_threshold,
98 void assignPorosity(
const Opm::EclipseState& eclState,
100 const int* global_cell);
102 std::vector<double> porosity_;
103 std::vector<double> permeability_;
104 std::vector<unsigned char> permfield_valid_;
112 #endif // OPM_ROCKFROMDECK_HEADER_INCLUDED
int numCells() const
Definition: RockFromDeck.hpp:61
void init(const Opm::EclipseState &eclState, int number_of_cells, const int *global_cell, const int *cart_dims)
Initialize from deck and cell mapping.
Definition: RockFromDeck.cpp:81
static void extractInterleavedPermeability(const Opm::EclipseState &eclState, const int number_of_cells, const int *global_cell, const int *cart_dims, const double perm_threshold, std::vector< double > &permeability)
Convert the permeabilites for the logically Cartesian grid in EclipseState to an array of size number...
Definition: RockFromDeck.cpp:110
const double * permeability() const
Definition: RockFromDeck.hpp:75
Definition: RockFromDeck.hpp:32
RockFromDeck()
Default constructor.
Definition: RockFromDeck.cpp:70
int numDimensions() const
Definition: RockFromDeck.hpp:55
const double * porosity() const
Definition: RockFromDeck.hpp:67