test_staggeredfvgeometry.cc 7.15 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
 *   (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
 * \brief Test for finite volume element geometry, sub control volume, and sub
          control volume faces
 */
#include <config.h>

#include <iostream>
#include <utility>
28
#include <iomanip>
29
30
31
32

#include <dune/common/test/iteratortest.hh>
#include <dune/grid/utility/structuredgridfactory.hh>
#include <dune/grid/yaspgrid.hh>
33
#include <dumux/common/parameters.hh>
34

35
36
37
38
39
40
41
42
43
44
45
46
#include <dumux/common/intersectionmapper.hh>
#include <dumux/common/defaultmappertraits.hh>
#include <dumux/discretization/cellcentered/subcontrolvolume.hh>
#include <dumux/discretization/staggered/fvelementgeometry.hh>
#include <dumux/discretization/staggered/fvgridgeometry.hh>
#include <dumux/discretization/staggered/subcontrolvolumeface.hh>

#ifndef DOXYGEN
namespace Dumux {
namespace Detail {
template<class T>
class NoopFunctor {
47
public:
48
49
  NoopFunctor() {}
  void operator()(const T& t){}
50
};
51
} // end namespace Detail
52

53
54
55
//! the fv grid geometry traits for this test
template<class GridView>
struct TestFVGGTraits : public DefaultMapperTraits<GridView>
56
{
57
58
59
60
    using SubControlVolume = CCSubControlVolume<GridView>;
    using SubControlVolumeFace = StaggeredSubControlVolumeFace<GridView>;
    using IntersectionMapper = ConformingGridIntersectionMapper<GridView>;
    using GeometryHelper = BaseStaggeredGeometryHelper<GridView>;
61

62
63
64
65
66
67
    struct DofTypeIndices
    {
        using CellCenterIdx = Dune::index_constant<0>;
        using FaceIdx = Dune::index_constant<1>;
    };

68
69
    //! Dummy connectivity map, required by GridGeometry
    template<class GridGeometry>
70
    struct MockConnectivityMap
71
    {
72
        void update(const GridGeometry& gridGeometry) {}
73
74
        void setStencilOrder(const int order) {}
    };
75

76
77
    template<class GridGeometry>
    using ConnectivityMap = MockConnectivityMap<GridGeometry>;
78

79
80
    template<class GridGeometry, bool enableCache>
    using LocalView = StaggeredFVElementGeometry<GridGeometry, enableCache>;
81
82
83
84
85
86
87
88
89

    struct PublicTraits
    {
        using CellSubControlVolume = SubControlVolume;
        using CellSubControlVolumeFace = SubControlVolumeFace;
        using FaceSubControlVolume = SubControlVolume;
        using FaceLateralSubControlVolumeFace = SubControlVolumeFace;
        using FaceFrontalSubControlVolumeFace = SubControlVolumeFace;
    };
90
91
};

92
93
94
} // end namespace Dumux
#endif

95
int main (int argc, char *argv[])
96
{
97
98
    using namespace Dumux;

99
100
101
    // maybe initialize mpi
    Dune::MPIHelper::instance(argc, argv);

102
103
104
    // parse command line arguments and input file
    Parameters::init(argc, argv);

105
106
    std::cout << "Checking the FVGeometries, SCVs and SCV faces" << std::endl;

107
    using Grid = Dune::YaspGrid<2>;
108

109
110
    constexpr int dim = Grid::dimension;
    constexpr int dimworld = Grid::dimensionworld;
111

112
    using GlobalPosition = typename Dune::FieldVector<Grid::ctype, dimworld>;
113
    using GridGeometry = StaggeredFVGridGeometry<typename Grid::LeafGridView, /*enable caching=*/ true,
114
                                                   TestFVGGTraits<typename Grid::LeafGridView> >;
115
    using FVElementGeometry = typename GridGeometry::LocalView;
116
117
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
118
119
120
121

    // make a grid
    GlobalPosition lower(0.0);
    GlobalPosition upper(1.0);
122
    std::array<unsigned int, dim> els{{2, 4}};
123
124
125
    std::shared_ptr<Grid> grid = Dune::StructuredGridFactory<Grid>::createCubeGrid(lower, upper, els);
    auto leafGridView = grid->leafGridView();

126
127
    GridGeometry gridGeometry(leafGridView);
    gridGeometry.update();
128
129
130
131

    // iterate over elements. For every element get fv geometry and loop over scvs and scvfaces
    for (const auto& element : elements(leafGridView))
    {
132
        auto eIdx = gridGeometry.elementMapper().index(element);
133
        std::cout << std::endl << "Checking fvGeometry of element " << eIdx << std::endl;
134
        auto fvGeometry = localView(gridGeometry);
135
136
137
138
139

        // bind the local view to the element
        if (fvGeometry.isBound())
            DUNE_THROW(Dune::Exception, "Local view should not be bound at this point");

140
141
        fvGeometry.bind(element);

142
143
144
145
146
147
148
149
        if (!fvGeometry.isBound())
            DUNE_THROW(Dune::Exception, "Local view should be bound at this point");

        // make sure the bound element fits
        auto eIdxBound = gridGeometry.elementMapper().index(fvGeometry.element());
        if (eIdx != eIdxBound)
            DUNE_THROW(Dune::Exception, "Bound element index does not match");

150
        auto range = scvs(fvGeometry);
151
        Detail::NoopFunctor<SubControlVolume> op;
152
        if(0 != testForwardIterator(range.begin(), range.end(), op))
153
154
155
156
            DUNE_THROW(Dune::Exception, "Iterator does not fulfill the forward iterator concept");

        for (auto&& scv : scvs(fvGeometry))
        {
157
            std::cout << "-- scv " << scv.localDofIndex() << " center at: " << scv.center() << " , volume: " << scv.volume()  << std::endl;
158
159
160
        }

        auto range2 = scvfs(fvGeometry);
161
        Detail::NoopFunctor<SubControlVolumeFace> op2;
162
        if(0 != testForwardIterator(range2.begin(), range2.end(), op2))
163
164
            DUNE_THROW(Dune::Exception, "Iterator does not fulfill the forward iterator concept");

165
        std::size_t boundaryCount = 0;
166
167
        for (auto&& scvf : scvfs(fvGeometry))
        {
168
            std::cout << "-- scvf " << scvf.index() << " ip at: " << scvf.ipGlobal() << " normal: " << scvf.unitOuterNormal();
169
170
171
172
173
            if (scvf.boundary())
            {
                ++boundaryCount;
                std::cout << " (on boundary).";
            }
174
175
            std::cout << std::endl;
        }
176
177
178
179

        if ((boundaryCount>0) != fvGeometry.hasBoundaryScvf())
            DUNE_THROW(Dune::InvalidStateException, "fvGeometry.hasBoundaryScvf() reports " << fvGeometry.hasBoundaryScvf()
                            << " but the number of boundary scvfs is " << boundaryCount);
180
181
    }
}