[firedrake] Mapping nodal values from parent to extruded mesh
Miklos Homolya
m.homolya14 at imperial.ac.uk
Fri Apr 10 12:04:24 BST 2015
Hi,
facet_support_dofs is the same infrastructure we use to support
geometric boundary conditions.
I have implemented geometric boundary conditions for (unstructured)
quadrilaterals, it wasn't too complicated:
https://bitbucket.org/mapdes/fiat/pull-request/29/geometric-bcs-for-quadrilaterals/diff
For general tensor product elements it would be a bit more complicated,
as horizontal and vertical facets are in general different. For example,
a prism has quadrilaterals as vertical facets, but triangles as
horizontal facets.
Regards,
Miklos
On 10/04/15 11:23, David Ham wrote:
> Hi Tuomas,
>
> Concerning the sets of spaces which this will work for, I think that
> the current implementation probably works for any horizontal space but
> only for vertical CG. One could just try relaxing the test in that case.
>
> The reason for the restriction to vertical CG is that fs.bt_masks uses
> the topological association of nodes with mesh entities in order to
> work out which nodes are on the top or bottom of the cell. In order to
> allow DG in the vertical, it would be necessary to support the
> geometric definition (ie which basis functions do not vanish on the
> top/bottom).
>
> If you look at functionspace.py:78 you can see where the bottom and
> top masks are generated. This uses entity_closure_dofs() from FIAT. In
> order to use the geometric definition of dofs one would need to
> support using facet_support_dofs() (which is how BC maps are set up at
> functionspace.py:375). Currently FIAT TensorFiniteElement objects do
> support entity_closure_dofs() but nobody has done the legwork to get
> them to support facet_support_dofs.
>
> Concerning documentation, COFFEE is deplorably underdocumented
> (although the author is on this list so maybe this will change ;). The
> actual available AST nodes are only "documented" by reading the
> source:
> https://github.com/coneoproject/COFFEE/blob/master/coffee/base.py
> However, you are fundamentally writing a PyOP2 kernel, and the C api
> for that *is* documented at: http://op2.github.io/PyOP2/kernels.html
> so that hopefully helps somewhat.
>
> Otherwise, keep asking here or in IRC (IRC tends to get faster
> responses, at least between about 10am and 8pm UK time).
>
> BTW, any chance of seeing you at FEniCS '15?
>
> Cheers,
>
> David
>
>
>
> On 10 April 2015 at 02:44, Tuomas Karna <tuomas.karna at gmail.com
> <mailto:tuomas.karna at gmail.com>> wrote:
>
> Hi all,
>
> A couple of questions regarding those extruded mesh -> parent mesh
> copy operations,
>
> A while back I got the reverse 2d->3d copy working with a
> hard-coded pyop2 kernel (like in mesh extrusion). I'm not familiar
> with COFFEE syntax, is there documentation/examples somewhere?
>
> The 3d->2d copy routine in extrusion_extraction branch checks that
> the function space is CG. Does this method easily generalize to
> other spaces? DGxCG prisms seem to work OK. I'm also interested in
> using DGxDG and RTxCG/DG spaces.
>
>
> Cheers,
>
> Tuomas
>
>
> On 03/05/2015 10:24 AM, Tuomas Karna wrote:
>> Thanks David,
>>
>> This is great, I'll try doing the reverse operation.
>>
>> - Tuomas
>>
>> On 03/05/2015 08:20 AM, David Ham wrote:
>>> Hi Tuomas,
>>>
>>> This wasn't there this morning but I've implemented one of the
>>> cases (pulling out the top and bottom maps). The result is in
>>> the extrusion_extraction branches of both PyOP2 and Firedrake.
>>> The 2d->3d operation would be rather similar except that you'd
>>> have to interrogate the fiat_element on the extruded Function in
>>> order to determine which extruded nodes are "on top of" the 2d
>>> nodes you have. Unfortunately I won't have time to do that one
>>> soon (huge amounts to do before SIAM CSE next week) but feel
>>> free to have a try and complain when it doesn't work!
>>>
>>> Cheers,
>>>
>>> David
>>>
>>> On 4 March 2015 at 22:23, Tuomas Karna <tuomas.karna at gmail.com
>>> <mailto:tuomas.karna at gmail.com>> wrote:
>>>
>>> Hi all,
>>>
>>> I'd need to copy nodal values between fields on parent and
>>> extruded
>>> meshes. For example, copy 2d->3d (constant over vertical) or
>>> 3d->2d
>>> (extract surface/bottom level). The horizontal function
>>> space is the
>>> same. Is there an easy way to do this? Is there a map of
>>> extruded nodes
>>> somewhere?
>>>
>>> Thanks,
>>>
>>> Tuomas
>>>
>>>
>>>
>>> _______________________________________________
>>> firedrake mailing list
>>> firedrake at imperial.ac.uk <mailto:firedrake at imperial.ac.uk>
>>> https://mailman.ic.ac.uk/mailman/listinfo/firedrake
>>>
>>>
>>>
>>>
>>> --
>>> Dr David Ham
>>> Departments of Mathematics and Computing
>>> Imperial College London
>>>
>>> http://www.imperial.ac.uk/people/david.ham
>>>
>>>
>>> _______________________________________________
>>> firedrake mailing list
>>> firedrake at imperial.ac.uk <mailto:firedrake at imperial.ac.uk>
>>> https://mailman.ic.ac.uk/mailman/listinfo/firedrake
>>
>
>
>
>
> --
> Dr David Ham
> Departments of Mathematics and Computing
> Imperial College London
>
> http://www.imperial.ac.uk/people/david.ham
>
>
> _______________________________________________
> firedrake mailing list
> firedrake at imperial.ac.uk
> https://mailman.ic.ac.uk/mailman/listinfo/firedrake
-------------- next part --------------
HTML attachment scrubbed and removed
More information about the firedrake
mailing list