BoundaryPeriodicity.hpp
Go to the documentation of this file.
1/*
2 Copyright 2011 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_BOUNDARYPERIODICITY_HEADER_INCLUDED
21#define OPM_BOUNDARYPERIODICITY_HEADER_INCLUDED
22
23#include <dune/common/fvector.hh>
24
25#include <opm/common/ErrorMacros.hpp>
26
27#include <algorithm>
28#include <array>
29#include <iostream>
30#include <numbers>
31#include <vector>
32
33namespace Opm
34{
35
39 {
43 int bid;
56 double area;
58 Dune::FieldVector<double,3> centroid;
59
63 bool operator<(const BoundaryFaceInfo& other) const
64 {
65 return cmpval() < other.cmpval();
66 }
67
68 private:
69 double cmpval() const
70 {
71 return centroid[(canon_pos/2 + 1)%3]
72 + std::numbers::pi*centroid[(canon_pos/2 + 2)%3];
73 }
74 };
75
76
83 bool match(std::vector<BoundaryFaceInfo>& bfaces, int face, int lower, int upper);
84
86 template <class GridView>
87 void findPeriodicPartners(std::vector<BoundaryFaceInfo>& bfinfo,
88 std::array<double, 6>& side_areas,
89 const GridView& g,
90 const std::array<bool, 2*GridView::dimension>& is_periodic,
91 double spatial_tolerance = 1e-6)
92 {
93 // Pick out all boundary faces, simultaneously find the
94 // bounding box of their centroids, and the max id.
95 const int dim = GridView::dimension;
96 typedef typename GridView::template Codim<0>::Iterator CI;
97 typedef typename GridView::IntersectionIterator FI;
98 typedef Dune::FieldVector<double, dim> Vector;
99 std::vector<FI> bface_iters;
100 Vector low(1e100);
101 Vector hi(-1e100);
102 int max_bid = 0;
103 for (CI c = g.template begin<0>(); c != g.template end<0>(); ++c) {
104 for (FI f = g.ibegin(*c); f != g.iend(*c); ++f) {
105 if (f->boundaryId()) {
106 bface_iters.push_back(f);
107 Vector fcent = f->geometry().center();
108 for (int dd = 0; dd < dim; ++dd) {
109 low[dd] = std::min(low[dd], fcent[dd]);
110 hi[dd] = std::max(hi[dd], fcent[dd]);
111 }
112 max_bid = std::max(max_bid, f->boundaryId());
113 }
114 }
115 }
116 int num_bdy = bface_iters.size();
117 if (max_bid != num_bdy) {
118 OPM_THROW(std::runtime_error, "createPeriodic() assumes that every boundary face has a unique boundary id. That seems to be violated.");
119 }
120
121 // Store boundary face info in a suitable structure. Also find side total volumes.
122 std::ranges::fill(side_areas, 0.0);
123 bfinfo.clear();
124 bfinfo.reserve(num_bdy);
125 for (int i = 0; i < num_bdy; ++i) {
127 bf.face_index = i;
128 bf.bid = bface_iters[i]->boundaryId();
129 bf.canon_pos = -1;
130 bf.partner_face_index = -1;
131 bf.partner_bid = 0;
132 bf.area = bface_iters[i]->geometry().volume();
133 bf.centroid = bface_iters[i]->geometry().center();
134 for (int dd = 0; dd < dim; ++dd) {
135 double coord = bf.centroid[dd];
136 if (fabs(coord - low[dd]) <= spatial_tolerance) {
137 bf.canon_pos = 2*dd;
138 break;
139 } else if (fabs(coord - hi[dd]) <= spatial_tolerance) {
140 bf.canon_pos = 2*dd + 1;
141 break;
142 }
143 }
144 if (bf.canon_pos == -1) {
145 std::cerr << "Centroid: " << bf.centroid << "\n";
146 std::cerr << "Bounding box min: " << low << "\n";
147 std::cerr << "Bounding box max: " << hi << "\n";
148 OPM_THROW(std::runtime_error, "Boundary face centroid not on bounding box. Maybe the grid is not an axis-aligned shoe-box?");
149 }
150 side_areas[bf.canon_pos] += bf.area;
151 bf.centroid[bf.canon_pos/2] = 0.0;
152 bfinfo.push_back(bf);
153 }
154 assert(bfinfo.size() == bface_iters.size());
155
156 // Sort the infos so that partners end up close.
157 std::sort(bfinfo.begin(), bfinfo.end());
158
159 // Identify partners.
160 for (int i = 0; i < num_bdy; ++i) {
161 if (bfinfo[i].partner_face_index != -1) {
162 continue;
163 }
164 if (!is_periodic[bfinfo[i].canon_pos]) {
165 continue;
166 }
167 int lower = std::max(0, i - 10);
168 int upper = std::min(num_bdy, i + 10);
169 bool ok = match(bfinfo, i, lower, upper);
170 if (!ok) {
171 // We have not found a partner.
172 ok = match(bfinfo, i, 0, num_bdy);
173 if (!ok) {
174 OPM_MESSAGE("Warning: No partner found for boundary id " << bfinfo[i].bid);
175 // OPM_THROW(std::runtime_error, "No partner found.");
176 }
177 }
178 }
179 }
180
181 /*
183 template <class GridInterface>
184 void findPeriodicPartners(std::vector<BoundaryFaceInfo>& bfinfo,
185 std::array<double, 6>& side_areas,
186 const GridInterface& g,
187 const std::array<bool, 2*GridInterface::Dimension>& is_periodic,
188 double spatial_tolerance = 1e-6)
189 {
190 // Pick out all boundary faces, simultaneously find the
191 // bounding box of their centroids, and the max id.
192 typedef typename GridInterface::CellIterator CI;
193 typedef typename CI::FaceIterator FI;
194 typedef typename GridInterface::Vector Vector;
195 std::vector<FI> bface_iters;
196 Vector low(1e100);
197 Vector hi(-1e100);
198 int max_bid = 0;
199 for (CI c = g.cellbegin(); c != g.cellend(); ++c) {
200 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
201 if (f->boundaryId()) {
202 bface_iters.push_back(f);
203 Vector fcent = f->centroid();
204 for (int dd = 0; dd < GridInterface::Dimension; ++dd) {
205 low[dd] = std::min(low[dd], fcent[dd]);
206 hi[dd] = std::max(hi[dd], fcent[dd]);
207 }
208 max_bid = std::max(max_bid, f->boundaryId());
209 }
210 }
211 }
212 int num_bdy = bface_iters.size();
213 if (max_bid != num_bdy) {
214 OPM_THROW(std::runtime_error, "createPeriodic() assumes that every boundary face has a unique boundary id. That seems to be violated.");
215 }
216
217 // Store boundary face info in a suitable structure. Also find side total volumes.
218 std::ranges::fill(side_areas, 0.0);
219 bfinfo.clear();
220 bfinfo.reserve(num_bdy);
221 for (int i = 0; i < num_bdy; ++i) {
222 BoundaryFaceInfo bf;
223 bf.face_index = i;
224 bf.bid = bface_iters[i]->boundaryId();
225 bf.canon_pos = -1;
226 bf.partner_face_index = -1;
227 bf.partner_bid = 0;
228 bf.area = bface_iters[i]->area();
229 bf.centroid = bface_iters[i]->centroid();
230 for (int dd = 0; dd < GridInterface::Dimension; ++dd) {
231 double coord = bf.centroid[dd];
232 if (fabs(coord - low[dd]) <= spatial_tolerance) {
233 bf.canon_pos = 2*dd;
234 break;
235 } else if (fabs(coord - hi[dd]) <= spatial_tolerance) {
236 bf.canon_pos = 2*dd + 1;
237 break;
238 }
239 }
240 if (bf.canon_pos == -1) {
241 std::cerr << "Centroid: " << bf.centroid << "\n";
242 std::cerr << "Bounding box min: " << low << "\n";
243 std::cerr << "Bounding box max: " << hi << "\n";
244 OPM_THROW(std::runtime_error, "Boundary face centroid not on bounding box. Maybe the grid is not an axis-aligned shoe-box?");
245 }
246 side_areas[bf.canon_pos] += bf.area;
247 bf.centroid[bf.canon_pos/2] = 0.0;
248 bfinfo.push_back(bf);
249 }
250 assert(bfinfo.size() == bface_iters.size());
251
252 // Sort the infos so that partners end up close.
253 std::sort(bfinfo.begin(), bfinfo.end());
254
255 // Identify partners.
256 for (int i = 0; i < num_bdy; ++i) {
257 if (bfinfo[i].partner_face_index != -1) {
258 continue;
259 }
260 if (!is_periodic[bfinfo[i].canon_pos]) {
261 continue;
262 }
263 int lower = std::max(0, i - 10);
264 int upper = std::min(num_bdy, i + 10);
265 bool ok = match(bfinfo, i, lower, upper);
266 if (!ok) {
267 // We have not found a partner.
268 ok = match(bfinfo, i, 0, num_bdy);
269 if (!ok) {
270 OPM_MESSAGE("Warning: No partner found for boundary id " << bfinfo[i].bid);
271 // OPM_THROW(std::runtime_error, "No partner found.");
272 }
273 }
274 }
275 }
276*/
277
278
279
280} // namespace Opm
281
282#endif // OPM_BOUNDARYPERIODICITY_HEADER_INCLUDED
min[0]
Definition: elasticity_upscale_impl.hpp:146
Dune::BlockVector< Dune::FieldVector< double, 1 > > Vector
A vector holding our RHS.
Definition: matrixops.hpp:33
Definition: ImplicitAssembly.hpp:43
bool match(std::vector< BoundaryFaceInfo > &bfaces, int face, int lower, int upper)
Find a match (periodic partner) for the given face.
void findPeriodicPartners(std::vector< BoundaryFaceInfo > &bfinfo, std::array< double, 6 > &side_areas, const GridView &g, const std::array< bool, 2 *GridView::dimension > &is_periodic, double spatial_tolerance=1e-6)
Common implementation for the various createPeriodic functions.
Definition: BoundaryPeriodicity.hpp:87
Definition: BoundaryPeriodicity.hpp:39
double area
Face area.
Definition: BoundaryPeriodicity.hpp:56
bool operator<(const BoundaryFaceInfo &other) const
Definition: BoundaryPeriodicity.hpp:63
int bid
Boundary id of this face.
Definition: BoundaryPeriodicity.hpp:43
int face_index
Face index in [0, ..., #faces - 1].
Definition: BoundaryPeriodicity.hpp:41
int partner_bid
Boundary id of periodic partner face, or 0 if no parner.
Definition: BoundaryPeriodicity.hpp:54
int partner_face_index
Face index of periodic partner, or -1 if no partner.
Definition: BoundaryPeriodicity.hpp:52
Dune::FieldVector< double, 3 > centroid
Face centroid.
Definition: BoundaryPeriodicity.hpp:58
int canon_pos
Definition: BoundaryPeriodicity.hpp:50