GridAdapter.hpp
Go to the documentation of this file.
1/*
2 Copyright 2010 SINTEF ICT, Applied Mathematics.
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
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_GRIDADAPTER_HEADER_INCLUDED
21#define OPM_GRIDADAPTER_HEADER_INCLUDED
22
23
25#include <opm/grid/CpGrid.hpp>
26#include <stdexcept>
27#include <vector>
28
30{
31public:
36 template <class Grid>
37 void init(const Grid& grid)
38 {
39 buildTopology(grid);
40 buildGeometry(grid);
41 buildGlobalCell(grid);
42 copyCartDims(grid);
43 }
44
46 {
47 return &g_;
48 }
50 const UnstructuredGrid* c_grid() const
51 {
52 return &g_;
53 }
54
55 // ------ Forwarding the same interface that init() expects ------
56 //
57 // This is only done in order to verify that init() works correctly.
58 //
59
60 enum { dimension = 3 }; // This is actually a hack used for testing (dim is a runtime parameter).
61
62 struct Vector
63 {
64 explicit Vector(const double* source)
65 {
66 for (int i = 0; i < dimension; ++i) {
67 data[i] = source[i];
68 }
69 }
70 double& operator[] (const int ix)
71 {
72 return data[ix];
73 }
74 double operator[] (const int ix) const
75 {
76 return data[ix];
77 }
78 double data[dimension];
79 };
80
81 // Topology
82 int numCells() const
83 {
84 return g_.number_of_cells;
85 }
86 int numFaces() const
87 {
88 return g_.number_of_faces;
89 }
90 int numVertices() const
91 {
92 return g_.number_of_nodes;
93 }
94
95 int numCellFaces(int cell) const
96 {
97 return cell_facepos_[cell + 1] - cell_facepos_[cell];
98 }
99 int cellFace(int cell, int local_index) const
100 {
101 return cell_faces_[cell_facepos_[cell] + local_index];
102 }
103 int faceCell(int face, int local_index) const
104 {
105 return face_cells_[2*face + local_index];
106 }
107 int numFaceVertices(int face) const
108 {
109 return face_nodepos_[face + 1] - face_nodepos_[face];
110 }
111 int faceVertex(int face, int local_index) const
112 {
113 return face_nodes_[face_nodepos_[face] + local_index];
114 }
115
116 // Geometry
117 Vector vertexPosition(int vertex) const
118 {
119 return Vector(&node_coordinates_[g_.dimensions*vertex]);
120 }
121 double faceArea(int face) const
122 {
123 return face_areas_[face];
124 }
125 Vector faceCentroid(int face) const
126 {
127 return Vector(&face_centroids_[g_.dimensions*face]);
128 }
129 Vector faceNormal(int face) const
130 {
131 Vector fn(&face_normals_[g_.dimensions*face]);
132 // We must renormalize since the stored normals are
133 // 'unit normal * face area'.
134 double invfa = 1.0 / faceArea(face);
135 for (int i = 0; i < dimension; ++i) {
136 fn[i] *= invfa;
137 }
138 return fn;
139 }
140 double cellVolume(int cell) const
141 {
142 return cell_volumes_[cell];
143 }
144 Vector cellCentroid(int cell) const
145 {
146 return Vector(&cell_centroids_[g_.dimensions*cell]);
147 }
148
149 bool operator==(const GridAdapter& other)
150 {
151 return face_nodes_ == other.face_nodes_
152 && face_nodepos_ == other.face_nodepos_
153 && face_cells_ == other.face_cells_
154 && cell_faces_ == other.cell_faces_
155 && cell_facepos_ == other.cell_facepos_
156 && node_coordinates_ == other.node_coordinates_
157 && face_centroids_ == other.face_centroids_
158 && face_areas_ == other.face_areas_
159 && face_normals_ == other.face_normals_
160 && cell_centroids_ == other.cell_centroids_
161 && cell_volumes_ == other.cell_volumes_;
162 }
163 // make a grid which looks periodic but do not have 2 half faces for each
164 // periodic boundary
165 void makeQPeriodic(const std::vector<int>& hf_ind,const std::vector<int>& periodic_cells){
166 for(int i=0; i<int(hf_ind.size());++i){
167 //std::array<int,2> cells;
168 int& cell0=face_cells_[2*cell_faces_[ hf_ind[i] ]+0];
169 int& cell1=face_cells_[2*cell_faces_[ hf_ind[i] ]+1];
170 assert(periodic_cells[2*i+1]>=0);
171 if(periodic_cells[2*i+0] == cell0){
172 assert(cell1==-1);
173 cell1=periodic_cells[2*i+1];
174 }else{
175 assert(periodic_cells[2*i+0] == cell1);
176 assert(cell0==-1);
177 cell0=periodic_cells[2*i+1];
178 }
179 }
180 }
181private:
183 // Topology storage.
184 std::vector<int> face_nodes_;
185 std::vector<int> face_nodepos_;
186 std::vector<int> face_cells_;
187 std::vector<int> cell_faces_;
188 std::vector<int> cell_facepos_;
189 // Geometry storage.
190 std::vector<double> node_coordinates_;
191 std::vector<double> face_centroids_;
192 std::vector<double> face_areas_;
193 std::vector<double> face_normals_;
194 std::vector<double> cell_centroids_;
195 std::vector<double> cell_volumes_;
196 // The global cell information
197 std::vector<int> global_cell_;
199 void buildGlobalCell(const Dune::CpGrid& grid)
200 {
201 bool all_active=true;
202 int old_cell=-1;
203 global_cell_.resize(grid.numCells());
204 for(int c=0; c<grid.numCells(); ++c)
205 {
206 int new_cell=global_cell_[c]=grid.globalCell()[c];
207 all_active = all_active && (new_cell==old_cell+1);
208 old_cell=new_cell;
209 }
210 if(all_active){
211 g_.global_cell=0;
212 // really release memory
213 std::vector<int>().swap(global_cell_);
214 }
215 else
216 g_.global_cell=&(global_cell_[0]);
217 }
219 template<class Grid>
220 void buildGlobalCell(const Grid& /*grid*/)
221 {
222 g_.global_cell=0;
223 }
225 void copyCartDims(const Dune::CpGrid& grid)
226 {
227 for(int i=0; i<3; ++i)
228 g_.cartdims[i] = grid.logicalCartesianSize()[i];
229 }
231 template <class Grid>
232 void copyCartDims(const Grid& /*grid*/)
233 {}
235 template <class Grid>
236 void buildTopology(const Grid& grid)
237 {
238 // Face topology.
239 int num_cells = grid.numCells();
240 int num_faces = grid.numFaces();
241 face_nodepos_.resize(num_faces + 1);
242 int facenodecount = 0;
243 for (int f = 0; f < num_faces; ++f) {
244 face_nodepos_[f] = facenodecount;
245 facenodecount += grid.numFaceVertices(f);
246 }
247 face_nodepos_.back() = facenodecount;
248 face_nodes_.resize(facenodecount);
249 for (int f = 0; f < num_faces; ++f) {
250 for (int local = 0; local < grid.numFaceVertices(f); ++local) {
251 face_nodes_[face_nodepos_[f] + local] = grid.faceVertex(f, local);
252 }
253 }
254 face_cells_.resize(2*num_faces);
255 for (int f = 0; f < num_faces; ++f) {
256 face_cells_[2*f] = grid.faceCell(f, 0);
257 face_cells_[2*f + 1] = grid.faceCell(f, 1);
258 }
259
260 // Cell topology.
261 int cellfacecount = 0;
262 cell_facepos_.resize(num_cells + 1);
263 for (int c = 0; c < num_cells; ++c) {
264 cell_facepos_[c] = cellfacecount;
265 cellfacecount += grid.numCellFaces(c);
266 }
267 cell_facepos_.back() = cellfacecount;
268 cell_faces_.resize(cellfacecount);
269 for (int c = 0; c < num_cells; ++c) {
270 for (int local = 0; local < grid.numCellFaces(c); ++local) {
271 cell_faces_[cell_facepos_[c] + local] = grid.cellFace(c, local);
272 }
273 }
274
275 // Set C grid members.
276 g_.dimensions = Grid::dimension;
277 g_.number_of_cells = grid.numCells();
278 g_.number_of_faces = grid.numFaces();
279 g_.number_of_nodes = grid.numVertices();
280 g_.face_nodes = &face_nodes_[0];
281 g_.face_nodepos = &face_nodepos_[0];
282 g_.face_cells = &face_cells_[0];
283 g_.cell_faces = &cell_faces_[0];
284 g_.cell_facepos = &cell_facepos_[0];
285 }
286
287
290 template <class Grid>
291 void buildGeometry(const Grid& grid)
292 {
293 // Node geometry.
294 int num_cells = grid.numCells();
295 int num_nodes = grid.numVertices();
296 int num_faces = grid.numFaces();
297 int dim = Grid::dimension;
298 node_coordinates_.resize(dim*num_nodes);
299 for (int n = 0; n < num_nodes; ++n) {
300 for (int dd = 0; dd < dim; ++dd) {
301 node_coordinates_[dim*n + dd] = grid.vertexPosition(n)[dd];
302 }
303 }
304
305 // Face geometry.
306 face_centroids_.resize(dim*num_faces);
307 face_areas_.resize(num_faces);
308 face_normals_.resize(dim*num_faces);
309 for (int f = 0; f < num_faces; ++f) {
310 face_areas_[f] = grid.faceArea(f);
311 for (int dd = 0; dd < dim; ++dd) {
312 face_centroids_[dim*f + dd] = grid.faceCentroid(f)[dd];
313 face_normals_[dim*f + dd] = grid.faceNormal(f)[dd]*face_areas_[f];
314 }
315 }
316
317 // Cell geometry.
318 cell_centroids_.resize(dim*num_cells);
319 cell_volumes_.resize(num_cells);
320 for (int c = 0; c < num_cells; ++c) {
321 cell_volumes_[c] = grid.cellVolume(c);
322 for (int dd = 0; dd < dim; ++dd) {
323 cell_centroids_[dim*c + dd] = grid.cellCentroid(c)[dd];
324 }
325 }
326
327 // Set C grid members.
328 g_.node_coordinates = &node_coordinates_[0];
329 g_.face_centroids = &face_centroids_[0];
330 g_.face_areas = &face_areas_[0];
331 g_.face_normals = &face_normals_[0];
332 g_.cell_centroids = &cell_centroids_[0];
333 g_.cell_volumes = &cell_volumes_[0];
334 }
335};
336
337
338#endif // OPM_GRIDADAPTER_HEADER_INCLUDED
[ provides Dune::Grid ]
Definition: CpGrid.hpp:238
const std::vector< int > & globalCell() const
const std::array< int, 3 > & logicalCartesianSize() const
int numCells() const
Get the number of cells.
Definition: GridAdapter.hpp:30
int faceCell(int face, int local_index) const
Definition: GridAdapter.hpp:103
int faceVertex(int face, int local_index) const
Definition: GridAdapter.hpp:111
int numCellFaces(int cell) const
Definition: GridAdapter.hpp:95
void makeQPeriodic(const std::vector< int > &hf_ind, const std::vector< int > &periodic_cells)
Definition: GridAdapter.hpp:165
int numFaceVertices(int face) const
Definition: GridAdapter.hpp:107
UnstructuredGrid * c_grid()
Definition: GridAdapter.hpp:45
void init(const Grid &grid)
Initialize the grid.
Definition: GridAdapter.hpp:37
const UnstructuredGrid * c_grid() const
Access the underlying C grid.
Definition: GridAdapter.hpp:50
double faceArea(int face) const
Definition: GridAdapter.hpp:121
@ dimension
Definition: GridAdapter.hpp:60
Vector cellCentroid(int cell) const
Definition: GridAdapter.hpp:144
int numVertices() const
Definition: GridAdapter.hpp:90
Vector faceNormal(int face) const
Definition: GridAdapter.hpp:129
int cellFace(int cell, int local_index) const
Definition: GridAdapter.hpp:99
int numCells() const
Definition: GridAdapter.hpp:82
int numFaces() const
Definition: GridAdapter.hpp:86
Vector faceCentroid(int face) const
Definition: GridAdapter.hpp:125
double cellVolume(int cell) const
Definition: GridAdapter.hpp:140
bool operator==(const GridAdapter &other)
Definition: GridAdapter.hpp:149
Vector vertexPosition(int vertex) const
Definition: GridAdapter.hpp:117
Dune::FieldVector< double, 3 > Vector
Definition: cpgrid/GridHelpers.hpp:385
Definition: GridAdapter.hpp:63
double data[dimension]
Definition: GridAdapter.hpp:78
Vector(const double *source)
Definition: GridAdapter.hpp:64
double & operator[](const int ix)
Definition: GridAdapter.hpp:70
Definition: UnstructuredGrid.h:99
int * face_nodes
Definition: UnstructuredGrid.h:121
int * face_nodepos
Definition: UnstructuredGrid.h:127
int number_of_cells
Definition: UnstructuredGrid.h:109
double * face_areas
Definition: UnstructuredGrid.h:173
int * cell_faces
Definition: UnstructuredGrid.h:146
int number_of_faces
Definition: UnstructuredGrid.h:111
double * cell_centroids
Definition: UnstructuredGrid.h:192
int * cell_facepos
Definition: UnstructuredGrid.h:152
double * cell_volumes
Definition: UnstructuredGrid.h:197
int * face_cells
Definition: UnstructuredGrid.h:138
double * face_centroids
Definition: UnstructuredGrid.h:168
int number_of_nodes
Definition: UnstructuredGrid.h:113
int dimensions
Definition: UnstructuredGrid.h:106
int cartdims[3]
Definition: UnstructuredGrid.h:227
int * global_cell
Definition: UnstructuredGrid.h:214
double * face_normals
Definition: UnstructuredGrid.h:184
double * node_coordinates
Definition: UnstructuredGrid.h:160