From: leicht Date: Tue, 18 Mar 2008 10:32:49 +0000 (+0000) Subject: Modify SolutionTransfer to store data in an internal map instead of using user_pointe... X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=00931495db2495997e285fcaa9fb43740a8fd59b;p=dealii-svn.git Modify SolutionTransfer to store data in an internal map instead of using user_pointers of the underlying Triangulation. This way we can use several instances for different DoFHandlers simultaneously. git-svn-id: https://svn.dealii.org/trunk@15901 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/deal.II/include/numerics/solution_transfer.h b/deal.II/deal.II/include/numerics/solution_transfer.h index a67b648e64..74aa73e26d 100644 --- a/deal.II/deal.II/include/numerics/solution_transfer.h +++ b/deal.II/deal.II/include/numerics/solution_transfer.h @@ -126,9 +126,6 @@ DEAL_II_NAMESPACE_OPEN * For deleting all stored data in @p SolutionTransfer and reinitializing it * use the clear() 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. * @@ -142,13 +139,12 @@ DEAL_II_NAMESPACE_OPEN * 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 refine_interpolate(in,out) 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 @@ -168,7 +164,7 @@ DEAL_II_NAMESPACE_OPEN * 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 @@ -180,10 +176,10 @@ DEAL_II_NAMESPACE_OPEN * As we need two different kinds of pointers (vector * for the Dof * indices and vector > * 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 - * void * user_pointer as vector * at one time and as + * hence we could use a + * void * pointer as vector * at one time and as * vector > * 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. @@ -194,7 +190,7 @@ DEAL_II_NAMESPACE_OPEN * 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 > *' 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 interpolate(all_in, all_out) only can be called with * the vector > all_in that previously was the parameter * of the prepare_for_coarsening_and_refinement(all_in) function. Hence @@ -450,6 +446,7 @@ class SolutionTransfer * 'multiplied' by this structure. */ struct Pointerstruct { + Pointerstruct(); unsigned int memory_consumption () const; std::vector *indices_ptr; @@ -457,15 +454,15 @@ class SolutionTransfer }; /** - * Vector of all @p Pointerstructs (cf. there). - * It makes it - * easier to delete all these structs - * (without going over all cell->user_pointer) - * 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 all_pointerstructs; + std::map, Pointerstruct> cell_map; /** * Is used for diff --git a/deal.II/deal.II/source/numerics/solution_transfer.cc b/deal.II/deal.II/source/numerics/solution_transfer.cc index f11630816a..1510619d88 100644 --- a/deal.II/deal.II/source/numerics/solution_transfer.cc +++ b/deal.II/deal.II/source/numerics/solution_transfer.cc @@ -41,15 +41,20 @@ SolutionTransfer::~SolutionTransfer() } +template +SolutionTransfer::Pointerstruct::Pointerstruct(): + indices_ptr(0), + dof_values_ptr(0) +{} + + + template void SolutionTransfer::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; } @@ -68,8 +73,10 @@ void SolutionTransfer::prepare_for_pure_refinement() 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 > (n_active_cells, - std::vector (dofs_per_cell)); + // efficient reallocation of indices_on_cell + std::vector > (n_active_cells, + std::vector (dofs_per_cell)) + .swap(indices_on_cell); typename DoFHandler::cell_iterator cell = dof_handler->begin(), endc = dof_handler->end(); @@ -86,12 +93,9 @@ void SolutionTransfer::prepare_for_pure_refinement() // 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; } @@ -116,11 +120,15 @@ SolutionTransfer::refine_interpolate(const Vector &in, typename DoFHandler::cell_iterator cell = dof_handler->begin(), endc = dof_handler->end(); - std::vector *indexptr; + typename std::map, 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 @@ -128,9 +136,8 @@ SolutionTransfer::refine_interpolate(const Vector &in, // which is both done by one // function { - indexptr=static_cast *>(cell->user_pointer()); for (unsigned int i=0; ioperator[](i)); + local_values(i)=in((*pointerstruct->second.indices_ptr)[i]); cell->set_dof_values_by_interpolation(local_values, out); } } @@ -201,32 +208,21 @@ prepare_for_coarsening_and_refinement(const std::vector > &all_in // allocate the needed memory. initialize // the following arrays in an efficient // way, without copying much - { - std::vector > - tmp(n_cells_to_stay_or_refine, - std::vector (dofs_per_cell)); - indices_on_cell.swap (tmp); - } - { - std::vector > > - tmp(n_coarsen_fathers, - std::vector > (in_size, - Vector (dofs_per_cell))); - dof_values_on_cell.swap (tmp); - } - { - std::vector - tmp(n_cells_to_stay_or_refine+n_coarsen_fathers); - all_pointerstructs.swap (tmp); - } - - + std::vector > + (n_cells_to_stay_or_refine, + std::vector (dofs_per_cell)) + .swap(indices_on_cell); + + std::vector > > + (n_coarsen_fathers, + std::vector > (in_size, + Vector (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::cell_iterator cell = dof_handler->begin(); for (; cell!=endc; ++cell) { @@ -237,11 +233,8 @@ prepare_for_coarsening_and_refinement(const std::vector > &all_in // 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()) { @@ -260,16 +253,9 @@ prepare_for_coarsening_and_refinement(const std::vector > &all_in 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()); @@ -311,26 +297,29 @@ interpolate (const std::vector > &all_in, Vector local_values(dofs_per_cell); std::vector dofs(dofs_per_cell); + typename std::map, Pointerstruct>::const_iterator + pointerstruct, + cell_map_end=cell_map.end(); + typename DoFHandler::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(cell->user_pointer()); - const std::vector * const indexptr - =structptr->indices_ptr; + =pointerstruct->second.indices_ptr; const std::vector > * 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 @@ -351,7 +340,7 @@ interpolate (const std::vector > &all_in, else if (valuesptr) { Assert (!cell->has_children(), ExcInternalError()); - Assert (structptr->indices_ptr == 0, + Assert (indexptr == 0, ExcInternalError()); // get the local @@ -398,11 +387,14 @@ template unsigned int SolutionTransfer::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)); } diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index a18014e31a..4f767829a8 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -283,6 +283,15 @@ constraints individually.

deal.II

    +
  1. Changed: From now on the SolutionTransfer class does not use any user + pointers of the underlying Triangulation object. Thus several SolutionTransfer + instances with different DoFHandler objects can be used simultaneously. It is + no longer necessary to assemble all solution vectors into a new one with a + combined DoFHandler. +
    + (Tobias Leicht, 2008/03/18) +

  2. +
  3. New: There is now a namespace DoFRenumbering::boost that contains the implementation of interfaces to three reordering strategies provided by the Boost Graph Library, namely DoFRenumbering::boost::Cuthill_McKee,