File: LagrangeInterpolator.py

package info (click to toggle)
dolfin 1.4.0%2Bdfsg-4
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 38,696 kB
  • ctags: 29,035
  • sloc: cpp: 409,305; python: 22,067; sh: 230; makefile: 215; xml: 126
file content (90 lines) | stat: -rw-r--r-- 2,922 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
"""Unit tests for interpolation using LagrangeInterpolator"""

# Copyright (C) 2014 Mikael Mortensen
#
# This file is part of DOLFIN.
#
# DOLFIN is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# DOLFIN 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 Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
#
# First added:  2014-02-18
# Last changed:

import unittest
import numpy
from dolfin import *

class Quadratic2D(Expression):
    def eval(self, values, x):
        values[0] = x[0]*x[0] + x[1]*x[1] + 1.0

class Quadratic3D(Expression):
    def eval(self, values, x):
        values[0] = x[0]*x[0] + x[1]*x[1] + x[2]*x[2] + 1.0

class LagrangeInterpolatorTest(unittest.TestCase):

    def test_functional2D(self):
        """Test integration of function interpolated in non-matching meshes"""

        f = Quadratic2D()
        
        ll = LagrangeInterpolator()

        # Interpolate quadratic function on course mesh
        mesh0 = UnitSquareMesh(8, 8)
        V0 = FunctionSpace(mesh0, "Lagrange", 2)
        u0 = Function(V0)
        ll.interpolate(u0, f)

        # Interpolate FE function on finer mesh
        mesh1 = UnitSquareMesh(31, 31)
        V1 = FunctionSpace(mesh1, "Lagrange", 2)
        u1 = Function(V1)
        ll.interpolate(u1, u0)
        self.assertAlmostEqual(assemble(u0*dx), assemble(u1*dx), 10)

        mesh1 = UnitSquareMesh(30, 30)
        V1 = FunctionSpace(mesh1, "Lagrange", 2)
        u1 = Function(V1)
        ll.interpolate(u1, u0)
        self.assertAlmostEqual(assemble(u0*dx), assemble(u1*dx), 10)

    def test_functional3D(self):
        """Test integration of function interpolated in non-matching meshes"""

        f = Quadratic3D()

        ll = LagrangeInterpolator()

        # Interpolate quadratic function on course mesh
        mesh0 = UnitCubeMesh(4, 4, 4)
        V0 = FunctionSpace(mesh0, "Lagrange", 2)
        u0 = Function(V0)
        ll.interpolate(u0, f)

        # Interpolate FE function on finer mesh
        mesh1 = UnitCubeMesh(11, 11, 11)
        V1 = FunctionSpace(mesh1, "Lagrange", 2)
        u1 = Function(V1)
        ll.interpolate(u1, u0)
        self.assertAlmostEqual(assemble(u0*dx), assemble(u1*dx), 10)

        mesh1 = UnitCubeMesh(10, 11, 10)
        V1 = FunctionSpace(mesh1, "Lagrange", 2)
        u1 = Function(V1)
        ll.interpolate(u1, u0)
        self.assertAlmostEqual(assemble(u0*dx), assemble(u1*dx), 10)

if __name__ == "__main__":
    unittest.main()