test_staggeredfvgeometry.cc 6.47 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
} // end namespace Dumux
#endif

86
87
int main (int argc, char *argv[]) try
{
88
89
    using namespace Dumux;

90
91
92
    // maybe initialize mpi
    Dune::MPIHelper::instance(argc, argv);

93
94
95
    // parse command line arguments and input file
    Parameters::init(argc, argv);

96
97
    std::cout << "Checking the FVGeometries, SCVs and SCV faces" << std::endl;

98
    using Grid = Dune::YaspGrid<2>;
99

100
101
    constexpr int dim = Grid::dimension;
    constexpr int dimworld = Grid::dimensionworld;
102

103
    using GlobalPosition = typename Dune::FieldVector<Grid::ctype, dimworld>;
104
    using GridGeometry = StaggeredFVGridGeometry<typename Grid::LeafGridView, /*enable caching=*/ true,
105
                                                   TestFVGGTraits<typename Grid::LeafGridView> >;
106
    using FVElementGeometry = typename GridGeometry::LocalView;
107
108
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
109
110
111
112

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

117
118
    GridGeometry gridGeometry(leafGridView);
    gridGeometry.update();
119
120
121
122

    // iterate over elements. For every element get fv geometry and loop over scvs and scvfaces
    for (const auto& element : elements(leafGridView))
    {
123
        auto eIdx = gridGeometry.elementMapper().index(element);
124
        std::cout << std::endl << "Checking fvGeometry of element " << eIdx << std::endl;
125
        auto fvGeometry = localView(gridGeometry);
126
127
128
        fvGeometry.bind(element);

        auto range = scvs(fvGeometry);
129
        Detail::NoopFunctor<SubControlVolume> op;
130
        if(0 != testForwardIterator(range.begin(), range.end(), op))
131
132
133
134
            DUNE_THROW(Dune::Exception, "Iterator does not fulfill the forward iterator concept");

        for (auto&& scv : scvs(fvGeometry))
        {
135
            std::cout << "-- scv " << scv.localDofIndex() << " center at: " << scv.center() << " , volume: " << scv.volume()  << std::endl;
136
137
138
        }

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

143
        std::size_t boundaryCount = 0;
144
145
        for (auto&& scvf : scvfs(fvGeometry))
        {
146
            std::cout << "-- scvf " << scvf.index() << " ip at: " << scvf.ipGlobal() << " normal: " << scvf.unitOuterNormal();
147
148
149
150
151
            if (scvf.boundary())
            {
                ++boundaryCount;
                std::cout << " (on boundary).";
            }
152
153
            std::cout << std::endl;
        }
154
155
156
157

        if ((boundaryCount>0) != fvGeometry.hasBoundaryScvf())
            DUNE_THROW(Dune::InvalidStateException, "fvGeometry.hasBoundaryScvf() reports " << fvGeometry.hasBoundaryScvf()
                            << " but the number of boundary scvfs is " << boundaryCount);
158
159
160
161
162
    }
}
// //////////////////////////////////
//   Error handler
// /////////////////////////////////
163
catch (Dune::Exception& e) {
164
165
166
167

    std::cout << e << std::endl;
    return 1;
}