BoundaryConditions.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: BoundaryConditions.hpp
4//
5// Created: Mon Jun 29 15:19:42 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// Bård Skaflestad <bard.skaflestad@sintef.no>
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
18 Copyright 2009, 2010 Statoil ASA.
19
20 This file is part of The Open Reservoir Simulator Project (OpenRS).
21
22 OpenRS is free software: you can redistribute it and/or modify
23 it under the terms of the GNU General Public License as published by
24 the Free Software Foundation, either version 3 of the License, or
25 (at your option) any later version.
26
27 OpenRS is distributed in the hope that it will be useful,
28 but WITHOUT ANY WARRANTY; without even the implied warranty of
29 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 GNU General Public License for more details.
31
32 You should have received a copy of the GNU General Public License
33 along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
34*/
35
36#ifndef OPENRS_BOUNDARYCONDITIONS_HEADER
37#define OPENRS_BOUNDARYCONDITIONS_HEADER
38
39
40#include <vector>
41#include <ostream>
42#include <boost/mpl/if.hpp>
43#include <opm/common/ErrorMacros.hpp>
44#include <dune/common/fvector.hh>
45
46namespace Opm
47{
48
50 template <typename T>
51 class BCBase
52 {
53 public:
59
64 template<typename charT, class traits>
65 void write(std::basic_ostream<charT,traits>& os) const
66 {
67 os << "Type: " << type_ << " Value: " << value_;
68 }
69
70 protected: // methods
73 : type_(Neumann), value_(0.0)
74 {
75 }
79 BCBase(BCType type, T value)
80 : type_(type), value_(value)
81 {
82 }
85 bool isDirichlet() const
86 {
87 return type_ == Dirichlet;
88 }
91 bool isNeumann() const
92 {
93 return type_ == Neumann;
94 }
97 bool isPeriodic() const
98 {
99 return type_ == Periodic;
100 }
101
102 protected: // data
105 };
106
108 template<typename charT, class traits, typename T>
109 std::basic_ostream<charT,traits>&
110 operator<<(std::basic_ostream<charT,traits>& os,
111 const BCBase<T>& bc)
112 {
113 bc.write(os);
114 return os;
115 }
116
117
118
119
121 class FlowBC : public BCBase<double>
122 {
123 public:
126 : BCBase<double>(Neumann, 0.0)
127 {
128 }
132 FlowBC(BCType type, double value)
133 : BCBase<double>(type, value)
134 {
135 assert(isNeumann() || isDirichlet() || isPeriodic());
136 }
137
139 using BCBase<double>::isDirichlet;
140 using BCBase<double>::isNeumann;
141 using BCBase<double>::isPeriodic;
142
145 double pressure() const
146 {
147#ifndef NDEBUG
148 OPM_ERROR_IF(!isDirichlet(), "Pressure boundary conditions are only valid for Dirichlet boundaries");
149#endif
150 return value_;
151 }
154 double outflux() const
155 {
156#ifndef NDEBUG
157 OPM_ERROR_IF(!isNeumann(), "Outflux boundary conditions are only valid for Neumann boundaries");
158#endif
159 return value_;
160 }
163 double pressureDifference() const
164 {
165#ifndef NDEBUG
166 OPM_ERROR_IF(!isPeriodic(), "Pressure difference boundary conditions are only valid for periodic boundaries");
167#endif
168 return value_;
169 }
170 };
171
172
173
175 class SatBC : public BCBase<double>
176 {
177 public:
180 : BCBase<double>(Dirichlet, 1.0)
181 {
182 }
186 SatBC(BCType type, double value)
187 : BCBase<double>(type, value)
188 {
189 assert(isDirichlet() || isPeriodic());
190 }
192 using BCBase<double>::isDirichlet;
193 using BCBase<double>::isPeriodic;
194
197 double saturation() const
198 {
199 assert (isDirichlet());
200 return value_;
201 }
204 double saturationDifference() const
205 {
206 assert (isPeriodic());
207 return value_;
208 }
209 };
210
211
213 template <int numComponents>
214 class SurfvolBC : public BCBase<Dune::FieldVector<double, numComponents> >
215 {
216 public:
217 typedef Dune::FieldVector<double, numComponents> CompVec;
221 : Base(Base::Dirichlet, CompVec(0.0))
222 {
223 }
226 explicit SurfvolBC(Dune::FieldVector<double, numComponents> value)
227 : Base(Base::Dirichlet, value)
228 {
229 assert(isDirichlet());
230 }
232 using Base::isDirichlet;
233
237 {
238 return Base::value_;
239 }
240 };
241
242
243
245 {
246 public:
248 {
249 }
250
251 PeriodicConditionHandler(int num_different_boundary_ids)
252 : periodic_partner_bid_(num_different_boundary_ids, 0),
253 canonical_bid_(num_different_boundary_ids, 0)
254 {
255 }
256
257 void resize(int new_size)
258 {
259 periodic_partner_bid_.resize(new_size, 0);
260 canonical_bid_.resize(new_size, 0);
261 }
262
263 bool empty() const
264 {
265 return periodic_partner_bid_.empty() && canonical_bid_.empty();
266 }
267
268 void clear()
269 {
270 periodic_partner_bid_.clear();
271 canonical_bid_.clear();
272 }
273
274 int size() const
275 {
276 return periodic_partner_bid_.size();
277 }
278
279 void setPeriodicPartners(int boundary_id_1, int boundary_id_2)
280 {
281 assert(boundary_id_1 >= 0 && boundary_id_1 < int(periodic_partner_bid_.size()));
282 assert(boundary_id_2 >= 0 && boundary_id_2 < int(periodic_partner_bid_.size()));
283 assert(periodic_partner_bid_[boundary_id_1] == 0);
284 assert(periodic_partner_bid_[boundary_id_2] == 0);
285 periodic_partner_bid_[boundary_id_1] = boundary_id_2;
286 periodic_partner_bid_[boundary_id_2] = boundary_id_1;
287 }
288
289 int getPeriodicPartner(int boundary_id) const
290 {
291 assert(boundary_id >= 0 && boundary_id < int(periodic_partner_bid_.size()));
292 return periodic_partner_bid_[boundary_id];
293 }
294
295 void setCanonicalBoundaryId(int boundary_id, int canonical_bid)
296 {
297 assert(boundary_id >= 0 && boundary_id < int(canonical_bid_.size()));
298 canonical_bid_[boundary_id] = canonical_bid;
299 }
300
301 int getCanonicalBoundaryId(int boundary_id) const
302 {
303 assert(boundary_id >= 0 && boundary_id < int(canonical_bid_.size()));
304 return canonical_bid_[boundary_id];
305 }
306
307 template<typename charT, class traits>
308 void write(std::basic_ostream<charT,traits>& os) const
309 {
310 for (int i = 0; i < int(periodic_partner_bid_.size()); ++i) {
311 os << "Partner of bid " << i << " is " << periodic_partner_bid_[i]
312 << " (canonical bid is " << canonical_bid_[i] << '\n';
313 }
314 os << std::endl;
315 }
316 private:
317 std::vector<int> periodic_partner_bid_;
318 std::vector<int> canonical_bid_;
319 };
320
321 template<typename charT, class traits>
322 std::basic_ostream<charT,traits>&
323 operator<<(std::basic_ostream<charT,traits>& os,
324 const PeriodicConditionHandler& pch)
325 {
326 pch.write(os);
327 return os;
328 }
329
330
331 template <typename T>
333 {
334 public:
336 DummyVec(int) {}
337 void resize(int) {}
338 void clear() {}
339 };
340
341 template <bool FC = false, bool SC = false, bool ZC = false, int numComponents = 3>
343 private boost::mpl::if_c<FC, std::vector<FlowBC>, DummyVec<FlowBC> >::type,
344 private boost::mpl::if_c<SC, std::vector<SatBC>, DummyVec<SatBC> >::type,
345 private boost::mpl::if_c<ZC, std::vector<SurfvolBC<numComponents> >, DummyVec<SurfvolBC<numComponents> > >::type
346 {
347 public:
348 typedef typename boost::mpl::if_c<FC, std::vector<FlowBC>, DummyVec<FlowBC> >::type FlowConds;
349 typedef typename boost::mpl::if_c<SC, std::vector<SatBC>, DummyVec<SatBC> >::type SatConds;
350 typedef typename boost::mpl::if_c<ZC, std::vector<SurfvolBC<numComponents> >, DummyVec<SurfvolBC<numComponents> > >::type SurfvolConds;
351 const static bool HasFlowConds = FC;
352 const static bool HasSatConds = SC;
353 const static bool HasSurfvolConds = SC;
354
356 {
357 }
358
359 BasicBoundaryConditions(int num_different_boundary_ids)
360 : PeriodicConditionHandler(num_different_boundary_ids),
361 FlowConds(num_different_boundary_ids),
362 SatConds(num_different_boundary_ids),
363 SurfvolConds(num_different_boundary_ids)
364 {
365 }
366
367 void resize(int new_size)
368 {
370 FlowConds::resize(new_size);
371 SatConds::resize(new_size);
372 SurfvolConds::resize(new_size);
373 }
374
375 bool empty() const
376 {
378 }
379
380 void clear()
381 {
383 FlowConds::clear();
384 SatConds::clear();
385 SurfvolConds::clear();
386 }
387
388 int size() const
389 {
391 }
392
394 {
395 return FlowConds::operator[](i);
396 }
397
398 const FlowBC& flowCond(int i) const
399 {
400 return FlowConds::operator[](i);
401 }
402
403 template <class BoundaryFace>
404 const FlowBC& flowCond(const BoundaryFace& bf) const
405 {
406 assert(bf.boundary());
407 return FlowConds::operator[](bf.boundaryId());
408 }
409
411 {
412 return SatConds::operator[](i);
413 }
414
415 const SatBC& satCond(int i) const
416 {
417 return SatConds::operator[](i);
418 }
419
420 template <class BoundaryFace>
421 const SatBC& satCond(const BoundaryFace& bf) const
422 {
423 assert(bf.boundary());
424 return SatConds::operator[](bf.boundaryId());
425 }
426
428 {
429 return SurfvolConds::operator[](i);
430 }
431
433 {
434 return SurfvolConds::operator[](i);
435 }
436
437 template <class BoundaryFace>
438 const SurfvolBC<numComponents>& surfvolCond(const BoundaryFace& bf) const
439 {
440 assert(bf.boundary());
441 return SurfvolConds::operator[](bf.boundaryId());
442 }
443
444 template<typename charT, class traits>
445 void write(std::basic_ostream<charT,traits>& os) const
446 {
448 }
449 };
450
451 template<typename charT, class traits, bool F, bool S> //, bool P>
452 std::basic_ostream<charT,traits>&
453 operator<<(std::basic_ostream<charT,traits>& os,
454 const BasicBoundaryConditions<F,S>& bcs)
455 {
456 bcs.write(os);
457 return os;
458 }
459
460
461} // namespace Opm
462
463
464#endif // OPENRS_BOUNDARYCONDITIONS_HEADER
A class for building boundary conditions in a uniform way.
Definition: BoundaryConditions.hpp:52
BCType type_
Definition: BoundaryConditions.hpp:103
BCType
Enum for the allowed boundary condition types. So far, we support Dirichlet, Neumann and Periodic con...
Definition: BoundaryConditions.hpp:58
@ Periodic
Definition: BoundaryConditions.hpp:58
@ Neumann
Definition: BoundaryConditions.hpp:58
@ Dirichlet
Definition: BoundaryConditions.hpp:58
bool isPeriodic() const
Type query.
Definition: BoundaryConditions.hpp:97
BCBase()
Default constructor, that makes a Neumann condition with value 0.0.
Definition: BoundaryConditions.hpp:72
T value_
Definition: BoundaryConditions.hpp:104
BCBase(BCType type, T value)
Constructor taking a type and value.
Definition: BoundaryConditions.hpp:79
bool isDirichlet() const
Type query.
Definition: BoundaryConditions.hpp:85
void write(std::basic_ostream< charT, traits > &os) const
Write type and value to an output stream.
Definition: BoundaryConditions.hpp:65
bool isNeumann() const
Type query.
Definition: BoundaryConditions.hpp:91
Definition: BoundaryConditions.hpp:346
boost::mpl::if_c< SC, std::vector< SatBC >, DummyVec< SatBC > >::type SatConds
Definition: BoundaryConditions.hpp:349
static const bool HasSatConds
Definition: BoundaryConditions.hpp:352
BasicBoundaryConditions(int num_different_boundary_ids)
Definition: BoundaryConditions.hpp:359
const SurfvolBC< numComponents > & surfvolCond(int i) const
Definition: BoundaryConditions.hpp:432
const FlowBC & flowCond(int i) const
Definition: BoundaryConditions.hpp:398
const SurfvolBC< numComponents > & surfvolCond(const BoundaryFace &bf) const
Definition: BoundaryConditions.hpp:438
void clear()
Definition: BoundaryConditions.hpp:380
static const bool HasFlowConds
Definition: BoundaryConditions.hpp:351
SurfvolBC< numComponents > & surfvolCond(int i)
Definition: BoundaryConditions.hpp:427
const SatBC & satCond(int i) const
Definition: BoundaryConditions.hpp:415
boost::mpl::if_c< ZC, std::vector< SurfvolBC< numComponents > >, DummyVec< SurfvolBC< numComponents > > >::type SurfvolConds
Definition: BoundaryConditions.hpp:350
const FlowBC & flowCond(const BoundaryFace &bf) const
Definition: BoundaryConditions.hpp:404
static const bool HasSurfvolConds
Definition: BoundaryConditions.hpp:353
int size() const
Definition: BoundaryConditions.hpp:388
bool empty() const
Definition: BoundaryConditions.hpp:375
void resize(int new_size)
Definition: BoundaryConditions.hpp:367
SatBC & satCond(int i)
Definition: BoundaryConditions.hpp:410
boost::mpl::if_c< FC, std::vector< FlowBC >, DummyVec< FlowBC > >::type FlowConds
Definition: BoundaryConditions.hpp:348
void write(std::basic_ostream< charT, traits > &os) const
Definition: BoundaryConditions.hpp:445
FlowBC & flowCond(int i)
Definition: BoundaryConditions.hpp:393
const SatBC & satCond(const BoundaryFace &bf) const
Definition: BoundaryConditions.hpp:421
BasicBoundaryConditions()
Definition: BoundaryConditions.hpp:355
Definition: BoundaryConditions.hpp:333
DummyVec(int)
Definition: BoundaryConditions.hpp:336
DummyVec()
Definition: BoundaryConditions.hpp:335
void clear()
Definition: BoundaryConditions.hpp:338
void resize(int)
Definition: BoundaryConditions.hpp:337
A class for representing a flow boundary condition.
Definition: BoundaryConditions.hpp:122
double pressure() const
Query a Dirichlet condition.
Definition: BoundaryConditions.hpp:145
double outflux() const
Query a Neumann condition.
Definition: BoundaryConditions.hpp:154
FlowBC(BCType type, double value)
Constructor taking a type and value.
Definition: BoundaryConditions.hpp:132
double pressureDifference() const
Query a Periodic condition.
Definition: BoundaryConditions.hpp:163
FlowBC()
Default constructor, that makes a noflow condition (Neumann, value 0.0).
Definition: BoundaryConditions.hpp:125
Definition: BoundaryConditions.hpp:245
void resize(int new_size)
Definition: BoundaryConditions.hpp:257
PeriodicConditionHandler()
Definition: BoundaryConditions.hpp:247
bool empty() const
Definition: BoundaryConditions.hpp:263
void clear()
Definition: BoundaryConditions.hpp:268
int getCanonicalBoundaryId(int boundary_id) const
Definition: BoundaryConditions.hpp:301
void setPeriodicPartners(int boundary_id_1, int boundary_id_2)
Definition: BoundaryConditions.hpp:279
PeriodicConditionHandler(int num_different_boundary_ids)
Definition: BoundaryConditions.hpp:251
void setCanonicalBoundaryId(int boundary_id, int canonical_bid)
Definition: BoundaryConditions.hpp:295
void write(std::basic_ostream< charT, traits > &os) const
Definition: BoundaryConditions.hpp:308
int getPeriodicPartner(int boundary_id) const
Definition: BoundaryConditions.hpp:289
int size() const
Definition: BoundaryConditions.hpp:274
A class for representing a saturation boundary condition.
Definition: BoundaryConditions.hpp:176
double saturation() const
Query a Dirichlet condition.
Definition: BoundaryConditions.hpp:197
double saturationDifference() const
Query a Periodic condition.
Definition: BoundaryConditions.hpp:204
SatBC()
Default constructor, that makes a Dirichlet condition with value 1.0.
Definition: BoundaryConditions.hpp:179
SatBC(BCType type, double value)
Constructor taking a type and value.
Definition: BoundaryConditions.hpp:186
A class for representing a surface volume boundary condition.
Definition: BoundaryConditions.hpp:215
Dune::FieldVector< double, numComponents > CompVec
Definition: BoundaryConditions.hpp:217
CompVec surfvol() const
Query a Dirichlet condition.
Definition: BoundaryConditions.hpp:236
SurfvolBC()
Default constructor, that makes a Dirichlet condition with all zero component values.
Definition: BoundaryConditions.hpp:220
BCBase< CompVec > Base
Definition: BoundaryConditions.hpp:218
bool isDirichlet() const
Forwarding the relevant type query.
Definition: BoundaryConditions.hpp:85
SurfvolBC(Dune::FieldVector< double, numComponents > value)
Constructor taking a value.
Definition: BoundaryConditions.hpp:226
Definition: BlackoilFluid.hpp:32
std::basic_ostream< charT, traits > & operator<<(std::basic_ostream< charT, traits > &os, const BCBase< T > &bc)
Stream insertion for BCBase.
Definition: BoundaryConditions.hpp:110