From: leicht Date: Thu, 3 Apr 2008 13:56:08 +0000 (+0000) Subject: SolutionTransfer now has an additional template parameter for the dof handler to... X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=7b3b6926decb30e6ec91151da2ed1b031318572a;p=dealii-svn.git SolutionTransfer now has an additional template parameter for the dof handler to be used. Instantiations are provided for DoFHandler (the default) and hp::DoFHandler. git-svn-id: https://svn.dealii.org/trunk@15980 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 74aa73e26d..cccef5b82b 100644 --- a/deal.II/deal.II/include/numerics/solution_transfer.h +++ b/deal.II/deal.II/include/numerics/solution_transfer.h @@ -201,7 +201,7 @@ DEAL_II_NAMESPACE_OPEN * @ingroup numerics * @author Ralf Hartmann, 1999 */ -template +template > class SolutionTransfer { public: @@ -210,7 +210,7 @@ class SolutionTransfer * Constructor, takes the current DoFHandler * as argument. */ - SolutionTransfer(const DoFHandler &dof); + SolutionTransfer(const DH &dof); /** * Destructor @@ -374,6 +374,11 @@ class SolutionTransfer */ DeclException0(ExcVectorsDifferFromInVectors); + /** + * Exception + */ + DeclException0(ExcNumberOfDoFsPerCellHasChanged); + /** * Exception */ @@ -388,7 +393,7 @@ class SolutionTransfer * Pointer to the degree of freedom handler * to work with. */ - SmartPointer > dof_handler; + SmartPointer dof_handler; /** * Stores the number of DoFs before the diff --git a/deal.II/deal.II/source/numerics/solution_transfer.cc b/deal.II/deal.II/source/numerics/solution_transfer.cc index 1510619d88..8ad91606f6 100644 --- a/deal.II/deal.II/source/numerics/solution_transfer.cc +++ b/deal.II/deal.II/source/numerics/solution_transfer.cc @@ -26,31 +26,31 @@ DEAL_II_NAMESPACE_OPEN -template -SolutionTransfer::SolutionTransfer(const DoFHandler &dof): +template +SolutionTransfer::SolutionTransfer(const DH &dof): dof_handler(&dof), n_dofs_old(0), prepared_for(none) {} -template -SolutionTransfer::~SolutionTransfer() +template +SolutionTransfer::~SolutionTransfer() { clear (); } -template -SolutionTransfer::Pointerstruct::Pointerstruct(): +template +SolutionTransfer::Pointerstruct::Pointerstruct(): indices_ptr(0), dof_values_ptr(0) {} -template -void SolutionTransfer::clear () +template +void SolutionTransfer::clear () { indices_on_cell.clear(); dof_values_on_cell.clear(); @@ -60,8 +60,8 @@ void SolutionTransfer::clear () } -template -void SolutionTransfer::prepare_for_pure_refinement() +template +void SolutionTransfer::prepare_for_pure_refinement() { Assert(prepared_for!=pure_refinement, ExcAlreadyPrepForRef()); Assert(prepared_for!=coarsening_and_refinement, @@ -70,40 +70,34 @@ void SolutionTransfer::prepare_for_pure_refinement() clear(); const unsigned int n_active_cells = dof_handler->get_tria().n_active_cells(); - const unsigned int dofs_per_cell = dof_handler->get_fe().dofs_per_cell; n_dofs_old=dof_handler->n_dofs(); // efficient reallocation of indices_on_cell - std::vector > (n_active_cells, - std::vector (dofs_per_cell)) + std::vector > (n_active_cells) .swap(indices_on_cell); - typename DoFHandler::cell_iterator cell = dof_handler->begin(), - endc = dof_handler->end(); + typename DH::active_cell_iterator cell = dof_handler->begin_active(), + endc = dof_handler->end(); - for (unsigned int i=0; cell!=endc; ++cell) + for (unsigned int i=0; cell!=endc; ++cell, ++i) { - if (cell->active()) - { - // on each cell store the indices - // of the dofs. after refining we - // get the values on the children - // by taking these indices, getting - // the respective values out of - // the data vectors and prolonging - // them to the children - cell->get_dof_indices(indices_on_cell[i]); - cell_map[std::make_pair(cell->level(),cell->index())].indices_ptr=&indices_on_cell[i]; - ++i; - } + indices_on_cell[i].resize(cell->get_fe().dofs_per_cell); + // on each cell store the indices of the + // dofs. after refining we get the values + // on the children by taking these + // indices, getting the respective values + // out of the data vectors and prolonging + // them to the children + cell->get_dof_indices(indices_on_cell[i]); + cell_map[std::make_pair(cell->level(),cell->index())].indices_ptr=&indices_on_cell[i]; } prepared_for=pure_refinement; } -template +template void -SolutionTransfer::refine_interpolate(const Vector &in, +SolutionTransfer::refine_interpolate(const Vector &in, Vector &out) const { Assert(prepared_for==pure_refinement, ExcNotPrepared()); @@ -114,11 +108,10 @@ SolutionTransfer::refine_interpolate(const Vector &in, ExcMessage ("Vectors cannot be used as input and output" " at the same time!")); - unsigned int dofs_per_cell=dof_handler->get_fe().dofs_per_cell; - Vector local_values(dofs_per_cell); + Vector local_values(0); - typename DoFHandler::cell_iterator cell = dof_handler->begin(), - endc = dof_handler->end(); + typename DH::cell_iterator cell = dof_handler->begin(), + endc = dof_handler->end(); typename std::map, Pointerstruct>::const_iterator pointerstruct, @@ -136,6 +129,18 @@ SolutionTransfer::refine_interpolate(const Vector &in, // which is both done by one // function { + const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell; + local_values.reinit(dofs_per_cell, true); // fast reinit, i.e. the + // entries are not set to 0. + // make sure that the size of the + // stored indices is the same as + // dofs_per_cell. this is kind of a + // test if we use the same fe in the + // hp case. to really do that test we + // would have to store the fe_index + // of all cells + Assert(dofs_per_cell==(*pointerstruct->second.indices_ptr).size(), + ExcNumberOfDoFsPerCellHasChanged()); for (unsigned int i=0; isecond.indices_ptr)[i]); cell->set_dof_values_by_interpolation(local_values, out); @@ -144,8 +149,8 @@ SolutionTransfer::refine_interpolate(const Vector &in, } -template -void SolutionTransfer::refine_interpolate (Vector &vec) const +template +void SolutionTransfer::refine_interpolate (Vector &vec) const { Assert(vec.size()==n_dofs_old, ExcWrongVectorSize(vec.size(),n_dofs_old)); @@ -156,9 +161,9 @@ void SolutionTransfer::refine_interpolate (Vector &vec) con } -template +template void -SolutionTransfer:: +SolutionTransfer:: prepare_for_coarsening_and_refinement(const std::vector > &all_in) { Assert(prepared_for!=pure_refinement, ExcAlreadyPrepForRef()); @@ -171,7 +176,6 @@ prepare_for_coarsening_and_refinement(const std::vector > &all_in clear(); const unsigned int n_active_cells = dof_handler->get_tria().n_active_cells(); - const unsigned int dofs_per_cell = dof_handler->get_fe().dofs_per_cell; n_dofs_old=dof_handler->n_dofs(); for (unsigned int i=0; i > &all_in // and that'll stay or be refined unsigned int n_cells_to_coarsen=0; unsigned int n_cells_to_stay_or_refine=0; - typename DoFHandler::active_cell_iterator + typename DH::active_cell_iterator act_cell = dof_handler->begin_active(), endc = dof_handler->end(); for (; act_cell!=endc; ++act_cell) @@ -209,25 +213,25 @@ prepare_for_coarsening_and_refinement(const std::vector > &all_in // the following arrays in an efficient // way, without copying much std::vector > - (n_cells_to_stay_or_refine, - std::vector (dofs_per_cell)) + (n_cells_to_stay_or_refine) .swap(indices_on_cell); std::vector > > (n_coarsen_fathers, - std::vector > (in_size, - Vector (dofs_per_cell))) + std::vector > (in_size)) .swap(dof_values_on_cell); // we need counters for // the 'to_stay_or_refine' cells 'n_sr' and // the 'coarsen_fathers' cells 'n_cf', unsigned int n_sr=0, n_cf=0; - typename DoFHandler::cell_iterator cell = dof_handler->begin(); + typename DH::cell_iterator cell = dof_handler->begin(); for (; cell!=endc; ++cell) { if (cell->active() && !cell->coarsen_flag_set()) { + const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell; + indices_on_cell[n_sr].resize(dofs_per_cell); // cell will not be coarsened, // so we get away by storing the // dof indices and later @@ -245,7 +249,10 @@ prepare_for_coarsening_and_refinement(const std::vector > &all_in for (unsigned int i=1; in_children(); ++i) Assert(cell->child(i)->coarsen_flag_set(), ExcTriaPrepCoarseningNotCalledBefore()); - + const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell; + std::vector >(in_size, Vector(dofs_per_cell)) + .swap(dof_values_on_cell[n_cf]); + for (unsigned int j=0; j > &all_in } -template +template void -SolutionTransfer::prepare_for_coarsening_and_refinement(const Vector &in) +SolutionTransfer::prepare_for_coarsening_and_refinement(const Vector &in) { std::vector > all_in=std::vector >(1, in); prepare_for_coarsening_and_refinement(all_in); } -template -void SolutionTransfer:: +template +void SolutionTransfer:: interpolate (const std::vector > &all_in, std::vector > &all_out) const { @@ -293,16 +300,15 @@ interpolate (const std::vector > &all_in, ExcMessage ("Vectors cannot be used as input and output" " at the same time!")); - const unsigned int dofs_per_cell=dof_handler->get_fe().dofs_per_cell; - Vector local_values(dofs_per_cell); - std::vector dofs(dofs_per_cell); + Vector local_values; + std::vector dofs; 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(); + typename DH::cell_iterator cell = dof_handler->begin(), + endc = dof_handler->end(); for (; cell!=endc; ++cell) { pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index())); @@ -315,22 +321,34 @@ interpolate (const std::vector > &all_in, const std::vector > * const valuesptr =pointerstruct->second.dof_values_ptr; + const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell; + // cell stayed or is // refined if (indexptr) { Assert (valuesptr == 0, ExcInternalError()); - + // make sure that the size of the + // stored indices is the same as + // dofs_per_cell. this is kind of + // a test if we use the same fe + // in the hp case. to really do + // that test we would have to + // store the fe_index of all + // cells + Assert(dofs_per_cell==(*indexptr).size(), + ExcNumberOfDoFsPerCellHasChanged()); // get the values of // each of the input // data vectors on this // cell and prolong it // to its children + local_values.reinit(dofs_per_cell, true); for (unsigned int j=0; joperator[](i)); + local_values(i)=all_in[j]((*indexptr)[i]); cell->set_dof_values_by_interpolation(local_values, all_out[j]); } @@ -342,7 +360,7 @@ interpolate (const std::vector > &all_in, Assert (!cell->has_children(), ExcInternalError()); Assert (indexptr == 0, ExcInternalError()); - + dofs.resize(dofs_per_cell); // get the local // indices cell->get_dof_indices(dofs); @@ -351,8 +369,23 @@ interpolate (const std::vector > &all_in, // stored data to the // new vectors for (unsigned int j=0; joperator[](j)(i); + { + // make sure that the size of + // the stored indices is the + // same as + // dofs_per_cell. this is + // kind of a test if we use + // the same fe in the hp + // case. to really do that + // test we would have to + // store the fe_index of all + // cells + Assert(dofs_per_cell==(*valuesptr)[j].size(), + ExcNumberOfDoFsPerCellHasChanged()); + + for (unsigned int i=0; i > &all_in, -template -void SolutionTransfer::interpolate(const Vector &in, +template +void SolutionTransfer::interpolate(const Vector &in, Vector &out) const { Assert (in.size()==n_dofs_old, @@ -383,9 +416,9 @@ void SolutionTransfer::interpolate(const Vector &in, -template +template unsigned int -SolutionTransfer::memory_consumption () const +SolutionTransfer::memory_consumption () const { // at the moment we do not include the memory // consumption of the cell_map as we have no @@ -400,17 +433,19 @@ SolutionTransfer::memory_consumption () const -template +template unsigned int -SolutionTransfer::Pointerstruct::memory_consumption () const +SolutionTransfer::Pointerstruct::memory_consumption () const { return sizeof(*this); } -template class SolutionTransfer; -template class SolutionTransfer; +template class SolutionTransfer >; +template class SolutionTransfer >; +template class SolutionTransfer >; +template class SolutionTransfer >; /*---------------------------- solution_transfer.cc ----------------------*/ diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 8566a72f00..ec1fe27f9a 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -290,6 +290,13 @@ constraints individually.
    +
  1. Extendend: SolutionTransfer has an additional template parameters DH + for the dof handler class to use, which defaults to DoFHandler. However, an + instantiation is also provided for hp::DoFHandler. +
    + (Tobias Leicht, 2008/04/03) +

  2. +
  3. Extended: DataOut has an extended interface now which enables a user decision, which cells shall be written as curved cells using the provided mapping: none, only cells at the boundary, or all cells. The last choice can