Reference documentation for deal.II version GIT relicensing-245-g36f19064f7 2024-03-29 07: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 | Typedefs | Enumerations | Functions
MeshWorker Namespace Reference

Namespaces

namespace  Assembler
 
namespace  internal
 

Classes

struct  CopyData
 
class  DoFInfo
 
class  DoFInfoBox
 
class  IntegrationInfo
 
class  IntegrationInfoBox
 
class  LocalIntegrator
 
class  LocalResults
 
class  LoopControl
 
class  MGVectorData
 
class  ScratchData
 
class  VectorData
 
class  VectorDataBase
 
class  VectorSelector
 

Typedefs

using CellWorkerFunctionType = std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)>
 
using CopierFunctionType = std::function< void(const CopyData &)>
 
using BoundaryWorkerFunctionType = std::function< void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>
 
using FaceWorkerFunctionType = std::function< void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>
 

Enumerations

enum  AssembleFlags {
  assemble_nothing = 0 , assemble_own_cells = 0x0001 , assemble_ghost_cells = 0x0002 , assemble_own_interior_faces_once = 0x0004 ,
  assemble_own_interior_faces_both = 0x0008 , assemble_ghost_faces_once = 0x0010 , assemble_ghost_faces_both = 0x0020 , assemble_boundary_faces = 0x0040 ,
  cells_after_faces = 0x0080 , work_on_cells = assemble_own_cells | assemble_ghost_cells , work_on_faces , work_on_boundary = assemble_boundary_faces
}
 

Functions

template<typename StreamType >
StreamType & operator<< (StreamType &s, AssembleFlags u)
 
AssembleFlags operator| (AssembleFlags f1, AssembleFlags f2)
 
AssembleFlagsoperator|= (AssembleFlags &f1, AssembleFlags f2)
 
AssembleFlags operator& (AssembleFlags f1, AssembleFlags f2)
 
AssembleFlagsoperator&= (AssembleFlags &f1, AssembleFlags f2)
 
template<class INFOBOX , class DOFINFO , int dim, int spacedim, typename IteratorType >
void 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 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 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 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 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 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 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 functions and classes for the mesh loops that are an ubiquitous part of each finite element program.

The workhorse of this namespace is the loop() function, which implements a completely generic loop over all mesh cells. Since the calls to loop() are error-prone due to its generality, for many applications it is advisable to derive a class from MeshWorker::LocalIntegrator and use the less general integration_loop() instead.

The loop() depends on certain objects handed to it as arguments. These objects are of two types, info objects like DoFInfo and IntegrationInfo and worker objects like LocalWorker and IntegrationWorker.

Worker objects usually do two different jobs: first, they compute the local contribution of a cell or face to the global operation. Second, they assemble this local contribution into the global result, whether a functional, a form or a bilinear form. While the first job is particular to the problem being solved, the second is generic and only depends on the data structures. Therefore, base classes for workers assembling into global data are provided in the namespace Assembler.

Template argument types

The functions loop() and cell_action() take some arguments which are template parameters. Let us list the minimum requirements for these classes here and describe their properties.

ITERATOR

Any object that has an operator++() and points to a TriaAccessor or derived class.

DOFINFO

For an example implementation, refer to the class template DoFInfo. In order to work with cell_action() and loop(), DOFINFO needs to follow the following interface.

class DOFINFO
{
private:
DOFINFO();
DOFINFO(const DOFINFO&);
DOFINFO& operator=(const DOFINFO&);
public:
template <class CellIt>
void reinit(const CellIt& c);
template <class CellIt, class FaceIt>
void reinit(const CellIt& c, const FaceIt& f, const unsigned int n);
template <class CellIt, class FaceIt>
void reinit(const CellIt& c, const FaceIt& f, const unsigned int n,
const unsigned int s);
friend template class DoFInfoBox<int dim, DOFINFO>;
};

The three private functions are called by DoFInfoBox and should not be needed elsewhere. Obviously, they can be made public and then the friend declaration at the end may be missing.

Additionally, you will need at least one public constructor. Furthermore DOFINFO is pretty useless yet: functions to interface with INTEGRATIONINFO and ASSEMBLER are needed.

DOFINFO objects are gathered in a DoFInfoBox. In those objects, we store the results of local operations on each cell and its faces. Once all this information has been gathered, an ASSEMBLER is used to assemble it into global data.

INFOBOX

This type is exemplified in IntegrationInfoBox. It collects the input data for actions on cells and faces in INFO objects (see below). It provides the following interface to loop() and cell_action():

class INFOBOX
{
public:
template <int dim, class DOFINFO>
void post_cell(const DoFInfoBox<dim, DOFINFO>&);
template <int dim, class DOFINFO>
void post_faces(const DoFInfoBox<dim, DOFINFO>&);
INFO cell;
INFO boundary;
INFO face;
INFO subface;
INFO neighbor;
};

The main purpose of this class is gathering the five INFO objects, which contain the temporary data used on each cell or face. The requirements on these objects are listed below. Here, we only note that there need to be these 5 objects with the names listed above.

The two function templates are call back functions called in cell_action(). The first is called before the faces are worked on, the second after the faces.

INFO

See IntegrationInfo for an example of these objects. They contain the temporary data needed on each cell or face to compute the result. The MeshWorker only uses the interface

class INFO
{
public:
void reinit(const DOFINFO& i);
};

Simplified interfaces

Since the loop() is fairly general, a specialization integration_loop() is available, which is a wrapper around loop() with a simplified interface.

The integration_loop() function loop takes most of the information that it needs to pass to loop() from an IntegrationInfoBox object. Its use is explained in step-12, but in short it requires functions that do the local integration on a cell, interior or boundary face, and it needs an object (called "assembler") that copies these local contributions into the global matrix and right hand side objects.

Before we can run the integration loop, we have to initialize several data structures in our IntegrationWorker and assembler objects. For instance, we have to decide on the quadrature rule or we may need more than the default update flags.

Typedef Documentation

◆ CellWorkerFunctionType

using MeshWorker::CellWorkerFunctionType = typedef std::function< void(const CellIteratorBaseType &, ScratchData &, CopyData &)>

This alias introduces a friendly and short name for the function type for the cell worker used in mesh_loop().

Definition at line 106 of file mesh_loop.h.

◆ CopierFunctionType

using MeshWorker::CopierFunctionType = typedef std::function<void(const CopyData &)>

This alias introduces a friendly and short name for the function type for the cell worker used in mesh_loop().

Definition at line 113 of file mesh_loop.h.

◆ BoundaryWorkerFunctionType

using MeshWorker::BoundaryWorkerFunctionType = typedef std::function<void(const CellIteratorBaseType &, const unsigned int, ScratchData &, CopyData &)>

This alias introduces a friendly and short name for the function type for the boundary worker used in mesh_loop().

Definition at line 119 of file mesh_loop.h.

◆ FaceWorkerFunctionType

using MeshWorker::FaceWorkerFunctionType = typedef std::function<void(const CellIteratorBaseType &, const unsigned int, const unsigned int, const CellIteratorBaseType &, const unsigned int, const unsigned int, ScratchData &, CopyData &)>

This alias introduces a friendly and short name for the function type for the face worker used in mesh_loop().

Definition at line 129 of file mesh_loop.h.

Enumeration Type Documentation

◆ AssembleFlags

The enum type given to the mesh_loop() function, telling that function which elements need to be assembled.

You can select more than one flag by concatenation using the bitwise or operator|(AssembleFlags,AssembleFlags).

Enumerator
assemble_nothing 

Do Nothing.

assemble_own_cells 

Assemble on locally owned cells.

assemble_ghost_cells 

Assemble on ghost cells.

assemble_own_interior_faces_once 

Assemble on interior faces between two locally owned cells, visiting each face only once.

assemble_own_interior_faces_both 

Assemble on interior faces between two locally owned cells, visiting each interior face twice, once from each of the two adjacent cells.

assemble_ghost_faces_once 

Assemble on faces between a locally owned cell and a ghost cell, making sure that only one of the processes will assemble these faces (from the finer side or the process with the lower mpi rank).

assemble_ghost_faces_both 

Assemble on faces between a locally owned cell and a ghost cell. Both processes will assemble these faces. Note that they are never assembled from both sides on a single process.

assemble_boundary_faces 

Assemble on boundary faces of the locally owned cells.

cells_after_faces 

By default we assemble cell integrals before face integrals. If this flag is specified, cells will be assembled after faces and boundaries.

work_on_cells 

Combination of flags to determine if any work on cells is done.

work_on_faces 

Combination of flags to determine if any work is done on faces.

work_on_boundary 

Combination of flags to determine if any work is done on the boundary faces.

Definition at line 40 of file assemble_flags.h.

Function Documentation

◆ operator<<()

template<typename StreamType >
StreamType & MeshWorker::operator<< ( StreamType &  s,
AssembleFlags  u 
)
inline

Output operator which outputs assemble flags as a set of or'd text values.

AssembleFlags

Definition at line 114 of file assemble_flags.h.

◆ operator|()

AssembleFlags MeshWorker::operator| ( AssembleFlags  f1,
AssembleFlags  f2 
)
inline

Global operator which returns an object in which all bits are set which are either set in the first or the second argument. This operator exists since if it did not then the result of the bit-or operator | would be an integer which would in turn trigger a compiler warning when we tried to assign it to an object of type AssembleFlags.

AssembleFlags

Definition at line 146 of file assemble_flags.h.

◆ operator|=()

AssembleFlags & MeshWorker::operator|= ( AssembleFlags f1,
AssembleFlags  f2 
)
inline

Global operator which sets the bits from the second argument also in the first one.

AssembleFlags

Definition at line 161 of file assemble_flags.h.

◆ operator&()

AssembleFlags MeshWorker::operator& ( AssembleFlags  f1,
AssembleFlags  f2 
)
inline

Global operator which returns an object in which all bits are set which are set in the first as well as the second argument. This operator exists since if it did not then the result of the bit-and operator & would be an integer which would in turn trigger a compiler warning when we tried to assign it to an object of type AssembleFlags.

AssembleFlags

Definition at line 178 of file assemble_flags.h.

◆ operator&=()

AssembleFlags & MeshWorker::operator&= ( AssembleFlags f1,
AssembleFlags  f2 
)
inline

Global operator which clears all the bits in the first argument if they are not also set in the second argument.

AssembleFlags

Definition at line 192 of file assemble_flags.h.