DOLFIN-X
DOLFIN-X C++ interface
|
48 template <
typename ScalarType>
50 const std::function<
int(std::int32_t,
const std::int32_t*, std::int32_t,
51 const std::int32_t*,
const ScalarType*)>&
53 const Form& a,
const std::vector<bool>& bc0,
const std::vector<bool>& bc1);
56 template <
typename ScalarType>
58 const std::function<
int(std::int32_t,
const std::int32_t*, std::int32_t,
59 const std::int32_t*,
const ScalarType*)>&
61 const mesh::Mesh& mesh,
const std::vector<std::int32_t>& active_cells,
64 const std::vector<bool>& bc0,
const std::vector<bool>& bc1,
65 const std::function<
void(ScalarType*,
const ScalarType*,
const ScalarType*,
66 const double*,
const int*,
const std::uint8_t*,
67 const std::uint32_t)>& kernel,
68 const Eigen::Array<ScalarType, Eigen::Dynamic, Eigen::Dynamic,
69 Eigen::RowMajor>& coeffs,
70 const Eigen::Array<ScalarType, Eigen::Dynamic, 1>& constant_values);
73 template <
typename ScalarType>
75 const std::function<
int(std::int32_t,
const std::int32_t*, std::int32_t,
76 const std::int32_t*,
const ScalarType*)>&
78 const mesh::Mesh& mesh,
const std::vector<std::int32_t>& active_facets,
79 const DofMap& dofmap0,
const DofMap& dofmap1,
const std::vector<bool>& bc0,
80 const std::vector<bool>& bc1,
81 const std::function<
void(ScalarType*,
const ScalarType*,
const ScalarType*,
82 const double*,
const int*,
const std::uint8_t*,
83 const std::uint32_t)>& fn,
84 const Eigen::Array<ScalarType, Eigen::Dynamic, Eigen::Dynamic,
85 Eigen::RowMajor>& coeffs,
86 const Eigen::Array<ScalarType, Eigen::Dynamic, 1> constant_values);
89 template <
typename ScalarType>
91 const std::function<
int(std::int32_t,
const std::int32_t*, std::int32_t,
92 const std::int32_t*,
const ScalarType*)>&
94 const mesh::Mesh& mesh,
const std::vector<std::int32_t>& active_facets,
95 const DofMap& dofmap0,
const DofMap& dofmap1,
const std::vector<bool>& bc0,
96 const std::vector<bool>& bc1,
97 const std::function<
void(ScalarType*,
const ScalarType*,
const ScalarType*,
98 const double*,
const int*,
const std::uint8_t*,
99 const std::uint32_t)>& kernel,
100 const Eigen::Array<ScalarType, Eigen::Dynamic, Eigen::Dynamic,
101 Eigen::RowMajor>& coeffs,
102 const std::vector<int>& offsets,
103 const Eigen::Array<ScalarType, Eigen::Dynamic, 1>& constant_values);
This class provides a static adjacency list data structure. It is commonly used to store directed gra...
Definition: assemble_matrix_impl.h:26
Degree-of-freedom map.
Definition: DofMap.h:40
void assemble_matrix(const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const ScalarType *)> &mat_set_values_local, const Form &a, const std::vector< bool > &bc0, const std::vector< bool > &bc1)
The matrix A must already be initialised. The matrix may be a proxy, i.e. a view into a larger matrix...
Definition: assemble_matrix_impl.cpp:24
void assemble_interior_facets(const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const ScalarType *)> &mat_set_values_local, const mesh::Mesh &mesh, const std::vector< std::int32_t > &active_facets, const DofMap &dofmap0, const DofMap &dofmap1, const std::vector< bool > &bc0, const std::vector< bool > &bc1, const std::function< void(ScalarType *, const ScalarType *, const ScalarType *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> &kernel, const Eigen::Array< ScalarType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &coeffs, const std::vector< int > &offsets, const Eigen::Array< ScalarType, Eigen::Dynamic, 1 > &constant_values)
Execute kernel over interior facets and accumulate result in Mat.
Definition: assemble_matrix_impl.cpp:274
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:46
void assemble_cells(const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const ScalarType *)> &mat_set_values_local, const mesh::Mesh &mesh, const std::vector< std::int32_t > &active_cells, const graph::AdjacencyList< std::int32_t > &dofmap0, int num_dofs_per_cell0, const graph::AdjacencyList< std::int32_t > &dofmap1, int num_dofs_per_cell1, const std::vector< bool > &bc0, const std::vector< bool > &bc1, const std::function< void(ScalarType *, const ScalarType *, const ScalarType *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> &kernel, const Eigen::Array< ScalarType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &coeffs, const Eigen::Array< ScalarType, Eigen::Dynamic, 1 > &constant_values)
Execute kernel over cells and accumulate result in Mat.
Definition: assemble_matrix_impl.cpp:101
This class represents a function in a finite element function space , given by.
Definition: Function.h:41
void assemble_exterior_facets(const std::function< int(std::int32_t, const std::int32_t *, std::int32_t, const std::int32_t *, const ScalarType *)> &mat_set_values_local, const mesh::Mesh &mesh, const std::vector< std::int32_t > &active_facets, const DofMap &dofmap0, const DofMap &dofmap1, const std::vector< bool > &bc0, const std::vector< bool > &bc1, const std::function< void(ScalarType *, const ScalarType *, const ScalarType *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> &fn, const Eigen::Array< ScalarType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &coeffs, const Eigen::Array< ScalarType, Eigen::Dynamic, 1 > constant_values)
Execute kernel over exterior facets and accumulate result in Mat.
Definition: assemble_matrix_impl.cpp:178