Reference documentation for deal.II version GIT relicensing-478-g3275795167 2024-04-24 07:10: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
Public Types | Public Member Functions | List of all members
MeshWorker::IntegrationInfoBox< dim, spacedim > Class Template Reference

#include <deal.II/meshworker/integration_info.h>

Public Types

using CellInfo = IntegrationInfo< dim, spacedim >
 

Public Member Functions

 IntegrationInfoBox ()
 
void initialize (const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const BlockInfo *block_info=nullptr)
 
template<typename VectorType >
void initialize (const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const AnyData &data, const VectorType &dummy, const BlockInfo *block_info=nullptr)
 
template<typename VectorType >
void initialize (const FiniteElement< dim, spacedim > &el, const Mapping< dim, spacedim > &mapping, const AnyData &data, const MGLevelObject< VectorType > &dummy, const BlockInfo *block_info=nullptr)
 
template<typename VectorType >
void initialize (const FiniteElement< dim, sdim > &el, const Mapping< dim, sdim > &mapping, const AnyData &data, const VectorType &, const BlockInfo *block_info)
 
template<typename VectorType >
void initialize (const FiniteElement< dim, sdim > &el, const Mapping< dim, sdim > &mapping, const AnyData &data, const MGLevelObject< VectorType > &, const BlockInfo *block_info)
 

Public Attributes

Data vectors
MeshWorker::VectorSelector cell_selector
 
MeshWorker::VectorSelector boundary_selector
 
MeshWorker::VectorSelector face_selector
 
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > cell_data
 
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > boundary_data
 
std::shared_ptr< MeshWorker::VectorDataBase< dim, spacedim > > face_data
 

FEValues setup

UpdateFlags cell_flags
 
UpdateFlags boundary_flags
 
UpdateFlags face_flags
 
UpdateFlags neighbor_flags
 
Quadrature< dim > cell_quadrature
 
Quadrature< dim - 1 > boundary_quadrature
 
Quadrature< dim - 1 > face_quadrature
 
void initialize_update_flags (bool neighbor_geometry=false)
 
void add_update_flags_all (const UpdateFlags flags)
 
void add_update_flags_cell (const UpdateFlags flags)
 
void add_update_flags_boundary (const UpdateFlags flags)
 
void add_update_flags_face (const UpdateFlags flags)
 
void add_update_flags (const UpdateFlags flags, const bool cell=true, const bool boundary=true, const bool face=true, const bool neighbor=true)
 
void initialize_gauss_quadrature (unsigned int n_cell_points, unsigned int n_boundary_points, unsigned int n_face_points, const bool force=true)
 
std::size_t memory_consumption () const
 

Interface for MeshWorker::loop()

CellInfo cell
 
CellInfo boundary
 
CellInfo face
 
CellInfo subface
 
CellInfo neighbor
 
template<class DOFINFO >
void post_cell (const DoFInfoBox< dim, DOFINFO > &)
 
template<class DOFINFO >
void post_faces (const DoFInfoBox< dim, DOFINFO > &)
 

Detailed Description

template<int dim, int spacedim = dim>
class MeshWorker::IntegrationInfoBox< dim, spacedim >

The object holding the scratch data for integrating over cells and faces. IntegrationInfoBox serves three main purposes:

  1. It provides the interface needed by MeshWorker::loop(), namely the two functions post_cell() and post_faces(), as well as the data members cell, boundary, face, subface, and neighbor.

  2. It contains all information needed to initialize the FEValues and FEFaceValues objects in the IntegrationInfo data members.

  3. It stores information on finite element vectors and whether their data should be used to compute values or derivatives of functions at quadrature points.

  4. It makes educated guesses on quadrature rules and update flags, so that minimal code has to be written when default parameters are sufficient.

In order to allow for sufficient generality, a few steps have to be undertaken to use this class.

First, you should consider if you need values from any vectors in a AnyData object. If so, fill the VectorSelector objects cell_selector, boundary_selector and face_selector with their names and the data type (value, gradient, Hessian) to be extracted.

Afterwards, you will need to consider UpdateFlags for FEValues objects. A good start is initialize_update_flags(), which looks at the selectors filled before and adds all the flags needed to get the selection. Additional flags can be set with add_update_flags().

Finally, we need to choose quadrature formulas. In the simplest case, you might be happy with the default settings, which are n-point Gauss formulas. If only derivatives of the shape functions are used (update_values is not set) n equals the highest polynomial degree in the FiniteElement, if update_values is set, n is one higher than this degree. If you choose to use Gauss formulas of other size, use initialize_gauss_quadrature() with appropriate values. Otherwise, you can fill the variables cell_quadrature, boundary_quadrature and face_quadrature directly.

In order to save time, you can set the variables boundary_fluxes and interior_fluxes of the base class to false, thus telling the Meshworker::loop() not to loop over those faces.

All the information in here is used to set up IntegrationInfo objects correctly, typically in an IntegrationInfoBox.

Definition at line 296 of file integration_info.h.

Member Typedef Documentation

◆ CellInfo

template<int dim, int spacedim = dim>
using MeshWorker::IntegrationInfoBox< dim, spacedim >::CellInfo = IntegrationInfo<dim, spacedim>

The type of the info object for cells.

Definition at line 302 of file integration_info.h.

Constructor & Destructor Documentation

◆ IntegrationInfoBox()

template<int dim, int spacedim = dim>
MeshWorker::IntegrationInfoBox< dim, spacedim >::IntegrationInfoBox ( )

Default constructor.

Member Function Documentation

◆ initialize() [1/5]

template<int dim, int spacedim = dim>
void MeshWorker::IntegrationInfoBox< dim, spacedim >::initialize ( const FiniteElement< dim, spacedim > &  el,
const Mapping< dim, spacedim > &  mapping,
const BlockInfo block_info = nullptr 
)

Initialize the IntegrationInfo objects contained.

Before doing so, add update flags necessary to produce the data needed and also set uninitialized quadrature rules to Gauss formulas, which integrate polynomial bilinear forms exactly.

◆ initialize() [2/5]

template<int dim, int spacedim = dim>
template<typename VectorType >
void MeshWorker::IntegrationInfoBox< dim, spacedim >::initialize ( const FiniteElement< dim, spacedim > &  el,
const Mapping< dim, spacedim > &  mapping,
const AnyData data,
const VectorType &  dummy,
const BlockInfo block_info = nullptr 
)

Initialize the IntegrationInfo objects contained.

Before doing so, add update flags necessary to produce the data needed and also set uninitialized quadrature rules to Gauss formulas, which integrate polynomial bilinear forms exactly.

◆ initialize() [3/5]

template<int dim, int spacedim = dim>
template<typename VectorType >
void MeshWorker::IntegrationInfoBox< dim, spacedim >::initialize ( const FiniteElement< dim, spacedim > &  el,
const Mapping< dim, spacedim > &  mapping,
const AnyData data,
const MGLevelObject< VectorType > &  dummy,
const BlockInfo block_info = nullptr 
)

Initialize the IntegrationInfo objects contained.

Before doing so, add update flags necessary to produce the data needed and also set uninitialized quadrature rules to Gauss formulas, which integrate polynomial bilinear forms exactly.

◆ initialize_update_flags()

template<int dim, int spacedim = dim>
void MeshWorker::IntegrationInfoBox< dim, spacedim >::initialize_update_flags ( bool  neighbor_geometry = false)

Call this function before initialize() in order to guess the update flags needed, based on the data selected.

When computing face fluxes, we normally can use the geometry (integration weights and normal vectors) from the original cell and thus can avoid updating these values on the neighboring cell. Set neighbor_geometry to true in order to initialize these values as well.

◆ add_update_flags_all()

template<int dim, int sdim>
void MeshWorker::IntegrationInfoBox< dim, sdim >::add_update_flags_all ( const UpdateFlags  flags)
inline

Add FEValues UpdateFlags for integration on all objects (cells, boundary faces and all interior faces).

Definition at line 751 of file integration_info.h.

◆ add_update_flags_cell()

template<int dim, int sdim>
void MeshWorker::IntegrationInfoBox< dim, sdim >::add_update_flags_cell ( const UpdateFlags  flags)
inline

Add FEValues UpdateFlags for integration on cells.

Definition at line 759 of file integration_info.h.

◆ add_update_flags_boundary()

template<int dim, int sdim>
void MeshWorker::IntegrationInfoBox< dim, sdim >::add_update_flags_boundary ( const UpdateFlags  flags)
inline

Add FEValues UpdateFlags for integration on boundary faces.

Definition at line 767 of file integration_info.h.

◆ add_update_flags_face()

template<int dim, int sdim>
void MeshWorker::IntegrationInfoBox< dim, sdim >::add_update_flags_face ( const UpdateFlags  flags)
inline

Add FEValues UpdateFlags for integration on interior faces.

Definition at line 776 of file integration_info.h.

◆ add_update_flags()

template<int dim, int spacedim = dim>
void MeshWorker::IntegrationInfoBox< dim, spacedim >::add_update_flags ( const UpdateFlags  flags,
const bool  cell = true,
const bool  boundary = true,
const bool  face = true,
const bool  neighbor = true 
)

Add additional update flags to the ones already set in this program. The four boolean flags indicate whether the additional flags should be set for cell, boundary, interelement face for the cell itself or neighbor cell, or any combination thereof.

◆ initialize_gauss_quadrature()

template<int dim, int sdim>
void MeshWorker::IntegrationInfoBox< dim, sdim >::initialize_gauss_quadrature ( unsigned int  n_cell_points,
unsigned int  n_boundary_points,
unsigned int  n_face_points,
const bool  force = true 
)
inline

Assign n-point Gauss quadratures to each of the quadrature rules. Here, a size of zero points means that no loop over these grid entities should be performed.

If the parameter force is true, then all quadrature sets are filled with new quadrature rules. If it is false, then only empty rules are changed.

Definition at line 735 of file integration_info.h.

◆ memory_consumption()

template<int dim, int spacedim = dim>
std::size_t MeshWorker::IntegrationInfoBox< dim, spacedim >::memory_consumption ( ) const

The memory used by this object.

◆ post_cell()

template<int dim, int sdim>
template<class DOFINFO >
void MeshWorker::IntegrationInfoBox< dim, sdim >::post_cell ( const DoFInfoBox< dim, DOFINFO > &  )

A callback function which is called in the loop over all cells, after the action on a cell has been performed and before the faces are dealt with.

In order for this function to have this effect, at least either of the arguments boundary_worker or face_worker arguments of loop() should be nonzero. Additionally, cells_first should be true. If cells_first is false, this function is called before any action on a cell is taken.

And empty function in this class, but can be replaced in other classes given to loop() instead.

See loop() and cell_action() for more details of how this function can be used.

Definition at line 891 of file integration_info.h.

◆ post_faces()

template<int dim, int sdim>
template<class DOFINFO >
void MeshWorker::IntegrationInfoBox< dim, sdim >::post_faces ( const DoFInfoBox< dim, DOFINFO > &  )

A callback function which is called in the loop over all cells, after the action on the faces of a cell has been performed and before the cell itself is dealt with (assumes cells_first is false).

In order for this function to have a reasonable effect, at least either of the arguments boundary_worker or face_worker arguments of loop() should be nonzero. Additionally, cells_first should be false.

And empty function in this class, but can be replaced in other classes given to loop() instead.

See loop() and cell_action() for more details of how this function can be used.

Definition at line 898 of file integration_info.h.

◆ initialize() [4/5]

template<int dim, int spacedim = dim>
template<typename VectorType >
void MeshWorker::IntegrationInfoBox< dim, spacedim >::initialize ( const FiniteElement< dim, sdim > &  el,
const Mapping< dim, sdim > &  mapping,
const AnyData data,
const VectorType &  ,
const BlockInfo block_info 
)

Definition at line 818 of file integration_info.h.

◆ initialize() [5/5]

template<int dim, int spacedim = dim>
template<typename VectorType >
void MeshWorker::IntegrationInfoBox< dim, spacedim >::initialize ( const FiniteElement< dim, sdim > &  el,
const Mapping< dim, sdim > &  mapping,
const AnyData data,
const MGLevelObject< VectorType > &  ,
const BlockInfo block_info 
)

Definition at line 854 of file integration_info.h.

Member Data Documentation

◆ cell_flags

template<int dim, int spacedim = dim>
UpdateFlags MeshWorker::IntegrationInfoBox< dim, spacedim >::cell_flags

The set of update flags for boundary cell integration.

Defaults to update_JxW_values.

Definition at line 431 of file integration_info.h.

◆ boundary_flags

template<int dim, int spacedim = dim>
UpdateFlags MeshWorker::IntegrationInfoBox< dim, spacedim >::boundary_flags

The set of update flags for boundary face integration.

Defaults to update_JxW_values and update_normal_vectors.

Definition at line 437 of file integration_info.h.

◆ face_flags

template<int dim, int spacedim = dim>
UpdateFlags MeshWorker::IntegrationInfoBox< dim, spacedim >::face_flags

The set of update flags for interior face integration.

Defaults to update_JxW_values and update_normal_vectors.

Definition at line 444 of file integration_info.h.

◆ neighbor_flags

template<int dim, int spacedim = dim>
UpdateFlags MeshWorker::IntegrationInfoBox< dim, spacedim >::neighbor_flags

The set of update flags for interior face integration.

Defaults to update_default, since quadrature weights are taken from the other cell.

Definition at line 452 of file integration_info.h.

◆ cell_quadrature

template<int dim, int spacedim = dim>
Quadrature<dim> MeshWorker::IntegrationInfoBox< dim, spacedim >::cell_quadrature

The quadrature rule used on cells.

Definition at line 457 of file integration_info.h.

◆ boundary_quadrature

template<int dim, int spacedim = dim>
Quadrature<dim - 1> MeshWorker::IntegrationInfoBox< dim, spacedim >::boundary_quadrature

The quadrature rule used on boundary faces.

Definition at line 462 of file integration_info.h.

◆ face_quadrature

template<int dim, int spacedim = dim>
Quadrature<dim - 1> MeshWorker::IntegrationInfoBox< dim, spacedim >::face_quadrature

The quadrature rule used on interior faces.

Definition at line 467 of file integration_info.h.

◆ cell_selector

template<int dim, int spacedim = dim>
MeshWorker::VectorSelector MeshWorker::IntegrationInfoBox< dim, spacedim >::cell_selector

Initialize the VectorSelector objects cell_selector, boundary_selector and face_selector in order to save computational effort. If no selectors are used, then values for all named vectors in DoFInfo::global_data will be computed in all quadrature points.

This function will also add UpdateFlags to the flags stored in this class. Select the vectors from DoFInfo::global_data that should be computed in the quadrature points on cells.

Definition at line 488 of file integration_info.h.

◆ boundary_selector

template<int dim, int spacedim = dim>
MeshWorker::VectorSelector MeshWorker::IntegrationInfoBox< dim, spacedim >::boundary_selector

Select the vectors from DoFInfo::global_data that should be computed in the quadrature points on boundary faces.

Definition at line 494 of file integration_info.h.

◆ face_selector

template<int dim, int spacedim = dim>
MeshWorker::VectorSelector MeshWorker::IntegrationInfoBox< dim, spacedim >::face_selector

Select the vectors from DoFInfo::global_data that should be computed in the quadrature points on interior faces.

Definition at line 500 of file integration_info.h.

◆ cell_data

template<int dim, int spacedim = dim>
std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > MeshWorker::IntegrationInfoBox< dim, spacedim >::cell_data

Definition at line 502 of file integration_info.h.

◆ boundary_data

template<int dim, int spacedim = dim>
std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > MeshWorker::IntegrationInfoBox< dim, spacedim >::boundary_data

Definition at line 503 of file integration_info.h.

◆ face_data

template<int dim, int spacedim = dim>
std::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > MeshWorker::IntegrationInfoBox< dim, spacedim >::face_data

Definition at line 504 of file integration_info.h.

◆ cell

template<int dim, int spacedim = dim>
CellInfo MeshWorker::IntegrationInfoBox< dim, spacedim >::cell

The info object for a cell.

Definition at line 555 of file integration_info.h.

◆ boundary

template<int dim, int spacedim = dim>
CellInfo MeshWorker::IntegrationInfoBox< dim, spacedim >::boundary

The info object for a boundary face.

Definition at line 559 of file integration_info.h.

◆ face

template<int dim, int spacedim = dim>
CellInfo MeshWorker::IntegrationInfoBox< dim, spacedim >::face

The info object for a regular interior face, seen from the first cell.

Definition at line 563 of file integration_info.h.

◆ subface

template<int dim, int spacedim = dim>
CellInfo MeshWorker::IntegrationInfoBox< dim, spacedim >::subface

The info object for the refined side of an interior face seen from the first cell.

Definition at line 568 of file integration_info.h.

◆ neighbor

template<int dim, int spacedim = dim>
CellInfo MeshWorker::IntegrationInfoBox< dim, spacedim >::neighbor

The info object for an interior face, seen from the other cell.

Definition at line 572 of file integration_info.h.


The documentation for this class was generated from the following file: