Facet integrals on RT element do not work on manifolds

Issue #22 resolved
Jan Blechta created an issue

Originally described by @chris_richardson and @dnarnold in this QA thread. This issue seems distinct from other RT issues #10 (which has not been properly resolved), #14, #21 (probably duplicate of #10). Example bellow fails with message error : 'cell_orientation' was not declared in this scope in foo_integral::tabulate_tensor (when using quadrature representation). With other representations it mostly fails gently during code generation or UFL analysis.

from dolfin import *

parameters['form_compiler']['representation'] = 'quadrature'
#parameters['form_compiler']['representation'] = 'uflacs'
#parameters['form_compiler']['representation'] = 'tensor'

# Create unit square embedded in 3D
mesh = Mesh()
me = MeshEditor()
me.open(mesh, 2, 3)
me.init_vertices(4)
me.init_cells(2)
me.add_vertex(0, 0., 0., 0.)
me.add_vertex(1, 1., 0., 0.)
me.add_vertex(2, 0., 1., 0.)
me.add_vertex(3, 1., 1., 0.)
me.add_cell(0, 0, 1, 3)
me.add_cell(1, 0, 2, 3)
me.close()
mesh.order()
global_normal = Expression(("0.", "0.", "1."))
mesh.init_cell_orientations(global_normal)

# Uncomment this to make it working (with quadrature repre)
#mesh = UnitSquareMesh(1, 1)

RT = FunctionSpace(mesh, "RT", 1)
n = FacetNormal(mesh)
u = TrialFunction(RT)
v = TestFunction(RT)
w = Function(RT)
assemble(dot(n, v)*ds) # Fails
assemble(dot(u, v)*ds) # Fails
assemble(dot(w, v)*ds) # Fails
assemble(dot(w, w)*ds) # Fails
assemble(dot(n('+'), v('+'))*dS) # Fails
assemble(dot(u('+'), v('+'))*dS) # Fails
assemble(dot(w('+'), v('+'))*dS) # Fails
assemble(dot(w('+'), w('+'))*dS) # Fails

Comments (5)

  1. Marie Elisabeth Rognes

    The problem is a missing argument in the tabulate_tensor signature. Will fix.

    Here is pure UFL to reproduce

    triangle = Cell("triangle", 3)
    RT = FiniteElement("RT", triangle, 1)
    n = FacetNormal(triangle)
    u = TrialFunction(RT)
    v = TestFunction(RT)
    w = Coefficient(RT)
    L = dot(n, v)*ds # Fails
    
  2. Log in to comment