Hello,

I came across this old thread, and successfully used the provided code in 
my own code (https://github.com/l-korous/mhdeal) - thank you Andreas - and 
I want to share how I used it:

- here are the (slightly changed) files provided by Andreas:
    https://github.com/l-korous/mhdeal/blob/master/DealiiExtensions.h
https://github.com/l-korous/mhdeal/blob/master/DealiiExtensions.cpp
* I changed the custom class PeriodicBoundary's functionality into simple 
arrays of boundary_ids and directions - 
https://github.com/l-korous/mhdeal/blob/master/DealiiExtensions.h#L54

- here is the usage of collect_periodic_faces & add_periodicity: 
https://github.com/l-korous/mhdeal/blob/master/parameters.cpp#L65
- and here is the usage of the custom functions from Andreas: 
https://github.com/l-korous/mhdeal/blob/master/problem.cpp#L50

Many thanks



On Monday, February 1, 2016 at 3:36:20 PM UTC+1, Nectarios wrote:
>
> Dear Andreas,
>
> Thank you very much for your help. I will look at it and let you know how 
> things go.
>
> Best Regards,
>
> Nectarios 
>
> On Mon, Feb 1, 2016 at 11:08 AM, Andreas Krämer <[email protected] 
> <javascript:>> wrote:
>
>> Hi Nectarios,
>>
>> the  way I did it works as follows (starting with an unrefined grid):
>>
>>    - Collect periodic faces to create a halo layer of ghost cells across 
>>    the periodic boundary:    
>>    // collect periodic faces and extend halo layer of ghost cells across 
>>    periodic boundaries
>>        // take care: collect_periodic_faces must not be called for a 
>>    refined grid
>>        if (m_triangulation->n_levels() > 1) {
>>            throw PeriodicBoundaryNotPossible(
>>                    "The periodic boundaries have to be set before the 
>>    refinement.");
>>        }
>>        std::vector<
>>                dealii::GridTools::PeriodicFacePair<
>>                        typename Mesh<dim>::cell_iterator> > matched_pairs
>>    ;
>>        dealii::GridTools::collect_periodic_faces(*m_triangulation,
>>                m_boundaryIndicator1, m_boundaryIndicator2, m_direction,
>>                matched_pairs);
>>        // add halo layer
>>        m_triangulation->add_periodicity(matched_pairs);
>>    
>>    
>>    
>>    - Refine grid (if you want to) and distribute DoFs
>>    - Create a cell map that identifies faces at the two sides of a 
>>    periodic boundary
>>    PeriodicCellMap<dim> m_cells;    
>>    DealIIExtensions::make_periodicity_map_dg(doFHandler, 
>>    m_boundaryIndicator1, m_boundaryIndicator2, m_direction, m_cells);
>>    
>>    - Make the sparsity pattern, incorporating the PeriodicCellMap (see 
>>    function make_sparser_flux_sparsity_pattern to get the idea). In general 
>>    you could just call deal.II's make_flux_sparsity_pattern and iterate over 
>>    the PeriodicCellMap to add the missing entries.
>>    - Assemble the system, incorporating the PeriodicCellMap to identify 
>>    neighbors across periodic faces. I did something like
>>    // loop over all cells
>>    ...
>>    // loop over all faces
>>        for (size_t j = 0; j < dealii::GeometryInfo<dim>::faces_per_cell; 
>>    j++) {
>>            //Faces at boundary
>>            if (cell->face(j)->at_boundary()) {
>>                size_t boundaryIndicator = cell->face(j)->boundary_id();
>>                if (m_boundaries->isPeriodic(boundaryIndicator)) {
>>                    // Apply periodic boundaries
>>                    const boost::shared_ptr<PeriodicBoundary<dim> >& 
>>    periodicBoundary =
>>                            m_boundaries->getPeriodicBoundary(
>>    boundaryIndicator);
>>                    assert(periodicBoundary->isFaceInBoundary(cell, j));
>>                    typename dealii::DoFHandler<dim>::cell_iterator 
>>    neighborCell;
>>                    size_t opposite_face =
>>                            periodicBoundary->
>>    getOppositeCellAtPeriodicBoundary(
>>                                    cell, neighborCell); // this function 
>>    extracts the neighbor cell from the PeriodicCellMap
>>                    // now assemble And Distribute Internal Face with the 
>>    periodic neighbor cell
>>    
>>    
>>    The names of the functions should give you an idea what actually 
>>    happens here.
>>    
>>
>> PeriodicCellMap, make_periodicity_map_dg and make_sparser_flux_sparsity 
>> pattern are defined in the attached file to give you a rough idea (this 
>> version works better than the one posted before).
>>
>> The code is part of a bigger framework, so you will have to change it 
>> quite a bit before using it. However, the most important function, 
>> make_periodicity_map_dg, should also work for your code.
>>
>>
>> If you have further questions, I'm happy to help.
>>
>> Best,
>>
>> Andreas
>>
>>
>>>
>>>> -- 
>> The deal.II project is located at http://www.dealii.org/
>> For mailing list/forum options, see 
>> https://groups.google.com/d/forum/dealii?hl=en
>> --- 
>> You received this message because you are subscribed to a topic in the 
>> Google Groups "deal.II User Group" group.
>> To unsubscribe from this topic, visit 
>> https://groups.google.com/d/topic/dealii/zdiXh-NRcYk/unsubscribe.
>> To unsubscribe from this group and all its topics, send an email to 
>> [email protected] <javascript:>.
>> For more options, visit https://groups.google.com/d/optout.
>>
>
>

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to [email protected].
For more options, visit https://groups.google.com/d/optout.

Reply via email to