FEniCS / dolfinx

Next generation FEniCS problem solving environment
https://fenicsproject.org
GNU Lesser General Public License v3.0
731 stars 177 forks source link

Problem assembling rectangular matrix in parallel #2177

Closed hermanmakhm closed 2 years ago

hermanmakhm commented 2 years ago

FEniCSx: Last pull from main on 2022-03-15 via spack

I attempted to assemble a rectangular matrix which worked in serial (mpirun -np 1 python *.py) but didn't work when I attempted to run the same code in parallel (mpirun -np 2 python *.py). Here is the sample code to reproduce the problem:

import gmsh
import ufl
import dolfinx.fem
import multiphenicsx.mesh

r = 3
lcar = 1. / 4.

gmsh.initialize()
gmsh.model.add("mesh")

p0 = gmsh.model.geo.addPoint(0.0, 0.0, 0.0, lcar)
p1 = gmsh.model.geo.addPoint(0.0, +r, 0.0, lcar)
p2 = gmsh.model.geo.addPoint(0.0, -r, 0.0, lcar)
c0 = gmsh.model.geo.addCircleArc(p1, p0, p2)
c1 = gmsh.model.geo.addCircleArc(p2, p0, p1)
boundary = gmsh.model.geo.addCurveLoop([c0, c1])
domain = gmsh.model.geo.addPlaneSurface([boundary])

gmsh.model.geo.synchronize()
gmsh.model.addPhysicalGroup(1, [c0, c1], 1)
gmsh.model.addPhysicalGroup(2, [boundary], 0)
gmsh.model.mesh.generate(2)

mesh, subdomains, boundaries = multiphenicsx.mesh.gmsh_to_fenicsx(gmsh.model, gdim=2)
gmsh.finalize()

# Define a function space
V = dolfinx.fem.FunctionSpace(mesh, ("Lagrange", 2))

u = ufl.TrialFunction(V)
v = ufl.TestFunction(V)

a = [[ufl.inner(u,v) * ufl.dx],
     [ufl.inner(u,v) * ufl.dx]]
a_cpp = dolfinx.fem.form(a)
A = dolfinx.fem.petsc.assemble_matrix_block(a_cpp, bcs=[])
A.assemble()

and below is the error message associated with running the code in parallel

Traceback (most recent call last):
  File "/stck/hmmak/fenicsx-scripts/10deg_naca12/test_rectangular_assembly.py", line 37, in <module>
    A = dolfinx.fem.petsc.assemble_matrix_block(a_cpp, bcs=[])
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/._view/c7ohj445pwtcbtlnhace2txzvkz7qrsg/lib/python3.9/functools.py", line 888, in wrapper
    return dispatch(args[0].__class__)(*args, **kw)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/view/lib/python3.9/site-packages/dolfinx/fem/petsc.py", line 383, in assemble_matrix_block
    return assemble_matrix_block(A, a, bcs, diagonal, constants, coeffs)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/._view/c7ohj445pwtcbtlnhace2txzvkz7qrsg/lib/python3.9/functools.py", line 888, in wrapper
    return dispatch(args[0].__class__)(*args, **kw)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/view/lib/python3.9/site-packages/dolfinx/fem/petsc.py", line 410, in _
    A.assemble(PETSc.Mat.AssemblyType.FLUSH)
  File "PETSc/Mat.pyx", line 1118, in petsc4py.PETSc.Mat.assemble
petsc4py.PETSc.Error: error code 63
[1] MatAssemblyEnd() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/interface/matrix.c:5671
[1] MatAssemblyEnd_MPIAIJ() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/impls/aij/mpi/mpiaij.c:707
[1] MatSetValues_MPIAIJ() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/impls/aij/mpi/mpiaij.c:482
[1] Argument out of range
[1] Column too large: col 2475 max 2268
Traceback (most recent call last):
  File "/stck/hmmak/fenicsx-scripts/10deg_naca12/test_rectangular_assembly.py", line 37, in <module>
    A = dolfinx.fem.petsc.assemble_matrix_block(a_cpp, bcs=[])
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/._view/c7ohj445pwtcbtlnhace2txzvkz7qrsg/lib/python3.9/functools.py", line 888, in wrapper
    return dispatch(args[0].__class__)(*args, **kw)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/view/lib/python3.9/site-packages/dolfinx/fem/petsc.py", line 383, in assemble_matrix_block
    return assemble_matrix_block(A, a, bcs, diagonal, constants, coeffs)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/._view/c7ohj445pwtcbtlnhace2txzvkz7qrsg/lib/python3.9/functools.py", line 888, in wrapper
    return dispatch(args[0].__class__)(*args, **kw)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/view/lib/python3.9/site-packages/dolfinx/fem/petsc.py", line 410, in _
    A.assemble(PETSc.Mat.AssemblyType.FLUSH)
  File "PETSc/Mat.pyx", line 1118, in petsc4py.PETSc.Mat.assemble
petsc4py.PETSc.Error: error code 63
[0] MatAssemblyEnd() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/interface/matrix.c:5671
[0] MatAssemblyEnd_MPIAIJ() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/impls/aij/mpi/mpiaij.c:707
[0] MatSetValues_MPIAIJ() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/impls/aij/mpi/mpiaij.c:482
[0] Argument out of range
[0] Column too large: col 3271 max 2268
--------------------------------------------------------------------------
Primary job  terminated normally, but 1 process returned
a non-zero exit code. Per user-direction, the job has been aborted.
garth-wells commented 2 years ago

Could you update the example so that it doesn't require multiphenicsx?

I can't test because I don't have multiphenicsx installed, but from inspection I think you need to change

V = dolfinx.fem.FunctionSpace(mesh, ("Lagrange", 2))
u = ufl.TrialFunction(V)
v = ufl.TestFunction(V)
a = [[ufl.inner(u, v)*ufl.dx],
     [ufl.inner(u, v)*ufl.dx]]

to

V0 = dolfinx.fem.FunctionSpace(mesh, ("Lagrange", 2))
V1 = V0.clone()
u = ufl.TrialFunction(V0)
v0, v1 = ufl.TrialFunction(V0), ufl.TestFunction(V1)
a = [[ufl.inner(u, v0)*ufl.dx],
     [ufl.inner(u, v1)*ufl.dx]]
hermanmakhm commented 2 years ago

I have updated the example (with the mesh from this tutorial) and also implemented your suggestions as follows

import mpi4py.MPI as MPI
import ufl
import dolfinx.fem
import dolfinx.mesh

mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 8, 8, dolfinx.mesh.CellType.quadrilateral)

# Define a function space
V0 = dolfinx.fem.FunctionSpace(mesh, ("Lagrange", 2))
V1 = V0.clone()

u = ufl.TrialFunction(V0)
v0, v1 = ufl.TestFunction(V0), ufl.TestFunction(V1)

a = [[ufl.inner(u,v0) * ufl.dx],
     [ufl.inner(u,v1) * ufl.dx]]
a_cpp = dolfinx.fem.form(a)
A = dolfinx.fem.petsc.assemble_matrix_block(a_cpp, bcs=[])
A.assemble()
print("Finished")

It again runs to completion with mpirun -np 1 and gives the following error with mpirun -np 2

Traceback (most recent call last):
  File "/stck/hmmak/fenicsx-scripts/10deg_naca12/test_rectangular_assembly.py", line 18, in <module>
    A = dolfinx.fem.petsc.assemble_matrix_block(a_cpp, bcs=[])
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/._view/c7ohj445pwtcbtlnhace2txzvkz7qrsg/lib/python3.9/functools.py", line 888, in wrapper
    return dispatch(args[0].__class__)(*args, **kw)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/view/lib/python3.9/site-packages/dolfinx/fem/petsc.py", line 383, in assemble_matrix_block
    return assemble_matrix_block(A, a, bcs, diagonal, constants, coeffs)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/._view/c7ohj445pwtcbtlnhace2txzvkz7qrsg/lib/python3.9/functools.py", line 888, in wrapper
    return dispatch(args[0].__class__)(*args, **kw)
  File "/stck/hmmak/tools/spack/var/spack/environments/2022-03-15-real-fenicsx/.spack-env/view/lib/python3.9/site-packages/dolfinx/fem/petsc.py", line 410, in _
    A.assemble(PETSc.Mat.AssemblyType.FLUSH)
  File "PETSc/Mat.pyx", line 1118, in petsc4py.PETSc.Mat.assemble
petsc4py.PETSc.Error: error code 63
[1] MatAssemblyEnd() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/interface/matrix.c:5671
[1] MatAssemblyEnd_MPIAIJ() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/impls/aij/mpi/mpiaij.c:707
[1] MatSetValues_MPIAIJ() at /tmp/hmmak/spack-stage/spack-stage-petsc-3.16.4-nhmtjugsebgrng3cgzdd5u5hiqdj7u2j/spack-src/src/mat/impls/aij/mpi/mpiaij.c:482
[1] Argument out of range
[1] Column too large: col 353 max 288
garth-wells commented 2 years ago

@hermanmakhm thanks for reporting - should be fixed now.