From: bangerth Date: Sun, 26 Jan 2014 19:36:10 +0000 (+0000) Subject: Disallow calling certain functions that have to do with active_fe_index on cells... X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=5ad6d0726c394a3299c8459ec38b3bcc7c56acf3;p=dealii-svn.git Disallow calling certain functions that have to do with active_fe_index on cells that are not in fact active. Fix the various places where we used to do so (and which probably produced something that made not much sense). git-svn-id: https://svn.dealii.org/trunk@32317 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 86aff38a9f..31489f7e85 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -50,6 +50,24 @@ inconvenience this causes.
    +
  1. Changed: It was possible to call DoFAccessor::set_active_fe_index() + on non-active cells. However, this made no sense: Since degrees of + freedoms only exist on active cells + for hp::DoFHandler (i.e., there is currently no implementation + of multilevel hp::DoFHandler objects), it does not make sense + to assign active FE indices to non-active cells since they + do not have finite element spaces associated with them without + having any degrees of freedom. +
    + The same of course is true for asking for the finite element active + on a non-active cell, i.e. using the functions + DoFAccessor::active_fe_index() and + DoFAccessor::get_fe(). All of these functions now produce exceptions on + non-active cells. +
    + (Wolfgang Bangerth, 2014/01/24) +
  2. +
  3. New: deal.II now links with the BOOST Iostreams library (at least if the libz and libbz2 libraries @@ -84,24 +102,18 @@ inconvenience this causes.

    Specific improvements

      -
    1. Fixed: A regression where a single whitespace accidentally added to - DEAL_II_LINKER_FLAGS internally prevented cmake-2.8.8 from configuring - sucessfully +
    2. Fixed: The SolutionTransfer class had all sorts of problems when + used with hp::DoFHandler that made its results at least questionable. + Several parts of this class have been rewritten to make the results + more predictable and, likely, more correct.
      - (Matthias Maier, Krysztof Bzowski 2014/01/26) + (Wolfgang Bangerth, 2014/01/26) -
    3. Changed: It was possible to call DoFAccessor::set_active_fe_index() - on non-active cells. However, this made no sense: Since degrees of - freedoms only exist on active cells - for hp::DoFHandler (i.e., there is currently no implementation - of multilevel hp::DoFHandler objects), it does not make sense - to assign active FE indices to non-active cells since they - do not have finite element spaces associated with them without - having any degrees of freedom. Consequently, this function will - now produce an exception when called on non-active cells. +
    4. Fixed: A regression where a single whitespace accidentally added to + DEAL_II_LINKER_FLAGS internally prevented cmake-2.8.8 from configuring + sucessfully.
      - (Wolfgang Bangerth, 2014/01/24) -
    5. + (Matthias Maier, Krysztof Bzowski, 2014/01/26)
    6. Fixed: SparsityPattern::max_entries_per_row() forgot to consider the last row of the matrix and consequently sometimes returned diff --git a/deal.II/include/deal.II/dofs/dof_accessor.h b/deal.II/include/deal.II/dofs/dof_accessor.h index d5c9cc46ee..038443436e 100644 --- a/deal.II/include/deal.II/dofs/dof_accessor.h +++ b/deal.II/include/deal.II/dofs/dof_accessor.h @@ -1660,6 +1660,14 @@ public: * this iterator. For non-hp DoF handlers, this is of course always * the same element, independent of the cell we are presently on, * but for hp DoF handlers, this may change from cell to cell. + * + * @note Since degrees of freedoms only exist on active cells + * for hp::DoFHandler (i.e., there is currently no implementation + * of multilevel hp::DoFHandler objects), it does not make sense + * to query the finite element on non-active cells since they + * do not have finite element spaces associated with them without + * having any degrees of freedom. Consequently, this function will + * produce an exception when called on non-active cells. */ const FiniteElement & get_fe () const; @@ -1669,6 +1677,14 @@ public: * FiniteElement used for this cell. This function is * only useful if the DoF handler object associated with * the current cell is an hp::DoFHandler. + * + * @note Since degrees of freedoms only exist on active cells + * for hp::DoFHandler (i.e., there is currently no implementation + * of multilevel hp::DoFHandler objects), it does not make sense + * to query active FE indices on non-active cells since they + * do not have finite element spaces associated with them without + * having any degrees of freedom. Consequently, this function will + * produce an exception when called on non-active cells. */ unsigned int active_fe_index () const; diff --git a/deal.II/include/deal.II/dofs/dof_accessor.templates.h b/deal.II/include/deal.II/dofs/dof_accessor.templates.h index d5a4ea6584..a0387765cc 100644 --- a/deal.II/include/deal.II/dofs/dof_accessor.templates.h +++ b/deal.II/include/deal.II/dofs/dof_accessor.templates.h @@ -3439,6 +3439,13 @@ inline const FiniteElement & DoFCellAccessor::get_fe () const { + Assert ((dynamic_cast*> + (this->dof_handler) != 0) + || + (this->has_children() == false), + ExcMessage ("In hp::DoFHandler objects, finite elements are only associated " + "with active cells. Consequently, you can not ask for the " + "active finite element on cells with children.")); return dealii::internal::DoFAccessor::get_fe (this->dof_handler->get_fe(), active_fe_index()); } @@ -3449,6 +3456,11 @@ inline unsigned int DoFCellAccessor::active_fe_index () const { + Assert (this->has_children() == false, + ExcMessage ("You can not ask for the active_fe_index on a cell that has " + "children because no degrees of freedom are assigned " + "to this cell and, consequently, no finite element " + "is associated with it.")); return dealii::internal::DoFCellAccessor::Implementation::active_fe_index (*this); } diff --git a/deal.II/source/dofs/dof_accessor_set.cc b/deal.II/source/dofs/dof_accessor_set.cc index 5a63df6a36..bc3a6434c7 100644 --- a/deal.II/source/dofs/dof_accessor_set.cc +++ b/deal.II/source/dofs/dof_accessor_set.cc @@ -101,7 +101,7 @@ set_dof_values_by_interpolation (const Vector &local_values, Assert (this->dof_handler != 0, typename BaseClass::ExcInvalidObject()); - Assert (&this->get_fe() != 0, + Assert (&this->get_dof_handler().get_fe() != 0, typename BaseClass::ExcInvalidObject()); Assert (local_values.size() == dofs_per_cell, typename BaseClass::ExcVectorDoesNotMatch()); diff --git a/deal.II/source/numerics/solution_transfer.cc b/deal.II/source/numerics/solution_transfer.cc index 2a6ba4b849..0f5ffc5873 100644 --- a/deal.II/source/numerics/solution_transfer.cc +++ b/deal.II/source/numerics/solution_transfer.cc @@ -135,10 +135,11 @@ 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 + const unsigned int this_fe_index = pointerstruct->second.active_fe_index; + const unsigned int dofs_per_cell=cell->get_dof_handler().get_fe()[this_fe_index].dofs_per_cell; + local_values.reinit(dofs_per_cell, true); + + // 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 @@ -150,7 +151,7 @@ SolutionTransfer::refine_interpolate(const VECTOR &in, for (unsigned int i=0; isecond.indices_ptr)[i]); cell->set_dof_values_by_interpolation(local_values, out, - cell->active_fe_index()); + this_fe_index); } } }