Andrew, Did you tell the triangulation that it should take periodic boundaries into account? My best bet (with the information given) would be that you didn't call parallel::distributed::Triangulation::add_periodicity() <https://www.dealii.org/current/doxygen/deal.II/classparallel_1_1distributed_1_1Triangulation.html#aa7b797070e5443a18f03a4a7f0267453> as described in https://www.dealii.org/current/doxygen/deal.II/step_45.html.
Best, Daniel Am Do., 27. Feb. 2020 um 17:45 Uhr schrieb Andrew Davis < andrew....@gmail.com>: > HI all, > > I'm trying to implement a time-dependent convection equation using DG > elements with periodic boundary conditions on and adaptive mesh. > > However, when I try to adapt the mesh using periodic using this code: > > // estimate the error in each cell > Vector<float> estimatedErrorPerCell(triangulation.n_active_cells()); > //KellyErrorEstimator<dim>::estimate(dofHandler, > QGauss<dim-1>(feSystem.degree + 1), {}, solution.block(0), > estimatedErrorPerCell); > KellyErrorEstimator<dim>::estimate(dofHandler, > QGauss<dim-1>(feSystem.degree + 1), {}, solution, estimatedErrorPerCell); > > // refine cells with the largest estimated error (80 per cent of the > error) and coarsens those cells with the smallest error (10 per cent of the > error) > GridRefinement::refine_and_coarsen_fixed_fraction(triangulation, > estimatedErrorPerCell, 0.8, 0.1); > > // to prevent the decrease in time step size limit the maximal > refinement depth of the meshlimit the maximal refinement depth of the mesh > if( triangulation.n_levels()>maxGridLevel ) { > for( auto& cell : > triangulation.active_cell_iterators_on_level(maxGridLevel) ) { > cell->clear_refine_flag(); } > } > > // store previous and current solution on the coarser mesh > std::vector<BlockVector<double> > tempSolution(2); > tempSolution[0] = oldSolution; tempSolution[1] = solution; > > // transfer the solution vectors from the old mesh to the new one > SolutionTransfer<dim, BlockVector<double> > transfer(dofHandler); > > // prepare for the refinement---we have to prepare the transfer as well > *triangulation.prepare_**coarsening_and_refinement();* // for some > reason, this doesn't tell the cells across a periodic boundary that they > need to be refined > transfer.prepare_for_coarsening_and_refinement(tempSolution); > > // do the refinement and reset the DoFs > triangulation.execute_coarsening_and_refinement(); > SetupSystem(); > > std::vector<BlockVector<double> > temp(2); > ReinitializeVector(temp[0]); ReinitializeVector(temp[1]); > transfer.interpolate(tempSolution, temp); > > oldSolution = temp[0]; solution = temp[1]; > > constraints.distribute(oldSolution); > constraints.distribute(solution); > > *I get this error:* > > -------------------------------------------------------- > An error occurred in line <992> of file > </home/add8536/Software/dealii/source/grid/tria.cc> in function > void dealii::{anonymous}::update_periodic_face_map_recursively(const > typename dealii::Triangulation<dim, spacedim>::cell_iterator&, const > typename dealii::Triangulation<dim, spacedim>::cell_iterator&, unsigned > int, unsigned int, const std::bitset<3>&, std::map<std::pair<typename > dealii::Triangulation<dim, spacedim>::cell_iterator, unsigned int>, > std::pair<std::pair<typename dealii::Triangulation<dim, > spacedim>::cell_iterator, unsigned int>, std::bitset<3> > >&) [with int dim > = 2; int spacedim = 2; typename dealii::Triangulation<dim, > spacedim>::cell_iterator = dealii::TriaIterator<dealii::CellAccessor<2, 2> > >] > The violated condition was: > dim == 1 || std::abs(cell_1->level() - cell_2->level()) < 2 > Additional information: > This exception -- which is used in many places in the library -- > usually indicates that some condition which the author of the code thought > must be satisfied at a certain point in an algorithm, is not fulfilled. An > example would be that the first part of an algorithm sorts elements of an > array in ascending order, and a second part of the algorithm later > encounters an element that is not larger than the previous one. > > There is usually not very much you can do if you encounter such an > exception since it indicates an error in deal.II, not in your own program. > Try to come up with the smallest possible program that still demonstrates > the error and contact the deal.II mailing lists with it to obtain help. > > For some reason the refinement does not seem to know that the boundaries > are periodic and that it should refine the cells on the other side of the > domain. I tried changing the refinement step so that instead of just > calling triangulation.prepare_coarsening_and_refinement();, I created my > own function that mimics the distributed triangulation (calling the > distributed version of enforce_mesh_balance_over_periodic_boundaries) > that does the following: > > while( true ) { > triangulation.prepare_coarsening_and_refinement(); > triangulation.update_periodic_face_map(); > bool changed = > enforce_mesh_balance_over_periodic_boundaries(triangulation); > if( !changed ) { break; } > } > > This fixed the error but I'm getting the strange behavior---it looks like > the cell face fluxes are no longer being computed correctly > > Does anyone have any idea about what the issue could be? I know I left out > a lot of code, but hopefully I included the important bits. > > Thanks, > Andy > ReplyForward > > <https://drive.google.com/u/0/settings/storage?hl=en&utm_medium=web&utm_source=gmail&utm_campaign=manage_storage> > <https://www.google.com/intl/en/policies/terms/> > <https://www.google.com/intl/en/policies/privacy/> > <https://www.google.com/gmail/about/policy/> > > -- > 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/fdcc9fe0-0a99-4db5-bf49-f447ed26aa4e%40googlegroups.com > <https://groups.google.com/d/msgid/dealii/fdcc9fe0-0a99-4db5-bf49-f447ed26aa4e%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/CAOYDWbJkuOy7qp2e4N03rQ1WYv%2BhsvyzqoLBqQD1F4-JbgsBgQ%40mail.gmail.com.