File: test_gf2.py

package info (click to toggle)
dune-grid 2.11.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 8,472 kB
  • sloc: cpp: 60,883; python: 1,438; perl: 191; makefile: 12; sh: 3
file content (145 lines) | stat: -rw-r--r-- 4,637 bytes parent folder | download
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
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
55
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
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
# SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
# SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception

from io import StringIO
import math
from dune.common import FieldVector
from dune.grid import gridFunction, structuredGrid

codeFunc = """
#include <cmath>
#include <dune/common/fvector.hh>
#include <dune/common/math.hh>

template <class GridView>
auto myFunction(double a)
{
  return [a](const auto& en,const auto& x) -> auto
  {
    auto y = en.geometry().global(x);
    return std::sin(a*Dune::MathematicalConstants<double>::pi()*(y[0]+y[1]));
  };
}
"""
codeVecFunc = """
#include <dune/common/fvector.hh>

template <class GridView, class GF>
auto myVecFunction(Dune::FieldVector<double,1> &a, const GF &gf)
{
  return [&a,lgf=localFunction(gf)](const auto& en,const auto& x) mutable -> auto
  {
    lgf.bind(en);
    auto v = lgf(x);
    auto y = en.geometry().global(x);
    return Dune::FieldVector<double,2>{v[0],(y[0]-0.5)*a[0]};
  };
}
"""

def testGF_second(gridView):
    gf1 = gridView.function(lambda e,x:\
              math.sin(math.pi*(e.geometry.toGlobal(x)[0]+e.geometry.toGlobal(x)[1])))

    if True:
        a = 2.
        gf1 = gridView.function(lambda e,x:\
                math.sin(a*math.pi*(e.geometry.toGlobal(x)[0]+e.geometry.toGlobal(x)[1])),
                name="gf1")
        lgf1 = gf1.localFunction()
        average1 = 0
        for e in gridView.elements:
            lgf1.bind(e)
            average1 += lgf1([0.5,0.5])*e.geometry.volume
        # print(average1)
        # gf1.plot()
        gf2 = gridView.function("myFunction",StringIO(codeFunc),a,name="gf2")
        lgf2 = gf2.localFunction()
        average2 = 0
        for e in gridView.elements:
           lgf2.bind(e)
           average2 += lgf2([0.5,0.5])*e.geometry.volume
        # print(average2)
        # gf2.plot()
        # assert abs(average1-average2)<1e-12
        diff = 0
        for e in gridView.elements:
           lgf1.bind(e)
           lgf2.bind(e)
           diff += abs(lgf1([0.5,0.5])-lgf2([0.5,0.5]))
        assert diff<1e-12

    if True:
        gf1 = gridView.function(lambda e,x:\
                [math.sin(2*math.pi*(e.geometry.toGlobal(x)[0]+e.geometry.toGlobal(x)[1])),\
                 (e.geometry.toGlobal(x)[0]-0.5)*2],
                name="gf1")
        lgf1 = gf1.localFunction()
        average1 = 0
        for e in gridView.elements:
            lgf1.bind(e)
            average1 += sum(lgf1([0.5,0.5]))*e.geometry.volume
        # print(average1)
        # gf1.plot()
        a = FieldVector([2])
        gf2 = gridView.function("myVecFunction",StringIO(codeVecFunc),a,gf1,name="gf2")
        lgf2 = gf2.localFunction()
        average2 = 0
        for e in gridView.elements:
            lgf2.bind(e)
            average2 += sum(lgf2([0.5,0.5]))*e.geometry.volume
        # print(average2)
        # gf2.plot()
        assert abs(average1-average2)<1e-12
        diff = 0
        for e in gridView.elements:
            lgf1.bind(e)
            lgf2.bind(e)
            diff += abs(lgf1([0.5,0.5]).two_norm-lgf2([0.5,0.5]).two_norm)
        assert diff<1e-12
        a[0] = 3
        diff = 0
        for e in gridView.elements:
            lgf1.bind(e)
            lgf2.bind(e)
            v = lgf1([0.5,0.5])
            v[1] *= 3./2.
            diff += abs(v.two_norm-lgf2([0.5,0.5]).two_norm)
        assert diff<1e-12

        gridView.writeVTK("test_gf",pointdata=[gf1,gf2])

    if False:
        a = 2.
        @gridFunction(gridView)
        def gf1(e,x):
            return math.sin(a*math.pi*(e.geometry.toGlobal(x)[0]+e.geometry.toGlobal(x)[1]))
        lgf1 = gf1.localFunction()
        average1 = 0
        for e in gridView.elements:
            lgf1.bind(e)
            average1 += lgf1([0.5,0.5])*e.geometry.volume
        # print(average1)
        # gf1.plot()
        @gridFunction(gridView)
        def gf2(x):
            return math.sin(a*math.pi*(x[0]+x[1]))
        gf2 = gridView.function("myFunction",StringIO(codeFunc),a,name="gf2")
        lgf2 = gf2.localFunction()
        average2 = 0
        for e in gridView.elements:
            lgf2.bind(e)
            average2 += lgf2([0.5,0.5])*e.geometry.volume
        # print(average2)
        # gf2.plot()
        assert abs(average1-average2)<1e-12
        diff = 0
        for e in gridView.elements:
            lgf1.bind(e)
            lgf2.bind(e)
            diff += abs(lgf1([0.5,0.5])-lgf2([0.5,0.5]))
        assert diff<1e-12

if __name__ == "__main__":
    gridView = structuredGrid([0,0],[1,1],[10,10])
    testGF_second(gridView)