test_2d2d_intersection.cc 10.5 KB
Newer Older
1
2
3
4
5
6
7
8
#include <config.h>

#include <iostream>

#include <dune/common/exceptions.hh>
#include <dune/common/fvector.hh>
#include <dune/geometry/multilineargeometry.hh>

9
#include <dumux/geometry/geometryintersection.hh>
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54

#ifndef DOXYGEN
template<int dimWorld>
Dune::FieldVector<double, dimWorld>
makePosition( double x, double y )
{
    if (dimWorld == 2)
        return {x, y};
    else if (dimWorld == 3)
        return {x, y, 1.0};
    else
        DUNE_THROW(Dune::InvalidStateException, "Invalid dimworld");
}

template<int dimWorld>
Dune::MultiLinearGeometry<double, 2, dimWorld>
makeTriangle( Dune::FieldVector<double, dimWorld>&& a,
              Dune::FieldVector<double, dimWorld>&& b,
              Dune::FieldVector<double, dimWorld>&& c )
{
    using CornerStorage = std::vector<Dune::FieldVector<double, dimWorld>>;
    return { Dune::GeometryTypes::triangle,
             CornerStorage({std::move(a), std::move(b), std::move(c)}) };
}

template<int dimWorld>
Dune::MultiLinearGeometry<double, 2, dimWorld>
makeQuadrilateral( Dune::FieldVector<double, dimWorld>&& a,
                   Dune::FieldVector<double, dimWorld>&& b,
                   Dune::FieldVector<double, dimWorld>&& c,
                   Dune::FieldVector<double, dimWorld>&& d )
{
    using CornerStorage = std::vector<Dune::FieldVector<double, dimWorld>>;
    return { Dune::GeometryTypes::quadrilateral,
             CornerStorage({std::move(a), std::move(b), std::move(c), std::move(d)}) };
}

template<int dimworld, class Polygon1, class Polygon2>
bool testPolygonIntersection(const Polygon1& pol1,
                             const Polygon2& pol2,
                             bool expectIntersection)
{
    using Policy = Dumux::IntersectionPolicy::PolygonPolicy<double, dimworld>;
    using Algorithm = Dumux::GeometryIntersection<Polygon1, Polygon2, Policy>;

55
    typename Algorithm::Intersection intersection;
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
100
101
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
127
128
129
130
    const bool found = Algorithm::intersection(pol1, pol2, intersection);
    if (found && !expectIntersection)
        std::cout << "Found false positive!" << std::endl;
    else if (!found && expectIntersection)
        std::cout << "Could not find intersection!" << std::endl;
    else if (found && expectIntersection)
        std::cout << "Found intersection" << std::endl;
    else
        std::cout << "No intersection" << std::endl;

    return (found == expectIntersection);
}

template<int dimWorld>
void testPolygonIntersections(std::vector<bool>& returns)
{
    for (auto scaling : {1.0, 1e3, 1e12, 1e-12})
    {
        std::cout << "Test with scaling " << scaling << std::endl;
        const auto tria1 = makeTriangle( makePosition<dimWorld>(0.0, 0.0),
                                         makePosition<dimWorld>(0.0, 1.0*scaling),
                                         makePosition<dimWorld>(1.0*scaling, 1.0*scaling) );
        const auto tria2 = makeTriangle( makePosition<dimWorld>(1.0*scaling, 0.0),
                                         makePosition<dimWorld>(0.0, 0.0),
                                         makePosition<dimWorld>(0.0, 1.0*scaling) );
        const auto tria3 = makeTriangle( makePosition<dimWorld>(1.0*scaling, 0.0),
                                         makePosition<dimWorld>(0.0, 0.5*scaling),
                                         makePosition<dimWorld>(1.0*scaling, 1.0*scaling) );
        const auto tria4 = makeTriangle( makePosition<dimWorld>(1.0*scaling, 0.0),
                                         makePosition<dimWorld>(1.0*scaling, 1.0*scaling),
                                         makePosition<dimWorld>(2.0*scaling, 1.0*scaling) );
        const auto tria5 = makeTriangle( makePosition<dimWorld>(1.0*scaling, 0.0),
                                         makePosition<dimWorld>(0.5*scaling, 0.5*scaling),
                                         makePosition<dimWorld>(2.0*scaling, 1.0*scaling) );
        const auto tria6 = makeTriangle( makePosition<dimWorld>(1.0*scaling, 0.0),
                                         makePosition<dimWorld>(0.501*scaling, 0.501*scaling),
                                         makePosition<dimWorld>(2.0*scaling, 1.0*scaling) );
        const auto tria7 = makeTriangle( makePosition<dimWorld>(1.0*scaling, 0.0),
                                         makePosition<dimWorld>(0.499*scaling, 0.501*scaling),
                                         makePosition<dimWorld>(2.0*scaling, 1.0*scaling) );

        returns.push_back(testPolygonIntersection<dimWorld>(tria1, tria2, true));
        returns.push_back(testPolygonIntersection<dimWorld>(tria1, tria3, true));
        returns.push_back(testPolygonIntersection<dimWorld>(tria1, tria4, false)); // touches in point
        returns.push_back(testPolygonIntersection<dimWorld>(tria1, tria5, false)); // touches in edge
        returns.push_back(testPolygonIntersection<dimWorld>(tria1, tria6, false)); // little bit outside triangle
        returns.push_back(testPolygonIntersection<dimWorld>(tria1, tria7, true));  // little bit inside triangle

        // test some quadrilaterals
        const auto quad1 = makeQuadrilateral( makePosition<dimWorld>(0.0, 0.0),
                                              makePosition<dimWorld>(1.0*scaling, 0.0),
                                              makePosition<dimWorld>(0.0, 1.0*scaling),
                                              makePosition<dimWorld>(1.0*scaling, 1.0*scaling) );
        const auto quad2 = makeQuadrilateral( makePosition<dimWorld>(1.0*scaling, 0.0),
                                              makePosition<dimWorld>(2.0*scaling, 0.0),
                                              makePosition<dimWorld>(1.0*scaling, 2.0*scaling),
                                              makePosition<dimWorld>(2.0*scaling, 2.0*scaling) );
        const auto quad3 = makeQuadrilateral( makePosition<dimWorld>(-1.0*scaling, 0.0),
                                              makePosition<dimWorld>(0.0, 0.0),
                                              makePosition<dimWorld>(-1.0*scaling, 1.0*scaling),
                                              makePosition<dimWorld>(0.0, 1.0*scaling) );
        const auto quad4 = makeQuadrilateral( makePosition<dimWorld>(0.5*scaling, 0.0),
                                              makePosition<dimWorld>(0.5*scaling, 0.501*scaling),
                                              makePosition<dimWorld>(1.0*scaling, 0.0),
                                              makePosition<dimWorld>(1.0*scaling, 1.0*scaling) );

        returns.push_back(testPolygonIntersection<dimWorld>(tria1, quad1, true));
        returns.push_back(testPolygonIntersection<dimWorld>(tria1, quad2, false)); // touches in point
        returns.push_back(testPolygonIntersection<dimWorld>(tria1, quad3, false)); // touches in edge
        returns.push_back(testPolygonIntersection<dimWorld>(tria1, quad4, true)); // small overlap region

        std::cout << std::endl;
    }
}

131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
void testParallelPolygons(std::vector<bool>& returns)
{
    using Point = Dune::FieldVector<double, 3>;

    for (auto scaling : {1.0, 1e3, 1e12, 1e-12})
    {
        const double unit = 1.0*scaling;
        const double offUnit = (1.0 + 1e-6)*scaling;

        std::cout << "Test with scaling " << scaling << std::endl;
        const auto tria1 = makeTriangle( Point{{0.0, 0.0, unit}},
                                         Point{{unit, 0.0, unit}},
                                         Point{{unit, unit, unit}} );
        const auto tria2 = makeTriangle( Point{{0.0, 0.0, offUnit}},
                                         Point{{0.0, unit, offUnit}},
                                         Point{{unit, 0.0, offUnit}} );
        returns.push_back(testPolygonIntersection<3>(tria1, tria2, false));
    }

    std::cout << std::endl;
}

void testNonParallelPolygons(std::vector<bool>& returns)
{
    using Point = Dune::FieldVector<double, 3>;

    for (auto scaling : {1.0, 1e3, 1e12, 1e-12})
    {
        const double unit = 1.0*scaling;
        const double offUnit = (1.0 + 1e-6)*scaling;

        std::cout << "Test with scaling " << scaling << std::endl;
        const auto tria1 = makeTriangle( Point{{0.0, 0.0, unit}},
                                         Point{{unit, 0.0, unit}},
                                         Point{{unit, unit, unit}} );
        const auto tria2 = makeTriangle( Point{{0.0, 0.0, unit}},
                                         Point{{0.0, unit, unit}},
                                         Point{{unit, 0.0, offUnit}} );
        returns.push_back(testPolygonIntersection<3>(tria1, tria2, false));
    }

    std::cout << std::endl;
}

175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
// unit test discovered in a simulation, that originally failed on commit 51fc9c1e3
void testSimulationCase(std::vector<bool>& returns)
{
    using Position = Dune::FieldVector<double, 3>;
    const auto quad1 = makeQuadrilateral(
        Position{{9.49264705847428835739e-01, 1.40073529420942810564e-01, 4.00000000000000077716e-01}},
        Position{{9.44669117627101373458e-01, 1.42095588240539000280e-01, 4.00000000000000077716e-01}},
        Position{{9.49632352923713174420e-01, 1.45036764710471721695e-01, 4.00000000000000022204e-01}},
        Position{{9.46446078418066827354e-01, 1.44730392160359544462e-01, 4.00000000000000077716e-01}}
    );
    const auto quad2 = makeQuadrilateral(
        Position{{9.39999999999999946709e-01, 1.39999999999999985567e-01, 4.00000000000000022204e-01}},
        Position{{9.59999999999999964473e-01, 1.39999999999999985567e-01, 4.00000000000000022204e-01}},
        Position{{9.39999999999999946709e-01, 1.60000000000000003331e-01, 4.00000000000000022204e-01}},
        Position{{9.59999999999999964473e-01, 1.60000000000000003331e-01, 4.00000000000000022204e-01}}
    );

    returns.push_back(testPolygonIntersection<3>(quad1, quad2, true));

    std::cout << std::endl;
}

197
198
#endif

199
int main(int argc, char* argv[])
200
201
{
    std::vector<bool> returns;
202
203

    std::cout << "Testing intersections in 2d space" << std::endl;
204
205
    testPolygonIntersections<2>(returns);

206
207
208
209
210
211
212
213
214
    std::cout << "Testing intersecions in 3d space" << std::endl;
    testPolygonIntersections<3>(returns);

    std::cout << "Testing parallel polygons in 3d space" << std::endl;
    testParallelPolygons(returns);

    std::cout << "Testing non-parallel polygons in 3d space" << std::endl;
    testNonParallelPolygons(returns);

215
216
217
    std::cout << "Testing 3d unit test from actual simulation" << std::endl;
    testSimulationCase(returns);

218
219
220
221
222
223
224
225
226
227
    // TODO: implement and test point and segment intersections

    // determine the exit code
    if (std::any_of(returns.begin(), returns.end(), [](bool i){ return !i; }))
        return 1;

    std::cout << "All tests passed!" << std::endl;

    return 0;
}