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

Move grouping of integral from build_integral_data to group_form_integrals #305

Conversation

jorgensd
Copy link
Member

Move group to resolve #304 by grouping prior to pull-back.
This will resolve FEniCS/ffcx#716 as well.

I think this is a better fix than what I came up with in https://github.com/FEniCS/ufl/pull/92/files
as we now group the integrals in the function that is supposed to do grouping, and it will reduce the load of the other functions (pull-back, geometry_lowering etc.) since it happens at the start of compute_form_data.

@jorgensd
Copy link
Member Author

jorgensd commented Sep 19, 2024

Seems like one should use the IndexRenumberingTransform:

class IndexRenumberingTransformer(VariableRenumberingTransformer):
"""Index renumbering transformer.
This is a poorly designed algorithm. It is used in some tests,
please do not use for anything else.

@wence- do you have any ideas on how to improve it?

ufl/algorithms/domain_analysis.py Outdated Show resolved Hide resolved
@wence-
Copy link
Collaborator

wence- commented Sep 26, 2024

Seems like one should use the IndexRenumberingTransform:

class IndexRenumberingTransformer(VariableRenumberingTransformer):
"""Index renumbering transformer.
This is a poorly designed algorithm. It is used in some tests,
please do not use for anything else.

@wence- do you have any ideas on how to improve it?

It doesn't look that bad. I would forward-port the model to the MultiFunction case and use map_expr_dags instead.

What you want to know is "are two dags equal to each other, up to index relabelling?". This is achievable by visiting the dag, child before parent, and replacing free indices with a unique new index.

class IndexRelabeller(MultiFunction):
    def __init__(self):
        super().__init__()
        count = itertools.count()
        self.index_cache = defaultdict(lambda: Index(next(count)))

    expr = MultiFunction.reuse_if_untouched
   
    def multi_index(self, o):
        return type(o)(tuple(self.index_cache[i] if isinstance(i, Index) else i for i in o.indices()))

    def zero(self, o):
        fi = o.ufl_free_indices
        fid = o.ufl_free_index_dimensions
        new_indices = [self.index_cache[Index(i)].count() for i in fi]
        new_fi, new_fid = zip(*sorted(zip(new_indices, fid), key=lambda x: x[0]))
        return type(o)(o.ufl_shape, tuple(new_fi), tuple(new_fid))

map_expr_dags(IndexRelabeller(), whatever)

Would probably work, untested.

@jorgensd
Copy link
Member Author

Seems like one should use the IndexRenumberingTransform:

class IndexRenumberingTransformer(VariableRenumberingTransformer):
"""Index renumbering transformer.
This is a poorly designed algorithm. It is used in some tests,
please do not use for anything else.

@wence- do you have any ideas on how to improve it?

It doesn't look that bad. I would forward-port the model to the MultiFunction case and use map_expr_dags instead.

What you want to know is "are two dags equal to each other, up to index relabelling?". This is achievable by visiting the dag, child before parent, and replacing free indices with a unique new index.

class IndexRelabeller(MultiFunction):
    def __init__(self):
        super().__init__()
        count = itertools.count()
        self.index_cache = defaultdict(lambda: Index(next(count)))

    expr = MultiFunction.reuse_if_untouched
   
    def multi_index(self, o):
        return type(o)(tuple(self.index_cache[i] if isinstance(i, Index) else i for i in o.indices()))

    def zero(self, o):
        fi = o.ufl_free_indices
        fid = o.ufl_free_index_dimensions
        new_indices = [self.index_cache[Index(i)].count() for i in fi]
        new_fi, new_fid = zip(*sorted(zip(new_indices, fid), key=lambda x: x[0]))
        return type(o)(o.ufl_shape, tuple(new_fi), tuple(new_fid))

map_expr_dags(IndexRelabeller(), whatever)

Would probably work, untested.

Thanks Lawrence! I will test this tomorrow!

@jorgensd
Copy link
Member Author

Thanks for the suggestions @wence- . I have now implemented this in the code, replacing the old transformer.
The test reduces the number of integrals computed in integral data from 3 on main to 1.
This should give a significant performance boost for both Firedrake and DOLFINx for problems with multiple subdomain ids, as there will be less forms to apply attach_estimated_degrees, apply_function_pullbacks, apply_integral_scaling etc. on, as well as fewer kernels generated by TSFC/FFCx

@jorgensd jorgensd requested a review from mscroggs September 27, 2024 08:14
@garth-wells
Copy link
Member

@jorgensd is this ready, or does it need more feedback?

@jorgensd
Copy link
Member Author

jorgensd commented Oct 2, 2024

@jorgensd is this ready, or does it need more feedback?

I believe this is ready for merging. @jpdean have you tested it for the use-cases that spawned revisiting this issue?

@jpdean
Copy link
Member

jpdean commented Oct 2, 2024

Yes, this has resolved FEniCS/ffcx#716!

@jorgensd jorgensd added this pull request to the merge queue Oct 2, 2024
Merged via the queue into main with commit 3b85665 Oct 2, 2024
21 checks passed
@jorgensd jorgensd deleted the dokken/move-grouping-of-integrals-prior-to-building-integral-data branch October 2, 2024 15:41
francesco-ballarin added a commit to RBniCS/ufl4rom that referenced this pull request Oct 6, 2024
…ing.renumber_indices in FEniCS/ufl#305

Regolding only plain ufl and dolfinx tests. firedrake tests will be regolded after its ufl fork gets updated to include the above PR.
francesco-ballarin added a commit to RBniCS/ufl4rom that referenced this pull request Nov 23, 2024
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Projects
None yet
Development

Successfully merging this pull request may close these issues.

Uniqueness of integrands after computing form data FFCx generates multiple copies of the same kernel
5 participants