Chucui,
> And I think what you means is the same as it above beforehand, but maybe not
> now. The dsp (BlockDynamicSparsityPattern) will be destroyed and cannot exist
> in all lifetime of BlockSparseMatrix, so we define : BlockSparsityPattern
> sparsity_pattern at the top, (this is what step-22 does) and let it copy
> dsp....
> But now, we know that BlockSparsityPattern sparsity_pattern also will be
> destroyed and cannot exist in all lifetime of BlockSparseMatrix. So I don't
> know how to correct it.
> And I cannot implement what you said "You'll first have to break that
> dependence before you destroy the sparsity pattern".
What I'm trying to say is this: Let's say you have code like this:
{
SparseMatrix m;
SparsityPattern sp;
...fill sp, for example by copying from a DynamicSparsityPattern...
m.reinit (sp); // m now points to sp
...some more code...
}
At the '}' brace, the compiler destroys objects in the reverse order as they
were created. So it first destroys 'sp'. But 'm' still points to it and so you
get an error: since 'sp' is still needed by the variable 'm', you can't
destroy 'sp' yet. If you want to destroy 'sp', you first have to tell 'm' to
not use 'sp' any more, for example by using code such as this:
{
SparseMatrix m;
SparsityPattern sp;
...fill sp, for example by copying from a DynamicSparsityPattern...
m.reinit (sp); // m now points to sp
...some more code...
m.reinit(); // ***
}
In the marked line, you tell 'm' to not use 'sp' any more. So at the '}'
brace, 'sp' can be destroyed now because no other variable is using it any
more and you will not get the error any more. Next 'm' is destroyed, but that
is no longer of concern to us.
> 2. As your suggestions, I change my code :
> |
> 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_active.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_active,
>
> estimated_error_per_cell,
> 0.2,0.1);
> std::cout << "triangulation.n_levels" << triangulation_active.n_levels() <<
> std::endl;
> if (triangulation_active.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 ();
> std::cout << "clear-refine" << triangulation_active.n_levels() << std::endl;
> if (triangulation_active.n_levels() < min_grid_level)
> for (typename Triangulation<dim>::active_cell_iterator
> cell = triangulation.end_active(min_grid_level);
> cell != triangulation.end(); ++cell)
> cell->clear_coarsen_flag ();
> std::cout << "clear-coarsen" << triangulation_active.n_levels() << std::endl;
> triangulation_active.execute_coarsening_and_refinement ();
> std::cout << "c-and-r" << triangulation_active.n_levels() << std::endl;
> }
> |
>
> Maybe I need to chage all"triangulation" in my originial
> StokesProblem<dim>::refine_mesh (const unsigned int max_grid_level, const
> unsigned int min_grid_level) into "triangulation_active", this function can
> run(but i don't know if it is right), and for implement periodic boundary
> condition, I use the "triangulation" and "triangulation_active" :
I don't know what the difference between triangulation and
triangulation_active is. It's really hard to debug a code you don't have or
know :-) But I will point out that something like this here:
if (triangulation_active.n_levels() > max_grid_level)
for (typename Triangulation<dim>::active_cell_iterator
cell = triangulation.begin_active(max_grid_level);
looks suspicious: In the 'if' condition, you are referencing
triangulation_active, but then when you are asking for a concrete cell, you
are referencing triangulation.begin_active(max_grid_level). This doesn't look
right.
> I don't know what I should use in whole code is "triangulation" or
> "triangulation_active" ?
A question like this suggests to me that you're not quite clear about what
these two objects are and where they should be used. That sounds like you
don't quite know what your code does or *should* do. Do you have colleagues
with whom you could talk through the design of your code by (i) explaining to
them what you want to do, and (ii) after that showing them what you have and
asking them for feedback?
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.