Refactor preprocessing and its use in ffc

Issue #36 resolved
Martin Sandve Alnæs created an issue

For several reasons the preprocessing needs to change significantly.

One reason is to reduce the overhead of symbolic computations in dolfin programs by minimizing the amount of processing done before a jit signature is computed and the information needed by dolfin assembly is available. Preprocessing in form.compute_form_data() should only include this minimal amount of symbolic processing. This symbolic overhead is particularly damaging in programs relying on large amounts of form generation, such as dolfin-adjoint programs, with non-trivial forms.

Another reason is that form compilers (read: uflacs) may want to do things slightly differently from what vanilla ffc does and some of the processing in preprocess() should therefore be moved to the form compiler analysis stage.

Comments (11)

  1. Martin Sandve Alnæs reporter

    Based on some source code searching, the minimal form data object needs:

    • Integration domain(s) in ordered list
    • Argument function spaces in ordered list
    • Original coefficient objects in ordered list
    • Subdomain data for each (domain x integral_type)
    • Signature

    Snippets of form_data usage in DOLFIN:

    domain = mesh.ufl_domain()
    form_data.subdomain_data[domain.label()]
    
    
    form_data.preprocessed_form.domains()
    
    function_spaces = [func.function_space() for func in form_data.original_arguments]
    
    assert all(isinstance(c, cpp.GenericFunction) for c in form_data.original_coefficients)
    coefficients = list(form_data.original_coefficients)
    
    domains = form_data.preprocessed_form.domains()
    assert len(domains) == 1, "This code assumes a single domain per form."
    domain, = domains
    subdomain_data = form_data.subdomain_data[domain.label()]
    

    Snippets of form_data usage in FFC if a pre-jitted module is found:

    jitcompiler.py:
            form_data = form.compute_form_data()
            jitobject = JITObject(form)
    jitobject.py:
        def signature(self):
            form_signature = self.form.form_data().signature
    
  2. Martin Sandve Alnæs reporter

    The best course of action seems to be

    1. Implement a new preprocess function doing the minimal amount of work
    2. Transition dolfin and ffc to use this function
    3. Implement a stage 2 preprocess function to produce the additional form data needed by ffc to compile a form
  3. Martin Sandve Alnæs reporter

    Here's an initial implementation of a new minimal preprocess function.

    from ufl import *
    from ufl.algorithms.analysis import extract_terminals
    from ufl.algorithms.signature import compute_form_signature
    
    def preprocess2(form):
    
        domains = form.domains()
    
        subdomain_data = {}
        for domain in domains:
            l = domain.label()
            subdomain_data[l] = {}
        for integral in form.integrals():
            l = integral.domain().label()
            t = integral.integral_type()
            data = subdomain_data[l].get(t)
            if data is None:
                data = integral.subdomain_data()
                subdomain_data[l][t] = data
            else:
                assert data.ufl_id() == integral.subdomain_data().ufl_id()
    
        terminals = extract_terminals(form)
        arguments = []
        coefficients = []
        for t in terminals:
            if isinstance(t, Argument):
                arguments.append(t)
            elif isinstance(t, Coefficient):
                coefficients.append(t)
    
        domains = sorted(domains, key=lambda x: x.label())
        arguments = sorted(arguments, key=lambda x: x.number())
        coefficients = sorted(coefficients, key=lambda x: x.count())
    
        domain_numbering = dict((d,i) for i,d in enumerate(domains))
        coefficient_numbering = dict((c,i) for i,c in enumerate(coefficients))
    
        # Temporary solution:
        function_replace_map = dict((c, c.reconstruct(count=i)) for c,i in coefficient_numbering.items())
        signature = compute_form_signature(form, function_replace_map)
        # Better solution:
        #signature = compute_form_signature2(form, domain_numbering, coefficient_numbering)
    
        return signature, domains, subdomain_data, arguments, coefficients
    

    I'd like to get the 1.4 release out first, then make dolfin.Mesh subclass ufl.Domain and see how that affects handling of domains here, then get this in place.

  4. Martin Sandve Alnæs reporter

    Preliminary work in

    https://bitbucket.org/fenics-project/ufl/branch/martinal/topic-redesign-form-preprocessing-mechanism

    This will be more merged into the Form class, making data available directly as form.coefficients(), form.arguments(), form.domains(), form.signature(), all computed when required and cached with the form object.

    Next steps:

    1. (DONE) Implement a modified ufl.preprocess to build a FormData object based on the data now available directly from Form.
    2. (DONE) Change ffc to call the new preprocess only if a cached module was not found in jit.
    3. (DONE) Change dolfin to use new form.* accessors instead of going via form_data.
    4. (DONE) Remove ufl.Form.[compute_]form_data().
  5. Martin Sandve Alnæs reporter

    Next: Test ufl/ffc branches in concert with fix to dolfin issue 306 when that is ready (after 1.4).

  6. Log in to comment