boundarytypes.hh 15.3 KB
Newer Older
1
2
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
3
/*****************************************************************************
4
 *   See the file COPYING for full copying permissions.                      *
5
 *                                                                           *
6
 *   This program 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
9
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (at your option) any later version.                                     *
10
 *                                                                           *
11
12
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
13
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
14
15
16
17
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
18
19
20
21
22
23
24
25
26
27
28
29
30
 *****************************************************************************/
/*!
 * \file
 * \brief Class to specify the type of a boundary.
 */
#ifndef BOUNDARY_TYPES_HH
#define BOUNDARY_TYPES_HH

#include <dumux/common/valgrind.hh>

namespace Dumux
{

31
/*!
32
 * \ingroup BC
33
34
 * \brief Class to specify the type of a boundary.
 */
35
36
37
38
39
40
41
42
43
44
45
template <int numEq>
class BoundaryTypes
{
public:
    BoundaryTypes()
    { reset(); }

    /*!
     * \brief Reset the boundary types.
     *
     * After this method no equations will be disabled and neither
46
47
     * Neumann nor Dirichlet conditions will be evaluated. This
     * corrosponds to a Neumann zero boundary.
48
49
50
51
52
53
54
     */
    void reset()
    {
        for (int i=0; i < numEq; ++i) {
            boundaryInfo_[i].visited = 0;

            boundaryInfo_[i].isDirichlet = 0;
Andreas Lauser's avatar
Andreas Lauser committed
55
            boundaryInfo_[i].isNeumann = 0;
56
            boundaryInfo_[i].isOutflow = 0;
57
58
            boundaryInfo_[i].isCouplingInflow = 0;
            boundaryInfo_[i].isCouplingOutflow = 0;
59
            boundaryInfo_[i].isMortarCoupling = 0;
60
61
62

            eq2pvIdx_[i] = i;
            pv2eqIdx_[i] = i;
63
        }
64
65
66
    }

    /*!
67
68
69
70
     * \brief Returns true if the boundary types for a given equation
     *        has been specified.
     *
     * \param eqIdx The index of the equation
71
72
     */
    bool isSet(int eqIdx) const
73
    { return boundaryInfo_[eqIdx].visited; }
74
75
76
77

    /*!
     * \brief Make sure the boundary conditions are well-posed.
     *
78
79
     * If they are not, an assertation fails and the program aborts!
     * (if the NDEBUG macro is not defined)
80
81
82
     */
    void checkWellPosed() const
    {
83
84
#ifndef NDEBUG
        for (int i=0; i < numEq; ++i)
85
86
            // if this fails, at least one condition is missing.
            assert(boundaryInfo_[i].visited);
87
#endif
88
    }
89
90

    /*!
91
     * \brief Set all boundary conditions to Neumann.
92
93
94
     */
    void setAllNeumann()
    {
95
96
97
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
        {
            setNeumann(eqIdx);
98
99
100
101
        }
    }

    /*!
102
     * \brief Set all boundary conditions to Dirichlet.
103
104
105
     */
    void setAllDirichlet()
    {
106
107
108
        for (int eqIdx = 0; eqIdx < numEq; ++ eqIdx)
        {
            setDirichlet(eqIdx);
109
110
111
        }
    }

112
    /*!
113
     * \brief Set all boundary conditions to Neumann.
114
115
116
     */
    void setAllOutflow()
    {
117
118
119
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
        {
            setOutflow(eqIdx);
120
121
122
123
124
125
        }
    }

    /*!
     * \brief Set all boundary conditions to coupling inflow.
     */
126
    DUNE_DEPRECATED_MSG("setAllCouplingInflow() is deprecated")
127
128
129
130
131
132
133
134
135
    void setAllCouplingInflow()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            boundaryInfo_[eqIdx].visited = 1;
            boundaryInfo_[eqIdx].isDirichlet = 0;
            boundaryInfo_[eqIdx].isNeumann = 0;
            boundaryInfo_[eqIdx].isOutflow = 0;
            boundaryInfo_[eqIdx].isCouplingInflow = 1;
            boundaryInfo_[eqIdx].isCouplingOutflow = 0;
136
            boundaryInfo_[eqIdx].isMortarCoupling = 0;
137
138
139
140
141
142
143
144

            Valgrind::SetDefined(boundaryInfo_[eqIdx]);
        }
    }

    /*!
     * \brief Set all boundary conditions to coupling outflow.
     */
145
    DUNE_DEPRECATED_MSG("setAllCouplingOutflow() is deprecated")
146
147
148
149
150
151
    void setAllCouplingOutflow()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            boundaryInfo_[eqIdx].visited = 1;
            boundaryInfo_[eqIdx].isDirichlet = 0;
            boundaryInfo_[eqIdx].isNeumann = 0;
152
            boundaryInfo_[eqIdx].isOutflow = 0;
153
154
            boundaryInfo_[eqIdx].isCouplingInflow = 0;
            boundaryInfo_[eqIdx].isCouplingOutflow = 1;
155
156
157
158
159
160
161
162
163
            boundaryInfo_[eqIdx].isMortarCoupling = 0;

            Valgrind::SetDefined(boundaryInfo_[eqIdx]);
        }
    }

    /*!
     * \brief Set all boundary conditions to mortar coupling.
     */
164
    DUNE_DEPRECATED_MSG("setAllMortarCoupling() is deprecated")
165
166
167
168
169
170
171
172
173
174
    void setAllMortarCoupling()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            boundaryInfo_[eqIdx].visited = 1;
            boundaryInfo_[eqIdx].isDirichlet = 0;
            boundaryInfo_[eqIdx].isNeumann = 0;
            boundaryInfo_[eqIdx].isOutflow = 0;
            boundaryInfo_[eqIdx].isCouplingInflow = 0;
            boundaryInfo_[eqIdx].isCouplingOutflow = 0;
            boundaryInfo_[eqIdx].isMortarCoupling = 1;
175
176
177
178
179

            Valgrind::SetDefined(boundaryInfo_[eqIdx]);
        }
    }

180
    /*!
181
     * \brief Set a Neumann boundary condition for a single a single
182
     *        equation.
183
184
     *
     * \param eqIdx The index of the equation
185
186
187
188
     */
    void setNeumann(int eqIdx)
    {
        boundaryInfo_[eqIdx].visited = 1;
189
190
191
        boundaryInfo_[eqIdx].isDirichlet = 0;
        boundaryInfo_[eqIdx].isNeumann = 1;
        boundaryInfo_[eqIdx].isOutflow = 0;
192
193
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 0;
194
        boundaryInfo_[eqIdx].isMortarCoupling = 0;
195
196
197
198
199

        Valgrind::SetDefined(boundaryInfo_[eqIdx]);
    }

    /*!
200
     * \brief Set a Dirichlet boundary condition for a single primary
201
202
     *        variable
     *
203
204
205
206
     * \param pvIdx The index of the primary variable for which the
     *              Dirichlet condition should apply.
     * \param eqIdx The index of the equation which should used to set
     *              the Dirichlet condition
207
     */
208
    void setDirichlet(int pvIdx, int eqIdx)
209
210
211
    {
        boundaryInfo_[eqIdx].visited = 1;
        boundaryInfo_[eqIdx].isDirichlet = 1;
Andreas Lauser's avatar
Andreas Lauser committed
212
        boundaryInfo_[eqIdx].isNeumann = 0;
213
        boundaryInfo_[eqIdx].isOutflow = 0;
214
215
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 0;
216
        boundaryInfo_[eqIdx].isMortarCoupling = 0;
217
218
219
220
221
222
223
224

        // update the equation <-> primary variable mapping
        eq2pvIdx_[eqIdx] = pvIdx;
        pv2eqIdx_[pvIdx] = eqIdx;

        Valgrind::SetDefined(boundaryInfo_[eqIdx]);
    }

225
    /*!
226
     * \brief Set a Neumann boundary condition for a single a single
227
     *        equation.
228
229
230
     *
     * \param eqIdx The index of the equation on which the outflow
     *              condition applies.
231
232
233
234
235
236
237
     */
    void setOutflow(int eqIdx)
    {
        boundaryInfo_[eqIdx].visited = 1;
        boundaryInfo_[eqIdx].isDirichlet = 0;
        boundaryInfo_[eqIdx].isNeumann = 0;
        boundaryInfo_[eqIdx].isOutflow = 1;
238
239
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 0;
240
        boundaryInfo_[eqIdx].isMortarCoupling = 0;
241
242
243
244

        Valgrind::SetDefined(boundaryInfo_[eqIdx]);
    }

245
246
247
    /*!
     * \brief Set a boundary condition for a single equation to coupling inflow.
     */
248
    DUNE_DEPRECATED_MSG("setCouplingInflow() is deprecated")
249
250
251
252
253
254
255
256
    void setCouplingInflow(int eqIdx)
    {
        boundaryInfo_[eqIdx].visited = 1;
        boundaryInfo_[eqIdx].isDirichlet = 0;
        boundaryInfo_[eqIdx].isNeumann = 0;
        boundaryInfo_[eqIdx].isOutflow = 0;
        boundaryInfo_[eqIdx].isCouplingInflow = 1;
        boundaryInfo_[eqIdx].isCouplingOutflow = 0;
257
        boundaryInfo_[eqIdx].isMortarCoupling = 0;
258
259
260
261
262
263
264

        Valgrind::SetDefined(boundaryInfo_[eqIdx]);
    }

    /*!
     * \brief Set a boundary condition for a single equation to coupling outflow.
     */
265
    DUNE_DEPRECATED_MSG("setCouplingOutflow() is deprecated")
266
267
268
269
270
271
    void setCouplingOutflow(int eqIdx)
    {
        boundaryInfo_[eqIdx].visited = 1;
        boundaryInfo_[eqIdx].isDirichlet = 0;
        boundaryInfo_[eqIdx].isNeumann = 0;
        boundaryInfo_[eqIdx].isOutflow = 0;
Klaus Mosthaf's avatar
Klaus Mosthaf committed
272
273
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 1;
274
275
276
277
278
279
280
281
        boundaryInfo_[eqIdx].isMortarCoupling = 0;

        Valgrind::SetDefined(boundaryInfo_[eqIdx]);
    }

    /*!
     * \brief Set a boundary condition for a single equation to mortar coupling.
     */
282
    DUNE_DEPRECATED_MSG("setMortarCoupling() is deprecated")
283
284
285
286
287
288
289
290
291
    void setMortarCoupling(int eqIdx)
    {
        boundaryInfo_[eqIdx].visited = 1;
        boundaryInfo_[eqIdx].isDirichlet = 0;
        boundaryInfo_[eqIdx].isNeumann = 0;
        boundaryInfo_[eqIdx].isOutflow = 0;
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 0;
        boundaryInfo_[eqIdx].isMortarCoupling = 1;
292
293
294
295

        Valgrind::SetDefined(boundaryInfo_[eqIdx]);
    }

296
    /*!
297
     * \brief Set a Dirichlet boundary condition for a single primary
298
     *        variable.
299
     *
300
301
     * Depending on the discretization, setting the Dirichlet condition
     * will replace the balance equation with index equal to pvIdx.
302
     *
303
304
     * \param pvIdx The index of the primary variable inside a
     *              PrimaryVariables object.
305
     */
306
    void setDirichlet(int pvIdx)
307
    {
308
        setDirichlet(pvIdx, pvIdx);
309
310
311
312
    }

    /*!
     * \brief Returns true if an equation is used to specify a
313
     *        Dirichlet condition.
314
315
     *
     * \param eqIdx The index of the equation
316
317
     */
    bool isDirichlet(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
318
    { return boundaryInfo_[eqIdx].isDirichlet; }
319
320
321

    /*!
     * \brief Returns true if some equation is used to specify a
322
     *        Dirichlet condition.
323
324
325
326
327
328
329
     */
    bool hasDirichlet() const
    {
        for (int i = 0; i < numEq; ++i)
            if (boundaryInfo_[i].isDirichlet)
                return true;
        return false;
Katherina Baber's avatar
Katherina Baber committed
330
    }
331
332
333

    /*!
     * \brief Returns true if an equation is used to specify a
334
     *        Neumann condition.
335
336
     *
     * \param eqIdx The index of the equation
337
338
     */
    bool isNeumann(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
339
    { return boundaryInfo_[eqIdx].isNeumann; }
340
341
342

    /*!
     * \brief Returns true if some equation is used to specify a
343
     *        Neumann condition.
344
345
346
347
     */
    bool hasNeumann() const
    {
        for (int i = 0; i < numEq; ++i)
Andreas Lauser's avatar
Andreas Lauser committed
348
            if (boundaryInfo_[i].isNeumann)
349
350
                return true;
        return false;
Katherina Baber's avatar
Katherina Baber committed
351
    }
352

353
354
355
    /*!
     * \brief Returns true if an equation is used to specify an
     *        outflow condition.
356
357
     *
     * \param eqIdx The index of the equation
358
359
     */
    bool isOutflow(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
360
    { return boundaryInfo_[eqIdx].isOutflow; }
361
362
363
364
365
366
367
368
369
370
371

    /*!
     * \brief Returns true if some equation is used to specify an
     *        outflow condition.
     */
    bool hasOutflow() const
    {
        for (int i = 0; i < numEq; ++i)
            if (boundaryInfo_[i].isOutflow)
                return true;
        return false;
Katherina Baber's avatar
Katherina Baber committed
372
    }
373

374
375
376
    /*!
     * \brief Returns true if an equation is used to specify an
     *        inflow coupling condition.
377
378
     *
     * \param eqIdx The index of the equation
379
     */
380
    DUNE_DEPRECATED_MSG("isCouplingInflow() is deprecated")
381
    bool isCouplingInflow(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
382
    { return boundaryInfo_[eqIdx].isCouplingInflow; }
383
384
385
386
387

    /*!
     * \brief Returns true if some equation is used to specify an
     *        inflow coupling condition.
     */
388
    DUNE_DEPRECATED_MSG("hasCouplingInflow() is deprecated")
389
390
391
392
393
394
    bool hasCouplingInflow() const
    {
        for (int i = 0; i < numEq; ++i)
            if (boundaryInfo_[i].isCouplingInflow)
                return true;
        return false;
Katherina Baber's avatar
Katherina Baber committed
395
    }
396
397
398
399

    /*!
     * \brief Returns true if an equation is used to specify an
     *        outflow coupling condition.
400
401
     *
     * \param eqIdx The index of the equation
402
     */
403
    DUNE_DEPRECATED_MSG("isCouplingOutflow() is deprecated")
404
    bool isCouplingOutflow(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
405
    { return boundaryInfo_[eqIdx].isCouplingOutflow; }
406
407
408
409
410

    /*!
     * \brief Returns true if some equation is used to specify an
     *        outflow coupling condition.
     */
411
    DUNE_DEPRECATED_MSG("hasCouplingOutflow() is deprecated")
412
413
414
415
416
417
    bool hasCouplingOutflow() const
    {
        for (int i = 0; i < numEq; ++i)
            if (boundaryInfo_[i].isCouplingOutflow)
                return true;
        return false;
Katherina Baber's avatar
Katherina Baber committed
418
    }
419

420
421
422
423
424
425
    /*!
     * \brief Returns true if an equation is used to specify a
     *        mortar coupling condition.
     *
     * \param eqIdx The index of the equation
     */
426
    DUNE_DEPRECATED_MSG("isMortarCoupling() is deprecated")
427
428
429
430
431
432
433
434
435
    bool isMortarCoupling(unsigned eqIdx) const
    {
        return boundaryInfo_[eqIdx].isMortarCoupling;
    }

    /*!
     * \brief Returns true if some equation is used to specify a
     *        mortar coupling condition.
     */
436
    DUNE_DEPRECATED_MSG("hasMortarCoupling() is deprecated")
437
438
439
440
441
442
443
444
    bool hasMortarCoupling() const
    {
        for (int i = 0; i < numEq; ++i)
            if (boundaryInfo_[i].isMortarCoupling)
                return true;
        return false;
    }

445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
    /*!
     * \brief Returns true if an equation is used to specify a
     *        coupling condition.
     *
     * \param eqIdx The index of the equation
     *
     * This only used for correcting the pressure in the Stokes equation.
     */
    bool isCoupling(unsigned eqIdx) const
    {
        return false;
    }

    /*!
     * \brief Returns true if some equation is used to specify a
     *        coupling condition.
     */
    bool hasCoupling() const
    {
        return false;
    }

467
468
    /*!
     * \brief Returns the index of the equation which should be used
469
     *        for the Dirichlet condition of the pvIdx's primary
470
     *        variable.
471
472
473
     *
     * \param pvIdx The index of the primary variable which is be set
     *              by the Dirichlet condition.
474
475
     */
    unsigned dirichletToEqIndex(unsigned pvIdx) const
Katherina Baber's avatar
Katherina Baber committed
476
    { return pv2eqIdx_[pvIdx]; }
477
478
479

    /*!
     * \brief Returns the index of the primary variable which should
480
     *        be used for the Dirichlet condition given an equation
481
     *        index.
482
483
484
     *
     * \param eqIdx The index of the equation which is used to set
     *              the Dirichlet condition.
485
486
     */
    unsigned eqToDirichletIndex(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
487
    { return eq2pvIdx_[eqIdx]; }
488
489
490

private:
    // this is a bitfield structure!
491
    struct __attribute__((__packed__)) {
492
493
        unsigned char visited : 1;
        unsigned char isDirichlet : 1;
Andreas Lauser's avatar
Andreas Lauser committed
494
        unsigned char isNeumann : 1;
495
        unsigned char isOutflow : 1;
496
497
        unsigned char isCouplingInflow : 1;
        unsigned char isCouplingOutflow : 1;
498
        unsigned char isMortarCoupling : 1;
499
500
501
502
503
504
505
506
507
    } boundaryInfo_[numEq];

    unsigned char eq2pvIdx_[numEq];
    unsigned char pv2eqIdx_[numEq];
};

}

#endif