Hi Rose,

the assert was definitely not positioned at the right place and checked on 
every face - even on non-active ones. I have moved it in the 
PR https://github.com/dealii/dealii/pull/14091. Now this assert is not 
triggered.

Generally, I am wondering if hp and PBC is working. At least, I cannot find 
a test in the hp test folder. Such a test should have triggered the assert.

Hope the patch works!

Peter

On Monday, 4 July 2022 at 18:09:41 UTC+2 jose.a...@gmail.com wrote:

> Dear dealii community,
>
> I am still trying to decipher what causes the error. It occurs in the 
> first assert of make_periodicity_constraints 
> <https://www.dealii.org/9.3.3/doxygen/deal.II/dof__tools__constraints_8cc_source.html#l02288>
> :
>
> AssertDimension(
> face_1->get_fe(face_1->nth_active_fe_index(0)).n_unique_faces(), 1);
>
> So I created vectors to store the active_fe_index and the 
> nth_active_fe_index(0), replicated the fe_index_is_active assert and the 
> get_fe() call:
>
> active_fe_index.reinit(triangulation.n_active_cells());
> nth_active_fe_index.reinit(triangulation.n_active_cells());
>
> active_fe_index = -1.0;
> nth_active_fe_index = -1.0;
>
> for (const auto &cell :
> dof_handler.active_cell_iterators())
> if (cell->is_locally_owned())
> {
> active_fe_index(cell->active_cell_index()) = cell->active_fe_index();
>
> if (cell->at_boundary())
> for (const auto &face_id : cell->face_indices())
> if (cell->face(face_id)->at_boundary())
> {
> AssertThrow(cell->face(face_id)->get_active_fe_indices().size() == 1, 
> dealii::ExcMessage("Error!"));
>
> nth_active_fe_index(cell->active_cell_index()) = cell->face(face_id)->
> nth_active_fe_index(0);
>
> AssertThrow(cell->face(face_id)->fe_index_is_active(cell->face(face_id)->
> nth_active_fe_index(0)), dealii::ExcMessage("Error!"));
>
> cell->face(face_id)->get_fe(cell->face(face_id)->nth_active_fe_index(0));
> }
> }
>
> I don't get any error from this loop and the from the visual inspection of 
> .vtu file containing both vectors one can see that they indeed coincide. 
> Maybe the problem is coming from GridTools::collect_periodic_faces as the 
> face iterators are created by it?
>
> Furthermore I also noticed that if no global refinement is performed the 
> executable runs without error in serial but still gets the error in 
> parallel.
>
> Attached is the updated MWE. A true/false argument can be passed to the 
> executable to perform one global refinement or not.
>
> Cheers,
> Jose
>
> On Wednesday, June 29, 2022 at 4:45:04 PM UTC+2 jose.a...@gmail.com wrote:
>
>> Dear dealii community,
>>
>> I am working on a gradient enhanced crystal plasticity code where the 
>> displacement field is continuous but the slip fields can be discontinuous. 
>> To this end, I am using a hp::FECollection<dim>. For example, in the 
>> case of two crystals with one slip system I have a hp::FECollection<2> 
>> with
>>
>> fe_collection[0] = [FE_Q<2> FE_Q<2> FE_Q<2>  FE_Nothing<2>]
>> fe_collection[1] = [FE_Q<2> FE_Q<2> FE_Nothing<2> FE_Q<2> ]
>>
>> where the first two FE_Q<2> correspond to the displacement field in 2D 
>> and the latter to the slip field.
>>
>> My reference problem consist of the unit square under simple shear with 
>> periodic boundary conditions on the x-direction. When I call the 
>> DoFTools::make_periodicity_constraints<dim, 
>> dim> method, I am getting the following error when I have subdomains 
>> with different values of fe_index
>>
>> An error occurred in line <1690> of file 
>> </calculate/temp/iwtm009/spack-stage-dealii-9.3.3-h7vrbckbqhxjllhyiwpmxjdp3242wjvv/spack-src/include/deal.II/dofs/dof_accessor.templates.h>
>>  
>> in function
>>     const dealii::FiniteElement<dim, spacedim>& 
>> dealii::DoFAccessor<structdim, dim, spacedim, lda>::get_fe(unsigned int) 
>> const [with int structdim = 1; int dim = 2; int spacedim = 2; bool 
>> level_dof_access = false]
>> The violated condition was: 
>>     fe_index_is_active(fe_index) == true
>> Additional information: 
>>     This function can only be called for active FE indices
>>
>> I wrote a minimum working example based on my code to showcase the error, 
>> which I have attached to this post. The executable takes a true or false 
>> argument to assign different material identifiers to different subdomains 
>> or not. If there are no subdomains and therefore, and more importantly, the 
>> size of the hp::FECollection<dim> is 1, the code runs with no problem.
>>
>> My question is if my algorithm is wrong or does the method 
>> DoFTools::make_periodicity_constraints<dim, 
>> dim> does not currently support hp::FECollection<dim> with sizes bigger 
>> than one.
>>
>> Cheers,
>> Jose
>>
>

-- 
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/fa901000-4471-4394-8acb-f759a9ab089fn%40googlegroups.com.

Reply via email to