Error in assembling a variable on interior facets

Issue #709 resolved
Corrado Maurini created an issue

It seems that there is a bug (maybe in UFL?) that prevents to assemble a UFL variable on the interior facets dS. Below a minimal code to reproduce the issue and the error I get. This is with dev, git commit 937568dce5b697aacac2df8e4addc203269f64a0.

from fenics import *
mesh = UnitSquareMesh(2,2)
element = FiniteElement('CG',triangle,1)
V = FunctionSpace(mesh,element)
u = Function(V)
# All of this works
print(assemble(variable(u)*dx))
print(assemble(variable(u)*ds))
print(assemble(u*dS))
# This gives an error
print(assemble(variable(u)*dS))

Output:

0.0
0.0
0.0
Calling FFC just-in-time (JIT) compiler, this may take some time.
Traceback (most recent call last):
  File "tmp.py", line 11, in <module>
    print(assemble(variable(u)*dS))
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/dolfin/fem/assembling.py", line 192, in assemble
    dolfin_form = _create_dolfin_form(form, form_compiler_parameters)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/dolfin/fem/assembling.py", line 66, in _create_dolfin_form
    function_spaces=function_spaces)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/dolfin/fem/form.py", line 94, in __init__
    mpi_comm=mesh.mpi_comm())
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/dolfin/compilemodules/jit.py", line 65, in mpi_jit
    return local_jit(*args, **kwargs)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/dolfin/compilemodules/jit.py", line 124, in jit
    result = ffc.jit(ufl_object, parameters=p)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/ffc/jitcompiler.py", line 201, in jit
    module = jit_build_with_instant(ufl_object, module_name, parameters)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/ffc/jitcompiler.py", line 98, in jit_build_with_instant
    code_h, code_c = jit_generate(ufl_object, module_name, parameters)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/ffc/jitcompiler.py", line 62, in jit_generate
    parameters=parameters, jit=True)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/ffc/compiler.py", line 154, in compile_form
    analysis = analyze_forms(forms, parameters)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/ffc/analysis.py", line 60, in analyze_forms
    parameters) for form in forms)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/ffc/analysis.py", line 60, in <genexpr>
    parameters) for form in forms)
  File "/opt/HPC/fenics-dev-petsc-stable/lib/python2.7/site-packages/ffc/analysis.py", line 142, in _analyze_form
    form_data = compute_form_data(form)
  File "build/bdist.macosx-10.11-intel/egg/ufl/algorithms/compute_form_data.py", line 285, in compute_form_data
  File "build/bdist.macosx-10.11-intel/egg/ufl/algorithms/apply_restrictions.py", line 201, in apply_restrictions
  File "build/bdist.macosx-10.11-intel/egg/ufl/algorithms/map_integrands.py", line 60, in map_integrand_dags
  File "build/bdist.macosx-10.11-intel/egg/ufl/algorithms/map_integrands.py", line 39, in map_integrands
  File "build/bdist.macosx-10.11-intel/egg/ufl/algorithms/map_integrands.py", line 46, in map_integrands
  File "build/bdist.macosx-10.11-intel/egg/ufl/algorithms/map_integrands.py", line 60, in <lambda>
  File "build/bdist.macosx-10.11-intel/egg/ufl/corealg/map_dag.py", line 35, in map_expr_dag
  File "build/bdist.macosx-10.11-intel/egg/ufl/corealg/map_dag.py", line 84, in map_expr_dags
TypeError: variable() takes exactly 3 arguments (4 given)
Abort trap: 6

Comments (3)

  1. Lawrence Mitchell

    This is a UFL bug in the restriction propagator. Try this patch?

    diff --git a/ufl/algorithms/apply_restrictions.py b/ufl/algorithms/apply_restrictions.py
    index 9bacdf2..712ac48 100644
    --- a/ufl/algorithms/apply_restrictions.py
    +++ b/ufl/algorithms/apply_restrictions.py
    @@ -89,7 +89,7 @@ class RestrictionPropagator(MultiFunction):
         cell_avg = _require_restriction
         facet_avg = _ignore_restriction
    
    -    def variable(self, o, op):
    +    def variable(self, o, op, label):
             "Strip variable."
             return op
    
  2. Log in to comment