Hi Michael

> Interesting. Just to clarify, since you are using Q_something every
> support point has a physical location, unlike say Raviart-Thomas
> elements where I believe they do not?

It's a mixed problem based on step-44 where:
            fe (FE_Q<dim>(parameters.poly_degree), dim,    // displacement
                FE_DGPMonomial<dim>(parameters.poly_degree-1), 1,  // pressure
                FE_DGPMonomial<dim>(parameters.poly_degree-1), 1), // dilatation
So only the displacement field has support points that coincide with physical 
points

All I'm actually trying to do is add periodic boundary conditions to code based 
on step-44.


> When you say the support points
> on the face are but not in the domain do you mean that you need
> information from support points in the cell that are not on the face
> of interest so they have a physical location that is not on the face?

I need information about the physical location of the displacement support 
points on faces of the domain. The pressure and dilitation are internal and 
thus do not appear as dof on the face.

> 
> I believe you can query a face's orientation from a CellAccessor
> (face_orientation(const unsigned int face) const) so maybe if not
> knowing that is the snag you can check the orientation and then have a
> switch statement to handle the various conditions. This would mean you
> can tie up the relevant dofs based on their support point number.

I already use the normal to the surface to sort faces. As it is in 3D and I 
need to impose periodic boundary conditions on the displacement for all 
opposite pairs of faces I already use this. I have numerous asserts that check 
this sort of logic as well.


> 
> In case it is relevant, you can get an FEValues object to provide
> support point locations rather than quadrature point locations. This
> is better than relying on methods that access vertices because it
> extends to higher shape functions. You use a call like:
> Quadrature<dim-1> dummy_quadrature
> ((fe.base_element(0)).get_unit_face_support_points());
> FEFaceValues<dim>   fe_face_values (fe.base_element (0),
> dummy_quadrature, update_quadrature_points);

thanks. this is of use. 
> 
> Not really answering your question but just pointing out some methods
> that seem like they could help. Hopefully someone who has tried what
> you are interested in doing will jump in here!

many thanks

Andrew
> 
> Cheers,
> Michael
> 
> On Tue, Mar 15, 2011 at 11:58 PM, Andrew McBride
> <[email protected]> wrote:
>> Thanks for the quick response
>> 
>> The problem is that identifying a dof with it's local and global dof index 
>> is not enough. You need to be able to have some way to compare support 
>> points on the face of elements on opposite sides of the domain. The 
>> comparison needs to be done on position. You don't know the orientation of 
>> the elements. (am I missing something here?)
>> 
>> This is fine if all support points are associated with physical points. In 
>> my case all support points on the face are; but not in the domain.
>> 
>> The hack would be to limit myself to trilinear interpolations for the 
>> displacement field and then use the link between vertices and dof but there 
>> must be a more elegant way to do this.
>> 
>> Cheers
>> Andrew
>> 
>> 
>> 
>> 
>> On 15.03.2011, at 17:42, Michael Rapson wrote:
>> 
>>> Hi there Andrew,
>>> 
>>> I think that fe.component_to_system_index could give you what you
>>> want. You would use it in a loop over the desired faces. It takes two
>>> unsigned ints if I remember right: one is the component you want and
>>> the other which support point for that component. I forget how it
>>> handles the situation where one component has fewer support points,
>>> like pressure often does. It gives you the local_index that can be
>>> used in the standard local_dof_indices (local_index) to get the global
>>> dof index (filled using cell->get_dof_indices (local_dof_indices). The
>>> method also has an inverse (system_to_component_index) that takes the
>>> local_dof_index (on the cell) and gives you the two variables I
>>> mentioned above. You can then test to see if it is the component you
>>> are interested in but this is probably more expensive. In either case
>>> you need to use it in a loop over the desired faces.
>>> 
>>> Cheers,
>>> Michael
>>> 
>>> On Tue, Mar 15, 2011 at 4:22 PM, Andrew McBride
>>> <[email protected]> wrote:
>>>> Hi all
>>>> I'm trying to implement periodic boundary conditions on pairs of opposite
>>>> faces of a cube in 3D. The problem is mixed in that each element has
>>>> vectorial (continuous) displacements,  and (discontinuous) pressures and
>>>> dilatation as unknowns. I wish to impose periodicity between the
>>>> displacement degrees of freedom only.
>>>> I've adapted a 2d periodic boundary conditions algorithm. I've sorted the
>>>> faces of the element on the opposite sides based on the position of the
>>>> centroid of the elements face. The problem is now to match degrees of
>>>> freedom on opposite faces. I can determine if they are associated with the
>>>> displacements and to which local support point they correspond but this 
>>>> does
>>>> not give me a geometric measure to order the dof.
>>>> I can't restrict the problem to only Q_1 elements for the displacement but
>>>> it will always be Q_something.
>>>> Is there any way to use a function
>>>> like DoFTools::map_dofs_to_support_points  but to restrict its operation to
>>>> only certain components? Or, if I know a dof index (local or global) can I
>>>> determine the position of the support point? In other words, is there a
>>>> function like DoFTools::map_dofs_to_support_point that only works with one
>>>> dof at a time.
>>>> I would appreciate any suggestions as I'm sure I'm not the first to do this
>>>> Thanks as ever
>>>> andrew
>>>> _______________________________________________
>>>> dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii
>>>> 
>>>> 
>> 
>> 
> _______________________________________________
> dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii

_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii

Reply via email to