boundarytypes.hh 13.1 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
template <int numEq>
class BoundaryTypes
{
public:
    BoundaryTypes()
    { reset(); }

    /*!
43
     * \brief Reset the boundary types for all equations.
44
45
     *
     * After this method no equations will be disabled and neither
46
     * Neumann nor Dirichlet conditions will be evaluated. This
47
     * corresponds to a Neumann zero boundary.
48
49
50
     */
    void reset()
    {
51
52
53
        for (int eqIdx=0; eqIdx < numEq; ++eqIdx)
        {
            resetEq(eqIdx);
54
        }
55
56
    }

57
58
59
60
61
    /*!
     * \brief Reset the boundary types for one equation.
     */
    void resetEq(int eqIdx)
    {
62
          boundaryInfo_[eqIdx].visited = false;
63

64
65
66
67
68
69
          boundaryInfo_[eqIdx].isDirichlet = false;
          boundaryInfo_[eqIdx].isNeumann = false;
          boundaryInfo_[eqIdx].isOutflow = false;
          boundaryInfo_[eqIdx].isCouplingDirichlet = false;
          boundaryInfo_[eqIdx].isCouplingNeumann = false;
          boundaryInfo_[eqIdx].isCouplingMortar = false;
70
71
72
73
74

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

75
    /*!
76
77
78
79
     * \brief Returns true if the boundary types for a given equation
     *        has been specified.
     *
     * \param eqIdx The index of the equation
80
81
     */
    bool isSet(int eqIdx) const
82
    { return boundaryInfo_[eqIdx].visited; }
83
84
85
86

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

    /*!
100
     * \brief Set all boundary conditions to Neumann.
101
102
103
     */
    void setAllNeumann()
    {
104
105
106
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
        {
            setNeumann(eqIdx);
107
108
109
110
        }
    }

    /*!
111
     * \brief Set all boundary conditions to Dirichlet.
112
113
114
     */
    void setAllDirichlet()
    {
115
116
117
        for (int eqIdx = 0; eqIdx < numEq; ++ eqIdx)
        {
            setDirichlet(eqIdx);
118
119
120
        }
    }

121
    /*!
122
     * \brief Set all boundary conditions to Neumann.
123
124
125
     */
    void setAllOutflow()
    {
126
127
128
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
        {
            setOutflow(eqIdx);
129
130
131
        }
    }

132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
    /*!
     * \brief Set all boundary conditions to Dirichlet-like coupling
     */
    void setAllCouplingDirichlet()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
        {
            setCouplingDirichlet(eqIdx);
        }
    }

    /*!
     * \brief Set all boundary conditions to Neumann-like coupling.
     */
    void setAllCouplingNeumann()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
        {
            setCouplingNeumann(eqIdx);
        }
    }
    /*!
     * \brief Set all boundary conditions to mortar coupling.
     */
    void setAllCouplingMortar()
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
        {
            setCouplingDirichlet(eqIdx);
        }
    }

164
    /*!
165
     * \brief Set a Neumann boundary condition for a single a single
166
     *        equation.
167
168
     *
     * \param eqIdx The index of the equation
169
170
171
     */
    void setNeumann(int eqIdx)
    {
172
        resetEq(eqIdx);
173
174
        boundaryInfo_[eqIdx].visited = true;
        boundaryInfo_[eqIdx].isNeumann = true;
175
176
177
178
179

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

    /*!
180
     * \brief Set a Dirichlet boundary condition for a single primary
181
182
     *        variable
     *
183
184
185
186
     * \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
187
     */
188
    void setDirichlet(int pvIdx, int eqIdx)
189
    {
190
        resetEq(eqIdx);
191
192
        boundaryInfo_[eqIdx].visited = true;
        boundaryInfo_[eqIdx].isDirichlet = true;
193
194
195
196
197
198
199
200

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

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

201
    /*!
202
     * \brief Set a Neumann boundary condition for a single a single
203
     *        equation.
204
205
206
     *
     * \param eqIdx The index of the equation on which the outflow
     *              condition applies.
207
208
209
     */
    void setOutflow(int eqIdx)
    {
210
        resetEq(eqIdx);
211
212
        boundaryInfo_[eqIdx].visited = true;
        boundaryInfo_[eqIdx].isOutflow = true;
213
214
215
216
217
218
219
220
221
222
223

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

    /*!
     * \brief Set a boundary condition for a single equation to
     *        a Dirichlet-like coupling condition.
     */
    void setCouplingDirichlet(int eqIdx)
    {
        resetEq(eqIdx);
224
225
        boundaryInfo_[eqIdx].visited = true;
        boundaryInfo_[eqIdx].isCouplingDirichlet = true;
226
227
228
229
230
231
232
233
234
235
236

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

    /*!
     * \brief Set a boundary condition for a single equation to
     *        a Neumann-like coupling condition.
     */
    void setCouplingNeumann(int eqIdx)
    {
        resetEq(eqIdx);
237
238
        boundaryInfo_[eqIdx].visited = true;
        boundaryInfo_[eqIdx].isCouplingNeumann = true;
239
240
241
242
243
244
245
246
247
248
249

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

    /*!
     * \brief Set a boundary condition for a single equation to
     *        a mortar coupling condition.
     */
    void setCouplingMortar(int eqIdx)
    {
        resetEq(eqIdx);
250
251
        boundaryInfo_[eqIdx].visited = true;
        boundaryInfo_[eqIdx].isCouplingMortar = true;
252
253
254
255

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

256
    /*!
257
     * \brief Set a Dirichlet boundary condition for a single primary
258
     *        variable.
259
     *
260
261
     * Depending on the discretization, setting the Dirichlet condition
     * will replace the balance equation with index equal to pvIdx.
262
     *
263
264
     * \param pvIdx The index of the primary variable inside a
     *              PrimaryVariables object.
265
     */
266
    void setDirichlet(int pvIdx)
267
    {
268
        setDirichlet(pvIdx, pvIdx);
269
270
271
272
    }

    /*!
     * \brief Returns true if an equation is used to specify a
273
     *        Dirichlet condition.
274
275
     *
     * \param eqIdx The index of the equation
276
277
     */
    bool isDirichlet(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
278
    { return boundaryInfo_[eqIdx].isDirichlet; }
279

280
281
282
283
284
285
286
287
288
289
290
291
292
293
    /*!
     * \brief Returns true if an equation is used to specify a
     *        Dirichlet condition.
     *
     * \param eqIdx The index of the equation
     */
    bool hasOnlyDirichlet() const
    {
        return std::all_of(boundaryInfo_.begin(),
                           boundaryInfo_.end(),
                           [](const BoundaryInfo& b){ return b.isDirichlet; }
                           );
    }

294
295
    /*!
     * \brief Returns true if some equation is used to specify a
296
     *        Dirichlet condition.
297
298
299
300
301
302
303
     */
    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
304
    }
305
306
307

    /*!
     * \brief Returns true if an equation is used to specify a
308
     *        Neumann condition.
309
310
     *
     * \param eqIdx The index of the equation
311
312
     */
    bool isNeumann(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
313
    { return boundaryInfo_[eqIdx].isNeumann; }
314
315
316

    /*!
     * \brief Returns true if some equation is used to specify a
317
     *        Neumann condition.
318
319
320
321
     */
    bool hasNeumann() const
    {
        for (int i = 0; i < numEq; ++i)
Andreas Lauser's avatar
Andreas Lauser committed
322
            if (boundaryInfo_[i].isNeumann)
323
324
                return true;
        return false;
Katherina Baber's avatar
Katherina Baber committed
325
    }
326

327
328
329
    /*!
     * \brief Returns true if an equation is used to specify an
     *        outflow condition.
330
331
     *
     * \param eqIdx The index of the equation
332
333
     */
    bool isOutflow(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
334
    { return boundaryInfo_[eqIdx].isOutflow; }
335
336
337
338
339
340
341
342
343
344
345

    /*!
     * \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
346
    }
347

348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
    /*!
     * \brief Returns true if an equation is used to specify an
     *        Dirichlet coupling condition.
     *
     * \param eqIdx The index of the equation
     */
    bool isCouplingDirichlet(unsigned eqIdx) const
    { return boundaryInfo_[eqIdx].isCouplingDirichlet; }

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

    /*!
     * \brief Returns true if an equation is used to specify an
     *        Neumann coupling condition.
     *
     * \param eqIdx The index of the equation
     */
    bool isCouplingNeumann(unsigned eqIdx) const
    { return boundaryInfo_[eqIdx].isCouplingNeumann; }

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

390
391
    /*!
     * \brief Returns true if an equation is used to specify a
392
     *        Mortar coupling condition.
393
394
     *
     * \param eqIdx The index of the equation
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
     */
    bool isCouplingMortar(unsigned eqIdx) const
    {
        return boundaryInfo_[eqIdx].isCouplingMortar;
    }

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

    /*!
     * \brief Returns true if an equation is used to specify a
     *        coupling condition.
416
     *
417
     * \param eqIdx The index of the equation
418
419
420
     */
    bool isCoupling(unsigned eqIdx) const
    {
421
422
423
        return boundaryInfo_[eqIdx].isCouplingDirichlet
               || boundaryInfo_[eqIdx].isCouplingNeumann
               || boundaryInfo_[eqIdx].isCouplingMortar;
424
425
426
427
428
429
430
431
    }

    /*!
     * \brief Returns true if some equation is used to specify a
     *        coupling condition.
     */
    bool hasCoupling() const
    {
432
433
434
        for (int i = 0; i < numEq; ++i)
            if (isCoupling(i))
                return true;
435
436
437
        return false;
    }

438
439
    /*!
     * \brief Returns the index of the equation which should be used
440
     *        for the Dirichlet condition of the pvIdx's primary
441
     *        variable.
442
443
444
     *
     * \param pvIdx The index of the primary variable which is be set
     *              by the Dirichlet condition.
445
446
     */
    unsigned dirichletToEqIndex(unsigned pvIdx) const
Katherina Baber's avatar
Katherina Baber committed
447
    { return pv2eqIdx_[pvIdx]; }
448
449
450

    /*!
     * \brief Returns the index of the primary variable which should
451
     *        be used for the Dirichlet condition given an equation
452
     *        index.
453
454
455
     *
     * \param eqIdx The index of the equation which is used to set
     *              the Dirichlet condition.
456
457
     */
    unsigned eqToDirichletIndex(unsigned eqIdx) const
Katherina Baber's avatar
Katherina Baber committed
458
    { return eq2pvIdx_[eqIdx]; }
459

460
protected:
461
462
463
464
465
466
467
468
469
470
471
472
    struct BoundaryInfo {
        bool visited;
        bool isDirichlet;
        bool isNeumann;
        bool isOutflow;
        bool isCouplingDirichlet;
        bool isCouplingNeumann;
        bool isCouplingMortar;
    };

    std::array<BoundaryInfo, numEq> boundaryInfo_;
    std::array<unsigned int, numEq> eq2pvIdx_, pv2eqIdx_;
473
474
475
476
477
};

}

#endif