Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add functionality to print patch statistics in ASMStarPC and friends #3875

Open
wants to merge 1 commit into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
26 changes: 25 additions & 1 deletion firedrake/preconditioners/asm.py
Original file line number Diff line number Diff line change
Expand Up @@ -33,7 +33,7 @@ def initialize(self, pc):
# Extract function space and mesh to obtain plex and indexing functions
V = get_function_space(dm)

# Obtain patches from user defined funtion
# Obtain patches from user defined function
ises = self.get_patches(V)
# PCASM expects at least one patch, so we define an empty one on idle processes
if len(ises) == 0:
Expand Down Expand Up @@ -103,6 +103,30 @@ def get_patches(self, V):

def view(self, pc, viewer=None):
self.asmpc.view(viewer=viewer)
self.prefix = pc.getOptionsPrefix() + self._prefix

if viewer is not None:
opts = PETSc.Options(self.prefix)
print_statistics = opts.getBool("print_patch_statistics", default=False)

if print_statistics:
from mpi4py import MPI

dm = pc.getDM()
V = get_function_space(dm)

# Obtain patches from user defined function
ises = self.get_patches(V)
Copy link
Contributor

Choose a reason for hiding this comment

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

It feels slightly extravagant to completely recompute the patch indices just for the sake of this monitor, but it is true that we don't stash these anywhere convenient for us to find. Is it worth stashing them in __init__ so that we don't need to recompute them here?

Copy link
Contributor Author

Choose a reason for hiding this comment

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

I agree it's a bit wasteful. I agree we could stash them in __init__ if the option is present. I'll make this change.

Copy link
Contributor

Choose a reason for hiding this comment

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

Isn't the correct approach here to compute the statistics during __init__ if the options are present? Then just stash the statistics to save memory.

# PCASM expects at least one patch, so we define an empty one on idle processes
if len(ises) == 0:
ises = [PETSc.IS().createGeneral(numpy.empty(0, dtype=IntType), comm=PETSc.COMM_SELF)]

max_local_patch = max(is_.getSize() for is_ in ises)
min_local_patch = min(is_.getSize() for is_ in ises)
max_global_patch = pc.comm.tompi4py().allreduce(max_local_patch, op=MPI.MAX)
min_global_patch = pc.comm.tompi4py().allreduce(min_local_patch, op=MPI.MIN)

viewer.printfASCII(f"Minimum / maximum patch sizes: {min_global_patch} / {max_global_patch}\n")

def update(self, pc):
# This is required to update an inplace ILU factorization
Expand Down
2 changes: 2 additions & 0 deletions tests/regression/test_star_pc.py
Original file line number Diff line number Diff line change
Expand Up @@ -41,6 +41,7 @@ def test_star_equivalence(problem_type, backend):

star_params = {"mat_type": "aij",
"snes_type": "ksponly",
"snes_view": None,
"ksp_type": "richardson",
"pc_type": "mg",
"pc_mg_type": "multiplicative",
Expand All @@ -50,6 +51,7 @@ def test_star_equivalence(problem_type, backend):
"mg_levels_ksp_max_it": 1,
"mg_levels_pc_type": "python",
"mg_levels_pc_python_type": "firedrake.ASMStarPC",
"mg_levels_pc_star_print_patch_statistics": None,
"mg_levels_pc_star_construct_dim": 0}

comp_params = {"mat_type": "aij",
Expand Down
Loading