jacobianpattern.hh 7.8 KB
Newer Older
1
2
3
4
5
6
7
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
 *   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
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
 *   (at your option) any later version.                                     *
 *                                                                           *
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
 *   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/>.   *
 *****************************************************************************/
/*!
 * \file
 * \ingroup Assembly
 * \brief Helper function to generate Jacobian pattern for different discretization methods
 */
#ifndef DUMUX_JACOBIAN_PATTERN_HH
#define DUMUX_JACOBIAN_PATTERN_HH

#include <type_traits>
#include <dune/istl/matrixindexset.hh>
29
#include <dumux/discretization/method.hh>
30
31
32
33
34
35
36
37

namespace Dumux {

/*!
 * \ingroup Assembly
 * \brief Helper function to generate Jacobian pattern for the box method
 */
template<bool isImplicit, class GridGeometry,
38
         typename std::enable_if_t<(GridGeometry::discMethod == DiscretizationMethod::box), int> = 0>
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
Dune::MatrixIndexSet getJacobianPattern(const GridGeometry& gridGeometry)
{
    const auto numDofs = gridGeometry.numDofs();
    Dune::MatrixIndexSet pattern;
    pattern.resize(numDofs, numDofs);

    // matrix pattern for implicit Jacobians
    if (isImplicit)
    {
        static constexpr int dim = std::decay_t<decltype(gridGeometry.gridView())>::dimension;
        for (const auto& element : elements(gridGeometry.gridView()))
        {
            for (unsigned int vIdx = 0; vIdx < element.subEntities(dim); ++vIdx)
            {
                const auto globalI = gridGeometry.vertexMapper().subIndex(element, vIdx, dim);
54
                for (unsigned int vIdx2 = 0; vIdx2 < element.subEntities(dim); ++vIdx2)
55
56
57
                {
                    const auto globalJ = gridGeometry.vertexMapper().subIndex(element, vIdx2, dim);
                    pattern.add(globalI, globalJ);
58
59
60
61
62
63
64
65
66

                    if (gridGeometry.dofOnPeriodicBoundary(globalI) && globalI != globalJ)
                    {
                        const auto globalIP = gridGeometry.periodicallyMappedDof(globalI);
                        pattern.add(globalIP, globalI);
                        pattern.add(globalI, globalIP);
                        if (globalI > globalIP)
                            pattern.add(globalIP, globalJ);
                    }
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
                }
            }
        }
    }

    // matrix pattern for explicit Jacobians -> diagonal matrix
    else
    {
        for (unsigned int globalI = 0; globalI < numDofs; ++globalI)
            pattern.add(globalI, globalI);
    }

    return pattern;
}

/*!
 * \ingroup Assembly
 * \brief Helper function to generate Jacobian pattern for cell-centered methods
 */
template<bool isImplicit, class GridGeometry,
87
88
         typename std::enable_if_t<( (GridGeometry::discMethod == DiscretizationMethod::cctpfa)
                                     || (GridGeometry::discMethod == DiscretizationMethod::ccmpfa) ), int> = 0>
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
Dune::MatrixIndexSet getJacobianPattern(const GridGeometry& gridGeometry)
{
    const auto numDofs = gridGeometry.numDofs();
    Dune::MatrixIndexSet pattern;
    pattern.resize(numDofs, numDofs);

    // matrix pattern for implicit Jacobians
    if (isImplicit)
    {
        for (unsigned int globalI = 0; globalI < numDofs; ++globalI)
        {
            pattern.add(globalI, globalI);
            for (const auto& dataJ : gridGeometry.connectivityMap()[globalI])
                pattern.add(dataJ.globalJ, globalI);
        }
    }

    // matrix pattern for explicit Jacobians -> diagonal matrix
    else
    {
        for (unsigned int globalI = 0; globalI < numDofs; ++globalI)
            pattern.add(globalI, globalI);
    }

    return pattern;
}

116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
/*!
 * \ingroup Assembly
 * \brief Helper function to generate Jacobian pattern for the staggered method
 */
template<bool isImplicit, class GridGeometry,
         typename std::enable_if_t<( (GridGeometry::discMethod == DiscretizationMethod::staggered) ), int> = 0>
auto getJacobianPattern(const GridGeometry& gridGeometry)
{
    // resize the jacobian and the residual
    const auto numDofs = gridGeometry.numDofs();
    Dune::MatrixIndexSet pattern(numDofs, numDofs);

    const auto& connectivityMap = gridGeometry.connectivityMap();

    // evaluate the acutal pattern
    for (const auto& element : elements(gridGeometry.gridView()))
    {
        if(gridGeometry.isCellCenter())
        {
            // the global index of the element at hand
            static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
            const auto ccGlobalI = gridGeometry.elementMapper().index(element);
138
            pattern.add(ccGlobalI, ccGlobalI);
139
140
141
142
143
144
145
146
147
148
149
150
151
            for (auto&& ccGlobalJ : connectivityMap(cellCenterIdx, cellCenterIdx, ccGlobalI))
                pattern.add(ccGlobalI, ccGlobalJ);
        }
        else
        {
            static constexpr auto faceIdx = GridGeometry::faceIdx();
            auto fvGeometry = localView(gridGeometry);
            fvGeometry.bindElement(element);

            // loop over sub control faces
            for (auto&& scvf : scvfs(fvGeometry))
            {
                const auto faceGlobalI = scvf.dofIndex();
152
                pattern.add(faceGlobalI, faceGlobalI);
153
154
155
156
157
158
159
160
161
                for (auto&& faceGlobalJ : connectivityMap(faceIdx, faceIdx, scvf.index()))
                    pattern.add(faceGlobalI, faceGlobalJ);
            }
        }
    }

    return pattern;
}

162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
/*!
 * \ingroup Assembly
 * \brief Helper function to generate Jacobian pattern for finite element scheme
 */
template<class FEBasis>
Dune::MatrixIndexSet getFEJacobianPattern(const FEBasis& feBasis)
{
    const auto numDofs = feBasis.size();

    Dune::MatrixIndexSet pattern;
    pattern.resize(numDofs, numDofs);

    // matrix pattern for implicit Jacobians
    for (const auto& element : elements(feBasis.gridView()))
    {
        auto localView = feBasis.localView();
        localView.bind(element);

        const auto& finiteElement = localView.tree().finiteElement();
        const auto numLocalDofs = finiteElement.localBasis().size();
182
        for (std::size_t i = 0; i < numLocalDofs; i++)
183
184
        {
            const auto dofIdxI = localView.index(i);
185
            for (std::size_t j = 0; j < numLocalDofs; j++)
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
            {
                const auto dofIdxJ = localView.index(j);
                pattern.add(dofIdxI, dofIdxJ);
            }
        }
    }

    return pattern;
}

/*!
 * \ingroup Assembly
 * \brief Helper function to generate Jacobian pattern for finite element scheme
 * \note This interface is for compatibility with the other schemes. The pattern
 *       in fem is the same independent of the time discretization scheme.
 */
template<bool isImplicit, class GridGeometry,
         typename std::enable_if_t<(GridGeometry::discMethod == DiscretizationMethod::fem), int> = 0>
Dune::MatrixIndexSet getJacobianPattern(const GridGeometry& gridGeometry)
{ return getFEJacobianPattern(gridGeometry.feBasis()); }

207
208
209
} // namespace Dumux

#endif