Does the test function(s) vanish when imposing Dirichlet BCs with mixed finite elements?

As far as I understand, in a FE problem the test function vanishes on the boundaries where Dirichlet BCs are imposed.

Now, given a sytem of PDEs on a domain \Omega for two functions u({\bf x}) and v({\bf x}) with Dirichlet boundary conditions for u and v on some parts of the boundary, \partial\Omega_A and \partial \Omega_B, respectively, which I solve with mixed FEs

P_V = VectorElement('P', triangle, 2)
P_Q = FiniteElement('P', triangle, 1)
element = MixedElement([P_V, P_Q])
VQ = FunctionSpace(mesh, element)

uv = Functions(VQ)
nu_u, nu_v = TestFunctions(VQ)
Fu = ...
Fv = ... 
F=Fu+Fv

bc_u = DirichletBC(VQ.sub(0), Constant(1.), partial_omega_A)
bc_v = DirichletBC(VQ.sub(1), Constant(2.), partial_omega_B)
bc = [bc_u, bc_v]

solve(F==0, uv, bc)

Do \nu_u, \nu_v vanish on some parts of \partial \Omega? If so, which \nu vanishes on which part?

Thanks

See for instance Application of Dirichlet boundary conditions — FEniCS Tutorial @ Sorbonne for a few ways of imposing Dirichlet BCs in practice.

Thank you, but I don’t see an answer to my question there.

The referenced text describes how the mathematical phrasing “vanishes” is handled for either \nu_u or \nu_v on an implementational level.
So yes, If you apply boundary conditions in a mixed space, DOLFINx will remove contributions from test-functions on those boundaries.

See for instance: Manually apply Dirichlet BC after assembly for mixed problem - #4 by dokken for an illustration

I don’t see this in the test, at least not explicitly stated.

In order to get a clear answer, may you please reply to the question in the original post?

Thank you

When the system is assembled with a DirichletCondition, either all columns corresponding to the dof of the DirichletBCs is set to the identity column. This is effectively the same as setting nu to zero when assembling over those dofs.

If you want to do this symmetrically, one uses the lifting procedure, which is also explained in the aforementioned document.

Yes, but my question in the original post was a

Let me clarify: there point here is that there may be different BCs for u and v. For example, if I set some Dirichlet BCs for u on \partial \Omega _A and some other BCs for v on \partial \Omega_B, where \partial \Omega_A \cap \partial \Omega_B = \emptyset, are these four statements true?

  1. \nu_u = 0 on \partial \Omega_A,
  2. \nu_v = 0 on \partial \Omega_B,
  3. In general, it is not true that \nu_v = 0 on \partial \Omega_A
  4. In general, it is not true that \nu_u = 0 on \partial \Omega_B

If these are true, are they true for both linear and nonlinear variational problems ? (in your answer you talk about ‘matrix’ so we would better clarify).
Is this a mathematical property independent of the implementation (DolfinX , legacy Dolfin) ?

Thank you

As I mentioned earlier, it is only the test function for the space that you apply the bc to that vanishes.
I.e. the statements you are mentioning are true.

This is then enforced in a linear or non-linear problem by modifying rows and columns of the matrix, as described earlier.

The property that a TestFunction disappears when a Dirichlet condition is applied is up to the implementation, but either implementation (Legacy Dolfin or DOLFINx) will enforce this in some way or another.

You can easily test this with:

import numpy as np
import dolfin


class Left(dolfin.SubDomain):
    def inside(self, x, on_boundary):
        return on_boundary and dolfin.near(x[0], 0.0)


class Right(dolfin.SubDomain):
    def inside(self, x, on_boundary):
        return on_boundary and dolfin.near(x[0], 1.0)


left = Left()
right = Right()


mesh = dolfin.UnitIntervalMesh(4)
facets = dolfin.MeshFunction("size_t", mesh, 0)
left.mark(facets, 1)
right.mark(facets, 2)


el = dolfin.FiniteElement("Lagrange", mesh.ufl_cell(), 1)
me = dolfin.MixedElement([el, el])
V = dolfin.FunctionSpace(mesh, me)
u_a, u_b = dolfin.TrialFunctions(V)
v_a, v_b = dolfin.TestFunctions(V)

dx = dolfin.dx
ds = dolfin.Measure("ds", domain=mesh, subdomain_data=facets)
a = dolfin.inner(u_a, v_a) * dx + dolfin.inner(u_b, v_b) * dx

bc_a = dolfin.DirichletBC(V.sub(0), dolfin.Constant(1), left)
bc_b = dolfin.DirichletBC(V.sub(1), dolfin.Constant(1), right)
bcs = [bc_a, bc_b]

A = dolfin.assemble(a)
[bc.apply(A) for bc in bcs]

# A2 should be equivalent to A
a2 = a + dolfin.inner(u_a, v_b) * ds(2) + dolfin.inner(u_b, v_a) * ds(1)
A2 = dolfin.assemble(a2)
[bc.apply(A2) for bc in bcs]

np.testing.assert_allclose(A.array(), A2.array())
print(A.array() - A2.array())
# A3 should be different from A
a3 = a + dolfin.inner(u_a, v_b) * ds(1) + dolfin.inner(u_b, v_a) * ds(2)
A3 = dolfin.assemble(a3)
[bc.apply(A3) for bc in bcs]
print(A.array() - A3.array())

yielding

[[0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]
 [0. 0. 0. 0. 0. 0. 0. 0. 0. 0.]]
[[ 0. -1.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0.  0.  0.]
 [ 0.  0.  0.  0.  0.  0.  0.  0. -1.  0.]]

That is crystal clear, thank you. The code also allowed me to learn some new syntax which will be useful. :slightly_smiling_face: