TSFC representation is not compatible with DOLFIN no-cell hack

Issue #144 resolved
Jan Blechta created an issue

demo_poisson1D-in-2D.py with parameters.form_compiler.representation = 'tsfc' fails with the message

[...]
  File "/home/fenics/local/lib/python3.5/site-packages/ffc/tsfc/tsfcrepresentation.py", line 45, in compute_integral_ir
    interface=ufc_interface)
  File "/home/fenics/local/lib/python3.5/site-packages/tsfc/driver.py", line 108, in compile_integral
    builder.set_coefficients(integral_data, form_data)
  File "/home/fenics/local/lib/python3.5/site-packages/tsfc/kernel_interface/ufc.py", line 98, in set_coefficients
    expressions = prepare_coefficients(coeffs, n, name, self.interior_facet)
  File "/home/fenics/local/lib/python3.5/site-packages/tsfc/kernel_interface/ufc.py", line 174, in prepare_coefficients
    elements = [create_element(coeff.ufl_element()) for coeff in coefficients]
  File "/home/fenics/local/lib/python3.5/site-packages/tsfc/kernel_interface/ufc.py", line 174, in <listcomp>
    elements = [create_element(coeff.ufl_element()) for coeff in coefficients]
  File "/home/fenics/local/lib/python3.5/site-packages/tsfc/finatinterface.py", line 177, in create_element
    raise ValueError("Don't know how to build element when cell is not given")
ValueError: Don't know how to build element when cell is not given

Comments (27)

  1. Jan Blechta reporter

    It's not just "some DOLFIN demos" - the same problem exists even in poisson demo.

    This is caused by DOLFIN no-cell hack in Expressions. Changing degree=2 to element=FiniteElement("P", mesh.ufl_cell(), 2) fixes the problem. Need to find a way how to make tsfc-repr compatible with the hack. One would expect that cell is already attached by UFL preprocessing. @martinal, any hints?

  2. Jan Blechta reporter

    @miklos1, any clue why this happens? Preprocessing goes through here which should take care of the problem but TSFC later complains. Any possibility that it works with unpreprocessed form?

  3. Miklós Homolya

    @blechta: I suspect this should be easy to solve if only I could inspect it in the debugger. But to do that, I first have to figure out how to install DOLFIN...

  4. Miklós Homolya

    OK, just reading the source code I suspect the problem is that the element map is never used in TSFC.

  5. Jan Blechta reporter

    But to do that, I first have to figure out how to install DOLFIN...

    Docker would be a good one-shot alternative.

  6. Jan Blechta reporter

    I can try fixing it. The question is

    • whether use just element_replace_map or function_replace_map,
    • whether to fix it just in UFC kernel interface and leave Firedrake pipeline intact,
    • whether coefficient numbering is correct now.

    Or should I leave it to you?

  7. Miklós Homolya

    I see UFLACS "just" (not considering that futzing about coefficient numbering) calls replace on the preprocessed integrand expression, and then everything just works later. That's why I didn't see anything related around FFC's fiatinterface.py.

    whether to fix it just in UFC kernel interface and leave Firedrake pipeline intact

    I don't think we have an API for kernel interface specific UFL munging.

    I'm trying to understand what happens in the UFLACS example before the ufl.replace call. I first thought it was there to undo some UFL renumbering, but looking at assert i == g.count() it looks more like a sanity check, since coefficient_numbering is completely redundant since coefficient_numbering[g] == g.count().

    If you have a patch that fixes this in DOLFIN, I can review it and test it in Firedrake. Meanwhile I'm trying to get a working DOLFIN installation...

  8. Miklós Homolya

    This patch seems to fix poisson and poisson1D-in-2D. @blechta: I don't know how to test all things in DOLFIN with parameters.form_compiler.representation = 'tsfc'.

  9. Jan Blechta reporter

    To be hones it is a mess. I found it difficult to get understand which members of FormData should be remapped. And UFLACS part demonstrates it too. Maybe there should be a clear mechanism in UFL to perform the mapping.

    @martinal, why was an application of replace maps removed from UFL and put to representations?

  10. Jan Blechta reporter

    To enforce optimize=False through env to test without optimizations (default is True).

  11. Miklós Homolya

    I'm happy to have optimisations on, that will test more things. I guess I can always disable that for debugging by editing some global default in a Python file, right?

  12. Jan Blechta reporter

    Sure, you can do parameters.form_compiler.optimize = False.

    I've taken a risk and merged the branch into master and next.

  13. Jan Blechta reporter
    cd test/unit
    DIJITSO_CACHE_DIR=$PWD FFC_FORCE_REPRESENTATION=tsfc py.test-3 python/multimesh/test_multimesh.py
    

    fails as expected but more informative message would be nice

    E         File "/vivid/home/jan/dev/fenics-master/src/ffc/ffc/representation.py", line 453, in _compute_integral_ir
    E           parameters)
    E         File "/vivid/home/jan/dev/fenics-master/src/ffc/ffc/tsfc/tsfcrepresentation.py", line 45, in compute_integral_ir
    E           interface=ufc_interface)
    E         File "/vivid/home/jan/dev/fenics-master/src/tsfc/tsfc/driver.py", line 93, in compile_integral
    E           integration_dim, entity_ids = lower_integral_type(fiat_cell, integral_type)
    E         File "/vivid/home/jan/dev/fenics-master/src/tsfc/tsfc/driver.py", line 352, in lower_integral_type
    E           basedim, extrdim = dim
    E       TypeError: 'int' object is not iterable
    
    DIJITSO_CACHE_DIR=$PWD FFC_FORCE_REPRESENTATION=tsfc py.test-3 python/fem/test_form.py -k test_derivative
    

    fails with

    E         File "/vivid/home/jan/dev/fenics-master/src/finat/finat/fiat_elements.py", line 24, in degree
    E           return self._element.degree()
    E       AttributeError: 'MixedElement' object has no attribute 'degree'
    

    This might be resolved with #143 (so far tested only with masters + this patch).

    Parallel tests

    DIJITSO_CACHE_DIR=$PWD FFC_FORCE_REPRESENTATION=tsfc mpirun -n 3 py.test-3
    

    are fine.

  14. Jan Blechta reporter

    I'd like to get a confidence that there's nothing wrong with coefficient numbering (after the map is applied). Do you understand the logic behind that, @miklos1?

    I guess that if it was wrong we might not see a crash in DOLFIN regression tests but FFC regression test should catch it.

  15. Miklós Homolya

    I'm confident that the TSFC does not change the ordering of coefficients, so it was correct so far, it shall stay that way. Yes, I also believe that the FFC regression tests ought to have caught such an issue.

  16. Log in to comment