Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Functions
deal.II and Modern C++ standards

Classes

class  IteratorRange< Iterator >
 

Functions

template<typename FunctionObjectType >
auto Threads::new_task (FunctionObjectType function_object) -> Task< decltype(function_object())>
 
template<typename BaseIterator , typename Predicate >
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators (IteratorRange< BaseIterator > i, const Predicate &p)
 
template<typename BaseIterator , typename Predicate , typename... Targs>
IteratorRange< typename internal::FilteredIteratorImplementation::NestFilteredIterators< BaseIterator, std::tuple< Predicate, Targs... > >::type > filter_iterators (IteratorRange< BaseIterator > i, const Predicate &p, const Targs... args)
 
template<typename BaseIterator , typename Predicate >
IteratorRange< FilteredIterator< BaseIterator > > operator| (IteratorRange< BaseIterator > i, const Predicate &p)
 
template<typename BaseIterator , typename Predicate >
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators (IteratorRange< BaseIterator > i, const Predicate &p)
 
template<typename BaseIterator , typename Predicate , typename... Targs>
IteratorRange< typename internal::FilteredIteratorImplementation::NestFilteredIterators< BaseIterator, std::tuple< Predicate, Targs... > >::type > filter_iterators (IteratorRange< BaseIterator > i, const Predicate &p, const Targs... args)
 
template<typename BaseIterator , typename Predicate >
IteratorRange< FilteredIterator< BaseIterator > > operator| (IteratorRange< BaseIterator > i, const Predicate &p)
 

Cell iterator functions returning ranges of iterators

IteratorRange< cell_iteratorDoFHandler< dim, spacedim >::cell_iterators () const
 
IteratorRange< active_cell_iteratorDoFHandler< dim, spacedim >::active_cell_iterators () const
 
IteratorRange< level_cell_iteratorDoFHandler< dim, spacedim >::mg_cell_iterators () const
 
IteratorRange< cell_iteratorDoFHandler< dim, spacedim >::cell_iterators_on_level (const unsigned int level) const
 
IteratorRange< active_cell_iteratorDoFHandler< dim, spacedim >::active_cell_iterators_on_level (const unsigned int level) const
 
IteratorRange< level_cell_iteratorDoFHandler< dim, spacedim >::mg_cell_iterators_on_level (const unsigned int level) const
 

Cell iterator functions returning ranges of iterators

IteratorRange< cell_iteratorTriangulation< dim, spacedim >::cell_iterators () const
 
IteratorRange< active_cell_iteratorTriangulation< dim, spacedim >::active_cell_iterators () const
 
IteratorRange< cell_iteratorTriangulation< dim, spacedim >::cell_iterators_on_level (const unsigned int level) const
 
IteratorRange< active_cell_iteratorTriangulation< dim, spacedim >::active_cell_iterators_on_level (const unsigned int level) const
 

Face iterator functions

IteratorRange< active_face_iteratorTriangulation< dim, spacedim >::active_face_iterators () const
 

Detailed Description

Since version 9.6, deal.II requires a compiler that supports at least C++17. Large parts of the library now depend on modern language constructs which are documented here.

One example is support for C++11 range-based for loops. deal.II-based codes often have many loops of the kind

...
cell = triangulation.begin_active(),
endc = triangulation.end();
for (; cell!=endc; ++cell)
cell->set_refine_flag();
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation

Using C++11's range-based for loops, you can now write this as follows:

...
for (auto &cell : triangulation.active_cell_iterators())
cell->set_refine_flag();

This works in the same way with Triangulation::active_cell_iterators() and DoFHandler::active_cell_iterators(). There are variants of these functions that provide iterator ranges for all cells (not just the active ones) and for cells on individual levels.

There are numerous other functions in the library that allow for the idiomatic use of range-based for loops. Examples are ReferenceCell::face_indices(), ReferenceCell::vertex_indices(), FEValuesBase::quadrature_point_indices(), among many others.

C++11 also introduces the concept of constexpr variables and function. The variables defined as constexpr are constant values that are computed during the compilation of the program and therefore have zero runtime cost associated with their initialization. Additionally, constexpr constants have properly defined lifetimes which prevent the so-called "static initialization order fiasco" completely. Functions can be marked as constexpr, indicating that they can produce compile-time constant return values if their input arguments are constant expressions. Additionally, classes with at least one constexpr constructor can be initialized as constexpr.

As an example, since the constructor Tensor::Tensor(const array_type &) is constexpr, we can initialize a tensor with an array during compile time as:

constexpr double[2][2] entries = {{1., 0.}, {0., 1.}};
constexpr Tensor<2, 2> A(entries);

Here, the contents of A are not stored on the stack. Rather, they are initialized during compile time and inserted into the .data portion of the executable program. The program can use these values at runtime without spending time for initialization. Initializing tensors can be simplified in one line.

constexpr Tensor<2, 2> A({{1., 0.}, {0., 1.}});

Some functions such as determinant() are specified as constexpr: these rely on the generalized constexpr support available in C++14. Some functions, such as unit_symmetric_tensor(), rely on further developments of constexpr only available in C++17 and newer. As such, this function is declared as

template <int dim, typename Number>
#define DEAL_II_CONSTEXPR
Definition config.h:236
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > unit_symmetric_tensor()

The macro DEAL_II_CONSTEXPR expands to constexpr if the compiler supports enough constexpr features (such as loops). If the compiler does not then this macro expands to nothing.

Functions declared as constexpr can be evaluated at compile time. Hence code like

constexpr double det_A = determinant(A);
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)

assuming A is declared with the constexpr specifier, will typically result in compile-time constants. This example shows the performance gains of using constexpr because here we performed an operation with \(O(\text{dim}^3)\) complexity during compile time, avoiding any runtime cost.

Function Documentation

◆ new_task()

template<typename FunctionObjectType >
auto Threads::new_task ( FunctionObjectType  function_object) -> Task<decltype(function_object())>
inline

Overload of the new_task function for objects that can be called like a function object without arguments. In particular, this function allows calling Threads::new_task() with either objects that result from using std::bind, or using lambda functions. For example, this function is called when writing code such as

task = Threads::new_task ( [] () {
do_this();
then_do_that();
return 42;
});
Task< RT > new_task(const std::function< RT()> &function)

Here, we schedule the call to the sequence of functions do_this() and then_do_that() on a separate task, by making the lambda function declared here the function to execute on the task. The lambda function then returns 42 (which is a bit pointless here, but it could of course be some computed number), and this is going to be the returned value you can later retrieve via task.return_value() once the task (i.e., the body of the lambda function) has completed.

Note
When MultithreadInfo::n_threads() returns 1, i.e., if the deal.II runtime system has been configured to only use one thread, then this function just executes the given function object immediately and stores the return value in the Task object returned by this function.
Every lambda function (or whatever else it is you pass to the new_task() function here, for example the result of a std::bind() expression) has a return type and consequently returns an object of this type. This type can be inferred using the C++11 decltype statement used in the declaration of this function, and it is then used as the template argument of the Threads::Task object returned by the current function. In the example above, because the lambda function returns 42 (which in C++ has data type int), the inferred type is int and the task object will have type Task<int>. In other words, it is not necessary to explicitly specify in user code what that return type of the lambda or std::bind expression will be, though it is possible to explicitly do so by (entirely equivalently) writing
task = Threads::new_task ( [] () -> int {
do_this();
then_do_that();
return 42;
});
In practice, the lambda functions you will pass to new_task() will of course typically be more complicated. In particular, they will likely capture variables from the surrounding context and use them within the lambda. See https://en.wikipedia.org/wiki/Anonymous_function#C.2B.2B_.28since_C.2B.2B11.29 for more on how lambda functions work.
If you pass a lambda function as an argument to the current function that captures a variable by reference, or if you use a std::bind that binds a function argument to a reference variable using std::ref() or std::cref(), then obviously you can only do this if the variables you reference or capture have a lifetime that extends at least until the time where the task finishes.
Threads::new_task() is, in essence, equivalent to calling std::async(std::launch::async, ...) in that it runs the given task in the background. See https://en.cppreference.com/w/cpp/thread/async for more information.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::invocable<FunctionObjectType>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 1177 of file thread_management.h.

◆ cell_iterators() [1/2]

template<int dim, int spacedim = dim>
IteratorRange< cell_iterator > DoFHandler< dim, spacedim >::cell_iterators ( ) const

Return an iterator range that contains all cells (active or not) that make up this DoFHandler. Such a range is useful to initialize range-based for loops as supported by C++11. See the example in the documentation of active_cell_iterators().

Returns
The half open range [this->begin(), this->end())

◆ active_cell_iterators() [1/2]

template<int dim, int spacedim = dim>
IteratorRange< active_cell_iterator > DoFHandler< dim, spacedim >::active_cell_iterators ( ) const

Return an iterator range that contains all active cells that make up this DoFHandler. Such a range is useful to initialize range-based for loops as supported by C++11, see also C++11 standard.

Range-based for loops are useful in that they require much less code than traditional loops (see here for a discussion of how they work). An example is that without range-based for loops, one often writes code such as the following:

DoFHandler<dim> dof_handler;
...
cell = dof_handler.begin_active(),
endc = dof_handler.end();
for (; cell!=endc; ++cell)
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
cell_iterator end() const
active_cell_iterator begin_active(const unsigned int level=0) const
typename ActiveSelector::active_cell_iterator active_cell_iterator

Using C++11's range-based for loops, this is now entirely equivalent to the following:

DoFHandler<dim> dof_handler;
...
for (const auto &cell : dof_handler.active_cell_iterators())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
void reinit(const Triangulation< dim, spacedim > &tria)
IteratorRange< active_cell_iterator > active_cell_iterators() const
Returns
The half open range [this->begin_active(), this->end())

◆ mg_cell_iterators()

template<int dim, int spacedim = dim>
IteratorRange< level_cell_iterator > DoFHandler< dim, spacedim >::mg_cell_iterators ( ) const

Return an iterator range that contains all cells (active or not) that make up this DoFHandler in their level-cell form. Such a range is useful to initialize range-based for loops as supported by C++11. See the example in the documentation of active_cell_iterators().

Returns
The half open range [this->begin_mg(), this->end_mg())

◆ cell_iterators_on_level() [1/2]

template<int dim, int spacedim = dim>
IteratorRange< cell_iterator > DoFHandler< dim, spacedim >::cell_iterators_on_level ( const unsigned int  level) const

Return an iterator range that contains all cells (active or not) that make up the given level of this DoFHandler. Such a range is useful to initialize range-based for loops as supported by C++11. See the example in the documentation of active_cell_iterators().

Parameters
[in]levelA given level in the refinement hierarchy of this triangulation.
Returns
The half open range [this->begin(level), this->end(level))
Precondition
level must be less than this->n_levels().

◆ active_cell_iterators_on_level() [1/2]

template<int dim, int spacedim = dim>
IteratorRange< active_cell_iterator > DoFHandler< dim, spacedim >::active_cell_iterators_on_level ( const unsigned int  level) const

Return an iterator range that contains all active cells that make up the given level of this DoFHandler. Such a range is useful to initialize range-based for loops as supported by C++11. See the example in the documentation of active_cell_iterators().

Parameters
[in]levelA given level in the refinement hierarchy of this triangulation.
Returns
The half open range [this->begin_active(level), this->end(level))
Precondition
level must be less than this->n_levels().

◆ mg_cell_iterators_on_level()

template<int dim, int spacedim = dim>
IteratorRange< level_cell_iterator > DoFHandler< dim, spacedim >::mg_cell_iterators_on_level ( const unsigned int  level) const

Return an iterator range that contains all cells (active or not) that make up the given level of this DoFHandler in their level-cell form. Such a range is useful to initialize range-based for loops as supported by C++11. See the example in the documentation of active_cell_iterators().

Parameters
[in]levelA given level in the refinement hierarchy of this triangulation.
Returns
The half open range [this->begin_mg(level), this->end_mg(level))
Precondition
level must be less than this->n_levels().

◆ filter_iterators() [1/4]

template<typename BaseIterator , typename Predicate >
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators ( IteratorRange< BaseIterator i,
const Predicate &  p 
)
inline

Filter the given range of iterators using a Predicate. This allows to replace:

DoFHandler<dim> dof_handler;
...
for (const auto &cell : dof_handler.active_cell_iterators())
{
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
}

by:

DoFHandler<dim> dof_handler;
...
const auto filtered_iterators_range =
for (const auto &cell : filtered_iterators_range)
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators(IteratorRange< BaseIterator > i, const Predicate &p)
Note
This function is obsolete and should be replaced by calling operator| instead. The latter is certainly more in the spirit of C++20 range adaptors and results in the following code instead of the one shown above:
DoFHandler<dim> dof_handler;
...
for (const auto &cell :
dof_handler.active_cell_iterators()
| IteratorFilters::LocallyOwnedCell())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::predicate<Predicate, BaseIterator>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 997 of file filtered_iterator.h.

◆ filter_iterators() [2/4]

template<typename BaseIterator , typename Predicate , typename... Targs>
IteratorRange< typename internal::FilteredIteratorImplementation::NestFilteredIterators< BaseIterator, std::tuple< Predicate, Targs... > >::type > filter_iterators ( IteratorRange< BaseIterator i,
const Predicate &  p,
const Targs...  args 
)

Filter the given range of iterators through an arbitrary number of Predicates. This allows to replace:

DoFHandler<dim> dof_handler;
...
for (const auto &cell : dof_handler.active_cell_iterators())
{
if (cell->is_locally_owned())
{
if (cell->at_boundary())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
}
}

by:

DoFHandler<dim> dof_handler;
...
const auto filtered_iterators_range =
for (const auto &cell : filter_iterators_range)
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
Note
This function is obsolete and should be replaced by calling operator| once or multiple times instead. The latter is certainly more in the spirit of C++20 range adaptors and results in the following code instead of the one shown above:
DoFHandler<dim> dof_handler;
...
for (const auto &cell :
dof_handler.active_cell_iterators()
| IteratorFilters::LocallyOwnedCell()
| IteratorFilters::AtBoundary())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::predicate<Predicate, BaseIterator>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 1071 of file filtered_iterator.h.

◆ operator|() [1/2]

template<typename BaseIterator , typename Predicate >
IteratorRange< FilteredIterator< BaseIterator > > operator| ( IteratorRange< BaseIterator i,
const Predicate &  p 
)
inline

Filter the given range of iterators using a predicate. This allows to replace:

DoFHandler<dim> dof_handler;
...
for (const auto &cell : dof_handler.active_cell_iterators())
{
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
}

by:

DoFHandler<dim> dof_handler;
...
const auto filtered_iterators_range =
for (const auto &cell :
dof_handler.active_cell_iterators() | IteratorFilters::LocallyOwnedCell())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}

Here, the operator| is to be interpreted in the same way as is done in the range adaptors feature that is part of C++20. It has the same meaning as the | symbol on the command line: It takes what is on its left as its inputs, and filters and transforms to produce some output. In the example above, it "filters" all of the active cell iterators and removes those that do not satisfy the predicate – that is, it produces a range of iterators that only contains those cells that are both active and locally owned.

Note
The | operator can be applied more than once. This results in code such as the following:
DoFHandler<dim> dof_handler;
...
for (const auto &cell :
dof_handler.active_cell_iterators()
| IteratorFilters::LocallyOwnedCell()
| IteratorFilters::AtBoundary())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
In this code, the loop executes over all active cells that are both locally owned and located at the boundary.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::predicate<Predicate, BaseIterator>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 1145 of file filtered_iterator.h.

◆ cell_iterators() [2/2]

template<int dim, int spacedim = dim>
IteratorRange< cell_iterator > Triangulation< dim, spacedim >::cell_iterators ( ) const

Return an iterator range that contains all cells (active or not) that make up this triangulation. Such a range is useful to initialize range-based for loops as supported by C++11. See the example in the documentation of active_cell_iterators().

Returns
The half open range [this->begin(), this->end())

◆ active_cell_iterators() [2/2]

template<int dim, int spacedim = dim>
IteratorRange< active_cell_iterator > Triangulation< dim, spacedim >::active_cell_iterators ( ) const

Return an iterator range that contains all active cells that make up this triangulation. Such a range is useful to initialize range-based for loops as supported by C++11, see also C++11 standard.

Range-based for loops are useful in that they require much less code than traditional loops (see here for a discussion of how they work). An example is that without range-based for loops, one often writes code such as the following (assuming for a moment that our goal is setting the user flag on every active cell):

...
cell = triangulation.begin_active(),
endc = triangulation.end();
for (; cell!=endc; ++cell)
cell->set_user_flag();

Using C++11's range-based for loops, this is now entirely equivalent to the following:

...
for (const auto &cell : triangulation.active_cell_iterators())
cell->set_user_flag();
IteratorRange< active_cell_iterator > active_cell_iterators() const
Returns
The half open range [this->begin_active(), this->end())

◆ cell_iterators_on_level() [2/2]

template<int dim, int spacedim = dim>
IteratorRange< cell_iterator > Triangulation< dim, spacedim >::cell_iterators_on_level ( const unsigned int  level) const

Return an iterator range that contains all cells (active or not) that make up the given level of this triangulation. Such a range is useful to initialize range-based for loops as supported by C++11. See the example in the documentation of active_cell_iterators().

Parameters
[in]levelA given level in the refinement hierarchy of this triangulation.
Returns
The half open range [this->begin(level), this->end(level))
Precondition
level must be less than this->n_levels().

◆ active_cell_iterators_on_level() [2/2]

template<int dim, int spacedim = dim>
IteratorRange< active_cell_iterator > Triangulation< dim, spacedim >::active_cell_iterators_on_level ( const unsigned int  level) const

Return an iterator range that contains all active cells that make up the given level of this triangulation. Such a range is useful to initialize range-based for loops as supported by C++11. See the example in the documentation of active_cell_iterators().

Parameters
[in]levelA given level in the refinement hierarchy of this triangulation.
Returns
The half open range [this->begin_active(level), this->end(level))
Precondition
level must be less than this->n_levels().

◆ active_face_iterators()

template<int dim, int spacedim = dim>
IteratorRange< active_face_iterator > Triangulation< dim, spacedim >::active_face_iterators ( ) const

Return an iterator range that contains all active faces that make up this triangulation. This function is the face version of Triangulation::active_cell_iterators(), and allows one to write code like, e.g.,

...
face->set_manifold_id(42);
IteratorRange< active_face_iterator > active_face_iterators() const
Returns
The half open range [this->begin_active_face(), this->end_face())

◆ filter_iterators() [3/4]

template<typename BaseIterator , typename Predicate >
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators ( IteratorRange< BaseIterator i,
const Predicate &  p 
)
related

Filter the given range of iterators using a Predicate. This allows to replace:

DoFHandler<dim> dof_handler;
...
for (const auto &cell : dof_handler.active_cell_iterators())
{
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
}

by:

DoFHandler<dim> dof_handler;
...
const auto filtered_iterators_range =
for (const auto &cell : filtered_iterators_range)
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
Note
This function is obsolete and should be replaced by calling operator| instead. The latter is certainly more in the spirit of C++20 range adaptors and results in the following code instead of the one shown above:
DoFHandler<dim> dof_handler;
...
for (const auto &cell :
dof_handler.active_cell_iterators()
| IteratorFilters::LocallyOwnedCell())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::predicate<Predicate, BaseIterator>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 997 of file filtered_iterator.h.

◆ filter_iterators() [4/4]

template<typename BaseIterator , typename Predicate , typename... Targs>
IteratorRange< typename internal::FilteredIteratorImplementation::NestFilteredIterators< BaseIterator, std::tuple< Predicate, Targs... > >::type > filter_iterators ( IteratorRange< BaseIterator i,
const Predicate &  p,
const Targs...  args 
)
related

Filter the given range of iterators through an arbitrary number of Predicates. This allows to replace:

DoFHandler<dim> dof_handler;
...
for (const auto &cell : dof_handler.active_cell_iterators())
{
if (cell->is_locally_owned())
{
if (cell->at_boundary())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
}
}

by:

DoFHandler<dim> dof_handler;
...
const auto filtered_iterators_range =
for (const auto &cell : filter_iterators_range)
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
Note
This function is obsolete and should be replaced by calling operator| once or multiple times instead. The latter is certainly more in the spirit of C++20 range adaptors and results in the following code instead of the one shown above:
DoFHandler<dim> dof_handler;
...
for (const auto &cell :
dof_handler.active_cell_iterators()
| IteratorFilters::LocallyOwnedCell()
| IteratorFilters::AtBoundary())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::predicate<Predicate, BaseIterator>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 1071 of file filtered_iterator.h.

◆ operator|() [2/2]

template<typename BaseIterator , typename Predicate >
IteratorRange< FilteredIterator< BaseIterator > > operator| ( IteratorRange< BaseIterator i,
const Predicate &  p 
)
related

Filter the given range of iterators using a predicate. This allows to replace:

DoFHandler<dim> dof_handler;
...
for (const auto &cell : dof_handler.active_cell_iterators())
{
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
}

by:

DoFHandler<dim> dof_handler;
...
const auto filtered_iterators_range =
for (const auto &cell :
dof_handler.active_cell_iterators() | IteratorFilters::LocallyOwnedCell())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}

Here, the operator| is to be interpreted in the same way as is done in the range adaptors feature that is part of C++20. It has the same meaning as the | symbol on the command line: It takes what is on its left as its inputs, and filters and transforms to produce some output. In the example above, it "filters" all of the active cell iterators and removes those that do not satisfy the predicate – that is, it produces a range of iterators that only contains those cells that are both active and locally owned.

Note
The | operator can be applied more than once. This results in code such as the following:
DoFHandler<dim> dof_handler;
...
for (const auto &cell :
dof_handler.active_cell_iterators()
| IteratorFilters::LocallyOwnedCell()
| IteratorFilters::AtBoundary())
{
fe_values.reinit (cell);
...do the local integration on 'cell'...;
}
In this code, the loop executes over all active cells that are both locally owned and located at the boundary.
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::predicate<Predicate, BaseIterator>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 1145 of file filtered_iterator.h.