EGrid.hpp
Go to the documentation of this file.
1/*
2 Copyright 2019 Equinor ASA.
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
9
10 OPM is distributed in the hope that it will be useful,
11 but WITHOUT ANY WARRANTY; without even the implied warranty of
12 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 GNU General Public License for more details.
14
15 You should have received a copy of the GNU General Public License
16 along with OPM. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19#ifndef OPM_IO_EGRID_HPP
20#define OPM_IO_EGRID_HPP
21
24
25#include <array>
26#include <iostream>
27#include <string>
28#include <fstream>
29#include <vector>
30#include <ctime>
31#include <map>
32
33namespace Opm { namespace EclIO {
34
35class EGrid : public EclFile
36{
37public:
38 explicit EGrid(const std::string& filename, std::string grid_name = "global");
39
40 int global_index(int i, int j, int k) const;
41 int active_index(int i, int j, int k) const;
42
43 const std::array<int, 3>& dimension() const { return nijk; }
44
45 std::array<int, 3> ijk_from_active_index(int actInd) const;
46 std::array<int, 3> ijk_from_global_index(int globInd) const;
47
48 void getCellCorners(int globindex, std::array<double,8>& X, std::array<double,8>& Y, std::array<double,8>& Z);
49 void getCellCorners(const std::array<int, 3>& ijk, std::array<double,8>& X, std::array<double,8>& Y, std::array<double,8>& Z);
50
51 int activeCells() const { return nactive; }
52 int totalNumberOfCells() const { return nijk[0] * nijk[1] * nijk[2]; }
53
56 bool is_radial() const { return m_radial; }
57
58 const std::vector<int>& hostCellsGlobalIndex() const { return host_cells; }
59 std::vector<std::array<int, 3>> hostCellsIJK();
60
61 // zero based: i1, j1,k1, i2,j2,k2, transmisibility
62 using NNCentry = std::tuple<int, int, int, int, int, int, float>;
63 std::vector<NNCentry> get_nnc_ijk();
64
65 const std::vector<std::string>& list_of_lgrs() const { return lgr_names; }
66
67private:
68 Opm::filesystem::path inputFileName, initFileName;
69 std::string m_grid_name;
70 bool m_radial;
71
72 std::array<int, 3> nijk;
73 std::array<int, 3> host_nijk;
74
75 int nactive;
76 mutable bool m_nncs_loaded;
77
78 std::vector<int> act_index;
79 std::vector<int> glob_index;
80
81 std::vector<float> coord_array;
82 std::vector<float> zcorn_array;
83
84 std::vector<int> nnc1_array;
85 std::vector<int> nnc2_array;
86 std::vector<float> transnnc_array;
87 std::vector<int> host_cells;
88
89 std::vector<std::string> lgr_names;
90
91 int zcorn_array_index;
92 int coord_array_index;
93 int actnum_array_index;
94 int nnc1_array_index;
95 int nnc2_array_index;
96};
97
98}} // namespace Opm::EclIO
99
100#endif // OPM_IO_EGRID_HPP
const char *const string
Definition: cJSON.h:170
Definition: EGrid.hpp:36
void getCellCorners(int globindex, std::array< double, 8 > &X, std::array< double, 8 > &Y, std::array< double, 8 > &Z)
const std::array< int, 3 > & dimension() const
Definition: EGrid.hpp:43
const std::vector< std::string > & list_of_lgrs() const
Definition: EGrid.hpp:65
std::tuple< int, int, int, int, int, int, float > NNCentry
Definition: EGrid.hpp:62
EGrid(const std::string &filename, std::string grid_name="global")
int global_index(int i, int j, int k) const
const std::vector< int > & hostCellsGlobalIndex() const
Definition: EGrid.hpp:58
void getCellCorners(const std::array< int, 3 > &ijk, std::array< double, 8 > &X, std::array< double, 8 > &Y, std::array< double, 8 > &Z)
int totalNumberOfCells() const
Definition: EGrid.hpp:52
int activeCells() const
Definition: EGrid.hpp:51
bool is_radial() const
Definition: EGrid.hpp:56
std::array< int, 3 > ijk_from_active_index(int actInd) const
std::vector< std::array< int, 3 > > hostCellsIJK()
std::vector< NNCentry > get_nnc_ijk()
int active_index(int i, int j, int k) const
std::array< int, 3 > ijk_from_global_index(int globInd) const
Definition: EclFile.hpp:36
Definition: A.hpp:4