FEniCS / dolfinx

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

[BUG]: Can't (easily) assemble a block vector without providing a left-hand side in Python #3301

Closed jpdean closed 3 days ago

jpdean commented 1 week ago

Summarize the issue

Calling assemble_vector_block in Python fails with RuntimeError: Unable to extract a mesh. if any rows of the bilinear forms provided are all None. This happens even when no boundary conditions need applying.

This problem occurs because assemble_vector_block calls apply_lifting on each row, which gets the mesh from the form (see here). Hence, if there forms are all nullptr, it can't extract a mesh.

There are a couple of solutions to this:

  1. If there is nothing to lift, apply_lifting should probably just return without doing anything. This would be done at the C++ level.
  2. Don't call apply_lifting when there is nothing to lift. This would mean adding a check in Python to assemble_vector_block.

I think I prefer the first option. Regardless, I think it would be good to separate the lifting and bc application out of assemble_vector_block so we'd have assemble_vector_block, apply_lifting_block, and set_bc_block. This is consistent with non-blocked problems and gives more control to the user.

How to reproduce the bug

Run the minimal example below.

Minimal Example (Python)

from mpi4py import MPI
from dolfinx import mesh, fem
import ufl
import dolfinx.fem.petsc

comm = MPI.COMM_WORLD
msh = mesh.create_unit_square(comm, 2, 2)

V = fem.functionspace(msh, ("Lagrange", 1))
W = fem.functionspace(msh, ("Lagrange", 2))

v = ufl.TestFunction(V)
w = ufl.TestFunction(W)

L = fem.form([v * ufl.dx, w * ufl.dx])

b = dolfinx.fem.petsc.assemble_vector_block(L, [[None, None], [None, None]])

Output (Python)

No response

Version

main branch

DOLFINx git commit

No response

Installation

No response

Additional information

No response

jorgensd commented 1 week ago

Does this also happen if you try to call dolfinx.fem.apply_lifting without a list of a's ?

jpdean commented 1 week ago

You can't call apply_lifting without providing a, it's a required argument.