> 1. I use the form to update my mesh:
> |
> [...]
> timestep_number =1;
> setup_dofs ();
> [...]
> for(timestep_number=2;timestep_number<=100;++timestep_number){
> [...]
>      refine_mesh (max_grid_level,min_grid_level);
>      setup_dofs ();
> [...]
> }
> |
> 
> But starting from timestep_number=2, the setup_dof () cannot work. the error 
> is:
> |
> --------------------------------------------------------
> An error occurred in line <128> of file 
> </home/yyyyy/boost/deal.ii-8.5.1/src/source/base/subscriptor.cc> in function
>      void dealii::Subscriptor::check_no_subscribers() const
> The violated condition was:
>      counter == 0
> Additional information:
>      Object of class N6dealii15SparsityPatternE is still used by 1 other 
> objects.
> 
> (Additional information:
>    from Subscriber SparseMatrix)
> 
> See the entry in the Frequently Asked Questions of deal.II (linked to from 
> http://www.dealii.org/) for a lot more information on what this error means 
> and how to fix programs in which it happens.

Have you read up on what this error message means? In your case, a 
SparsityPattern is still used by some other object but is being destroyed. 
You'll have to find out what other object that is -- likely a SparseMatrix 
object. You'll first have to break that dependence before you destroy the 
sparsity pattern.


> 2. My refine_mesh() is:
> |
>    template <int dim>
>    void
>    StokesProblem<dim>::refine_mesh (const unsigned int max_grid_level,
>                                     const unsigned int min_grid_level)
>    {
>      Vector<float> estimated_error_per_cell (triangulation.n_active_cells());
>      FEValuesExtractors::Scalar phase(0);
>      KellyErrorEstimator<dim>::estimate (dof_handler,
>                                          QGauss<dim-1>(degree+1),
>                                          typename FunctionMap<dim>::type(),
>                                          old_solution,
>                                          estimated_error_per_cell,
>                                          fe.component_mask(phase));
>      GridRefinement::refine_and_coarsen_fixed_number (triangulation,
>                                                       
> estimated_error_per_cell,
>                                                       0.2,0.1);
>     if (triangulation.n_levels() > max_grid_level)
>        for (typename Triangulation<dim>::active_cell_iterator
>             cell = triangulation.begin_active(max_grid_level);
>             cell != triangulation.end(); ++cell)
>          cell->clear_refine_flag ();
> /*  this part for clear_coarsen_flag() cannot work */
>     if (triangulation.n_levels() < min_grid_level)
>        for (typename Triangulation<dim>::active_cell_iterator
>             cell = triangulation.begin_active(min_grid_level);
>             cell != triangulation.end(); ++cell)

I suspect you want to use
   cell != triangulation.end_active(min_grid_level)
or some such. It will not be correct otherwise.


> |
> triangulation_active.refine_global (8)
> |
> then I use
> |
> refine_mesh (9,6);
> |
> the error is:
> |
> --------------------------------------------------------
> An error occurred in line <11031> of file 
> </home/yyyyy/boost/deal.ii-8.5.1/src/source/grid/tria.cc> in function
>      dealii::Triangulation<dim, spacedim>::raw_quad_iterator 
> dealii::Triangulation<dim, spacedim>::begin_raw_quad(unsigned int) const 
> [with 
> int dim = 2; int spacedim = 2; dealii::Triangulation<dim, 
> spacedim>::raw_quad_iterator = 
> dealii::TriaRawIterator<dealii::CellAccessor<2, 
> 2> >]
> The violated condition was:
>      level<n_global_levels() || level<levels.size()
> Additional information:
>      The given level 6 is not in the valid range!

I suspect that your mesh does not at this point have 6 refinement levels. You 
will have to find out why this is so. Are you calling the function before you 
execute the refine_global(8) call?

Best
  W.

-- 
------------------------------------------------------------------------
Wolfgang Bangerth          email:                 [email protected]
                            www: http://www.math.colostate.edu/~bangerth/

-- 
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