File: test_petsc4py.py

package info (click to toggle)
dolfin 2019.2.0~legacy20240219.1c52e83-18
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 31,700 kB
  • sloc: xml: 104,040; cpp: 102,227; python: 24,356; sh: 460; makefile: 330; javascript: 226
file content (101 lines) | stat: -rw-r--r-- 2,968 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
"""Unit tests for the pybind11 wrapping of PETSc / petsc4py"""

# Copyright (C) 2017 Tormod Landet
#
# 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/>.

import gc
from dolfin import (UnitSquareMesh, TrialFunction, TestFunction,
                    FunctionSpace, Function, assemble, dx,
                    parameters, as_backend_type, has_petsc)

if has_petsc() :
    from dolfin import (PETScVector, PETScMatrix)

from dolfin_utils.test import skip_if_not_petsc4py, pushpop_parameters


@skip_if_not_petsc4py
def test_petsc4py_vector(pushpop_parameters):
    "Test PETScVector <-> petsc4py.PETSc.Vec conversions"
    parameters["linear_algebra_backend"] = "PETSc"

    # Assemble a test matrix
    mesh = UnitSquareMesh(4, 4)
    V = FunctionSpace(mesh, "Lagrange", 1)
    v = TestFunction(V)
    a = v*dx
    b1 = assemble(a)

    # Test conversion dolfin.PETScVector -> petsc4py.PETSc.Vec
    b1 = as_backend_type(b1)
    v1 = b1.vec()

    # Copy and scale vector with petsc4py
    v2 = v1.copy()
    v2.scale(2.0)

    # Test conversion petsc4py.PETSc.Vec  -> PETScVector
    b2 = PETScVector(v2)

    assert (b1.get_local()*2.0 == b2.get_local()).all()


@skip_if_not_petsc4py
def test_petsc4py_matrix(pushpop_parameters):
    "Test PETScMatrix <-> petsc4py.PETSc.Mat conversions"
    parameters["linear_algebra_backend"] = "PETSc"

    # Assemble a test matrix
    mesh = UnitSquareMesh(4, 4)
    V = FunctionSpace(mesh, "Lagrange", 1)
    u, v = TrialFunction(V), TestFunction(V)
    a = u*v*dx
    A1 = assemble(a)

    # Test conversion dolfin.PETScMatrix -> petsc4py.PETSc.Mat
    A1 = as_backend_type(A1)
    M1 = A1.mat()

    # Copy and scale matrix with petsc4py
    M2 = M1.copy()
    M2.scale(2.0)

    # Test conversion petsc4py.PETSc.Mat  -> PETScMatrix
    A2 = PETScMatrix(M2)

    assert (A1.array()*2.0 == A2.array()).all()

@skip_if_not_petsc4py
def test_ref_count(pushpop_parameters):
    "Test petsc4py reference counting"
    parameters["linear_algebra_backend"] = "PETSc"

    mesh = UnitSquareMesh(3, 3)
    V = FunctionSpace(mesh, "P", 1)

    # Check u and x own the vector
    u = Function(V)
    x = as_backend_type(u.vector()).vec()
    assert x.refcount == 2

    # Check decref
    del u; gc.collect()  # destroy u
    assert x.refcount == 1

    # Check incref
    vec = PETScVector(x)
    assert x.refcount == 2