Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d9/dd0/assemble__scalar__impl_8h_source.html
DOLFINx 0.8.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{
23
25template <dolfinx::scalar T>
26T assemble_cells(mdspan2_t x_dofmap, std::span<const scalar_value_type_t<T>> x,
27 std::span<const std::int32_t> cells, FEkernel<T> auto fn,
28 std::span<const T> constants, std::span<const T> coeffs,
29 int cstride)
30{
31 T value(0);
32 if (cells.empty())
33 return value;
34
35 // Create data structures used in assembly
36 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
37
38 // Iterate over all cells
39 for (std::size_t index = 0; index < cells.size(); ++index)
40 {
41 std::int32_t c = cells[index];
42
43 // Get cell coordinates/geometry
44 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
45 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
46 for (std::size_t i = 0; i < x_dofs.size(); ++i)
47 {
48 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
49 std::next(coordinate_dofs.begin(), 3 * i));
50 }
51
52 const T* coeff_cell = coeffs.data() + index * cstride;
53 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(), nullptr,
54 nullptr);
55 }
56
57 return value;
58}
59
61template <dolfinx::scalar T>
62T assemble_exterior_facets(mdspan2_t x_dofmap,
63 std::span<const scalar_value_type_t<T>> x,
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{
68 T value(0);
69 if (facets.empty())
70 return value;
71
72 // Create data structures used in assembly
73 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
74
75 // Iterate over all facets
76 assert(facets.size() % 2 == 0);
77 for (std::size_t index = 0; index < facets.size(); index += 2)
78 {
79 std::int32_t cell = facets[index];
80 std::int32_t local_facet = facets[index + 1];
81
82 // Get cell coordinates/geometry
83 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
84 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
85 for (std::size_t i = 0; i < x_dofs.size(); ++i)
86 {
87 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
88 std::next(coordinate_dofs.begin(), 3 * i));
89 }
90
91 const T* coeff_cell = coeffs.data() + index / 2 * cstride;
92 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
93 &local_facet, nullptr);
94 }
95
96 return value;
97}
98
100template <dolfinx::scalar T>
101T assemble_interior_facets(mdspan2_t x_dofmap,
102 std::span<const scalar_value_type_t<T>> x,
103 int num_cell_facets,
104 std::span<const std::int32_t> facets,
105 FEkernel<T> auto fn, std::span<const T> constants,
106 std::span<const T> coeffs, int cstride,
107 std::span<const int> offsets,
108 std::span<const std::uint8_t> perms)
109{
110 T value(0);
111 if (facets.empty())
112 return value;
113
114 // Create data structures used in assembly
115 using X = scalar_value_type_t<T>;
116 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
117 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
118 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
119 x_dofmap.extent(1) * 3);
120
121 std::vector<T> coeff_array(2 * offsets.back());
122 assert(offsets.back() == cstride);
123
124 // Iterate over all facets
125 assert(facets.size() % 4 == 0);
126 for (std::size_t index = 0; index < facets.size(); index += 4)
127 {
128 std::array<std::int32_t, 2> cells = {facets[index], facets[index + 2]};
129 std::array<std::int32_t, 2> local_facet
130 = {facets[index + 1], facets[index + 3]};
131
132 // Get cell geometry
133 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
134 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
135 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
136 {
137 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
138 std::next(cdofs0.begin(), 3 * i));
139 }
140 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
141 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
142 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
143 {
144 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
145 std::next(cdofs1.begin(), 3 * i));
146 }
147
148 const std::array perm{perms[cells[0] * num_cell_facets + local_facet[0]],
149 perms[cells[1] * num_cell_facets + local_facet[1]]};
150 fn(&value, coeffs.data() + index / 2 * cstride, constants.data(),
151 coordinate_dofs.data(), local_facet.data(), perm.data());
152 }
153
154 return value;
155}
156
158template <dolfinx::scalar T, std::floating_point U>
159T assemble_scalar(
160 const fem::Form<T, U>& M, mdspan2_t x_dofmap,
161 std::span<const scalar_value_type_t<T>> x, std::span<const T> constants,
162 const std::map<std::pair<IntegralType, int>,
163 std::pair<std::span<const T>, int>>& coefficients)
164{
165 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
166 assert(mesh);
167
168 T value = 0;
169 for (int i : M.integral_ids(IntegralType::cell))
170 {
171 auto fn = M.kernel(IntegralType::cell, i);
172 assert(fn);
173 auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
174 std::span<const std::int32_t> cells = M.domain(IntegralType::cell, i);
175 value += impl::assemble_cells(x_dofmap, x, cells, fn, constants, coeffs,
176 cstride);
177 }
178
179 for (int i : M.integral_ids(IntegralType::exterior_facet))
180 {
181 auto fn = M.kernel(IntegralType::exterior_facet, i);
182 assert(fn);
183 auto& [coeffs, cstride]
184 = coefficients.at({IntegralType::exterior_facet, i});
185 value += impl::assemble_exterior_facets(
186 x_dofmap, x, M.domain(IntegralType::exterior_facet, i), fn, constants,
187 coeffs, cstride);
188 }
189
190 if (M.num_integrals(IntegralType::interior_facet) > 0)
191 {
192 mesh->topology_mutable()->create_entity_permutations();
193 const std::vector<std::uint8_t>& perms
194 = mesh->topology()->get_facet_permutations();
195 mesh::CellType cell_type = mesh->topology()->cell_type();
196 int num_cell_facets
197 = mesh::cell_num_entities(cell_type, mesh->topology()->dim() - 1);
198 const std::vector<int> c_offsets = M.coefficient_offsets();
199 for (int i : M.integral_ids(IntegralType::interior_facet))
200 {
201 auto fn = M.kernel(IntegralType::interior_facet, i);
202 assert(fn);
203 auto& [coeffs, cstride]
204 = coefficients.at({IntegralType::interior_facet, i});
205 value += impl::assemble_interior_facets(
206 x_dofmap, x, num_cell_facets,
207 M.domain(IntegralType::interior_facet, i), fn, constants, coeffs,
208 cstride, c_offsets, perms);
209 }
210 }
211
212 return value;
213}
214
215} // namespace dolfinx::fem::impl
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:15
IntegralType
Type of integral.
Definition Form.h:34
@ 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