test_staggeredfvgeometry.cc 6.4 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
34

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

35
36
#include <dumux/common/basicproperties.hh>
#include <dumux/discretization/staggered/properties.hh>
37
#include <dumux/discretization/staggered/fvgridgeometry.hh>
38
39
40
41
42
43
44
#include <dumux/discretization/staggered/fvelementgeometry.hh>

namespace Dumux
{

namespace Properties
{
45
NEW_TYPE_TAG(TestFVGeometry, INHERITS_FROM(StaggeredModel, NumericModel));
46
47
48

SET_TYPE_PROP(TestFVGeometry, Grid, Dune::YaspGrid<2>);

49
SET_BOOL_PROP(TestFVGeometry, EnableFVGridGeometryCache, ENABLE_CACHING);
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
}

}

template<class T>
class NoopFunctor {
public:
  NoopFunctor() {}
  void operator()(const T& t){}
};

int main (int argc, char *argv[]) try
{
    // maybe initialize mpi
    Dune::MPIHelper::instance(argc, argv);

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

    // aliases
    using TypeTag = TTAG(TestFVGeometry);
    using Grid = typename GET_PROP_TYPE(TypeTag, Grid);
    using GridView = typename Grid::LeafGridView;

    constexpr int dim = GridView::dimension;
    constexpr int dimworld = GridView::dimensionworld;

76
    using GlobalPosition = typename Dune::FieldVector<Grid::ctype, dimworld>;
77
78
    using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume);
    using SubControlVolumeFace = typename GET_PROP_TYPE(TypeTag, SubControlVolumeFace);
79
    using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
80
81
82
83

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

88
89
    FVGridGeometry fvGridGeometry(leafGridView);
    fvGridGeometry.update();
90

91
    std::cout << "Abbreviatons:\n"
92
93
94
95
              << "ip - postition of face center\n"
              << "face - face index\n"
              << "self/oppo - dofIdx on intersection (self/opposite)\n"
              << "norm in/out - dofIdx on side normal to intersection (within own element / in adjacent element)" << std::endl;
96

97
98
99
    // iterate over elements. For every element get fv geometry and loop over scvs and scvfaces
    for (const auto& element : elements(leafGridView))
    {
100
        auto eIdx = fvGridGeometry.elementMapper().index(element);
101
        std::cout << std::endl << "Checking fvGeometry of element " << eIdx << std::endl;
102
        auto fvGeometry = localView(fvGridGeometry);
103
104
105
106
        fvGeometry.bind(element);

        auto range = scvs(fvGeometry);
        NoopFunctor<SubControlVolume> op;
107
        if(0 != testForwardIterator(range.begin(), range.end(), [](const SubControlVolumeFace& scvf){}))
108
109
110
111
            DUNE_THROW(Dune::Exception, "Iterator does not fulfill the forward iterator concept");

        for (auto&& scv : scvs(fvGeometry))
        {
112
            std::cout << "-- scv " << scv.dofIndex() << " center at: " << scv.center() << std::endl;
113
114
115
116
        }

        auto range2 = scvfs(fvGeometry);
        NoopFunctor<SubControlVolumeFace> op2;
117
        if(0 != testForwardIterator(range2.begin(), range2.end(), [](const SubControlVolumeFace& scvf){}))
118
119
            DUNE_THROW(Dune::Exception, "Iterator does not fulfill the forward iterator concept");

120

121
122
        for (auto&& scvf : scvfs(fvGeometry))
        {
123
124
125
            std::cout <<  std::fixed << std::left << std::setprecision(2)
            << "ip "<< scvf.ipGlobal()
            << "; face "  << std::setw(3)  << scvf.index()
126
            << "; self/oppo " << std::setw(3) << scvf.dofIndex() << "/" << std::setw(3) <<scvf.dofIndexOpposingFace()
127
            << "; dist self/oppo " << std::setw(3) << scvf.selfToOppositeDistance()
128
129
            << ", norm1 in/out " << std::setw(3) << scvf.pairData(0).normalPair.first << "/" << std::setw(3) << scvf.pairData(0).normalPair.second
            << ", norm2 in/out " << std::setw(3) << scvf.pairData(1).normalPair.first << "/" << std::setw(3) << scvf.pairData(1).normalPair.second
130
131
            << ", par1 in/out " << std::setw(3) << scvf.dofIndex() << "/" << std::setw(3) << scvf.pairData(0).outerParallelFaceDofIdx
            << ", par2 in/out " << std::setw(3) << scvf.dofIndex() << "/" << std::setw(3) << scvf.pairData(1).outerParallelFaceDofIdx
132
133
134
135
            << ", normDist1 " << std::setw(3) << scvf.pairData(0).normalDistance
            << ", normDist2 " << std::setw(3) << scvf.pairData(1).normalDistance
            << ", parDist1 " << std::setw(3) << scvf.pairData(0).parallelDistance
            << ", parDist2 " << std::setw(3) << scvf.pairData(1).parallelDistance;
136
            if (scvf.boundary()) std::cout << " (on boundary)";
137
138
139
140
141
142
143
144
145
146
147
148
            std::cout << std::endl;
        }
    }
}
// //////////////////////////////////
//   Error handler
// /////////////////////////////////
catch (Dune::Exception e) {

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