Reference documentation for deal.II version 9.1.1
|
Namespaces | |
MeshWorker | |
MeshWorker::Assembler | |
Functions | |
template<class INFOBOX , class DOFINFO , int dim, int spacedim, class ITERATOR > | |
void | MeshWorker::cell_action (ITERATOR cell, DoFInfoBox< dim, DOFINFO > &dof_info, INFOBOX &info, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, const LoopControl &loop_control) |
template<int dim, int spacedim, class DOFINFO , class INFOBOX , class ASSEMBLER , class ITERATOR > | |
void | MeshWorker::loop (ITERATOR begin, typename identity< ITERATOR >::type end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl()) |
template<int dim, int spacedim, class ITERATOR , class ASSEMBLER > | |
void | MeshWorker::integration_loop (ITERATOR begin, typename identity< ITERATOR >::type end, DoFInfo< dim, spacedim > &dof_info, IntegrationInfoBox< dim, spacedim > &box, const LocalIntegrator< dim, spacedim > &integrator, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl()) |
template<class CellIteratorType , class ScratchData , class CopyData , class CellIteratorBaseType = typename internal::CellIteratorBaseType<CellIteratorType>::type> | |
void | MeshWorker::mesh_loop (const CellIteratorType &begin, const typename identity< CellIteratorType >::type &end, const typename identity< std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>::type &cell_worker, const typename identity< std::function< void(const CopyData &)>>::type &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>>::type &boundary_worker=std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>(), const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>>::type &face_worker=std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>(), const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8) |
template<class CellIteratorType , class ScratchData , class CopyData , class CellIteratorBaseType = typename internal::CellIteratorBaseType<CellIteratorType>::type> | |
void | MeshWorker::mesh_loop (IteratorRange< CellIteratorType > iterator_range, const typename identity< std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>::type &cell_worker, const typename identity< std::function< void(const CopyData &)>>::type &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>>::type &boundary_worker=std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>(), const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>>::type &face_worker=std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>(), const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8) |
template<class CellIteratorType , class ScratchData , class CopyData , class MainClass > | |
void | MeshWorker::mesh_loop (const CellIteratorType &begin, const typename identity< CellIteratorType >::type &end, MainClass &main_class, void(MainClass::*cell_worker)(const CellIteratorType &, ScratchData &, CopyData &), void(MainClass::*copier)(const CopyData &), const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, void(MainClass::*boundary_worker)(const CellIteratorType &, const unsigned int, ScratchData &, CopyData &)=nullptr, void(MainClass::*face_worker)(const CellIteratorType &, const unsigned int, const unsigned int, const CellIteratorType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)=nullptr, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8) |
template<class CellIteratorType , class ScratchData , class CopyData , class MainClass , class CellIteratorBaseType = typename internal::CellIteratorBaseType<CellIteratorType>::type> | |
void | MeshWorker::mesh_loop (IteratorRange< CellIteratorType > iterator_range, MainClass &main_class, void(MainClass::*cell_worker)(const CellIteratorBaseType &, ScratchData &, CopyData &), void(MainClass::*copier)(const CopyData &), const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, void(MainClass::*boundary_worker)(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)=nullptr, void(MainClass::*face_worker)(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)=nullptr, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8) |
A collection of classes and functions simplifying the coding of loops over all cells and faces. All classes and functions of this module are in the MeshWorker namespace, which also contains documentation on the usage.
void MeshWorker::cell_action | ( | ITERATOR | cell, |
DoFInfoBox< dim, DOFINFO > & | dof_info, | ||
INFOBOX & | info, | ||
const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> & | cell_worker, | ||
const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> & | boundary_worker, | ||
const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> & | face_worker, | ||
const LoopControl & | loop_control | ||
) |
The function called by loop() to perform the required actions on a cell and its faces. The three functions cell_worker
, boundary_worker
and face_worker
are the same ones handed to loop(). While there we only run the loop over all cells, here, we do a single cell and, if necessary, its faces, interior and boundary.
Upon return, the DoFInfo objects in the DoFInfoBox are filled with the data computed on the cell and each of the faces. Thus, after the execution of this function, we are ready to call DoFInfoBox::assemble() to distribute the local data into global data.
cell | is the cell we work on |
dof_info | is the object into which local results are entered. It is expected to have been set up for the right types of data. |
info | is the object containing additional data only needed for internal processing. |
cell_worker | defines the local action on each cell. |
boundary_worker | defines the local action on boundary faces |
face_worker | defines the local action on interior faces. |
loop_control | control structure to specify what actions should be performed. |
void MeshWorker::loop | ( | ITERATOR | begin, |
typename identity< ITERATOR >::type | end, | ||
DOFINFO & | dinfo, | ||
INFOBOX & | info, | ||
const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> & | cell_worker, | ||
const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> & | boundary_worker, | ||
const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> & | face_worker, | ||
ASSEMBLER & | assembler, | ||
const LoopControl & | lctrl = LoopControl() |
||
) |
The main work function of this namespace. It is a loop over all cells in an iterator range, in which cell_action() is called for each cell. Unilaterally refined interior faces are handled automatically by the loop. Most of the work in this loop is done in cell_action(), which also receives most of the parameters of this function. See the documentation there for more details.
If you don't want anything to be done on cells, interior or boundary faces to happen, simply pass the Null pointer to one of the function arguments.
void MeshWorker::integration_loop | ( | ITERATOR | begin, |
typename identity< ITERATOR >::type | end, | ||
DoFInfo< dim, spacedim > & | dof_info, | ||
IntegrationInfoBox< dim, spacedim > & | box, | ||
const LocalIntegrator< dim, spacedim > & | integrator, | ||
ASSEMBLER & | assembler, | ||
const LoopControl & | lctrl = LoopControl() |
||
) |
Simplified interface for loop() if specialized for integration, using the virtual functions in LocalIntegrator.
void MeshWorker::mesh_loop | ( | const CellIteratorType & | begin, |
const typename identity< CellIteratorType >::type & | end, | ||
const typename identity< std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>::type & | cell_worker, | ||
const typename identity< std::function< void(const CopyData &)>>::type & | copier, | ||
const ScratchData & | sample_scratch_data, | ||
const CopyData & | sample_copy_data, | ||
const AssembleFlags | flags = assemble_own_cells , |
||
const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>>::type & | boundary_worker = std::function<void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>() , |
||
const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>>::type & | face_worker = std::function<void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>() , |
||
const unsigned int | queue_length = 2 * MultithreadInfo::n_threads() , |
||
const unsigned int | chunk_size = 8 |
||
) |
This function extends the WorkStream concept to work on meshes (cells and/or faces) and handles the complicated logic for work on adaptively refined faces and parallel computation (work on faces to ghost neighbors for example). The mesh_loop
can be used to simplify operations on cells (for example assembly), on boundaries (Neumann type boundary conditions), or on interior faces (for example in discontinuous Galerkin methods).
For uniformly refined meshes, it would be relatively easy to use WorkStream::run() with a cell_worker
that also loops over faces, and takes care of assembling face terms depending on the current and neighbor cell. All user codes that do these loops would then need to insert manually the logic that identifies, for every face of the current cell, the neighboring cell, and the face index on the neighboring cell that corresponds to the current face.
This is more complicated if local refinement is enabled and the current or neighbor cells have hanging nodes. In this case it is also necessary to identify the corresponding subface on either the current or the neighbor faces.
This method externalises that logic (which is independent from user codes) and separates the assembly of face terms (internal faces, boundary faces, or faces between different subdomain ids on parallel computations) from the assembling on cells, allowing the user to specify two additional workers (a cell_worker
, a boundary_worker
, and a face_worker
) that are called automatically in each cell
, according to the specific AssembleFlags flags
that are passed. The cell_worker
is passed the cell identifier, a ScratchData object, and a CopyData object, following the same principles of WorkStream::run. Internally the function passes to boundary_worker
, in addition to the above, also a face_no
parameter that identifies the face on which the integration should be performed. The face_worker
instead needs to identify the current face unambiguously both on the cell and on the neighboring cell, and it is therefore called with six arguments (three for each cell: the actual cell, the face index, and the subface_index. If no subface integration is needed, then the subface_index is numbers::invalid_unsigned_int) in addition to the usual ScratchData and CopyData objects.
If the flag AssembleFlags::assemble_own_cells is passed, then the default behavior is to first loop over faces and do the work there, and then compute the actual work on the cell. It is possible to perform the integration on the cells after working on faces, by adding the flag AssembleFlags::cells_after_faces.
If the flag AssembleFlags::assemble_own_interior_faces_once is specified, then each interior face is visited only once, and the face_worker
is assumed to integrate all face terms at once (and add contributions to both sides of the face in a discontinuous Galerkin setting).
This method is equivalent to the WorkStream::run() method when AssembleFlags contains only assemble_own_cells
, and can be used as a drop-in replacement for that method.
The two data types ScratchData and CopyData need to have a working copy constructor. ScratchData is only used in the worker function, while CopyData is the object passed from the worker to the copier.
The queue_length argument indicates the number of items that can be live at any given time. Each item consists of chunk_size elements of the input stream that will be worked on by the worker and copier functions one after the other on the same thread.
If your data objects are large, or their constructors are expensive, it is helpful to keep in mind that queue_length copies of the ScratchData object and queue_length*chunk_size copies of the CopyData object are generated.
queue_length
and chunk_size
can be found in the documentation of the WorkStream namespace and its members.Definition at line 179 of file mesh_loop.h.
void MeshWorker::mesh_loop | ( | IteratorRange< CellIteratorType > | iterator_range, |
const typename identity< std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)>>::type & | cell_worker, | ||
const typename identity< std::function< void(const CopyData &)>>::type & | copier, | ||
const ScratchData & | sample_scratch_data, | ||
const CopyData & | sample_copy_data, | ||
const AssembleFlags | flags = assemble_own_cells , |
||
const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>>::type & | boundary_worker = std::function<void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>() , |
||
const typename identity< std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>>::type & | face_worker = std::function<void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>() , |
||
const unsigned int | queue_length = 2 * MultithreadInfo::n_threads() , |
||
const unsigned int | chunk_size = 8 |
||
) |
Same as the function above, but for iterator ranges (and, therefore, filtered iterators).
An example usage of the function for the serial case is given by
and an example usage of the function for the parallel distributed case, where the copier is only to be called on locally owned cells, is given by
Definition at line 531 of file mesh_loop.h.
void MeshWorker::mesh_loop | ( | const CellIteratorType & | begin, |
const typename identity< CellIteratorType >::type & | end, | ||
MainClass & | main_class, | ||
void(MainClass::*)(const CellIteratorType &, ScratchData &, CopyData &) | cell_worker, | ||
void(MainClass::*)(const CopyData &) | copier, | ||
const ScratchData & | sample_scratch_data, | ||
const CopyData & | sample_copy_data, | ||
const AssembleFlags | flags = assemble_own_cells , |
||
void(MainClass::*)(const CellIteratorType &, const unsigned int, ScratchData &, CopyData &) | boundary_worker = nullptr , |
||
void(MainClass::*)(const CellIteratorType &, const unsigned int, const unsigned int, const CellIteratorType &, const unsigned int, const unsigned int, ScratchData &, CopyData &) | face_worker = nullptr , |
||
const unsigned int | queue_length = 2 * MultithreadInfo::n_threads() , |
||
const unsigned int | chunk_size = 8 |
||
) |
This is a variant of the mesh_loop() function, that can be used for worker and copier functions that are member functions of a class.
The argument passed as end
must be convertible to the same type as begin
, but doesn't have to be of the same type itself. This allows to write code like mesh_loop(dof_handler.begin_active(), dof_handler.end(), ...)
where the first is of type DoFHandler::active_cell_iterator whereas the second is of type DoFHandler::raw_cell_iterator.
The queue_length
argument indicates the number of items that can be live at any given time. Each item consists of chunk_size
elements of the input stream that will be worked on by the worker and copier functions one after the other on the same thread.
queue_length
copies of the ScratchData
object and queue_length*chunk_size
copies of the CopyData
object are generated.An example usage of the function is given by
Definition at line 655 of file mesh_loop.h.
void MeshWorker::mesh_loop | ( | IteratorRange< CellIteratorType > | iterator_range, |
MainClass & | main_class, | ||
void(MainClass::*)(const CellIteratorBaseType &, ScratchData &, CopyData &) | cell_worker, | ||
void(MainClass::*)(const CopyData &) | copier, | ||
const ScratchData & | sample_scratch_data, | ||
const CopyData & | sample_copy_data, | ||
const AssembleFlags | flags = assemble_own_cells , |
||
void(MainClass::*)(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &) | boundary_worker = nullptr , |
||
void(MainClass::*)(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &) | face_worker = nullptr , |
||
const unsigned int | queue_length = 2 * MultithreadInfo::n_threads() , |
||
const unsigned int | chunk_size = 8 |
||
) |
Same as the function above, but for iterator ranges (and, therefore, filtered iterators).
An example usage of the function for the serial case is given by
and an example usage of the function for the parallel distributed case, where the copier is only to be called on locally owned cells, is given by
Definition at line 823 of file mesh_loop.h.