interior facet integrals do not work on RT element

Issue #14 resolved
Jan Blechta created an issue
from dolfin import *

mesh = UnitSquareMesh(4, 4)
RT1 = FunctionSpace(mesh, 'RT', 1)

u = Function(RT1)
v, w = TrialFunction(RT1), TestFunction(RT1)

assemble(dot(u, u)('+')*dS) # Pass
assemble(dot(u, w)('+')*dS) # Fail
assemble(dot(v, w)('+')*dS) # Fail

Both fails with

error: ‘detJ’ was not declared in this scope

Comments (9)

  1. Colin Cotter

    The following (which arises from reconstructing fluxes for DG advection) also fails:

    n = FacetNormal(mesh)
    
    V1 = FunctionSpace(mesh,"RT",1)
    v = TestFunction(V1)
    u = TrialFunction(V1)
    
    mat = assemble(inner(v('+'),n('+'))*inner(u('+'),n('+'))*dS)
    
  2. Marie Elisabeth Rognes

    This is a bug in tensor representation, in other words it runs with quadrature representation:

    form_compiler_parameters={"representation": "quadrature"}
    assemble(..., form_compiler_parameters=form_compiler_parameters)
    

    Also, here is the pure UFL example to help with debugging

    RT1 = FiniteElement('RT', triangle, 1)
    u = Coefficient(RT1)
    w = TestFunction(RT1)
    L = dot(u, w)('+')*dS # Fail
    
  3. Marie Elisabeth Rognes

    Good point. Fixed in this:

    commit f09815da164cd3741bfd00f3ef5c36c4563b5286
    Author: Marie E. Rognes <meg@simula.no>
    Date:   Mon Mar 10 21:28:12 2014 +0100
    
        Implement restrictions for monomial determinants.
    
  4. Prof Garth Wells

    Thanks - it you just put the commit ID in the message, Bitbucket is clever enough to automatically make it a link to the change.

  5. Log in to comment