[firedrake] Concatenation of function spaces?
Justin Chang
jychang48 at gmail.com
Fri Jul 24 12:13:52 BST 2015
Hi Lawrence,
Thanks for the responses (and same to you Fabio). I am attempting to
experiment with slip condition dirichlet boundary conditions for velocity
with these lines:
bcs1 = DirichletBC(W.sub(0).sub(0),
Expression(("2*pi*cos(2*pi*x[0])*sin(2*pi*x[1])*sin(2*pi*x[2])")), (1,2))
bcs2 = DirichletBC(W.sub(0).sub(1),
Expression(("2*pi*sin(2*pi*x[0])*cos(2*pi*x[1])*sin(2*pi*x[2])")), (3,4))
bcs3 = DirichletBC(W.sub(0).sub(2),
Expression(("2*pi*sin(2*pi*x[0])*sin(2*pi*x[1])*cos(2*pi*x[2])")), (5,6))
but I am getting these errors:
Traceback (most recent call last):
File "Compare_P2P0.py", line 66, in <module>
solver.solve()
File "<string>", line 2, in solve
File "/home/justin/Software/firedrake-deps/PyOP2/pyop2/profiling.py",
line 203, in wrapper
return f(*args, **kwargs)
File
"/home/justin/Software/firedrake-deps/firedrake/firedrake/variational_solver.py",
line 165, in solve
bc.apply(self._problem.u)
File "<string>", line 2, in apply
File "/home/justin/Software/firedrake-deps/PyOP2/pyop2/profiling.py",
line 203, in wrapper
return f(*args, **kwargs)
File "/home/justin/Software/firedrake-deps/firedrake/firedrake/bcs.py",
line 199, in apply
raise RuntimeError("%r defined on incompatible FunctionSpace!" % r)
RuntimeError: Coefficient(MixedElement(*[VectorElement('Lagrange',
Domain(Coefficient(VectorElement('Lagrange', Domain(Cell('tetrahedron', 3),
label=None, data='<data with id 140068500682768>'), 1, dim=3,
quad_scheme=None), 1)), 2, dim=3, quad_scheme=None),
FiniteElement('Discontinuous Lagrange',
Domain(Coefficient(VectorElement('Lagrange', Domain(Cell('tetrahedron', 3),
label=None, data='<data with id 140068500682768>'), 1, dim=3,
quad_scheme=None), 1)), 0, quad_scheme=None)], **{'value_shape': (4,) }),
7) defined on incompatible FunctionSpace!
Attached is the code. Any ideas on what this one may be?
Thanks,
Justin
On Fri, Jul 24, 2015 at 5:42 AM, Lawrence Mitchell <
lawrence.mitchell at imperial.ac.uk> wrote:
>
> > On 24 Jul 2015, at 11:25, Justin Chang <jychang48 at gmail.com> wrote:
> >
> > Thanks for the responses guys. However, I am running into issues
> regarding
> > Dirichlet bcs for P2/P0 elements. Attached is the code:
> >
> > I am trying to apply homogeneous boundary conditions in the P0 space, and
> > when I looked at the documentation online it said to use DirichletBC(Q,
> > <value>, <boundary_id>, method="geometric")
> > but I am getting these kinds of errors:
>
>
> Aha, this is a simple problem. When you define a MixedFunctionSpace and
> then want to apply boundary conditions on a part of the space, you must
> define the boundary condition on an indexed subspace, rather than the
> original space.
>
> That is:
>
> V = ...
> Q = ...
>
> W = V*Q
>
> bc = DirichletBC(W.sub(1), ...)
>
> rather than
>
> bc = DirichletBC(Q, ...)
>
> This is in the documentation (
> http://firedrakeproject.org/variational-problems.html#incorporating-boundary-conditions)
> but maybe not clearly enough.
>
> Making this change, and your code runs fine.
>
> I note in passing, that the DirichletBC constructor can take a list of
> boundary ids as well as just a single id. So, since all the boundary
> condition values are the same you can replace your six separate BC objects
> by:
>
> bcs = DirichletBC(W.sub(1), 0, (1, 2, 3, 4, 5, 6), method="geometric")
>
> Cheers,
>
> Lawrence
>
> _______________________________________________
> firedrake mailing list
> firedrake at imperial.ac.uk
> https://mailman.ic.ac.uk/mailman/listinfo/firedrake
>
>
-------------- next part --------------
HTML attachment scrubbed and removed
-------------- next part --------------
A non-text attachment was scrubbed...
Name: Compare_P2P0.py
Type: text/x-python
Size: 3872 bytes
Desc: not available
URL: <http://mailman.ic.ac.uk/pipermail/firedrake/attachments/20150724/19aa226f/attachment.py>
More information about the firedrake
mailing list