![]() |
Reference documentation for deal.II version Git c9feb145b3 2021-03-01 21:22:52 +0100
|
Classes | |
class | IteratorRange< Iterator > |
Functions | |
template<typename FunctionObjectType > | |
auto | Threads::new_thread (FunctionObjectType function_object) -> Thread< decltype(function_object())> |
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 > > | 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) |
Cell iterator functions returning ranges of iterators | |
IteratorRange< cell_iterator > | Triangulation< dim, spacedim >::cell_iterators () const |
IteratorRange< active_cell_iterator > | Triangulation< dim, spacedim >::active_cell_iterators () const |
IteratorRange< cell_iterator > | Triangulation< dim, spacedim >::cell_iterators_on_level (const unsigned int level) const |
IteratorRange< active_cell_iterator > | Triangulation< dim, spacedim >::active_cell_iterators_on_level (const unsigned int level) const |
Face iterator functions | |
IteratorRange< active_face_iterator > | Triangulation< dim, spacedim >::active_face_iterators () const |
Since version 9.0, deal.II requires a compiler that supports at least C++11. As part of this, many places in the internal implementation of deal.II are now using features that were only introduced in C++11. That said, deal.II also has functions and classes that make using it with C++11 features easier.
One example is support for C++11 range-based for loops. deal.II-based codes often have many loops of the kind
Using C++11's range-based for loops, you can now write this as follows:
This relies on functions such as Triangulation::active_cell_iterators(), and equivalents in the DoF handler classes, DoFHandler::active_cell_iterators(), hp::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 GeometryInfo::face_indices(), GeometryInfo::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:
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.
Some functions such as determinant() are specified as constexpr
but they require a compiler with C++14 capability. As such, this function is internally declared as:
The macro DEAL_II_CONSTEXPR simplifies to constexpr
if a C++14-capable compiler is available. Otherwise, for old compilers, it ignores DEAL_II_CONSTEXPR altogether. Therefore, with newer compilers, the user can write
assuming A
is declared with the constexpr
specifier. 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.
|
inline |
Overload of the new_thread() function for objects that can be called like a function object without arguments. In particular, this function allows calling Threads::new_thread() with either objects that result from using std::bind, or using lambda functions. For example, this function is called when writing code such as
Here, we run the sequence of functions do_this()
and then_do_that()
on a separate thread, by making the lambda function declared here the function to execute on the thread. 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 thread.return_value()
once the thread (i.e., the body of the lambda function) has completed.
decltype
statement used in the declaration of this function, and it is then used as the template argument of the Threads::Thread 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 Definition at line 822 of file thread_management.h.
|
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
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.
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 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.) The only difference is if you configured deal.II with DEAL_II_WITH_THREADS=OFF
, then the operation described by the arguments of this function are executed immediately and the returned value is placed in the Task object returned here. This is useful for cases where one wants to run a program in a way where deal.II does not internally create parallel tasks, for example because one is already using one MPI process per core in a parallel computation. Definition at line 1414 of file thread_management.h.
IteratorRange< typename DoFHandler< dim, spacedim >::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().
Definition at line 2186 of file dof_handler.cc.
IteratorRange< typename DoFHandler< dim, spacedim >::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:
Using C++11's range-based for loops, this is now entirely equivalent to the following:
[this->begin_active(), this->end())
Definition at line 2196 of file dof_handler.cc.
IteratorRange< typename DoFHandler< dim, spacedim >::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().
[this->begin_mg(), this->end_mg())
Definition at line 2207 of file dof_handler.cc.
IteratorRange< typename DoFHandler< dim, spacedim >::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().
[in] | level | A given level in the refinement hierarchy of this triangulation. |
[this->begin(level), this->end(level))
Definition at line 2217 of file dof_handler.cc.
IteratorRange< typename DoFHandler< dim, spacedim >::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().
[in] | level | A given level in the refinement hierarchy of this triangulation. |
[this->begin_active(level), this->end(level))
Definition at line 2228 of file dof_handler.cc.
IteratorRange< typename DoFHandler< dim, spacedim >::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().
[in] | level | A given level in the refinement hierarchy of this triangulation. |
[this->begin_mg(level), this->end_mg(level))
Definition at line 2240 of file dof_handler.cc.
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators | ( | IteratorRange< BaseIterator > | i, |
const Predicate & | p | ||
) |
Filter the given range of iterators using a Predicate. This allows to replace:
by:
Definition at line 860 of file filtered_iterator.h.
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:
by:
Definition at line 910 of file filtered_iterator.h.
IteratorRange< typename Triangulation< dim, spacedim >::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().
IteratorRange< typename Triangulation< dim, spacedim >::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):
Using C++11's range-based for loops, this is now entirely equivalent to the following:
[this->begin_active(), this->end())
IteratorRange< typename Triangulation< dim, spacedim >::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().
[in] | level | A given level in the refinement hierarchy of this triangulation. |
[this->begin(level), this->end(level))
IteratorRange< typename Triangulation< dim, spacedim >::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().
[in] | level | A given level in the refinement hierarchy of this triangulation. |
[this->begin_active(level), this->end(level))
IteratorRange< typename Triangulation< dim, spacedim >::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.,
[this->begin_active_face(), this->end_face())
|
related |
Filter the given range of iterators using a Predicate. This allows to replace:
by:
Definition at line 860 of file filtered_iterator.h.
|
related |
Filter the given range of iterators through an arbitrary number of Predicates. This allows to replace:
by:
Definition at line 910 of file filtered_iterator.h.
types::global_dof_index DoFHandler< dim, spacedim >::n_dofs | ( | ) | const |
Return the global number of degrees of freedom. If the current object handles all degrees of freedom itself (even if you may intend to solve your linear system in parallel, such as in step-17 or step-18), then this number equals the number of locally owned degrees of freedom since this object doesn't know anything about what you want to do with it and believes that it owns every degree of freedom it knows about.
On the other hand, if this object operates on a parallel::distributed::Triangulation object, then this function returns the global number of degrees of freedom, accumulated over all processors.
In either case, included in the returned number are those DoFs which are constrained by hanging nodes, see Constraints on degrees of freedom.
Mathematically speaking, the number returned by this function equals the dimension of the finite element space (without taking into account constraints) that corresponds to (i) the mesh on which it is defined, and (ii) the finite element that is used by the current object. It also, of course, equals the number of shape functions that span this space.
types::global_dof_index DoFHandler< dim, spacedim >::n_dofs | ( | const unsigned int | level | ) | const |
The (global) number of multilevel degrees of freedom on a given level.
If no level degrees of freedom have been assigned to this level, returns numbers::invalid_dof_index. Else returns the number of degrees of freedom on this level.
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs | ( | ) | const |
Return the number of locally owned degrees of freedom located on the boundary.
Definition at line 2255 of file dof_handler.cc.
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs | ( | const std::map< types::boundary_id, const Function< spacedim, number > *> & | boundary_ids | ) | const |
Return the number of locally owned degrees of freedom located on those parts of the boundary which have a boundary indicator listed in the given set. The reason that a map
rather than a set
is used is the same as described in the documentation of that variant of DoFTools::make_boundary_sparsity_pattern() that takes a map.
There is, however, another overload of this function that takes a set
argument (see below).
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs | ( | const std::set< types::boundary_id > & | boundary_ids | ) | const |
Return the number of degrees of freedom located on those parts of the boundary which have a boundary indicator listed in the given set. The
Definition at line 2304 of file dof_handler.cc.
const BlockInfo& DoFHandler< dim, spacedim >::block_info | ( | ) | const |
Access to an object informing of the block structure of the dof handler.
If an FESystem is used in distribute_dofs(), degrees of freedom naturally split into several blocks. For each base element as many blocks appear as its multiplicity.
At the end of distribute_dofs(), the number of degrees of freedom in each block is counted, and stored in a BlockInfo object, which can be accessed here. If you have previously called distribute_mg_dofs(), the same is done on each level of the multigrid hierarchy. Additionally, the block structure on each cell can be generated in this object by calling initialize_local_block_info().
types::global_dof_index DoFHandler< dim, spacedim >::n_locally_owned_dofs | ( | ) | const |
Return the number of degrees of freedom that belong to this process.
If this is a sequential DoFHandler, then the result equals that produced by n_dofs(). (Here, "sequential" means that either the whole program does not use MPI, or that it uses MPI but only uses a single MPI process, or that there are multiple MPI processes but the Triangulation on which this DoFHandler builds works only on one MPI process.) On the other hand, if we are operating on a parallel::distributed::Triangulation or parallel::shared::Triangulation, then it includes only the degrees of freedom that the current processor owns. Note that in this case this does not include all degrees of freedom that have been distributed on the current processor's image of the mesh: in particular, some of the degrees of freedom on the interface between the cells owned by this processor and cells owned by other processors may be theirs, and degrees of freedom on ghost cells are also not necessarily included.
const IndexSet& DoFHandler< dim, spacedim >::locally_owned_dofs | ( | ) | const |
Return an IndexSet describing the set of locally owned DoFs as a subset of 0..n_dofs(). The number of elements of this set equals n_locally_owned_dofs().
const IndexSet& DoFHandler< dim, spacedim >::locally_owned_mg_dofs | ( | const unsigned int | level | ) | const |
Return an IndexSet describing the set of locally owned DoFs used for the given multigrid level as a subset of 0..n_dofs(level).
const std::vector<IndexSet>& DoFHandler< dim, spacedim >::locally_owned_dofs_per_processor | ( | ) | const |
Return a vector that stores the locally owned DoFs of each processor.
Utilities::MPI::all_gather(comm, locally_owned_dofs())
upon the first invocation, including global communication. Use Utilities::MPI::all_gather(comm, dof_handler.locally_owned_dofs())
instead if using up to a few thousands of MPI ranks or some variant involving local communication with more processors. const std::vector<types::global_dof_index>& DoFHandler< dim, spacedim >::n_locally_owned_dofs_per_processor | ( | ) | const |
Return a vector that stores the number of degrees of freedom each processor that participates in this triangulation owns locally. The sum of all these numbers equals the number of degrees of freedom that exist globally, i.e. what n_dofs() returns.
Utilities::MPI::all_gather(comm, n_locally_owned_dofs()
upon the first invocation, including global communication. Use Utilities::MPI::all_gather(comm, dof_handler.n_locally_owned_dofs()
instead if using up to a few thousands of MPI ranks or some variant involving local communication with more processors. const std::vector<IndexSet>& DoFHandler< dim, spacedim >::locally_owned_mg_dofs_per_processor | ( | const unsigned int | level | ) | const |
Return a vector that stores the locally owned DoFs of each processor on the given level level
.
Utilities::MPI::all_gather(comm, locally_owned_dofs_mg())
upon the first invocation, including global communication. Use Utilities::MPI::all_gather(comm, dof_handler.locally_owned_dofs_mg())
instead if using up to a few thousands of MPI ranks or some variant involving local communication with more processors. const FiniteElement<dim, spacedim>& DoFHandler< dim, spacedim >::get_fe | ( | const unsigned int | index = 0 | ) | const |
Return a constant reference to the selected finite element object. Since there is only one FiniteElement index
must be equal to zero which is also the default value.
const hp::FECollection<dim, spacedim>& DoFHandler< dim, spacedim >::get_fe_collection | ( | ) | const |
Return a constant reference to the set of finite element objects that are used by this DoFHandler
. Since this object only contains one FiniteElement, only this one object is returned wrapped in a hp::FECollection.
const Triangulation<dim, spacedim>& DoFHandler< dim, spacedim >::get_triangulation | ( | ) | const |
Return a constant reference to the triangulation underlying this object.
MPI_Comm DoFHandler< dim, spacedim >::get_communicator | ( | ) | const |
Return MPI communicator used by the underlying triangulation.
void DoFHandler< dim, spacedim >::prepare_for_serialization_of_active_fe_indices | ( | ) |
Whenever serialization with a parallel::distributed::Triangulation as the underlying triangulation is considered, we also need to consider storing the active FE indices on all active cells as well.
This function registers that these indices are to be stored whenever the parallel::distributed::Triangulation::save() function is called on the underlying triangulation.
Definition at line 3327 of file dof_handler.cc.
void DoFHandler< dim, spacedim >::deserialize_active_fe_indices | ( | ) |
Whenever serialization with a parallel::distributed::Triangulation as the underlying triangulation is considered, we also need to consider storing the active FE indices on all active cells as well.
This function deserializes and distributes the previously stored active FE indices on all active cells.
Definition at line 3385 of file dof_handler.cc.
|
virtual |
Determine an estimate for the memory consumption (in bytes) of this object.
This function is made virtual, since a dof handler object might be accessed through a pointers to this base class, although the actual object might be a derived class.
Definition at line 2356 of file dof_handler.cc.
void DoFHandler< dim, spacedim >::save | ( | Archive & | ar, |
const unsigned int | version | ||
) | const |
Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.
void DoFHandler< dim, spacedim >::load | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.
void DoFHandler< dim, spacedim >::serialize | ( | Archive & | archive, |
const unsigned int | version | ||
) |
Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.
|
static |
Exception
|
static |
Exception
|
static |
Exception
|
static |
Exception
|
static |
Exception
|
static |
Exception used when a certain feature doesn't make sense when DoFHandler does not hp-capabilities.
|
static |
Exception used when a certain feature is not implemented when the DoFHandler has hp-capabilities.