* For deleting all stored data in @p SolutionTransfer and reinitializing it
* use the <tt>clear()</tt> function.
*
- * Note that the @p user_pointer of some cells are used.
- * Be sure that you don't need them otherwise.
- *
* The template argument @p number denotes the data type of the vectors you want
* to transfer.
*
* DoFs of the discretisation. If we now refine the grid then the calling of
* DoFHandler::distribute_dofs() will change at least some of the
* DoF indices. Hence we need to store the DoF indices of all active cells
- * before the refinement. The @p user_pointer of each active cell
- * is used to point to the vector of these DoF indices of that cell, all other
- * @p user_pointers are cleared. All this is
- * done by prepare_for_pure_refinement().
+ * before the refinement. A pointer for each active cell
+ * is used to point to the vector of these DoF indices of that cell.
+ * This is done by prepare_for_pure_refinement().
*
* In the function <tt>refine_interpolate(in,out)</tt> and on each cell where the
- * @p user_pointer is set (i.e. the cells that were active in the original grid)
+ * pointer is set (i.e. the cells that were active in the original grid)
* we can now access the local values of the solution vector @p in
* on that cell by using the stored DoF indices. These local values are
* interpolated and set into the vector @p out that is at the end the
* the children cells are needed. Hence this interpolation
* and the storing of the interpolated values of each of the discrete functions
* that we want to interpolate needs to take place before these children cells
- * are coarsened (and deleted!!). Again the @p user_pointers of the cells are
+ * are coarsened (and deleted!!). Again a pointers for the relevant cells is
* set to point to these values (see below).
* Additionally the DoF indices of the cells
* that will not be coarsened need to be stored according to the solution
* As we need two different kinds of pointers (<tt>vector<unsigned int> *</tt> for the Dof
* indices and <tt>vector<Vector<number> > *</tt> for the interpolated DoF values)
* we use the @p Pointerstruct that includes both of these pointers and
- * the @p user_pointer of each cell points to these @p Pointerstructs.
+ * the pointer for each cell points to these @p Pointerstructs.
* On each cell only one of the two different pointers is used at one time
- * hence we could use the
- * <tt>void * user_pointer</tt> as <tt>vector<unsigned int> *</tt> at one time and as
+ * hence we could use a
+ * <tt>void * pointer</tt> as <tt>vector<unsigned int> *</tt> at one time and as
* <tt>vector<Vector<number> > *</tt> at the other but using this @p Pointerstruct
* in between makes the use of these pointers more safe and gives better
* possibility to expand their usage.
* the values of the discrete
* functions in @p all_out are set to the stored local interpolated values
* that are accessible due to the 'vector<Vector<number> > *' pointer in
- * @p Pointerstruct that is pointed to by the @p user_pointer of that cell.
+ * @p Pointerstruct that is pointed to by the pointer of that cell.
* It is clear that <tt>interpolate(all_in, all_out)</tt> only can be called with
* the <tt>vector<Vector<number> > all_in</tt> that previously was the parameter
* of the <tt>prepare_for_coarsening_and_refinement(all_in)</tt> function. Hence
* 'multiplied' by this structure.
*/
struct Pointerstruct {
+ Pointerstruct();
unsigned int memory_consumption () const;
std::vector<unsigned int> *indices_ptr;
};
/**
- * Vector of all @p Pointerstructs (cf. there).
- * It makes it
- * easier to delete all these structs
- * (without going over all <tt>cell->user_pointer</tt>)
- * after they are not used any more, and
- * collecting all these structures in a vector
- * helps avoiding fraqmentation of the memory.
+ * Map mapping from level and index of cell
+ * to the @p Pointerstructs (cf. there).
+ * This map makes it possible to keep all
+ * the information needed to transfer the
+ * solution inside this object rather than
+ * using user pointers of the Triangulation
+ * for this purpose.
*/
- std::vector<Pointerstruct> all_pointerstructs;
+ std::map<std::pair<unsigned int, unsigned int>, Pointerstruct> cell_map;
/**
* Is used for
}
+template<int dim, typename number>
+SolutionTransfer<dim, number>::Pointerstruct::Pointerstruct():
+ indices_ptr(0),
+ dof_values_ptr(0)
+{}
+
+
+
template<int dim, typename number>
void SolutionTransfer<dim, number>::clear ()
{
- if (indices_on_cell.size())
- indices_on_cell.erase(indices_on_cell.begin(), indices_on_cell.end());
- if (all_pointerstructs.size())
- all_pointerstructs.erase(all_pointerstructs.begin(), all_pointerstructs.end());
- if (dof_values_on_cell.size())
- dof_values_on_cell.erase(dof_values_on_cell.begin(), dof_values_on_cell.end());
+ indices_on_cell.clear();
+ dof_values_on_cell.clear();
+ cell_map.clear();
prepared_for=none;
}
const unsigned int dofs_per_cell = dof_handler->get_fe().dofs_per_cell;
n_dofs_old=dof_handler->n_dofs();
- indices_on_cell=std::vector<std::vector<unsigned int> > (n_active_cells,
- std::vector<unsigned int> (dofs_per_cell));
+ // efficient reallocation of indices_on_cell
+ std::vector<std::vector<unsigned int> > (n_active_cells,
+ std::vector<unsigned int> (dofs_per_cell))
+ .swap(indices_on_cell);
typename DoFHandler<dim>::cell_iterator cell = dof_handler->begin(),
endc = dof_handler->end();
// the data vectors and prolonging
// them to the children
cell->get_dof_indices(indices_on_cell[i]);
- cell->set_user_pointer(&indices_on_cell[i]);
-
+ cell_map[std::make_pair(cell->level(),cell->index())].indices_ptr=&indices_on_cell[i];
++i;
}
- else
- cell->clear_user_pointer();
}
prepared_for=pure_refinement;
}
typename DoFHandler<dim>::cell_iterator cell = dof_handler->begin(),
endc = dof_handler->end();
- std::vector<unsigned int> *indexptr;
+ typename std::map<std::pair<unsigned int, unsigned int>, Pointerstruct>::const_iterator
+ pointerstruct,
+ cell_map_end=cell_map.end();
for (; cell!=endc; ++cell)
{
- if (cell->user_pointer())
+ pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index()));
+
+ if (pointerstruct!=cell_map_end)
// this cell was refined or not
// touched at all, so we can get
// the new values by just setting
// which is both done by one
// function
{
- indexptr=static_cast<std::vector<unsigned int> *>(cell->user_pointer());
for (unsigned int i=0; i<dofs_per_cell; ++i)
- local_values(i)=in(indexptr->operator[](i));
+ local_values(i)=in((*pointerstruct->second.indices_ptr)[i]);
cell->set_dof_values_by_interpolation(local_values, out);
}
}
// allocate the needed memory. initialize
// the following arrays in an efficient
// way, without copying much
- {
- std::vector<std::vector<unsigned int> >
- tmp(n_cells_to_stay_or_refine,
- std::vector<unsigned int> (dofs_per_cell));
- indices_on_cell.swap (tmp);
- }
- {
- std::vector<std::vector<Vector<number> > >
- tmp(n_coarsen_fathers,
- std::vector<Vector<number> > (in_size,
- Vector<number> (dofs_per_cell)));
- dof_values_on_cell.swap (tmp);
- }
- {
- std::vector<Pointerstruct>
- tmp(n_cells_to_stay_or_refine+n_coarsen_fathers);
- all_pointerstructs.swap (tmp);
- }
-
-
+ std::vector<std::vector<unsigned int> >
+ (n_cells_to_stay_or_refine,
+ std::vector<unsigned int> (dofs_per_cell))
+ .swap(indices_on_cell);
+
+ std::vector<std::vector<Vector<number> > >
+ (n_coarsen_fathers,
+ std::vector<Vector<number> > (in_size,
+ Vector<number> (dofs_per_cell)))
+ .swap(dof_values_on_cell);
+
// we need counters for
- // the 'to_stay_or_refine' cells 'n_sr',
+ // the 'to_stay_or_refine' cells 'n_sr' and
// the 'coarsen_fathers' cells 'n_cf',
- // and all the cells where a
- // @p{Pointerstruct} is needed 'n'
- unsigned int n_sr=0, n_cf=0, n=0;
+ unsigned int n_sr=0, n_cf=0;
typename DoFHandler<dim>::cell_iterator cell = dof_handler->begin();
for (; cell!=endc; ++cell)
{
// dof indices and later
// interpolating to the children
cell->get_dof_indices(indices_on_cell[n_sr]);
- all_pointerstructs[n].indices_ptr=&indices_on_cell[n_sr];
- all_pointerstructs[n].dof_values_ptr=0;
- cell->set_user_pointer(&all_pointerstructs[n]);
+ cell_map[std::make_pair(cell->level(), cell->index())].indices_ptr=&indices_on_cell[n_sr];
++n_sr;
- ++n;
}
else if (cell->has_children() && cell->child(0)->coarsen_flag_set())
{
cell->get_interpolated_dof_values(all_in[j],
dof_values_on_cell[n_cf][j]);
}
- all_pointerstructs[n].indices_ptr=0;
- all_pointerstructs[n].dof_values_ptr=&dof_values_on_cell[n_cf];
- cell->set_user_pointer(&all_pointerstructs[n]);
+ cell_map[std::make_pair(cell->level(), cell->index())].dof_values_ptr=&dof_values_on_cell[n_cf];
++n_cf;
- ++n;
}
- else
- // some cell on the lower levels to
- // which nothing will happen
- cell->clear_user_pointer();
}
Assert(n_sr==n_cells_to_stay_or_refine, ExcInternalError());
Assert(n_cf==n_coarsen_fathers, ExcInternalError());
Vector<number> local_values(dofs_per_cell);
std::vector<unsigned int> dofs(dofs_per_cell);
+ typename std::map<std::pair<unsigned int, unsigned int>, Pointerstruct>::const_iterator
+ pointerstruct,
+ cell_map_end=cell_map.end();
+
typename DoFHandler<dim>::cell_iterator cell = dof_handler->begin(),
endc = dof_handler->end();
for (; cell!=endc; ++cell)
{
- if (cell->user_pointer())
+ pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index()));
+
+ if (pointerstruct!=cell_map_end)
{
- const Pointerstruct * const structptr
- =static_cast<Pointerstruct *>(cell->user_pointer());
-
const std::vector<unsigned int> * const indexptr
- =structptr->indices_ptr;
+ =pointerstruct->second.indices_ptr;
const std::vector<Vector<number> > * const valuesptr
- =structptr->dof_values_ptr;
+ =pointerstruct->second.dof_values_ptr;
// cell stayed or is
// refined
if (indexptr)
{
- Assert (structptr->dof_values_ptr == 0,
+ Assert (valuesptr == 0,
ExcInternalError());
// get the values of
else if (valuesptr)
{
Assert (!cell->has_children(), ExcInternalError());
- Assert (structptr->indices_ptr == 0,
+ Assert (indexptr == 0,
ExcInternalError());
// get the local
unsigned int
SolutionTransfer<dim, number>::memory_consumption () const
{
+ // at the moment we do not include the memory
+ // consumption of the cell_map as we have no
+ // real idea about memory consumption of a
+ // std::map
return (MemoryConsumption::memory_consumption (dof_handler) +
MemoryConsumption::memory_consumption (n_dofs_old) +
sizeof (prepared_for) +
MemoryConsumption::memory_consumption (indices_on_cell) +
- MemoryConsumption::memory_consumption (all_pointerstructs) +
MemoryConsumption::memory_consumption (dof_values_on_cell));
}