DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
assemble_scalar_impl.h
1// Copyright (C) 2019-2020 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 "Constant.h"
10#include "Form.h"
11#include "FunctionSpace.h"
12#include "utils.h"
13#include <algorithm>
14#include <dolfinx/common/IndexMap.h>
15#include <dolfinx/mesh/Geometry.h>
16#include <dolfinx/mesh/Mesh.h>
17#include <dolfinx/mesh/Topology.h>
18#include <memory>
19#include <vector>
20
21namespace dolfinx::fem::impl
22{
24template <dolfinx::scalar T>
25T assemble_cells(mdspan2_t x_dofmap, std::span<const scalar_value_type_t<T>> x,
26 std::span<const std::int32_t> cells, FEkernel<T> auto fn,
27 std::span<const T> constants, std::span<const T> coeffs,
28 int cstride)
29{
30 T value(0);
31 if (cells.empty())
32 return value;
33
34 // Create data structures used in assembly
35 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
36
37 // Iterate over all cells
38 for (std::size_t index = 0; index < cells.size(); ++index)
39 {
40 std::int32_t c = cells[index];
41
42 // Get cell coordinates/geometry
43 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
44 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
45 for (std::size_t i = 0; i < x_dofs.size(); ++i)
46 {
47 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
48 std::next(coordinate_dofs.begin(), 3 * i));
49 }
50
51 const T* coeff_cell = coeffs.data() + index * cstride;
52 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(), nullptr,
53 nullptr);
54 }
55
56 return value;
57}
58
60template <dolfinx::scalar T>
61T assemble_exterior_facets(mdspan2_t x_dofmap,
62 std::span<const scalar_value_type_t<T>> x,
63 int num_facets_per_cell,
64 std::span<const std::int32_t> facets,
65 FEkernel<T> auto fn, std::span<const T> constants,
66 std::span<const T> coeffs, int cstride,
67 std::span<const std::uint8_t> perms)
68{
69 T value(0);
70 if (facets.empty())
71 return value;
72
73 // Create data structures used in assembly
74 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
75
76 // Iterate over all facets
77 assert(facets.size() % 2 == 0);
78 for (std::size_t index = 0; index < facets.size(); index += 2)
79 {
80 std::int32_t cell = facets[index];
81 std::int32_t local_facet = facets[index + 1];
82
83 // Get cell coordinates/geometry
84 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
85 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
86 for (std::size_t i = 0; i < x_dofs.size(); ++i)
87 {
88 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
89 std::next(coordinate_dofs.begin(), 3 * i));
90 }
91
92 // Permutations
93 std::uint8_t perm
94 = perms.empty() ? 0 : perms[cell * num_facets_per_cell + local_facet];
95 const T* coeff_cell = coeffs.data() + index / 2 * cstride;
96 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
97 &local_facet, &perm);
98 }
99
100 return value;
101}
102
104template <dolfinx::scalar T>
105T assemble_interior_facets(mdspan2_t x_dofmap,
106 std::span<const scalar_value_type_t<T>> x,
107 int num_facets_per_cell,
108 std::span<const std::int32_t> facets,
109 FEkernel<T> auto fn, std::span<const T> constants,
110 std::span<const T> coeffs, int cstride,
111 std::span<const int> offsets,
112 std::span<const std::uint8_t> perms)
113{
114 T value(0);
115 if (facets.empty())
116 return value;
117
118 // Create data structures used in assembly
119 using X = scalar_value_type_t<T>;
120 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
121 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
122 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
123 x_dofmap.extent(1) * 3);
124
125 std::vector<T> coeff_array(2 * offsets.back());
126 assert(offsets.back() == cstride);
127
128 // Iterate over all facets
129 assert(facets.size() % 4 == 0);
130 for (std::size_t index = 0; index < facets.size(); index += 4)
131 {
132 std::array<std::int32_t, 2> cells = {facets[index], facets[index + 2]};
133 std::array<std::int32_t, 2> local_facet
134 = {facets[index + 1], facets[index + 3]};
135
136 // Get cell geometry
137 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
138 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
139 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
140 {
141 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
142 std::next(cdofs0.begin(), 3 * i));
143 }
144 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
145 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
146 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
147 {
148 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
149 std::next(cdofs1.begin(), 3 * i));
150 }
151
152 std::array perm
153 = perms.empty()
154 ? std::array<std::uint8_t, 2>{0, 0}
155 : std::array{
156 perms[cells[0] * num_facets_per_cell + local_facet[0]],
157 perms[cells[1] * num_facets_per_cell + local_facet[1]]};
158 fn(&value, coeffs.data() + index / 2 * cstride, constants.data(),
159 coordinate_dofs.data(), local_facet.data(), perm.data());
160 }
161
162 return value;
163}
164
166template <dolfinx::scalar T, std::floating_point U>
167T assemble_scalar(
168 const fem::Form<T, U>& M, mdspan2_t x_dofmap,
169 std::span<const scalar_value_type_t<T>> x, std::span<const T> constants,
170 const std::map<std::pair<IntegralType, int>,
171 std::pair<std::span<const T>, int>>& coefficients)
172{
173 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
174 assert(mesh);
175
176 T value = 0;
177 for (int i : M.integral_ids(IntegralType::cell))
178 {
179 auto fn = M.kernel(IntegralType::cell, i);
180 assert(fn);
181 auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
182 std::span<const std::int32_t> cells = M.domain(IntegralType::cell, i);
183 value += impl::assemble_cells(x_dofmap, x, cells, fn, constants, coeffs,
184 cstride);
185 }
186
187 std::span<const std::uint8_t> perms;
188 if (M.needs_facet_permutations())
189 {
190 mesh->topology_mutable()->create_entity_permutations();
191 perms = std::span(mesh->topology()->get_facet_permutations());
192 }
193
194 mesh::CellType cell_type = mesh->topology()->cell_type();
195 int num_facets_per_cell
196 = mesh::cell_num_entities(cell_type, mesh->topology()->dim() - 1);
197 for (int i : M.integral_ids(IntegralType::exterior_facet))
198 {
199 auto fn = M.kernel(IntegralType::exterior_facet, i);
200 assert(fn);
201 auto& [coeffs, cstride]
202 = coefficients.at({IntegralType::exterior_facet, i});
203 value += impl::assemble_exterior_facets(
204 x_dofmap, x, num_facets_per_cell,
205 M.domain(IntegralType::exterior_facet, i), fn, constants, coeffs,
206 cstride, perms);
207 }
208
209 for (int i : M.integral_ids(IntegralType::interior_facet))
210 {
211 const std::vector<int> c_offsets = M.coefficient_offsets();
212 auto fn = M.kernel(IntegralType::interior_facet, i);
213 assert(fn);
214 auto& [coeffs, cstride]
215 = coefficients.at({IntegralType::interior_facet, i});
216 value += impl::assemble_interior_facets(
217 x_dofmap, x, num_facets_per_cell,
218 M.domain(IntegralType::interior_facet, i), fn, constants, coeffs,
219 cstride, c_offsets, perms);
220 }
221
222 return value;
223}
224
225} // namespace dolfinx::fem::impl
A representation of finite element variational forms.
Definition Form.h:139
Functions supporting finite element method operations.
void cells(la::SparsityPattern &pattern, std::array< std::span< const std::int32_t >, 2 > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:16
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
int cell_num_entities(CellType type, int dim)
Definition cell_types.cpp:139
CellType
Cell type identifier.
Definition cell_types.h:22