Integrals of a Form, including those defined over cells, interior and exterior facets, and vertices.
More...
#include <FormIntegrals.h>
|
| FormIntegrals () |
| Construct empty object.
|
|
const std::function< void(T *, const T *, const T *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> & | get_tabulate_tensor (IntegralType type, int i) const |
| Get the function for 'tabulate_tensor' for integral i of given type. More...
|
|
void | set_tabulate_tensor (IntegralType type, int i, std::function< void(T *, const T *, const T *, 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...
|
|
std::set< IntegralType > | types () const |
| Get types of integrals in the form. More...
|
|
int | num_integrals (IntegralType type) const |
| Number of integrals of given type. More...
|
|
std::vector< int > | integral_ids (IntegralType type) 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 (IntegralType type, 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 (IntegralType 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...
|
|
template<typename T>
class dolfinx::fem::FormIntegrals< T >
Integrals of a Form, including those defined over cells, interior and exterior facets, and vertices.
◆ get_tabulate_tensor()
template<typename T >
const std::function<void(T*, const T*, const T*, const double*, const int*, const std::uint8_t*, const std::uint32_t)>& dolfinx::fem::FormIntegrals< T >::get_tabulate_tensor |
( |
IntegralType |
type, |
|
|
int |
i |
|
) |
| const |
|
inline |
Get the function for 'tabulate_tensor' for integral i of given type.
- Parameters
-
[in] | type | Integral type |
[in] | i | Integral number |
- Returns
- Function to call for tabulate_tensor
◆ integral_domains()
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] | type | Integral type |
[in] | i | Integral number |
- Returns
- List of active entities for this integral
◆ integral_ids()
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
-
- Returns
- List of IDs for this integral
◆ num_integrals()
Number of integrals of given type.
- Parameters
-
- Returns
- Number of integrals
◆ set_default_domains()
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
-
◆ set_domains()
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] | type | Integral type |
[in] | marker | MeshTags mapping entities to integrals |
◆ set_tabulate_tensor()
template<typename T >
void dolfinx::fem::FormIntegrals< T >::set_tabulate_tensor |
( |
IntegralType |
type, |
|
|
int |
i, |
|
|
std::function< void(T *, const T *, const T *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> |
fn |
|
) |
| |
|
inline |
Set the function for 'tabulate_tensor' for integral i of given type.
- Parameters
-
[in] | type | Integral type |
[in] | i | Integral number |
[in] | fn | tabulate function |
◆ types()
Get types of integrals in the form.
- Returns
- Integrals types
The documentation for this class was generated from the following file:
- /build/dolfinx-g9rZs9/dolfinx-2019.2.0~git20200723.696fbc0/cpp/dolfinx/fem/FormIntegrals.h