DOLFIN-X
DOLFIN-X C++ interface
assemble_scalar_impl.h
1 // Copyright (C) 2019 Garth N. Wells
2 //
3 // This file is part of DOLFINX (https://www.fenicsproject.org)
4 //
5 // SPDX-License-Identifier: LGPL-3.0-or-later
6 
7 #pragma once
8 
9 #include <Eigen/Dense>
10 #include <dolfinx/common/types.h>
11 #include <memory>
12 #include <petscsys.h>
13 #include <vector>
14 
15 namespace dolfinx
16 {
17 
18 namespace function
19 {
20 class Function;
21 }
22 
23 namespace mesh
24 {
25 class Mesh;
26 }
27 
28 namespace fem
29 {
30 class DirichletBC;
31 class Form;
32 class DofMap;
33 
34 namespace impl
35 {
36 
38 PetscScalar assemble_scalar(const fem::Form& M);
39 
41 PetscScalar assemble_cells(
42  const mesh::Mesh& mesh, const std::vector<std::int32_t>& active_cells,
43  const std::function<void(PetscScalar*, const PetscScalar*,
44  const PetscScalar*, const double*, const int*,
45  const std::uint8_t*, const std::uint32_t)>& fn,
46  const Eigen::Array<PetscScalar, Eigen::Dynamic, Eigen::Dynamic,
47  Eigen::RowMajor>& coeffs,
48  const std::vector<PetscScalar>& constant_values);
49 
51 PetscScalar assemble_exterior_facets(
52  const mesh::Mesh& mesh, const std::vector<std::int32_t>& active_cells,
53  const std::function<void(PetscScalar*, const PetscScalar*,
54  const PetscScalar*, const double*, const int*,
55  const std::uint8_t*, const std::uint32_t)>& fn,
56  const Eigen::Array<PetscScalar, Eigen::Dynamic, Eigen::Dynamic,
57  Eigen::RowMajor>& coeffs,
58  const std::vector<PetscScalar>& constant_values);
59 
61 PetscScalar assemble_interior_facets(
62  const mesh::Mesh& mesh, const std::vector<std::int32_t>& active_cells,
63  const std::function<void(PetscScalar*, const PetscScalar*,
64  const PetscScalar*, const double*, const int*,
65  const std::uint8_t*, const std::uint32_t)>& fn,
66  const Eigen::Array<PetscScalar, Eigen::Dynamic, Eigen::Dynamic,
67  Eigen::RowMajor>& coeffs,
68  const std::vector<int>& offsets,
69  const std::vector<PetscScalar>& constant_values);
70 
71 } // namespace impl
72 } // namespace fem
73 } // namespace dolfinx
dolfinx::fem::impl::assemble_scalar
PetscScalar assemble_scalar(const fem::Form &M)
Assemble functional into an scalar.
Definition: assemble_scalar_impl.cpp:25
dolfinx::fem::impl::assemble_interior_facets
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
dolfinx::fem::impl::assemble_cells
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
dolfinx::fem::impl::assemble_exterior_facets
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