DOLFIN-X
DOLFIN-X C++ interface
Classes | Public Types | Public Member Functions | List of all members
dolfinx::fem::FormIntegrals Class Reference

Integrals of a Form, including those defined over cells, interior and exterior facets, and vertices. More...

#include <FormIntegrals.h>

Public Types

enum  Type : std::int8_t { cell = 0, exterior_facet = 1, interior_facet = 2, vertex = 3 }
 Type of integral.
 

Public Member Functions

 FormIntegrals ()
 Construct empty object.
 
const std::function< void(PetscScalar *, const PetscScalar *, const PetscScalar *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> & get_tabulate_tensor (FormIntegrals::Type type, int i) const
 Get the function for 'tabulate_tensor' for integral i of given type. More...
 
void set_tabulate_tensor (FormIntegrals::Type type, int i, std::function< void(PetscScalar *, const PetscScalar *, const PetscScalar *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> fn)
 Set the function for 'tabulate_tensor' for integral i of given type. More...
 
int num_integrals (FormIntegrals::Type t) const
 Number of integrals of given type. More...
 
std::vector< int > integral_ids (FormIntegrals::Type t) const
 Get the integer IDs of integrals of type t. The IDs correspond to the domains which the integrals are defined for in the form, except ID -1, which denotes the default integral. More...
 
const std::vector< std::int32_t > & integral_domains (FormIntegrals::Type t, int i) const
 Get the list of active entities for the ith integral of type t. Note, these are not retrieved by ID, but stored in order. The IDs can be obtained with "FormIntegrals::integral_ids()". For cell integrals, a list of cells. For facet integrals, a list of facets etc. More...
 
void set_domains (FormIntegrals::Type type, const mesh::MeshTags< int > &marker)
 Set the valid domains for the integrals of a given type from a MeshTags "marker". Note the MeshTags is not stored, so if there any changes to the integration domain this must be called again. More...
 
void set_default_domains (const mesh::Mesh &mesh)
 If there exists a default integral of any type, set the list of entities for those integrals from the mesh topology. For cell integrals, this is all cells. For facet integrals, it is either all interior or all exterior facets. More...
 

Detailed Description

Integrals of a Form, including those defined over cells, interior and exterior facets, and vertices.

Member Function Documentation

◆ get_tabulate_tensor()

const std::function< void(PetscScalar *, const PetscScalar *, const PetscScalar *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> & FormIntegrals::get_tabulate_tensor ( FormIntegrals::Type  type,
int  i 
) const

Get the function for 'tabulate_tensor' for integral i of given type.

Parameters
[in]typeIntegral type
[in]iIntegral number
Returns
Function to call for tabulate_tensor

◆ integral_domains()

const std::vector< std::int32_t > & FormIntegrals::integral_domains ( FormIntegrals::Type  t,
int  i 
) const

Get the list of active entities for the ith integral of type t. Note, these are not retrieved by ID, but stored in order. The IDs can be obtained with "FormIntegrals::integral_ids()". For cell integrals, a list of cells. For facet integrals, a list of facets etc.

Parameters
[in]tIntegral type
[in]iIntegral number
Returns
List of active entities for this integral

◆ integral_ids()

std::vector< int > FormIntegrals::integral_ids ( FormIntegrals::Type  t) const

Get the integer IDs of integrals of type t. The IDs correspond to the domains which the integrals are defined for in the form, except ID -1, which denotes the default integral.

Parameters
[in]tIntegral type
Returns
List of IDs for this integral

◆ num_integrals()

int FormIntegrals::num_integrals ( FormIntegrals::Type  t) const

Number of integrals of given type.

Parameters
[in]tIntegral type
Returns
Number of integrals

◆ set_default_domains()

void FormIntegrals::set_default_domains ( const mesh::Mesh mesh)

If there exists a default integral of any type, set the list of entities for those integrals from the mesh topology. For cell integrals, this is all cells. For facet integrals, it is either all interior or all exterior facets.

Parameters
[in]meshMesh

◆ set_domains()

void FormIntegrals::set_domains ( FormIntegrals::Type  type,
const mesh::MeshTags< int > &  marker 
)

Set the valid domains for the integrals of a given type from a MeshTags "marker". Note the MeshTags is not stored, so if there any changes to the integration domain this must be called again.

Parameters
[in]typeIntegral type
[in]markerMeshTags mapping entities to integrals

◆ set_tabulate_tensor()

void FormIntegrals::set_tabulate_tensor ( FormIntegrals::Type  type,
int  i,
std::function< void(PetscScalar *, const PetscScalar *, const PetscScalar *, const double *, const int *, const std::uint8_t *, const std::uint32_t)>  fn 
)

Set the function for 'tabulate_tensor' for integral i of given type.

Parameters
[in]typeIntegral type
[in]iIntegral number
[in]fntabulate function

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