Issues when trying to set different BCs along different boundaries

Hello,

I’m trying to solve this PDE:

\nabla^2 u = u - b, \,\,\, u\in \bf{R} \, ;

in a 2D domain \Omega, with Robin boundary conditions on \partial \Omega:

\left( \nabla u \cdot \hat{n} + c u \right)|_{\partial \Omega} = 0 .

The variational form is thus:

\int_\Omega \left( \nabla u \cdot \nabla v + u v \right) \, \text{d}^2 {\bf{x}} + c \int_{\partial \Omega} u v \, \text{d} {\bf{s}} = b \int_\Omega v \, \text{d}^2 {\bf{x}} .

I know how to solve this for a constant c over the boundary, i.e., using the pre-defined measures:

a = (dot(grad(u), grad(v)) + u*v)*dx + Constant(c)*u*v*ds
L = Constant(b)*v*dx

However, I’m getting nonsense when I try to set different c values for different boundaries. Here’s a MWE:

import dolfin as fem

# simple mesh
mesh = fem.UnitSquareMesh(20, 20)

# horizontal boundary
class HBoundary(fem.SubDomain):
    def __init__(self, Y0):
        self.Y0 = Y0
        super().__init__()
    def inside(self, x, on_boundary):
        return on_boundary and fem.near(x[1], self.Y0, eps= 1.e-12)

# vertical boundary
class VBoundary(fem.SubDomain):
    def __init__(self, X0):
        self.X0 = X0
        super().__init__()
    def inside(self, x, on_boundary):
        return on_boundary and fem.near(x[0], self.X0, eps= 1.e-12)

# a MeshFunction over facets
boundfn = fem.MeshFunction('size_t', mesh, mesh.topology().dim()-1)
boundfn.set_all(999) # dummy value for inner facets

# mark boundaries
left, right = VBoundary(0.), VBoundary(1.)
bottom, top = HBoundary(0.), HBoundary(1.)

bottom.mark(boundfn, 0)
right.mark(boundfn,  1)
top.mark(boundfn,    2)
left.mark(boundfn,   3)

# space, trial and test functions
V = fem.FunctionSpace(mesh, 'CG', 1)

u = fem.TrialFunction(V)
v = fem.TestFunction(V)

# redefine the ds measure
ds = fem.Measure('ds', domain= mesh, subdomain_data= boundfn)

# the variational form, volume terms
a = (fem.dot(fem.grad(u), fem.grad(v)) + u*v)*fem.dx
b = fem.Constant(1.); L = b*v*fem.dx

# the boundary terms
c = fem.Constant(1.) # for simplicity; what matters is the following
for i in sorted(list(set(boundfn.array())))[:-1]: # not to include inner facets
    a += c*u*v*ds(i)

# if instead of using the previous loop, we use this:
#a += c*u*v*ds
# it works; however this will not allow different c values

# further, setting the boundary terms for all
# but one of the boundaries produces the expected solution
# for c = 0 along the unspecified boundary term (3 here, left boundary)
#a += c*u*v*ds(0)
#a += c*u*v*ds(1)
#a += c*u*v*ds(2)

# compute solution
m = fem.Function(V)
fem.solve(a == L, m)

I am lost here, I will appreciate your input

Cheers,
Miguel

This looks like a bug in UFLACS. If I switch the form representation to quadrature (deprecated) or TSFC (which you need to manually install on top of a standard FEniCS distribution), I don’t get this strange behavior:

import dolfin as fem
fem.parameters["form_compiler"]["representation"] = "quadrature"
#fem.parameters["form_compiler"]["representation"] = "tsfc"

Hello,
I tried that but I still get a nonsensical solution, apart from the deprecation warnings:

    Ignoring precision in integral metadata compiled using quadrature representation. Not implemented.

/usr/lib/python3/dist-packages/ffc/quadrature/quadraturerepresentation.py:58: QuadratureRepresentationDeprecationWarning: 
*** ===================================================== ***
*** FFC: quadrature representation is deprecated! It will ***
*** likely be removed in 2018.2.0 release. Use uflacs     ***
*** representation instead.                               ***
*** ===================================================== ***
  issue_deprecation_warning()

    Ignoring precision in integral metadata compiled using quadrature representation. Not implemented.
    Ignoring precision in integral metadata compiled using quadrature representation. Not implemented.
    Ignoring precision in integral metadata compiled using quadrature representation. Not implemented.
    Ignoring precision in integral metadata compiled using quadrature representation. Not implemented.
    Ignoring precision in integral metadata compiled using quadrature representation. Not implemented.

/usr/lib/python3/dist-packages/ffc/jitcompiler.py:234: QuadratureRepresentationDeprecationWarning: 
*** ===================================================== ***
*** FFC: quadrature representation is deprecated! It will ***
*** likely be removed in 2018.2.0 release. Use uflacs     ***
*** representation instead.                               ***
*** ===================================================== ***
  issue_deprecation_warning()
/usr/lib/python3/dist-packages/ffc/quadrature/quadraturerepresentation.py:58: QuadratureRepresentationDeprecationWarning: 
*** ===================================================== ***
*** FFC: quadrature representation is deprecated! It will ***
*** likely be removed in 2018.2.0 release. Use uflacs     ***
*** representation instead.                               ***
*** ===================================================== ***
  issue_deprecation_warning()
/usr/lib/python3/dist-packages/ffc/quadrature/quadraturerepresentation.py:58: QuadratureRepresentationDeprecationWarning: 
*** ===================================================== ***
*** FFC: quadrature representation is deprecated! It will ***
*** likely be removed in 2018.2.0 release. Use uflacs     ***
*** representation instead.                               ***
*** ===================================================== ***
  issue_deprecation_warning()
/usr/lib/python3/dist-packages/ffc/quadrature/quadraturerepresentation.py:58: QuadratureRepresentationDeprecationWarning: 
*** ===================================================== ***
*** FFC: quadrature representation is deprecated! It will ***
*** likely be removed in 2018.2.0 release. Use uflacs     ***
*** representation instead.                               ***
*** ===================================================== ***
  issue_deprecation_warning()
/usr/lib/python3/dist-packages/ffc/quadrature/quadraturerepresentation.py:58: QuadratureRepresentationDeprecationWarning: 
*** ===================================================== ***
*** FFC: quadrature representation is deprecated! It will ***
*** likely be removed in 2018.2.0 release. Use uflacs     ***
*** representation instead.                               ***
*** ===================================================== ***
  issue_deprecation_warning()
/usr/lib/python3/dist-packages/ffc/quadrature/quadraturerepresentation.py:58: QuadratureRepresentationDeprecationWarning: 
*** ===================================================== ***
*** FFC: quadrature representation is deprecated! It will ***
*** likely be removed in 2018.2.0 release. Use uflacs     ***
*** representation instead.                               ***
*** ===================================================== ***
  issue_deprecation_warning()
/usr/lib/python3/dist-packages/ffc/jitcompiler.py:234: QuadratureRepresentationDeprecationWarning: 
*** ===================================================== ***
*** FFC: quadrature representation is deprecated! It will ***
*** likely be removed in 2018.2.0 release. Use uflacs     ***
*** representation instead.                               ***
*** ===================================================== ***
  issue_deprecation_warning()

I’m trying to avoid building over FEniCS if possible as I don’t enjoy dealing with installations and I’m trying to keep my program simple for a future release. This forced me to try a different solution. This is what I came up with:

import dolfin as FEM

# simple mesh
mesh = FEM.UnitSquareMesh(52, 52)

# horizontal boundary
class HBoundary(FEM.SubDomain):
    def __init__(self, Y0):
        self.Y0 = Y0
        super().__init__()
    def inside(self, x, on_boundary):
        return on_boundary and FEM.near(x[1], self.Y0, eps= 1.e-12)

# vertical boundary
class VBoundary(FEM.SubDomain):
    def __init__(self, X0):
        self.X0 = X0
        super().__init__()
    def inside(self, x, on_boundary):
        return on_boundary and FEM.near(x[0], self.X0, eps= 1.e-12)

# a MeshFunction over facets
boundfn = FEM.MeshFunction('size_t', mesh, mesh.topology().dim()-1)
boundfn.set_all(999) # dummy value that will remain assigned to the inner facets

# mark boundaries
left, right = VBoundary(0.), VBoundary(1.)
bottom, top = HBoundary(0.), HBoundary(1.)

bottom.mark(boundfn, 0)
right.mark(boundfn,  1)
top.mark(boundfn,    2)
left.mark(boundfn,   3)

# redefine ds
ds = FEM.Measure('ds', domain= mesh, subdomain_data= boundfn)

# a UserExpression subclass for the values on the boundary
class ValuesOnBoundary(FEM.UserExpression):
    def __init__(self, mesh, boundfn, list_of_values, **kwargs):
        self._cells = tuple(FEM.cells(mesh))
        self._boundfn = boundfn.array()
        self._values = list_of_values
        super().__init__(**kwargs)
    def eval_cell(self, values, x, cell):
        values[0] = self._values[self._boundfn[self._cells[cell.index].entities(1)[cell.local_facet]]]
    def value_shape(self):
        return ()

c_values = [1., 1., 1., 1.]
c = ValuesOnBoundary(mesh, boundfn, c_values, degree= 0)

# space, trial and test functions
V = FEM.FunctionSpace(mesh, 'CG', 1)

u = FEM.TrialFunction(V)
v = FEM.TestFunction(V)

# the variational form
a = (FEM.dot(FEM.grad(u), FEM.grad(v)) + u*v)*FEM.dx + c*u*v*ds
b = FEM.Constant(1.); L = b*v*FEM.dx

# compute solution
m = FEM.Function(V)
FEM.solve(a == L, m)

# visualization
FEM.plot(m)

Now the values on the boundary for my Robin conditions are handled by a subclass of UserExpression which I can instantiate with whatever c_values I may need and use that as a coefficient for the variational form like c*u*v*ds. This works. I tried it for meshes with over 250 different boundaries and it works perfectly. One thing I don’t know if this is efficient or if it will work in parallel (I’m developing my code on colab.google so I haven’t tried) so I would appreciate any thoughts on that

Cheers,
Miguel

You can shut off the deprecation warnings with the following:

import warnings
from ffc.quadrature.deprecation \
    import QuadratureRepresentationDeprecationWarning
warnings.simplefilter("ignore", QuadratureRepresentationDeprecationWarning)

TSFC is actually pretty easy to use, and doesn’t require re-building FEniCS. You just need to install a few extra things with pip3 for it to work:

pip3 install git+https://github.com/blechta/tsfc.git@2018.1.0
pip3 install git+https://github.com/blechta/COFFEE.git@2018.1.0
pip3 install git+https://github.com/blechta/FInAT.git@2018.1.0
pip3 install singledispatch networkx pulp

As for optimizing UserExpression, you can code your own C++ subclasses of Expression and use them in Python scripts, as in this example

although I’m not sure that counts as “trying to keep my program simple”.

I’ll check out tsfc, thanks for the tips on how to build. I know writing directly c++ code helps with performance, so I’ll give that a try if this approach turns awful for larger or 3D meshes.
Now this is a bit off-topic but, do you know if there’s anything I’d need to change for my solution to work in parallel? (I can’t check right now because the company where I work has a ridiculous proxy that blocks everything, including the anaconda servers so no local FEniCS installation, thus I run everything in colab.google). I know that the dofs are distributed to different processes in parallel so, for example I have this function that assigns the solution to the equation above as initial condition for the timestepping part, and for that I need to access the dofs and use m.vector().get_local(), do the cell and facet info also get distributed so I’d need something similar?

Cheers,
Miguel

For some purposes, you may need to manually update the ghost DoF values after modifying local DoF arrays, e.g.,

fem.as_backend_type(m.vector()).vec().ghostUpdate()

using the (default) PETSc backend. In parallel, meshes are partitioned by cells, and I believe these are just indexed separately on each MPI task. Anecdotally, the last version of your code with the UserExpression ran fine for me in parallel. For some cases, you may need to tinker with the ghost mode, e.g.,

fem.parameters["ghost_mode"]="shared_facet"

is necessary when assembling over interior facets with the dS measure in parallel. Facets on the external boundary are necessarily each associated to a unique cell, though, and should always be local.