20#ifndef OPM_ROCK_IMPL_HEADER_INCLUDED 
   21#define OPM_ROCK_IMPL_HEADER_INCLUDED 
   26#include <opm/core/io/eclipse/EclipseGridInspector.hpp> 
   45                         const std::vector<int>& global_cell,
 
   46                         const double perm_threshold)
 
   49        static_assert(dim == 3, 
"");
 
   51        permfield_valid_.assign(global_cell.size(), 
false);
 
   53        assignPorosity    (deck, global_cell);
 
   54        assignPermeability(deck, global_cell, perm_threshold);
 
   60                         const double uniform_poro,
 
   61                         const double uniform_perm)
 
   63        permfield_valid_.assign(num_cells, 
true);
 
   64        porosity_.assign(num_cells, uniform_poro);
 
   65        permeability_.assign(dim*dim*num_cells, 0.0);
 
   66        for (
int i = 0; i < num_cells; ++i) {
 
   68            for (
int dd = 0; dd < dim; ++dd) {
 
   69                K(dd, dd) = uniform_perm;
 
   78        return porosity_[cell_index];
 
   86        assert (permfield_valid_[cell_index]);
 
   88        const PermTensor K(dim, dim, &permeability_[dim*dim*cell_index]);
 
  101        permfield_valid_[cell_index] = std::vector<unsigned char>::value_type(1);
 
  116                                   const std::vector<int>& global_cell)
 
  118        porosity_.assign(global_cell.size(), 1.0);
 
  120        if (deck->hasKeyword(
"PORO")) {
 
  121            const std::vector<double>& poro = deck->getKeyword(
"PORO").getSIDoubleData();
 
  123            for (
int c = 0; c < int(porosity_.size()); ++c) {
 
  124                porosity_[c] = poro[global_cell[c]];
 
  133                                       const std::vector<int>& global_cell,
 
  134                                       double perm_threshold)
 
  136        Opm::EclipseGridInspector insp(deck);
 
  137        std::array<int, 3> dims = insp.gridSize();
 
  138        int num_global_cells = dims[0]*dims[1]*dims[2];
 
  139        assert (num_global_cells > 0);
 
  141        permeability_.assign(dim * dim * global_cell.size(), 0.0);
 
  143        std::vector<const std::vector<double>*> tensor;
 
  146        const std::vector<double> 
zero(num_global_cells, 0.0);
 
  147        tensor.push_back(&
zero);
 
  149        static_assert(dim == 3, 
"");
 
  150        std::array<int,9> kmap;
 
  151        permeability_kind_ = fillTensor(deck, tensor, kmap);
 
  160        if (tensor.size() > 1) {
 
  161            const int nc  = global_cell.size();
 
  164            for (
int c = 0; c < nc; ++c, off += dim*dim) {
 
  167                const int glob = global_cell[c];
 
  169                for (
int i = 0; i < dim; ++i) {
 
  170                    for (
int j = 0; j < dim; ++j, ++kix) {
 
  171                        K(i,j) = (*tensor[kmap[kix]])[glob];
 
  173                    K(i,i) = std::max(K(i,i), perm_threshold);
 
  176                permfield_valid_[c] = std::vector<unsigned char>::value_type(1);
 
PermTensor permeability(int cell_index) const
Read-access to permeability.
Definition: Rock_impl.hpp:84
 
SharedPermTensor permeabilityModifiable(int cell_index)
Read- and write-access to permeability. Use with caution.
Definition: Rock_impl.hpp:95
 
SharedCMatrix SharedPermTensor
Tensor type for read and write access to permeability.
Definition: Rock.hpp:43
 
ImmutableCMatrix PermTensor
Tensor type for read-only access to permeability.
Definition: Rock.hpp:39
 
Rock()
Default constructor.
Definition: Rock_impl.hpp:37
 
void init(Opm::DeckConstPtr deck, const std::vector< int > &global_cell, const double perm_threshold=0.0)
Initialize from a grdecl file.
Definition: Rock_impl.hpp:44
 
void assignPermeability(Opm::DeckConstPtr deck, const std::vector< int > &global_cell, const double perm_threshold)
Definition: Rock_impl.hpp:132
 
void assignPorosity(Opm::DeckConstPtr deck, const std::vector< int > &global_cell)
Definition: Rock_impl.hpp:115
 
double porosity(int cell_index) const
Read-access to porosity.
Definition: Rock_impl.hpp:76
 
Definition: BlackoilFluid.hpp:32
 
@ Invalid
Definition: ReservoirPropertyCommon.hpp:50
 
void zero(Matrix &A)
Zero-fill a.
Definition: Matrix.hpp:603