Invalid code is emmited

Issue #35 duplicate
Jan Blechta created an issue

Invalid code is generated for following DOLFIN code

from dolfin import *

parameters['form_compiler']['representation'] = 'uflacs'

mesh = UnitSquareMesh(3, 3)

V = VectorFunctionSpace(mesh, "Lagrange", 2)
P = FunctionSpace(mesh, "Discontinuous Lagrange", 0)
W = V*P

w  = Function(W)
w0 = Function(W)
v,  h  = split(w)
v0, h0 = split(w0)
v_, h_ = TestFunctions(W)

width = Function(P)


x = SpatialCoordinate(mesh)
F = as_matrix([[width, +x[0]*width.dx(1)],
               [ 0.0,                1.0]])
A = as_matrix([[1.0/width, -x[0]*width.dx(1)/width],
               [     0.0,                    1.0]])
detF = width
detA = 1.0/width
cofF = as_matrix([[             1.0,  0.0],
                  [-x[0]*width.dx(1), width]])
cofA = as_matrix([[                  1.0,      0.0],
                  [+x[0]*width.dx(1)/width, 1.0/width]])
n = FacetNormal(mesh)

Grad = lambda u: dot(grad(u), A)
Div  = lambda u: Grad(u)[..., Index(0), Index(0)]
class dX(object):
    """This is volumetric integral with built-in transformation."""
    def __init__(self, *args, **kwargs):
        self.dx = Measure('dx', *args, **kwargs)
    def __rmul__(self, integrand):
        return integrand*detF*self.dx
class dsX(object):
    """This is surface integral of vector field with built-in
    transformation and outter normal."""
    def __init__(self, *args, **kwargs):
        self.ds = Measure('ds', *args, **kwargs)
    def __rmul__(self, integrand):
        return inner(integrand, dot(cofF, n))*self.ds

mu    = Expression('1.0', domain=mesh)
rho_g = Expression('1.0', domain=mesh)
rho_t = Expression('1.0', domain=mesh)
g = Constant(9.80251)
idt = Constant(1.0)
gamma = Constant(0.2)

D = sym(Grad(v))
P = rho_g*(1.0-rho_g/rho_t)*g*h*h*0.5
I = Identity(2)
T_gs = -P*I + 2.0*mu*h*(Div(rho_g*v)/rho_g*I + D)
F = (   Div(rho_g*h*v)*h_
    + inner(idt*rho_g*(h*v-h0*v0) + Div(rho_g*h*outer(v, v)), v_)
    + inner(T_gs, Grad(v_))
    ) * dX() \
    + gamma*v_*dsX()

v, h = w.split()

problem = NonlinearVariationalProblem(F, w, [], derivative(F, w))
solver = NonlinearVariationalSolver(problem)
solver.solve()

C++ compiler stops with

In member function 'virtual void ffc_form_8567291d470ce8e0cbc2552effa459be5e98ee19_cell_integral_0_otherwise::tabulate_tensor(double*, const double* const*, const double*, int) const':
error: 'FE1_C2_D01' was not declared in this scope
error: 'FE1_C2_D01' was not declared in this scope

I will try to simplify the example later. Note that the Stokes problem on this element does not fail. Also note that without UFLACS it seems to work.

Comments (3)

  1. Jan Blechta reporter

    There are no second derivatives although it seems related. Problem is generating a derivative of piecewise constant:

    from dolfin import *
    parameters['form_compiler']['representation'] = 'uflacs'
    
    mesh = UnitSquareMesh(3, 3)
    
    P = FunctionSpace(mesh, "Discontinuous Lagrange", 0)
    p = TestFunction(P)
    assemble(p.dx(0)*dx)
    
  2. Log in to comment