#include <deal.II/base/config.h>
+#include <deal.II/grid/filtered_iterator.h>
+
#include <deal.II/numerics/data_out_dof_data.h>
#include <memory>
* <h3>Extensions</h3>
*
* By default, this class produces patches for all active cells. Sometimes,
- * this is not what you want, maybe because they are simply too many (and too
+ * this is not what you want, maybe because there are simply too many (and too
* small to be seen individually) or because you only want to see a certain
- * region of the domain (for example in parallel programs such as the step-18
- * example program), or for some other reason.
+ * region of the domain (for example only in the fluid part of the domain in
+ * step-46), or for some other reason.
*
* For this, internally build_patches() does not generate the sequence of
- * cells to be converted into patches itself, but relies on the two functions
- * first_cell() and next_cell(). By default, they return the first active
- * cell, and the next active cell, respectively. Since they are @p virtual
- * functions, you can write your own class derived from DataOut in which you
- * overload these two functions to select other cells for output. This may,
+ * cells to be converted into patches itself, but relies on the two function
+ * that we'll call first_cell() and next_cell(). By default, they return the
+ * first active cell, and the next active cell, respectively. But this can
+ * be changed using the set_cell_selection() function that allows you to
+ * replace this behavior. What set_cell_selection() wants to know is how
+ * you want to pick out the first cell on which output should be generated,
+ * and how given one cell on which output is generated you want to pick the
+ * next cell.
+ *
+ * This may,
* for example, include only cells that are in parts of a domain (e.g., if you
* don't care about the solution elsewhere, think for example a buffer region
* in which you attenuate outgoing waves in the Perfectly Matched Layer
* method) or if you don't want output to be generated at all levels of an
* adaptively refined mesh because this creates too much data (in this case,
- * the set of cells returned by your implementations of first_cell() and
- * next_cell() will include non-active cells, and DataOut::build_patches()
+ * the set of cells returned by your implementations of the `first_cell` and
+ * `next_cell` arguments to set_cell_selection() will include
+ * non-active cells, and DataOut::build_patches()
* will simply take interpolated values of the solution instead of the exact
- * values on these cells children for output). Once you derive your own class,
- * you would just create an object of this type instead of an object of type
- * DataOut, and everything else will remain the same.
- *
- * The two functions are not constant, so you may store information within
- * your derived class about the last accessed cell. This is useful if the
- * information of the last cell which was accessed is not sufficient to
- * determine the next one.
- *
- * There is one caveat, however: if you have cell data (in contrast to nodal,
- * or dof, data) such as error indicators, then you must make sure that
- * first_cell() and next_cell() only walk over active cells, since cell data
- * cannot be interpolated to a coarser cell. If you do have cell data and use
- * this pair of functions and they return a non-active cell, then an exception
- * will be thrown.
+ * values on these cells children for output).
*
* @ingroup output
* @author Wolfgang Bangerth, 1999
"The dimension given explicitly as a template argument to "
"this class must match the dimension of the DoFHandler "
"template argument");
+ static constexpr unsigned int spacedim = DoFHandlerType::space_dimension;
+
/**
* Typedef to the iterator type of the dof handler class under
* consideration.
curved_inner_cells
};
+ /**
+ * Constructor.
+ */
+ DataOut();
+
/**
* This is the central function of this class since it builds the list of
* patches to be written by the low-level functions of the base class. A
const unsigned int n_subdivisions = 0,
const CurvedCellRegion curved_region = curved_boundary);
+ /**
+ * A function that allows selecting for which cells output should be
+ * generated. This function takes two arguments, both `std::function`
+ * objects that can be used what the first cell on which output is
+ * generated is supposed to be, and what given one cell the next
+ * function is supposed to be. Through these function objects,
+ * it is possible to select a subset of cells on which output should
+ * be produced (e.g., only selecting those cells that belong to a
+ * part of the domain -- say, the fluid domain in a code such as step-46),
+ * or to completely change *where* output is produced (e.g., to produce
+ * output on non-active cells of a multigrid hierarchy or if the finest
+ * level of a mesh is so fine that generating graphical output would lead
+ * to an overwhelming amount of data).
+ *
+ * @param[in] first_cell A function object that takes as argument the
+ * triangulation this class works on and that should return the first cell
+ * on which output should be generated.
+ * @param[in] next_cell A function object that takes as arguments the
+ * triangulation as well as the last cell
+ * on which output was generated, and that should return the next
+ * cell on which output should be generated. If there is no next
+ * cell, i.e., if the input argument to the `next_cell` function object
+ * is the last cell on which output is to be generated, then `next_cell`
+ * must return `triangulation.end()`.
+ *
+ * These function objects are not difficult to write, but also not immediately
+ * obvious. As a consequence, there is a second variation of this function
+ * that takes a IteratorFilter argument and generates the corresponding
+ * functions itself.
+ *
+ * @note This function is also called in the constructor of this class,
+ * where the default behavior is set. By default, this class will select all
+ * @ref GlossLocallyOwnedCell "locally owned"
+ * and
+ * @ref GlossActive "active"
+ * cells for output.
+ *
+ * @note If you have cell data (in contrast to nodal, or dof, data) such as
+ * error indicators, then you must make sure that the `first_cell` and
+ * `next_cell` function objects only walk over active cells, since cell data
+ * cannot be interpolated to a coarser cell. If you do have cell data and
+ * use this pair of functions and they return a non-active cell, then an
+ * exception will be thrown.
+ */
+ void
+ set_cell_selection(
+ const std::function<cell_iterator(const Triangulation<dim, spacedim> &)>
+ & first_cell,
+ const std::function<cell_iterator(const Triangulation<dim, spacedim> &,
+ const cell_iterator &)> &next_cell);
+
+ /**
+ * A variation of the previous function that selects a subset of all
+ * cells for output based on the filter encoded in the FilteredIterator
+ * object given as argument. A typical way to generate the argument
+ * is via the make_filtered_iterator() function.
+ *
+ * @note Not all filters will result in subsets of cells for which
+ * output can actually be generated. For example, if you are working
+ * on parallel meshes where data is only available on some cells,
+ * then you better make sure that your `filtered_iterator` only
+ * loops over the
+ * @ref GlossLocallyOwnedCell "locally owned"
+ * cells; likewise, in most cases you will probably only want to work on
+ * @ref GlossActive "active"
+ * cells since this is where the solution actually lives. In particular,
+ * if you have added vectors that represent data defined on cells
+ * (instead of nodal data), then you can not generate output on non-active
+ * cells and your iterator filter should reflect this.
+ */
+ void
+ set_cell_selection(const FilteredIterator<cell_iterator> &filtered_iterator);
+
/**
* Return the first cell which we want output for. The default
* implementation returns the first active cell, but you might want to
* return other cells in a derived class.
+ *
+ * @deprecated Use the set_cell_selection() function instead.
*/
+ DEAL_II_DEPRECATED
virtual cell_iterator
first_cell();
/**
* Return the next cell after @p cell which we want output for. If there
- * are no more cells, <tt>#dofs->end()</tt> shall be returned.
+ * are no more cells, any implementation of this function should return
+ * <tt>dof_handler->end()</tt>.
*
* The default implementation returns the next active cell, but you might
* want to return other cells in a derived class. Note that the default
* guaranteed as long as first_cell() is also used from the default
* implementation. Overloading only one of the two functions might not be a
* good idea.
+ *
+ * @deprecated Use the set_cell_selection() function instead.
*/
+ DEAL_II_DEPRECATED
virtual cell_iterator
next_cell(const cell_iterator &cell);
private:
+ /**
+ * A function object that is used to select what the first cell is going to
+ * be on which to generate graphical output. See the set_cell_selection()
+ * function for more information.
+ */
+ std::function<cell_iterator(const Triangulation<dim, spacedim> &)>
+ first_cell_function;
+
+ /**
+ * A function object that is used to select what the next cell is going to
+ * be on which to generate graphical output, given a previous cell. See
+ * the set_cell_selection() function for more information.
+ */
+ std::function<cell_iterator(const Triangulation<dim, spacedim> &,
+ const cell_iterator &)>
+ next_cell_function;
+
/**
* Return the first cell produced by the first_cell()/next_cell() function
* pair that is locally owned. If this object operates on a non-distributed
* triangulation, the result equals what first_cell() returns.
+ *
+ * @deprecated Use the set_cell_selection() function instead.
*/
+ DEAL_II_DEPRECATED
virtual cell_iterator
first_locally_owned_cell();
* Return the next cell produced by the next_cell() function that is locally
* owned. If this object operates on a non-distributed triangulation, the
* result equals what first_cell() returns.
+ *
+ * @deprecated Use the set_cell_selection() function instead.
*/
+ DEAL_II_DEPRECATED
virtual cell_iterator
next_locally_owned_cell(const cell_iterator &cell);
+template <int dim, typename DoFHandlerType>
+DataOut<dim, DoFHandlerType>::DataOut()
+{
+ // For the moment, just call the existing virtual functions. This
+ // preserves backward compatibility. When these deprecated functions are
+ // removed, we'll just inline their functionality into the lambda
+ // functions created here.
+ set_cell_selection(
+ [this](const Triangulation<dim, spacedim> &) {
+ return this->first_locally_owned_cell();
+ },
+ [this](const Triangulation<dim, spacedim> &, const cell_iterator &cell) {
+ return this->next_locally_owned_cell(cell);
+ });
+}
+
+
+
template <int dim, typename DoFHandlerType>
void
DataOut<dim, DoFHandlerType>::build_one_patch(
// - active_index: index for a cell when counting from begin_active() using
// ++cell (identical to cell->active_cell_index())
// - cell_index: unique index of a cell counted using
- // next_locally_owned_cell() starting from first_locally_owned_cell()
+ // next_cell_function() starting from first_cell_function()
//
// It turns out that we create one patch for each selected cell, so
// patch_index==cell_index.
{
// max_index is the largest cell->index on level l
unsigned int max_index = 0;
- for (cell_iterator cell = first_locally_owned_cell();
+ for (cell_iterator cell = first_cell_function(*this->triangulation);
cell != this->triangulation->end();
- cell = next_locally_owned_cell(cell))
+ cell = next_cell_function(*this->triangulation, cell))
if (static_cast<unsigned int>(cell->level()) == l)
max_index =
std::max(max_index, static_cast<unsigned int>(cell->index()));
// important: we need to compute the active_index of the cell in the range
// 0..n_active_cells() because this is where we need to look up cell
// data from (cell data vectors do not have the length distance computed by
- // first_locally_owned_cell/next_locally_owned_cell because this might skip
+ // first_cell_function/next_cell_function because this might skip
// some values (FilteredIterator).
auto active_cell = this->triangulation->begin_active();
unsigned int active_index = 0;
- cell_iterator cell = first_locally_owned_cell();
+ cell_iterator cell = first_cell_function(*this->triangulation);
for (; cell != this->triangulation->end();
- cell = next_locally_owned_cell(cell))
+ cell = next_cell_function(*this->triangulation, cell))
{
// move forward until active_cell points at the cell (cell) we are
// looking at to compute the current active_index
+template <int dim, typename DoFHandlerType>
+void
+DataOut<dim, DoFHandlerType>::set_cell_selection(
+ const std::function<cell_iterator(const Triangulation<dim, spacedim> &)>
+ & first_cell,
+ const std::function<cell_iterator(const Triangulation<dim, spacedim> &,
+ const cell_iterator &)> &next_cell)
+{
+ first_cell_function = first_cell;
+ next_cell_function = next_cell;
+}
+
+
+
+template <int dim, typename DoFHandlerType>
+void
+DataOut<dim, DoFHandlerType>::set_cell_selection(
+ const FilteredIterator<cell_iterator> &filtered_iterator)
+{
+ const auto first_cell =
+ [filtered_iterator](const Triangulation<dim, spacedim> &triangulation) {
+ // Create a copy of the filtered iterator so that we can
+ // call a non-const function -- though we are really only
+ // interested in the return value of that function, not the
+ // state of the object
+ FilteredIterator<cell_iterator> x = filtered_iterator;
+ return x.set_to_next_positive(triangulation.begin());
+ };
+
+
+ const auto next_cell =
+ [filtered_iterator](const Triangulation<dim, spacedim> &,
+ const cell_iterator &cell) {
+ // Create a copy of the filtered iterator so that we can
+ // call a non-const function -- though we are really only
+ // interested in the return value of that function, not the
+ // state of the object
+ FilteredIterator<cell_iterator> x = filtered_iterator;
+
+ // Set the iterator to 'cell'. Since 'cell' must satisfy the
+ // predicate (that's how it was created), set_to_next_positive
+ // simply sets the iterator to 'cell'.
+ x.set_to_next_positive(cell);
+
+ // Advance by one:
+ ++x;
+
+ return x;
+ };
+
+ set_cell_selection(first_cell, next_cell);
+}
+
+
+
template <int dim, typename DoFHandlerType>
typename DataOut<dim, DoFHandlerType>::cell_iterator
DataOut<dim, DoFHandlerType>::first_cell()