File: test_lu_solver.py

package info (click to toggle)
dolfin 2018.1.0.post1-16
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 28,764 kB
  • sloc: xml: 104,040; cpp: 98,856; python: 22,511; makefile: 204; sh: 182
file content (119 lines) | stat: -rwxr-xr-x 3,670 bytes parent folder | download | duplicates (5)
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
"""Unit tests for the LUSolver interface"""

# Copyright (C) 2014 Garth N. Wells
#
# 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/>.

from dolfin import *
import pytest
from dolfin_utils.test import skip_if_not_PETSc, skip_in_parallel

backends = ["PETSc", pytest.param(("Eigen"), marks=skip_in_parallel)]

@pytest.mark.parametrize('backend', backends)
def test_lu_solver(backend):

    # Check whether backend is available
    if not has_linear_algebra_backend(backend):
        pytest.skip('Need %s as backend to run this test' % backend)

    # Set linear algebra backend
    prev_backend = parameters["linear_algebra_backend"]
    parameters["linear_algebra_backend"] = backend

    mesh = UnitSquareMesh(12, 12)
    V = FunctionSpace(mesh, "Lagrange", 1)
    u, v = TrialFunction(V), TestFunction(V)
    A = assemble(Constant(1.0)*u*v*dx)
    b = assemble(Constant(1.0)*v*dx)

    norm = 13.0

    solver = LUSolver()
    x = Vector()
    solver.solve(A, x, b)
    assert round(x.norm("l2") - norm, 10) == 0

    solver = LUSolver(A)
    x = Vector()
    solver.solve(x, b)
    assert round(x.norm("l2") - norm, 10) == 0

    solver = LUSolver()
    x = Vector()
    solver.set_operator(A)
    solver.solve(x, b)
    assert round(x.norm("l2") - norm, 10) == 0

    solver = LUSolver()
    x = Vector()
    solver.solve(A, x, b)
    assert round(x.norm("l2") - norm, 10) == 0

    # Reset backend
    parameters["linear_algebra_backend"] = prev_backend


@pytest.mark.parametrize('backend', backends)
def test_lu_solver_reuse(backend):
    """Test that LU re-factorisation is only performed after
    set_operator(A) is called"""

    # Test requires PETSc version 3.5 or later. Use petsc4py to check
    # version number.
    try:
        from petsc4py import PETSc
    except ImportError:
        pytest.skip("petsc4py required to check PETSc version")
    else:
        if not PETSc.Sys.getVersion() >= (3, 5, 0):
            pytest.skip("PETSc version must be 3.5  of higher")


    # Check whether backend is available
    if not has_linear_algebra_backend(backend):
        pytest.skip('Need %s as backend to run this test' % backend)

    # Set linear algebra backend
    prev_backend = parameters["linear_algebra_backend"]
    parameters["linear_algebra_backend"] = backend

    mesh = UnitSquareMesh(12, 12)
    V = FunctionSpace(mesh, "Lagrange", 1)
    u, v = TrialFunction(V), TestFunction(V)
    b = assemble(Constant(1.0)*v*dx)

    A = assemble(Constant(1.0)*u*v*dx)
    norm = 13.0

    solver = LUSolver(A)
    x = Vector()
    solver.solve(x, b)
    assert round(x.norm("l2") - norm, 10) == 0

    assemble(Constant(0.5)*u*v*dx, tensor=A)
    x = Vector()
    solver.solve(x, b)
    if backend != 'Eigen':
        # The eigen backend only recomputes the factorization once set_operator is called
        assert round(x.norm("l2") - 2.0*norm, 10) == 0

    solver.set_operator(A)
    solver.solve(x, b)
    assert round(x.norm("l2") - 2.0*norm, 10) == 0

    # Reset backend
    parameters["linear_algebra_backend"] = prev_backend