Aaditya,

 the documentation says
"[...] More precisely, faces with coordinates only differing in the
direction component are identified."
Hence, it is crucial to choose this correctly. For a cube with standard
coloring faces 2 and 3 share the x-coordinates so direction must be 1.

Best,
Daniel

Am Di., 8. Dez. 2020 um 20:14 Uhr schrieb Aaditya Lakshmanan <
aaditya8...@gmail.com>:

> Hi Mark,
>    Quick update. I tried running a simulation with an implementation using
> the 'Triangulation' object instead and I face the same issue as earlier. It
> returns the same error message as in the parallel case. Also, making the
> change in the *direction* argument from 0 to 1 for the constraints
> between faces with boundary indicators 2 and 3 resulted in the code running
> just fine and yielding periodic variation of the fields.
>
> Best,
> Aaditya
>
> On Tuesday, December 8, 2020 at 7:20:56 PM UTC-5 Aaditya Lakshmanan wrote:
>
>> Hi Mark,
>>     Thank you for responding. I haven't tried using a standard
>> 'Triangulation' object instead of its parallel counterpart. I will try
>> another implementation with the appropriate substitutions.
>>  Based on some trials I ran since I posted the question, my code runs
>> without issues(even on multiple processors) when I make the following
>> changes(enforcing periodicity on the triangulation and constraints) :
>>
>> GridTools::collect_periodic_faces(triangulation,
>>                                                             0,
>>                                                             1,
>>                                                             0,
>>
>>  periodicity_vector);
>> GridTools::collect_periodic_faces(triangulation,
>>                                                              2,
>>                                                              3,
>>                                                             * 1*,
>>                         // The earlier 0 changed to 1 does not yield errors
>> and results are periodic
>>
>>  periodicity_vector);
>>
>>    Although, as you mention, based on the documentation the direction
>> argument shouldn't play a role in the problem I am trying to simulate. I
>> have attached my code here. I am running my simulations with deal.II
>> version 9.2.0 compiled with Trilinos.
>>
>> Best,
>> Aaditya
>>
>> On Tuesday, December 8, 2020 at 7:01:07 PM UTC-5 mafe...@gmail.com wrote:
>>
>>> Hi Aaditya,
>>>
>>> on first look your implementation looks good to me. Does the same error
>>> occur when you are using a standard `Triangulation` object instead of a
>>> `parallel::distributed::Triangulation`?
>>>
>>> As far as I know, the direction parameter does not matter for scalar
>>> fields (see also step-45).
>>>
>>> Would you mind sharing your source code?
>>>
>>> Best,
>>> Marc
>>>
>>> On Saturday, December 5, 2020 at 10:12:00 PM UTC-7 aadit...@gmail.com
>>> wrote:
>>>
>>>> Hi,
>>>>    I am trying to simulate a reaction-diffusion system containing two
>>>> species on a square domain(structured mesh) with periodic boundary
>>>> conditions enforced on the concentration fields on opposite edges. After
>>>> testing my implementation on a single processor I obtain the following
>>>> error message :
>>>>
>>>> *---------------------------------------------------------
>>>>
>>>>                  TimerOutput objects finalize timed values printed to the
>>>>
>>>>                               screen by communicating over MPI in their
>>>> destructors.
>>>>                                               Since an exception is
>>>> currently uncaught, this
>>>>
>>>> synchronization (and subsequent output) will be skipped
>>>>
>>>>              to avoid a possible deadlock.
>>>>
>>>>
>>>>  ---------------------------------------------------------
>>>>
>>>>
>>>>
>>>>
>>>>
>>>>
>>>>  ----------------------------------------------------
>>>>
>>>>               Exception on processing:
>>>>
>>>>
>>>>
>>>>
>>>> --------------------------------------------------------
>>>>
>>>>             An error occurred in line <2107> of file
>>>> </home/aadi/dealii-9.2.0/source/grid/grid_tools_dof_handlers.cc> in
>>>> function                                                   void
>>>> dealii::GridTools::match_periodic_face_pairs(std::set<std::pair<CellIterator,
>>>> unsigned int> >&, std::set<std::pair<typename
>>>> dealii::identity<RangeType>::type, unsigned int> >&, int,
>>>> std::vector<dealii::GridTools::PeriodicFacePair<CellIterator> >&, const
>>>> dealii::Tensor<1, typename FaceIterator::AccessorType:: space_dimension>&,
>>>> const dealii::FullMatrix<double>&) [with CellIterator =
>>>> dealii::TriaIterator<dealii::CellAccessor<2, 2> >; typename
>>>> dealii::identity<RangeType>::type =
>>>> dealii::TriaIterator<dealii::CellAccessor<2, 2> >; typename
>>>> FaceIterator::AccessorType = dealii::CellAccessor<2, 2>]
>>>>                               The violated condition was:
>>>>
>>>>                                                n_matches == pairs1.size()
>>>> && pairs2.size() == 0
>>>>                                                           Additional
>>>> information:
>>>>
>>>>      Unmatched faces on periodic boundaries
>>>>
>>>>               --------------------------------------------------------
>>>>
>>>>
>>>>
>>>>                                       Aborting!  *
>>>>
>>>>
>>>> *Additional Details* : The error message is associated with the
>>>> *create_mesh()* method of the problem class whose implementation I
>>>> have included below. The part highlighted in red is the cause of the error
>>>> message :
>>>>
>>>> template <int dim>
>>>>   void Schnakenberg<dim>::create_mesh()
>>>>   {
>>>>       TimerOutput::Scope t(computing_timer, "setup");
>>>>
>>>>       GridGenerator::hyper_cube(triangulation, min_coord, max_coord,
>>>> true);
>>>>
>>>>       std::vector<GridTools::PeriodicFacePair<
>>>> typename parallel::distributed::Triangulation<dim>::cell_iterator>>
>>>> periodicity_vector;
>>>>
>>>>       GridTools::collect_periodic_faces(triangulation,
>>>> 0,
>>>> 1,
>>>> 0,
>>>> periodicity_vector);
>>>>       GridTools::collect_periodic_faces(triangulation,
>>>> 2,
>>>> 3,
>>>> 0,
>>>> periodicity_vector);
>>>>      triangulation.add_periodicity(periodicity_vector);
>>>>
>>>>       triangulation.refine_global(refine_factor);
>>>>
>>>>   }
>>>>
>>>> Also, I am curious as to what the '*direction*' argument in the
>>>> *GridTools::collect_periodic_faces* function should be for a scalar
>>>> solution field. I would appreciate some insight on this. Thank you.
>>>>
>>>> Best,
>>>> Aaditya
>>>>
>>> --
> 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 dealii+unsubscr...@googlegroups.com.
> To view this discussion on the web visit
> https://groups.google.com/d/msgid/dealii/02884162-bc73-4685-8ad8-4633b7ae2647n%40googlegroups.com
> <https://groups.google.com/d/msgid/dealii/02884162-bc73-4685-8ad8-4633b7ae2647n%40googlegroups.com?utm_medium=email&utm_source=footer>
> .
>

-- 
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 dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/CAOYDWbKZswfZ18OE-n1YGkh9rYnoJtG3oX9HKrpk4x9%3DWX-zWw%40mail.gmail.com.

Reply via email to