gmshgriddatahandle.hh 8.12 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 InputOutput
 * \brief A data handle for commucating grid data for gmsh grids
 */
#ifndef DUMUX_GMSH_GRID_DATA_HANDLE_HH
#define DUMUX_GMSH_GRID_DATA_HANDLE_HH

#include <memory>
#include <algorithm>
Timo Koch's avatar
Timo Koch committed
29
#include <map>
30

31
#include <dune/common/parallel/communication.hh>
Timo Koch's avatar
Timo Koch committed
32
#include <dune/geometry/dimension.hh>
33
34
35
36
#include <dune/grid/common/partitionset.hh>
#include <dune/grid/common/datahandleif.hh>

// UGGrid specific includes
37
#if HAVE_DUNE_UGGRID
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
#include <dune/grid/uggrid.hh>
#endif

namespace Dumux {

/*!
 * \ingroup InputOutput
 * \brief A data handle for commucating grid data for gmsh grids
 */
template<class Grid, class GridFactory, class Data>
struct GmshGridDataHandle : public Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>
{
    using GridView = typename Grid::LevelGridView;

    GmshGridDataHandle(const Grid& grid, const GridFactory& gridFactory, Data& elementMarkers, Data& boundaryMarkers, Data& faceMarkers)
    : gridView_(grid.levelGridView(0))
    , idSet_(grid.localIdSet())
    , elementMarkers_(elementMarkers)
    , boundaryMarkers_(boundaryMarkers)
    , faceMarkers_(faceMarkers)
    {
        const auto& indexSet = gridView_.indexSet();

        for (const auto& element : elements(gridView_, Dune::Partitions::interior))
           std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);

        for (const auto& face : entities(gridView_, Dune::Codim<1>()))
           std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);
    }

    ~GmshGridDataHandle()
    {
        const auto& indexSet = gridView_.indexSet();

        elementMarkers_.resize(indexSet.size(0));
        for (const auto& element : elements(gridView_))
           std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);

        faceMarkers_.resize(indexSet.size(1));
        for (const auto& face : entities(gridView_, Dune::Codim<1>()))
           std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);

        boundaryMarkers_.resize(gridView_.grid().numBoundarySegments(), 0);
        for (const auto& element : elements(gridView_.grid().leafGridView()))
        {
            for (const auto& intersection : intersections(gridView_.grid().leafGridView(), element))
            {
                if (intersection.boundary())
                {
                    const auto marker = faceMarkers_[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))];
                    boundaryMarkers_[intersection.boundarySegmentIndex()] = marker;
                }
            }
       }
    }

    Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>& interface()
    { return *this; }

    bool contains (int dim, int codim) const
    { return codim == 0 || codim == 1; }

100
101
    //! returns true if size per entity of given dim and codim is a constant
    bool fixedSize(int dim, int codim) const
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
    { return true; }

    template<class EntityType>
    std::size_t size (const EntityType& e) const
    { return 1; }

    template<class MessageBufferImp, class EntityType>
    void gather (MessageBufferImp& buff, const EntityType& e) const
    { buff.write(data_[idSet_.id(e)]); }

    template<class MessageBufferImp, class EntityType>
    void scatter (MessageBufferImp& buff, const EntityType& e, std::size_t n)
    { buff.read(data_[idSet_.id(e)]); }

private:
    using IdSet = typename Grid::LocalIdSet;

    const GridView gridView_;
    const IdSet &idSet_;
    Data& elementMarkers_;
    Data& boundaryMarkers_;
    Data& faceMarkers_;
    mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
};

127
#if HAVE_DUNE_UGGRID
128
129
130
131
132
133
134
135
136
137
138
139

/*!
 * \ingroup InputOutput
 * \brief A data handle for commucating grid data for gmsh grids (specialization for UGGrid)
 */
template<class GridFactory, class Data, int dimgrid>
struct GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>
: public Dune::CommDataHandleIF<GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>, typename Data::value_type>
{
    using Grid = Dune::UGGrid<dimgrid>;
    using GridView = typename Grid::LevelGridView;

140
    GmshGridDataHandle(const Grid& grid, const GridFactory& gridFactory, Data& elementMarkers, Data& boundaryMarkers)
141
142
143
    : gridView_(grid.levelGridView(0))
    , idSet_(grid.localIdSet())
    , elementMarkers_(elementMarkers)
144
    , boundaryMarkers_(boundaryMarkers)
145
146
147
    {
        for (const auto& element : elements(gridView_, Dune::Partitions::interior))
           std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162

        // Depending on the Dune version, the boundary markers are present on
        // all processes (<= 2.6) or on the root process only (>= 2.7). Try to
        // handle this in a flexible way: determine if the minimum size over
        // all processes of the boundary markers vector is zero. If yes, assume
        // that the root process contains all markers and broadcast them.
        auto bmSizeMin = boundaryMarkers_.size();
        Dune::MPIHelper::getCollectiveCommunication().min(&bmSizeMin, 1);
        if (bmSizeMin == 0)
        {
            auto bmSize = boundaryMarkers_.size();
            Dune::MPIHelper::getCollectiveCommunication().broadcast(&bmSize, 1, 0);
            boundaryMarkers_.resize(bmSize);
            Dune::MPIHelper::getCollectiveCommunication().broadcast(&boundaryMarkers_.front(), bmSize, 0);
        }
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
    }

    ~GmshGridDataHandle()
    {
        const auto& indexSet = gridView_.indexSet();
        elementMarkers_.resize(indexSet.size(0));
        for (const auto& element : elements(gridView_))
           std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);
    }

    Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>& interface()
    { return *this; }

    bool contains (int dim, int codim) const
    { return codim == 0 || codim == 1; }

179
180
    //! returns true if size per entity of given dim and codim is a constant
    bool fixedSize(int dim, int codim) const
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
    { return true; }

    template<class EntityType>
    std::size_t size (const EntityType& e) const
    { return 1; }

    template<class MessageBufferImp, class EntityType>
    void gather (MessageBufferImp& buff, const EntityType& e) const
    { buff.write(data_[idSet_.id(e)]); }

    template<class MessageBufferImp, class EntityType>
    void scatter (MessageBufferImp& buff, const EntityType& e, std::size_t n)
    { buff.read(data_[idSet_.id(e)]); }

private:
    using IdSet = typename Grid::LocalIdSet;

    const GridView gridView_;
    const IdSet &idSet_;
    Data& elementMarkers_;
201
    Data& boundaryMarkers_;
202
203
204
    mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
};

205
#endif // HAVE_DUNE_UGGRID
206
207
208
209

} // namespace Dumux

#endif