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
Namespaces | Classes | Functions
The MeshWorker interface
Collaboration diagram for The MeshWorker interface:

Namespaces

namespace  MeshWorker
 
namespace  MeshWorker::Assembler
 

Classes

class  MeshWorker::Assembler::ResidualLocalBlocksToGlobalBlocks< VectorType >
 
class  MeshWorker::Assembler::MatrixLocalBlocksToGlobalBlocks< MatrixType, number >
 
class  MeshWorker::Assembler::MGMatrixLocalBlocksToGlobalBlocks< MatrixType, number >
 
class  MeshWorker::DoFInfo< dim, spacedim, number >
 
class  MeshWorker::DoFInfoBox< dim, DOFINFO >
 
class  MeshWorker::Assembler::Functional< number >
 
class  MeshWorker::Assembler::CellsAndFaces< number >
 
class  MeshWorker::IntegrationInfo< dim, spacedim >
 
class  MeshWorker::IntegrationInfoBox< dim, spacedim >
 
class  MeshWorker::LocalIntegrator< dim, spacedim, number >
 
class  MeshWorker::LocalResults< number >
 
class  MeshWorker::Assembler::ResidualSimple< VectorType >
 
class  MeshWorker::Assembler::MatrixSimple< MatrixType >
 
class  MeshWorker::Assembler::MGMatrixSimple< MatrixType >
 
class  MeshWorker::Assembler::SystemSimple< MatrixType, VectorType >
 
class  MeshWorker::VectorSelector
 
class  MeshWorker::VectorDataBase< dim, spacedim, Number >
 
class  MeshWorker::VectorData< VectorType, dim, spacedim >
 
class  MeshWorker::MGVectorData< VectorType, dim, spacedim >
 

Functions

template<class INFOBOX , class DOFINFO , int dim, int spacedim, typename IteratorType >
void MeshWorker::cell_action (IteratorType 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 , typename AssemblerType , typename IteratorType >
void MeshWorker::loop (IteratorType begin, std_cxx20::type_identity_t< IteratorType > 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, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
 
template<int dim, int spacedim, typename IteratorType , typename AssemblerType >
void MeshWorker::integration_loop (IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DoFInfo< dim, spacedim > &dof_info, IntegrationInfoBox< dim, spacedim > &box, const LocalIntegrator< dim, spacedim > &integrator, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
 
template<typename CellIteratorType , class ScratchData , class CopyData , typename CellIteratorBaseType = typename internal::CellIteratorBaseType<CellIteratorType>::type>
void MeshWorker::mesh_loop (const CellIteratorType &begin, const CellIteratorType &end, const CellWorkerFunctionType &cell_worker, const CopierType &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const BoundaryWorkerFunctionType &boundary_worker=BoundaryWorkerFunctionType(), const FaceWorkerFunctionType &face_worker=FaceWorkerFunctionType(), const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
 
template<typename CellIteratorType , class ScratchData , class CopyData , typename CellIteratorBaseType = typename internal::CellIteratorBaseType<CellIteratorType>::type>
void MeshWorker::mesh_loop (IteratorRange< CellIteratorType > iterator_range, const std_cxx20::type_identity_t< std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)> > &cell_worker, const std_cxx20::type_identity_t< std::function< void(const CopyData &)> > &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const std_cxx20::type_identity_t< std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)> > &boundary_worker=std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>(), const std_cxx20::type_identity_t< std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)> > &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<typename CellIteratorType , class ScratchData , class CopyData , class MainClass >
void MeshWorker::mesh_loop (const CellIteratorType &begin, const std_cxx20::type_identity_t< CellIteratorType > &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<typename CellIteratorType , class ScratchData , class CopyData , class MainClass , typename 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)
 

Detailed Description

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.

Function Documentation

◆ cell_action()

template<class INFOBOX , class DOFINFO , int dim, int spacedim, typename IteratorType >
void MeshWorker::cell_action ( IteratorType  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.

Parameters
cellis the cell we work on
dof_infois the object into which local results are entered. It is expected to have been set up for the right types of data.
infois the object containing additional data only needed for internal processing.
cell_workerdefines the local action on each cell.
boundary_workerdefines the local action on boundary faces
face_workerdefines the local action on interior faces.
loop_controlcontrol structure to specify what actions should be performed.

Definition at line 194 of file loop.h.

◆ loop()

template<int dim, int spacedim, class DOFINFO , class INFOBOX , typename AssemblerType , typename IteratorType >
void MeshWorker::loop ( IteratorType  begin,
std_cxx20::type_identity_t< IteratorType >  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,
AssemblerType &  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.

Definition at line 442 of file loop.h.

◆ integration_loop()

template<int dim, int spacedim, typename IteratorType , typename AssemblerType >
void MeshWorker::integration_loop ( IteratorType  begin,
std_cxx20::type_identity_t< IteratorType >  end,
DoFInfo< dim, spacedim > &  dof_info,
IntegrationInfoBox< dim, spacedim > &  box,
const LocalIntegrator< dim, spacedim > &  integrator,
AssemblerType &  assembler,
const LoopControl lctrl = LoopControl() 
)

Simplified interface for loop() if specialized for integration, using the virtual functions in LocalIntegrator.

Definition at line 501 of file loop.h.

◆ mesh_loop() [1/4]

template<typename CellIteratorType , class ScratchData , class CopyData , typename CellIteratorBaseType = typename internal::CellIteratorBaseType<CellIteratorType>::type>
void MeshWorker::mesh_loop ( const CellIteratorType &  begin,
const CellIteratorType &  end,
const CellWorkerFunctionType cell_worker,
const CopierType &  copier,
const ScratchData sample_scratch_data,
const CopyData sample_copy_data,
const AssembleFlags  flags = assemble_own_cells,
const BoundaryWorkerFunctionType boundary_worker = BoundaryWorkerFunctionType(),
const FaceWorkerFunctionType face_worker = FaceWorkerFunctionType(),
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). The function is used in a number of tutorials, including step-12, step-16, and step-47, to name just a few.

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 externalizes 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 CopyData object is reset to the value provided to this function every time this function visits a new cell (where it then calls the cell and face workers). In other words, no state carries over between calling the copier on one cell and the cell_worker/face_worker/boundary_worker functions on the next cell, and user code needs not reset the copy object either at the beginning of the cell integration or end of the copy operation. Resetting the state of the copier inside of a face_worker or boundary_worker constitutes a bug, and may lead to some unexpected results. The following example shows what is not permissible, as the copier is potentially shared among numerous faces on a cell:

using CellIteratorType = decltype(dof_handler.begin_active());
ScratchData scratch(...);
CopyData copy(...);
std::function<void(const CellIteratorType &, ScratchData &, CopyData &)>
empty_cell_worker;
auto boundary_worker = [...] (
const CellIteratorType &cell,
const unsigned int face,
ScratchData &scratch_data,
CopyData &copy_data)
{
const auto &fe_face_values = scratch_data.reinit(cell, face);
copy_data = CopyData(...); // This is an error, as we lose the
// accumulation that has been performed on
// other boundary faces of the same cell.
for (unsigned int q_point = 0;
q_point < fe_face_values.n_quadrature_points;
++q_point)
{
copy_data.vectors[0][0] += 1.0 * fe_face_values.JxW(q_point);
}
};
double value = 0;
auto copier = [...](const CopyData &copy_data)
{
value += copy_data.vectors[0][0]; // Contributions from some faces may
// be missing.
};
MeshWorker::mesh_loop(dof_handler.active_cell_iterators(),
empty_cell_worker, copier,
scratch, copy,
boundary_worker);
void mesh_loop(const CellIteratorType &begin, const CellIteratorType &end, const CellWorkerFunctionType &cell_worker, const CopierType &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const BoundaryWorkerFunctionType &boundary_worker=BoundaryWorkerFunctionType(), const FaceWorkerFunctionType &face_worker=FaceWorkerFunctionType(), const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
Definition mesh_loop.h:281
@ assemble_boundary_faces

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.

Note
The types of the function arguments and the default values (empty worker functions) displayed in the Doxygen documentation here are slightly simplified compared to the real types.
More information about requirements on template types and meaning of queue_length and chunk_size can be found in the documentation of the WorkStream namespace and its members.

Definition at line 281 of file mesh_loop.h.

◆ mesh_loop() [2/4]

template<typename CellIteratorType , class ScratchData , class CopyData , typename CellIteratorBaseType = typename internal::CellIteratorBaseType<CellIteratorType>::type>
void MeshWorker::mesh_loop ( IteratorRange< CellIteratorType >  iterator_range,
const std_cxx20::type_identity_t< std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)> > &  cell_worker,
const std_cxx20::type_identity_t< std::function< void(const CopyData &)> > &  copier,
const ScratchData sample_scratch_data,
const CopyData sample_copy_data,
const AssembleFlags  flags = assemble_own_cells,
const std_cxx20::type_identity_t< std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)> > &  boundary_worker = std::function<void(const CellIteratorBaseType &,                         const unsigned int,                         ScratchData &,                         CopyData &)>(),
const std_cxx20::type_identity_t< std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)> > &  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

using CellIteratorType = decltype(dof_handler.begin_active());
ScratchData scratch(...);
CopyData copy(...);
auto cell_worker = [...] (
const CellIteratorType &cell,
ScratchData &scratch_data,
CopyData &copy_data)
{
...
};
auto copier = [...](const CopyData &copy_data)
{
...
};
MeshWorker::mesh_loop(dof_handler.active_cell_iterators(),
cell_worker, copier,
scratch, copy,

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

using CellIteratorType = decltype(dof_handler.begin_active());
ScratchData scratch(...);
CopyData copy(...);
auto cell_worker = [...] (
const CellIteratorType &cell,
ScratchData &scratch_data,
CopyData &copy_data)
{
...
};
auto copier = [...](const CopyData &copy_data)
{
...
};
const auto filtered_iterator_range =
filter_iterators(dof_handler.active_cell_iterators(),
MeshWorker::mesh_loop(filtered_iterator_range,
cell_worker, copier,
scratch, copy,
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators(IteratorRange< BaseIterator > i, const Predicate &p)

Definition at line 715 of file mesh_loop.h.

◆ mesh_loop() [3/4]

template<typename CellIteratorType , class ScratchData , class CopyData , class MainClass >
void MeshWorker::mesh_loop ( const CellIteratorType &  begin,
const std_cxx20::type_identity_t< CellIteratorType > &  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.

Note
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.

An example usage of the function is given by

struct ScratchData;
struct CopyData;
template <int dim, int spacedim>
class MyClass
{
public:
void
cell_worker(const CellIteratorType &cell, ScratchData &, CopyData &);
void
copier(const CopyData &);
...
};
...
MyClass<dim, spacedim> my_class;
ScratchData scratch;
CopyData copy;
tria.end(),
my_class,
&MyClass<dim, spacedim>::cell_worker,
&MyClass<dim, spacedim>::copier,
scratch,
copy,
cell_iterator end() const
active_cell_iterator begin_active(const unsigned int level=0) const

Definition at line 840 of file mesh_loop.h.

◆ mesh_loop() [4/4]

template<typename CellIteratorType , class ScratchData , class CopyData , class MainClass , typename CellIteratorBaseType = typename internal::CellIteratorBaseType<CellIteratorType>::type>
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

struct ScratchData;
struct CopyData;
template <int dim, int spacedim>
class MyClass
{
public:
void
cell_worker(const CellIteratorType &cell, ScratchData &, CopyData &);
void
copier(const CopyData &);
...
};
...
MyClass<dim, spacedim> my_class;
ScratchData scratch;
CopyData copy;
my_class,
&MyClass<dim, spacedim>::cell_worker,
&MyClass<dim, spacedim>::copier,
scratch,
copy,
IteratorRange< active_cell_iterator > active_cell_iterators() const

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

struct ScratchData;
struct CopyData;
template <int dim, int spacedim>
class MyClass
{
public:
void
cell_worker(const CellIteratorType &cell, ScratchData &, CopyData &);
void
copier(const CopyData &);
...
};
...
MyClass<dim, spacedim> my_class;
ScratchData scratch;
CopyData copy;
const auto filtered_iterator_range =
filter_iterators(distributed_tria.active_cell_iterators(),
mesh_loop(filtered_iterator_range,
my_class,
&MyClass<dim, spacedim>::cell_worker,
&MyClass<dim, spacedim>::copier,
scratch,
copy,

Definition at line 1022 of file mesh_loop.h.