- changed title to FunctionAssigner does not always call apply('insert')
FunctionAssigner does not always call apply('insert')
Issue #587
resolved
The following code shows the bug
from dolfin import *
mesh = UnitSquareMesh(4, 4)
V = FunctionSpace(mesh, 'CG', 1)
W = MixedFunctionSpace([V, V])
v1 = Function(V)
v2 = Function(V)
w = Function(W)
w.vector()[:] = 1
fa = FunctionAssigner([V, V], W)
fa.assign([v1, v2], w)
print v1.vector()[:].array()
# Removing this line gives PETSc error
v2.vector().apply('insert')
print v2.vector()[:].array()
Comments (7)
-
reporter -
- removed milestone
Removing milestone: 1.7 (automated comment)
-
reporter - changed status to resolved
I can no longer replicate this
-
reporter - changed status to open
I can replicate the bug with the following code (updated to make it work with current master), but only when running in parallel. In serial there is no problem.
from dolfin import * mesh = UnitSquareMesh(4, 4) V = FunctionSpace(mesh, 'CG', 1) e = MixedElement([V.ufl_element(), V.ufl_element()]) W = FunctionSpace(mesh, e) v1 = Function(V) v2 = Function(V) w = Function(W) w.vector()[:] = 1 fa = FunctionAssigner([V, V], W) fa.assign([v1, v2], w) print v1.vector()[:].array() # Removing this line gives PETSc error #v2.vector().apply('insert') print v2.vector()[:].array()
In parallel I get an error
*** Error: Unable to successfully call PETSc function 'VecCopy'. *** Reason: PETSc error code is: 73. *** Where: This error was encountered inside [...]/dolfin/dolfin/la/PETScVector.cpp.
-
-
-
assigned issue to
-
assigned issue to
-
- changed status to resolved
Should be fixed by 5d210681f6620b5ef00ef486ed13d283b733da06
- Log in to comment