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
31
32
 *****************************************************************************/
/*!
 * \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
{

33
/*!
34
 * \ingroup BC
35
36
 * \brief Class to specify the type of a boundary.
 */
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
template <int numEq>
class BoundaryTypes
{
public:
    BoundaryTypes()
    { reset(); }

    /*!
     * \brief Reset the boundary types.
     *
     * After this method no equations will be disabled and neither
     * neumann nor dirichlet conditions will be evaluated. This
     * corrosponds to a neumann zero boundary.
     */
    void reset()
    {
        for (int i=0; i < numEq; ++i) {
            boundaryInfo_[i].visited = 0;

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

            eq2pvIdx_[i] = i;
            pv2eqIdx_[i] = i;
65
        }
66
67
68
    }

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

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

    /*!
     * \brief Set all boundary conditions to neuman.
     */
    void setAllNeumann()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            boundaryInfo_[eqIdx].visited = 1;
            boundaryInfo_[eqIdx].isDirichlet = 0;
Andreas Lauser's avatar
Andreas Lauser committed
100
            boundaryInfo_[eqIdx].isNeumann = 1;
101
            boundaryInfo_[eqIdx].isOutflow = 0;
102
103
            boundaryInfo_[eqIdx].isCouplingInflow = 0;
            boundaryInfo_[eqIdx].isCouplingOutflow = 0;
104
            boundaryInfo_[eqIdx].isMortarCoupling = 0;
105
106
107
108
109
110
111
112
113
114
115
116
117

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

    /*!
     * \brief Set all boundary conditions to dirichlet.
     */
    void setAllDirichlet()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
            boundaryInfo_[eqIdx].visited = 1;
            boundaryInfo_[eqIdx].isDirichlet = 1;
Andreas Lauser's avatar
Andreas Lauser committed
118
            boundaryInfo_[eqIdx].isNeumann = 0;
119
            boundaryInfo_[eqIdx].isOutflow = 0;
120
121
            boundaryInfo_[eqIdx].isCouplingInflow = 0;
            boundaryInfo_[eqIdx].isCouplingOutflow = 0;
122
            boundaryInfo_[eqIdx].isMortarCoupling = 0;
123
124
125
126
127
128
129
130

            eq2pvIdx_[eqIdx] = eqIdx;
            pv2eqIdx_[eqIdx] = eqIdx;

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

131
132
133
134
135
136
137
138
139
140
    /*!
     * \brief Set all boundary conditions to neuman.
     */
    void setAllOutflow()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            boundaryInfo_[eqIdx].visited = 1;
            boundaryInfo_[eqIdx].isDirichlet = 0;
            boundaryInfo_[eqIdx].isNeumann = 0;
            boundaryInfo_[eqIdx].isOutflow = 1;
141
142
            boundaryInfo_[eqIdx].isCouplingInflow = 0;
            boundaryInfo_[eqIdx].isCouplingOutflow = 0;
143
            boundaryInfo_[eqIdx].isMortarCoupling = 0;
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160

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

    /*!
     * \brief Set all boundary conditions to coupling inflow.
     */
    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;
161
            boundaryInfo_[eqIdx].isMortarCoupling = 0;
162
163
164
165
166
167
168
169
170
171
172
173
174
175

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

    /*!
     * \brief Set all boundary conditions to coupling outflow.
     */
    void setAllCouplingOutflow()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            boundaryInfo_[eqIdx].visited = 1;
            boundaryInfo_[eqIdx].isDirichlet = 0;
            boundaryInfo_[eqIdx].isNeumann = 0;
176
            boundaryInfo_[eqIdx].isOutflow = 0;
177
178
            boundaryInfo_[eqIdx].isCouplingInflow = 0;
            boundaryInfo_[eqIdx].isCouplingOutflow = 1;
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
            boundaryInfo_[eqIdx].isMortarCoupling = 0;

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

    /*!
     * \brief Set all boundary conditions to mortar coupling.
     */
    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;
198
199
200
201
202

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

203
204
205
    /*!
     * \brief Set a neumann boundary condition for a single a single
     *        equation.
206
207
     *
     * \param eqIdx The index of the equation
208
209
210
211
     */
    void setNeumann(int eqIdx)
    {
        boundaryInfo_[eqIdx].visited = 1;
212
213
214
        boundaryInfo_[eqIdx].isDirichlet = 0;
        boundaryInfo_[eqIdx].isNeumann = 1;
        boundaryInfo_[eqIdx].isOutflow = 0;
215
216
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 0;
217
        boundaryInfo_[eqIdx].isMortarCoupling = 0;
218
219
220
221
222
223
224
225

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

    /*!
     * \brief Set a dirichlet boundary condition for a single primary
     *        variable
     *
226
227
228
229
     * \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
230
     */
231
    void setDirichlet(int pvIdx, int eqIdx)
232
233
234
    {
        boundaryInfo_[eqIdx].visited = 1;
        boundaryInfo_[eqIdx].isDirichlet = 1;
Andreas Lauser's avatar
Andreas Lauser committed
235
        boundaryInfo_[eqIdx].isNeumann = 0;
236
        boundaryInfo_[eqIdx].isOutflow = 0;
237
238
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 0;
239
        boundaryInfo_[eqIdx].isMortarCoupling = 0;
240
241
242
243
244
245
246
247

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

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

248
249
250
    /*!
     * \brief Set a neumann boundary condition for a single a single
     *        equation.
251
252
253
     *
     * \param eqIdx The index of the equation on which the outflow
     *              condition applies.
254
255
256
257
258
259
260
     */
    void setOutflow(int eqIdx)
    {
        boundaryInfo_[eqIdx].visited = 1;
        boundaryInfo_[eqIdx].isDirichlet = 0;
        boundaryInfo_[eqIdx].isNeumann = 0;
        boundaryInfo_[eqIdx].isOutflow = 1;
261
262
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 0;
263
        boundaryInfo_[eqIdx].isMortarCoupling = 0;
264
265
266
267

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

268
269
270
271
272
273
274
275
276
277
278
    /*!
     * \brief Set a boundary condition for a single equation to coupling inflow.
     */
    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;
279
        boundaryInfo_[eqIdx].isMortarCoupling = 0;
280
281
282
283
284
285
286
287
288
289
290
291
292

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

    /*!
     * \brief Set a boundary condition for a single equation to coupling outflow.
     */
    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
293
294
        boundaryInfo_[eqIdx].isCouplingInflow = 0;
        boundaryInfo_[eqIdx].isCouplingOutflow = 1;
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
        boundaryInfo_[eqIdx].isMortarCoupling = 0;

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

    /*!
     * \brief Set a boundary condition for a single equation to mortar coupling.
     */
    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;
312
313
314
315

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

316
317
    /*!
     * \brief Set a dirichlet boundary condition for a single primary
318
     *        variable.
319
     *
320
321
     * Depending on the discretization, setting the Dirichlet condition
     * will replace the balance equation with index equal to pvIdx.
322
     *
323
324
     * \param pvIdx The index of the primary variable inside a
     *              PrimaryVariables object.
325
     */
326
    void setDirichlet(int pvIdx)
327
    {
328
        setDirichlet(pvIdx, pvIdx);
329
330
331
332
333
    }

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

    /*!
     * \brief Returns true if some equation is used to specify a
     *        dirichlet condition.
     */
    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
350
    }
351
352
353
354

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

    /*!
     * \brief Returns true if some equation is used to specify a
     *        neumann condition.
     */
    bool hasNeumann() const
    {
        for (int i = 0; i < numEq; ++i)
Andreas Lauser's avatar
Andreas Lauser committed
368
            if (boundaryInfo_[i].isNeumann)
369
370
                return true;
        return false;
Katherina Baber's avatar
Katherina Baber committed
371
    }
372

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

    /*!
     * \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
392
    }
393

394
395
396
    /*!
     * \brief Returns true if an equation is used to specify an
     *        inflow coupling condition.
397
398
     *
     * \param eqIdx The index of the equation
399
400
     */
    bool isCouplingInflow(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
401
    { return boundaryInfo_[eqIdx].isCouplingInflow; }
402
403
404
405
406
407
408
409
410
411
412

    /*!
     * \brief Returns true if some equation is used to specify an
     *        inflow coupling condition.
     */
    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
413
    }
414
415
416
417

    /*!
     * \brief Returns true if an equation is used to specify an
     *        outflow coupling condition.
418
419
     *
     * \param eqIdx The index of the equation
420
421
     */
    bool isCouplingOutflow(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
422
    { return boundaryInfo_[eqIdx].isCouplingOutflow; }
423
424
425
426
427
428
429
430
431
432
433

    /*!
     * \brief Returns true if some equation is used to specify an
     *        outflow coupling condition.
     */
    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
434
    }
435

436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
    /*!
     * \brief Returns true if an equation is used to specify a
     *        mortar coupling condition.
     *
     * \param eqIdx The index of the equation
     */
    bool isMortarCoupling(unsigned eqIdx) const
    {
        return boundaryInfo_[eqIdx].isMortarCoupling;
    }

    /*!
     * \brief Returns true if some equation is used to specify a
     *        mortar coupling condition.
     */
    bool hasMortarCoupling() const
    {
        for (int i = 0; i < numEq; ++i)
            if (boundaryInfo_[i].isMortarCoupling)
                return true;
        return false;
    }

459
460
    /*!
     * \brief Returns the index of the equation which should be used
461
     *        for the Dirichlet condition of the pvIdx's primary
462
     *        variable.
463
464
465
     *
     * \param pvIdx The index of the primary variable which is be set
     *              by the Dirichlet condition.
466
467
     */
    unsigned dirichletToEqIndex(unsigned pvIdx) const
Katherina Baber's avatar
Katherina Baber committed
468
    { return pv2eqIdx_[pvIdx]; }
469
470
471

    /*!
     * \brief Returns the index of the primary variable which should
472
     *        be used for the Dirichlet condition given an equation
473
     *        index.
474
475
476
     *
     * \param eqIdx The index of the equation which is used to set
     *              the Dirichlet condition.
477
478
     */
    unsigned eqToDirichletIndex(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
479
    { return eq2pvIdx_[eqIdx]; }
480
481
482

private:
    // this is a bitfield structure!
483
    struct __attribute__((__packed__)) {
484
485
        unsigned char visited : 1;
        unsigned char isDirichlet : 1;
Andreas Lauser's avatar
Andreas Lauser committed
486
        unsigned char isNeumann : 1;
487
        unsigned char isOutflow : 1;
488
489
        unsigned char isCouplingInflow : 1;
        unsigned char isCouplingOutflow : 1;
490
        unsigned char isMortarCoupling : 1;
491
492
493
494
495
496
497
498
499
    } boundaryInfo_[numEq];

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

}

#endif