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

LiteralFloat with value Zero not supported in compilation of integral #699

Closed
jorgensd opened this issue May 27, 2024 · 7 comments · Fixed by #701
Closed

LiteralFloat with value Zero not supported in compilation of integral #699

jorgensd opened this issue May 27, 2024 · 7 comments · Fixed by #701

Comments

@jorgensd
Copy link
Member

jorgensd commented May 27, 2024

MWE:

import basix.ufl
import ufl

c_el = basix.ufl.element("Lagrange", "interval", 1, shape=(1,))
mesh = ufl.Mesh(c_el)
C = ufl.as_tensor([[0, 0], [0, 1]])

integrand = ufl.transpose(C)[0, 0]

dx = ufl.Measure("dx", domain=mesh)
J = integrand * dx
forms = [J]

Traceback

Traceback (most recent call last):
  File "/usr/lib/python3.10/runpy.py", line 196, in _run_module_as_main
    return _run_code(code, main_globals, None,
  File "/usr/lib/python3.10/runpy.py", line 86, in _run_code
    exec(code, run_globals)
  File "/root/shared/ffcx/ffcx/__main__.py", line 14, in <module>
    sys.exit(main())
  File "/root/shared/ffcx/ffcx/main.py", line 75, in main
    code_h, code_c = compiler.compile_ufl_objects(
  File "/root/shared/ffcx/ffcx/compiler.py", line 118, in compile_ufl_objects
    code = generate_code(ir, options)
  File "/root/shared/ffcx/ffcx/codegeneration/codegeneration.py", line 48, in generate_code
    code_integrals = [integral_generator(integral_ir, options) for integral_ir in ir.integrals]
  File "/root/shared/ffcx/ffcx/codegeneration/codegeneration.py", line 48, in <listcomp>
    code_integrals = [integral_generator(integral_ir, options) for integral_ir in ir.integrals]
  File "/root/shared/ffcx/ffcx/codegeneration/C/integrals.py", line 42, in generator
    parts = ig.generate()
  File "/root/shared/ffcx/ffcx/codegeneration/integral_generator.py", line 168, in generate
    all_quadparts += self.generate_quadrature_loop(rule)
  File "/root/shared/ffcx/ffcx/codegeneration/integral_generator.py", line 273, in generate_quadrature_loop
    tensor_comp, intermediates_fw = self.generate_dofblock_partition(quadrature_rule)
  File "/root/shared/ffcx/ffcx/codegeneration/integral_generator.py", line 387, in generate_dofblock_partition
    block_quadparts, intermediate = self.generate_block_parts(
  File "/root/shared/ffcx/ffcx/codegeneration/integral_generator.py", line 506, in generate_block_parts
    var = fw if isinstance(fw, L.Symbol) else fw.array
AttributeError: 'LiteralFloat' object has no attribute 'array'

Issue first reported at
https://fenicsproject.discourse.group/t/asemble-scalar-error-literalfloat-no-array-attribute/14807?u=dokken

@jorgensd
Copy link
Member Author

The issue might be:

block_is_transposed = False # FIXME: Handle transposes for these block types

@jorgensd
Copy link
Member Author

Or it seems like some optimization of zeros is not applied:

import basix.ufl
import ufl

c_el = basix.ufl.element("Lagrange", "triangle", 1, shape=(2,))
mesh = ufl.Mesh(c_el)
C = ufl.as_tensor([[0, 0], [0, 1]])

int = ufl.transpose(C)[0, 0]
# int = C[0, 0]
dx = ufl.Measure("dx", domain=mesh)
J = int * dx

forms = [J]

This should be reduced to 0s, but it doesn't seem like it is.

@jorgensd jorgensd changed the title Cannot access index of transposed tensor. LiteralFloat is not supported in language May 27, 2024
@jorgensd jorgensd changed the title LiteralFloat is not supported in language LiteralFloat is untested and not fully supported May 27, 2024
@chrisrichardson
Copy link
Contributor

Is this a problem only because the value is zero. The kernel works with non-zero values in C[0,0]?

@jorgensd
Copy link
Member Author

Yes, it is only an issue when C[0,0] is 0.

@jorgensd jorgensd changed the title LiteralFloat is untested and not fully supported LiteralFloat with value Zero not supported in compilation of integral May 28, 2024
@michalhabera
Copy link
Contributor

michalhabera commented May 28, 2024

An optimization won't be applied in this case, that is expected (UFL does not zero parts of expression tree at the granularity of each component). For expression you should see the 0.0 multiplication in the generated code?

I think what happens is https://github.com/FEniCS/ffcx/blob/main/ffcx/codegeneration/lnodes.py#L98 triggers and fw becomes LiteralFloat(0.0) which is neither a temporary symbol nor it has .array. A quick fix might be to check for literals and access .value in that case. Or in float_product return 0.0 for zero literal?

@chrisrichardson
Copy link
Contributor

I don't think this will work, as the value 0.0 is still neither a Symbol, nor have an array property, and it will fail in the same place. We could cut out those lines in https://github.com/FEniCS/ffcx/blob/main/ffcx/codegeneration/lnodes.py#L98 - then it will just return Product() which will actually multiply by zero.

@michalhabera
Copy link
Contributor

I don't think this will work, as the value 0.0 is still neither a Symbol, nor have an array property, and it will fail in the same place. We could cut out those lines in https://github.com/FEniCS/ffcx/blob/main/ffcx/codegeneration/lnodes.py#L98 - then it will just return Product() which will actually multiply by zero.

Thats sounds good. It will avoid extra work to check for possibly unused quadrature tables.

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 a pull request may close this issue.

3 participants