Skip to content

Fix IndexError in _pic_swarm_in_mesh on MPI ranks with empty cell partition#5078

Open
hardik-corintis wants to merge 6 commits intofiredrakeproject:mainfrom
Corintis:fix/pic-swarm-empty-rank
Open

Fix IndexError in _pic_swarm_in_mesh on MPI ranks with empty cell partition#5078
hardik-corintis wants to merge 6 commits intofiredrakeproject:mainfrom
Corintis:fix/pic-swarm-empty-rank

Conversation

@hardik-corintis
Copy link
Copy Markdown

@hardik-corintis hardik-corintis commented May 4, 2026

Description

Fixes #5077.

In _pic_swarm_in_mesh, parent_cell_nums_local contains -1 sentinels for points not
owned by the current rank. On ranks with an empty cell partition, cell_closure has shape
(0, k) and NumPy rejects -1 as an axis-0 index, causing an IndexError when
interpolating from a Submesh onto its parent mesh in parallel.

Fix: index cell_closure only on visible_idxs (the non-sentinel entries), pre-filling
the rest with -1. The same issue exists in the extruded branch: _parent_extrusion_numbering
floor-divides parent_cell_nums_local, and -1 // N == -1 in Python, so sentinels pass
through into base_parent_cell_nums and cause the same crash. Both branches are fixed.

Two regression tests added in tests/firedrake/submesh/test_submesh_interpolate.py using
8 MPI ranks:

  • test_submesh_interpolate_3Dcell_2Dfacet_empty_rank_8_processes: interpolates from a
    Submesh of the x=1 face of UnitCubeMesh(2, 2, 2) onto the parent mesh.
  • test_submesh_interpolate_3Dcell_extruded_empty_rank_8_processes: interpolates from
    ExtrudedMesh(UnitSquareMesh(1, 1), layers=3) (6 cells) onto a UnitCubeMesh(2, 2, 2),
    guaranteeing at least two ranks own zero extruded cells.

@connorjward connorjward requested a review from leo-collins May 5, 2026 11:00


@pytest.mark.parallel(nprocs=2)
def test_submesh_interpolate_3Dcell_2Dfacet_empty_rank_2_processes():
Copy link
Copy Markdown
Contributor

@leo-collins leo-collins May 5, 2026

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I've just had a look at this now (on the main branch). I don't get the IndexError you described. However, it does fail the assertion at the bottom, but this is because subm is an immersed manifold. We handle these by expanding the bounding boxes into an N-hypercube, where N is the geometric dimension of the mesh, and points not on x=1 were being found and evaluated on subm because they are within the default tolerance of the mesh. The test passes after setting subm.tolerance=0.1.

I'm also not sure that the reasoning for your fix is correct; see my comment here #5077

Copy link
Copy Markdown
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Actually, if I run this on 8 ranks I do get the IndexError you described. The tests in test_swarm.py I mentioned don't actually input any points onto the rank which has zero cells, which is what is happening here. In that case I think this fix is good.

Copy link
Copy Markdown
Contributor

@leo-collins leo-collins left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think this is good. The same thing probably happens for the extruded case too - can you add a test for this (and a fix if necessary)?



@pytest.mark.parallel(nprocs=2)
def test_submesh_interpolate_3Dcell_2Dfacet_empty_rank_2_processes():
Copy link
Copy Markdown
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Actually, if I run this on 8 ranks I do get the IndexError you described. The tests in test_swarm.py I mentioned don't actually input any points onto the rank which has zero cells, which is what is happening here. In that case I think this fix is good.

Comment thread tests/firedrake/submesh/test_submesh_interpolate.py Outdated
Comment thread tests/firedrake/submesh/test_submesh_interpolate.py Outdated
Comment thread tests/firedrake/submesh/test_submesh_interpolate.py
hardik-corintis and others added 4 commits May 5, 2026 16:05
Co-authored-by: Leo Collins <leocollins511@gmail.com>
Co-authored-by: Leo Collins <leocollins511@gmail.com>
Co-authored-by: Leo Collins <leocollins511@gmail.com>
@hardik-corintis
Copy link
Copy Markdown
Author

Thanks @leo-collins for the quick review. I have applied all the inline suggestions and added the test and fix for the extruded meshes.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

BUG: IndexError in _pic_swarm_in_mesh on MPI ranks with empty cell partition

2 participants