test_staggeredfvgeometry.cc 5.93 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
// -*- 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    *
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (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
33

#include <dune/common/test/iteratortest.hh>
#include <dune/grid/utility/structuredgridfactory.hh>
#include <dune/grid/yaspgrid.hh>

34
35
36
37
38
39
40
41
42
43
44
45
#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 {
46
public:
47
48
  NoopFunctor() {}
  void operator()(const T& t){}
49
};
50
} // end namespace Detail
51

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

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

67
68
69
70
    //! Dummy connectivity map, required by FVGridGeometry
    template<class FVGridGeometry>
    struct MockConnectivityMap
    { void update(const FVGridGeometry& fvGridGeometry) {} };
71

72
73
    template<class FVGridGeometry>
    using ConnectivityMap = MockConnectivityMap<FVGridGeometry>;
74

75
76
    template<class FVGridGeometry, bool enableCache>
    using LocalView = StaggeredFVElementGeometry<FVGridGeometry, enableCache>;
77
78
};

79
80
81
} // end namespace Dumux
#endif

82
83
int main (int argc, char *argv[]) try
{
84
85
    using namespace Dumux;

86
87
88
89
90
    // maybe initialize mpi
    Dune::MPIHelper::instance(argc, argv);

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

91
    using Grid = Dune::YaspGrid<2>;
92

93
94
    constexpr int dim = Grid::dimension;
    constexpr int dimworld = Grid::dimensionworld;
95

96
    using GlobalPosition = typename Dune::FieldVector<Grid::ctype, dimworld>;
97
98
99
    using FVGridGeometry = StaggeredFVGridGeometry<typename Grid::LeafGridView, /*enable caching=*/ true,
                                                   TestFVGGTraits<typename Grid::LeafGridView> >;
    using FVElementGeometry = typename FVGridGeometry::LocalView;
100
101
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
102
103
104
105

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

110
111
    FVGridGeometry fvGridGeometry(leafGridView);
    fvGridGeometry.update();
112
113
114
115

    // iterate over elements. For every element get fv geometry and loop over scvs and scvfaces
    for (const auto& element : elements(leafGridView))
    {
116
        auto eIdx = fvGridGeometry.elementMapper().index(element);
117
        std::cout << std::endl << "Checking fvGeometry of element " << eIdx << std::endl;
118
        auto fvGeometry = localView(fvGridGeometry);
119
120
121
        fvGeometry.bind(element);

        auto range = scvs(fvGeometry);
122
        Detail::NoopFunctor<SubControlVolume> op;
123
        if(0 != testForwardIterator(range.begin(), range.end(), op))
124
125
126
127
            DUNE_THROW(Dune::Exception, "Iterator does not fulfill the forward iterator concept");

        for (auto&& scv : scvs(fvGeometry))
        {
128
            std::cout << "-- scv " << scv.indexInElement() << " center at: " << scv.center() << " , volume: " << scv.volume()  << std::endl;
129
130
131
        }

        auto range2 = scvfs(fvGeometry);
132
        Detail::NoopFunctor<SubControlVolumeFace> op2;
133
        if(0 != testForwardIterator(range2.begin(), range2.end(), op2))
134
135
136
137
            DUNE_THROW(Dune::Exception, "Iterator does not fulfill the forward iterator concept");

        for (auto&& scvf : scvfs(fvGeometry))
        {
138
139
            std::cout << "-- scvf " << scvf.index() << " ip at: " << scvf.ipGlobal() << " normal: " << scvf.unitOuterNormal();
            if (scvf.boundary()) std::cout << " (on boundary).";
140
141
142
143
144
145
146
147
148
149
150
151
            std::cout << std::endl;
        }
    }
}
// //////////////////////////////////
//   Error handler
// /////////////////////////////////
catch (Dune::Exception e) {

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