Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d9/d40/assemble__vector__impl_8h_source.html
DOLFINx 0.6.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
assemble_vector_impl.h
1// Copyright (C) 2018-2021 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 "DirichletBC.h"
11#include "DofMap.h"
12#include "Form.h"
13#include "FunctionSpace.h"
14#include "utils.h"
15#include <algorithm>
16#include <dolfinx/common/IndexMap.h>
17#include <dolfinx/graph/AdjacencyList.h>
18#include <dolfinx/mesh/Geometry.h>
19#include <dolfinx/mesh/Mesh.h>
20#include <dolfinx/mesh/Topology.h>
21#include <functional>
22#include <memory>
23#include <span>
24#include <vector>
25
26namespace dolfinx::fem::impl
27{
28
30
38template <typename T, int _bs0 = -1, int _bs1 = -1>
39void _lift_bc_cells(
40 std::span<T> b, const mesh::Geometry& geometry, FEkernel<T> auto kernel,
41 std::span<const std::int32_t> cells,
42 const std::function<void(const std::span<T>&,
43 const std::span<const std::uint32_t>&,
44 std::int32_t, int)>& dof_transform,
45 const graph::AdjacencyList<std::int32_t>& dofmap0, int bs0,
46 const std::function<void(const std::span<T>&,
47 const std::span<const std::uint32_t>&,
48 std::int32_t, int)>& dof_transform_to_transpose,
49 const graph::AdjacencyList<std::int32_t>& dofmap1, int bs1,
50 std::span<const T> constants, std::span<const T> coeffs, int cstride,
51 std::span<const std::uint32_t> cell_info, std::span<const T> bc_values1,
52 std::span<const std::int8_t> bc_markers1, std::span<const T> x0,
53 double scale)
54{
55 assert(_bs0 < 0 or _bs0 == bs0);
56 assert(_bs1 < 0 or _bs1 == bs1);
57
58 if (cells.empty())
59 return;
60
61 // Prepare cell geometry
62 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
63 const std::size_t num_dofs_g = geometry.cmap().dim();
64 std::span<const double> x = geometry.x();
65
66 // Data structures used in bc application
67 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * num_dofs_g);
68 std::vector<T> Ae, be;
69 const scalar_value_type_t<T> _scale
70 = static_cast<scalar_value_type_t<T>>(scale);
71 for (std::size_t index = 0; index < cells.size(); ++index)
72 {
73 std::int32_t c = cells[index];
74
75 // Get dof maps for cell
76 auto dmap1 = dofmap1.links(c);
77
78 // Check if bc is applied to cell
79 bool has_bc = false;
80 for (std::size_t j = 0; j < dmap1.size(); ++j)
81 {
82 if constexpr (_bs1 > 0)
83 {
84 for (int k = 0; k < _bs1; ++k)
85 {
86 assert(_bs1 * dmap1[j] + k < (int)bc_markers1.size());
87 if (bc_markers1[_bs1 * dmap1[j] + k])
88 {
89 has_bc = true;
90 break;
91 }
92 }
93 }
94 else
95 {
96 for (int k = 0; k < bs1; ++k)
97 {
98 assert(bs1 * dmap1[j] + k < (int)bc_markers1.size());
99 if (bc_markers1[bs1 * dmap1[j] + k])
100 {
101 has_bc = true;
102 break;
103 }
104 }
105 }
106 }
107
108 if (!has_bc)
109 continue;
110
111 // Get cell coordinates/geometry
112 auto x_dofs = x_dofmap.links(c);
113 for (std::size_t i = 0; i < x_dofs.size(); ++i)
114 {
115 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
116 std::next(coordinate_dofs.begin(), 3 * i));
117 }
118
119 // Size data structure for assembly
120 auto dmap0 = dofmap0.links(c);
121
122 const int num_rows = bs0 * dmap0.size();
123 const int num_cols = bs1 * dmap1.size();
124
125 const T* coeff_array = coeffs.data() + index * cstride;
126 Ae.resize(num_rows * num_cols);
127 std::fill(Ae.begin(), Ae.end(), 0);
128 kernel(Ae.data(), coeff_array, constants.data(), coordinate_dofs.data(),
129 nullptr, nullptr);
130 dof_transform(Ae, cell_info, c, num_cols);
131 dof_transform_to_transpose(Ae, cell_info, c, num_rows);
132
133 // Size data structure for assembly
134 be.resize(num_rows);
135 std::fill(be.begin(), be.end(), 0);
136 for (std::size_t j = 0; j < dmap1.size(); ++j)
137 {
138 if constexpr (_bs1 > 0)
139 {
140 for (int k = 0; k < _bs1; ++k)
141 {
142 const std::int32_t jj = _bs1 * dmap1[j] + k;
143 assert(jj < (int)bc_markers1.size());
144 if (bc_markers1[jj])
145 {
146 const T bc = bc_values1[jj];
147 const T _x0 = x0.empty() ? 0.0 : x0[jj];
148 // const T _x0 = 0.0;
149 // be -= Ae.col(bs1 * j + k) * scale * (bc - _x0);
150 for (int m = 0; m < num_rows; ++m)
151 be[m] -= Ae[m * num_cols + _bs1 * j + k] * _scale * (bc - _x0);
152 }
153 }
154 }
155 else
156 {
157 for (int k = 0; k < bs1; ++k)
158 {
159 const std::int32_t jj = bs1 * dmap1[j] + k;
160 assert(jj < (int)bc_markers1.size());
161 if (bc_markers1[jj])
162 {
163 const T bc = bc_values1[jj];
164 const T _x0 = x0.empty() ? 0.0 : x0[jj];
165 // be -= Ae.col(bs1 * j + k) * scale * (bc - _x0);
166 for (int m = 0; m < num_rows; ++m)
167 be[m] -= Ae[m * num_cols + bs1 * j + k] * _scale * (bc - _x0);
168 }
169 }
170 }
171 }
172
173 for (std::size_t i = 0; i < dmap0.size(); ++i)
174 {
175 if constexpr (_bs0 > 0)
176 {
177 for (int k = 0; k < _bs0; ++k)
178 b[_bs0 * dmap0[i] + k] += be[_bs0 * i + k];
179 }
180 else
181 {
182 for (int k = 0; k < bs0; ++k)
183 b[bs0 * dmap0[i] + k] += be[bs0 * i + k];
184 }
185 }
186 }
187}
188
195template <typename T, int _bs = -1>
196void _lift_bc_exterior_facets(
197 std::span<T> b, const mesh::Geometry& geometry, FEkernel<T> auto kernel,
198 std::span<const std::int32_t> facets,
199 const std::function<void(const std::span<T>&,
200 const std::span<const std::uint32_t>&,
201 std::int32_t, int)>& dof_transform,
202 const graph::AdjacencyList<std::int32_t>& dofmap0, int bs0,
203 const std::function<void(const std::span<T>&,
204 const std::span<const std::uint32_t>&,
205 std::int32_t, int)>& dof_transform_to_transpose,
206 const graph::AdjacencyList<std::int32_t>& dofmap1, int bs1,
207 std::span<const T> constants, std::span<const T> coeffs, int cstride,
208 std::span<const std::uint32_t> cell_info, std::span<const T> bc_values1,
209 std::span<const std::int8_t> bc_markers1, std::span<const T> x0,
210 double scale)
211{
212 if (facets.empty())
213 return;
214
215 // Prepare cell geometry
216 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
217 const std::size_t num_dofs_g = geometry.cmap().dim();
218 std::span<const double> x = geometry.x();
219
220 // Data structures used in bc application
221 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * num_dofs_g);
222 std::vector<T> Ae, be;
223 assert(facets.size() % 2 == 0);
224 const scalar_value_type_t<T> _scale
225 = static_cast<scalar_value_type_t<T>>(scale);
226 for (std::size_t index = 0; index < facets.size(); index += 2)
227 {
228 std::int32_t cell = facets[index];
229 std::int32_t local_facet = facets[index + 1];
230
231 // Get dof maps for cell
232 auto dmap1 = dofmap1.links(cell);
233
234 // Check if bc is applied to cell
235 bool has_bc = false;
236 for (std::size_t j = 0; j < dmap1.size(); ++j)
237 {
238 for (int k = 0; k < bs1; ++k)
239 {
240 if (bc_markers1[bs1 * dmap1[j] + k])
241 {
242 has_bc = true;
243 break;
244 }
245 }
246 }
247
248 if (!has_bc)
249 continue;
250
251 // Get cell coordinates/geometry
252 auto x_dofs = x_dofmap.links(cell);
253 for (std::size_t i = 0; i < x_dofs.size(); ++i)
254 {
255 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
256 std::next(coordinate_dofs.begin(), 3 * i));
257 }
258
259 // Size data structure for assembly
260 auto dmap0 = dofmap0.links(cell);
261 const int num_rows = bs0 * dmap0.size();
262 const int num_cols = bs1 * dmap1.size();
263
264 const T* coeff_array = coeffs.data() + index / 2 * cstride;
265 Ae.resize(num_rows * num_cols);
266 std::fill(Ae.begin(), Ae.end(), 0);
267 kernel(Ae.data(), coeff_array, constants.data(), coordinate_dofs.data(),
268 &local_facet, nullptr);
269 dof_transform(Ae, cell_info, cell, num_cols);
270 dof_transform_to_transpose(Ae, cell_info, cell, num_rows);
271
272 // Size data structure for assembly
273 be.resize(num_rows);
274 std::fill(be.begin(), be.end(), 0);
275 for (std::size_t j = 0; j < dmap1.size(); ++j)
276 {
277 for (int k = 0; k < bs1; ++k)
278 {
279 const std::int32_t jj = bs1 * dmap1[j] + k;
280 if (bc_markers1[jj])
281 {
282 const T bc = bc_values1[jj];
283 const T _x0 = x0.empty() ? 0.0 : x0[jj];
284 // be -= Ae.col(bs1 * j + k) * scale * (bc - _x0);
285 for (int m = 0; m < num_rows; ++m)
286 be[m] -= Ae[m * num_cols + bs1 * j + k] * _scale * (bc - _x0);
287 }
288 }
289 }
290
291 for (std::size_t i = 0; i < dmap0.size(); ++i)
292 for (int k = 0; k < bs0; ++k)
293 b[bs0 * dmap0[i] + k] += be[bs0 * i + k];
294 }
295}
296
303template <typename T, int _bs = -1>
304void _lift_bc_interior_facets(
305 std::span<T> b, const mesh::Geometry& geometry, int num_cell_facets,
306 FEkernel<T> auto kernel, std::span<const std::int32_t> facets,
307 const std::function<void(const std::span<T>&,
308 const std::span<const std::uint32_t>&,
309 std::int32_t, int)>& dof_transform,
310 const graph::AdjacencyList<std::int32_t>& dofmap0, int bs0,
311 const std::function<void(const std::span<T>&,
312 const std::span<const std::uint32_t>&,
313 std::int32_t, int)>& dof_transform_to_transpose,
314 const graph::AdjacencyList<std::int32_t>& dofmap1, int bs1,
315 std::span<const T> constants, std::span<const T> coeffs, int cstride,
316 std::span<const std::uint32_t> cell_info,
317 const std::function<std::uint8_t(std::size_t)>& get_perm,
318 std::span<const T> bc_values1, std::span<const std::int8_t> bc_markers1,
319 std::span<const T> x0, double scale)
320{
321 if (facets.empty())
322 return;
323
324 // Prepare cell geometry
325 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
326 const std::size_t num_dofs_g = geometry.cmap().dim();
327 std::span<const double> x = geometry.x();
328
329 // Data structures used in assembly
330 using X = scalar_value_type_t<T>;
331 std::vector<X> coordinate_dofs(2 * num_dofs_g * 3);
332 std::span<X> cdofs0(coordinate_dofs.data(), num_dofs_g * 3);
333 std::span<X> cdofs1(coordinate_dofs.data() + num_dofs_g * 3, num_dofs_g * 3);
334 std::vector<T> Ae, be;
335
336 // Temporaries for joint dofmaps
337 std::vector<std::int32_t> dmapjoint0, dmapjoint1;
338 assert(facets.size() % 4 == 0);
339
340 const scalar_value_type_t<T> _scale
341 = static_cast<scalar_value_type_t<T>>(scale);
342 for (std::size_t index = 0; index < facets.size(); index += 4)
343 {
344 std::array<std::int32_t, 2> cells = {facets[index], facets[index + 2]};
345 std::array<std::int32_t, 2> local_facet
346 = {facets[index + 1], facets[index + 3]};
347
348 // Get cell geometry
349 auto x_dofs0 = x_dofmap.links(cells[0]);
350 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
351 {
352 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
353 std::next(cdofs0.begin(), 3 * i));
354 }
355 auto x_dofs1 = x_dofmap.links(cells[1]);
356 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
357 {
358 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
359 std::next(cdofs1.begin(), 3 * i));
360 }
361
362 // Get dof maps for cells and pack
363 std::span<const std::int32_t> dmap0_cell0 = dofmap0.links(cells[0]);
364 std::span<const std::int32_t> dmap0_cell1 = dofmap0.links(cells[1]);
365 dmapjoint0.resize(dmap0_cell0.size() + dmap0_cell1.size());
366 std::copy(dmap0_cell0.begin(), dmap0_cell0.end(), dmapjoint0.begin());
367 std::copy(dmap0_cell1.begin(), dmap0_cell1.end(),
368 std::next(dmapjoint0.begin(), dmap0_cell0.size()));
369
370 std::span<const std::int32_t> dmap1_cell0 = dofmap1.links(cells[0]);
371 std::span<const std::int32_t> dmap1_cell1 = dofmap1.links(cells[1]);
372 dmapjoint1.resize(dmap1_cell0.size() + dmap1_cell1.size());
373 std::copy(dmap1_cell0.begin(), dmap1_cell0.end(), dmapjoint1.begin());
374 std::copy(dmap1_cell1.begin(), dmap1_cell1.end(),
375 std::next(dmapjoint1.begin(), dmap1_cell0.size()));
376
377 // Check if bc is applied to cell0
378 bool has_bc = false;
379 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
380 {
381 for (int k = 0; k < bs1; ++k)
382 {
383 if (bc_markers1[bs1 * dmap1_cell0[j] + k])
384 {
385 has_bc = true;
386 break;
387 }
388 }
389 }
390
391 // Check if bc is applied to cell1
392 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
393 {
394 for (int k = 0; k < bs1; ++k)
395 {
396 if (bc_markers1[bs1 * dmap1_cell1[j] + k])
397 {
398 has_bc = true;
399 break;
400 }
401 }
402 }
403
404 if (!has_bc)
405 continue;
406
407 const int num_rows = bs0 * dmapjoint0.size();
408 const int num_cols = bs1 * dmapjoint1.size();
409
410 // Tabulate tensor
411 Ae.resize(num_rows * num_cols);
412 std::fill(Ae.begin(), Ae.end(), 0);
413 const std::array perm{
414 get_perm(cells[0] * num_cell_facets + local_facet[0]),
415 get_perm(cells[1] * num_cell_facets + local_facet[1])};
416 kernel(Ae.data(), coeffs.data() + index / 2 * cstride, constants.data(),
417 coordinate_dofs.data(), local_facet.data(), perm.data());
418
419 std::span<T> _Ae(Ae);
420 std::span<T> sub_Ae0 = _Ae.subspan(bs0 * dmap0_cell0.size() * num_cols,
421 bs0 * dmap0_cell1.size() * num_cols);
422 std::span<T> sub_Ae1
423 = _Ae.subspan(bs1 * dmap1_cell0.size(),
424 num_rows * num_cols - bs1 * dmap1_cell0.size());
425
426 dof_transform(_Ae, cell_info, cells[0], num_cols);
427 dof_transform(sub_Ae0, cell_info, cells[1], num_cols);
428 dof_transform_to_transpose(_Ae, cell_info, cells[0], num_rows);
429 dof_transform_to_transpose(sub_Ae1, cell_info, cells[1], num_rows);
430
431 be.resize(num_rows);
432 std::fill(be.begin(), be.end(), 0);
433
434 // Compute b = b - A*b for cell0
435 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
436 {
437 for (int k = 0; k < bs1; ++k)
438 {
439 const std::int32_t jj = bs1 * dmap1_cell0[j] + k;
440 if (bc_markers1[jj])
441 {
442 const T bc = bc_values1[jj];
443 const T _x0 = x0.empty() ? 0.0 : x0[jj];
444 // be -= Ae.col(bs1 * j + k) * scale * (bc - _x0);
445 for (int m = 0; m < num_rows; ++m)
446 be[m] -= Ae[m * num_cols + bs1 * j + k] * _scale * (bc - _x0);
447 }
448 }
449 }
450
451 // Compute b = b - A*b for cell1
452 const int offset = bs1 * dmap1_cell0.size();
453 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
454 {
455 for (int k = 0; k < bs1; ++k)
456 {
457 const std::int32_t jj = bs1 * dmap1_cell1[j] + k;
458 if (bc_markers1[jj])
459 {
460 const T bc = bc_values1[jj];
461 const T _x0 = x0.empty() ? 0.0 : x0[jj];
462 // be -= Ae.col(offset + bs1 * j + k) * scale * (bc - x0[jj]);
463 for (int m = 0; m < num_rows; ++m)
464 {
465 be[m] -= Ae[m * num_cols + offset + bs1 * j + k] * _scale
466 * (bc - _x0);
467 }
468 }
469 }
470 }
471
472 for (std::size_t i = 0; i < dmap0_cell0.size(); ++i)
473 for (int k = 0; k < bs0; ++k)
474 b[bs0 * dmap0_cell0[i] + k] += be[bs0 * i + k];
475
476 const int offset_be = bs0 * dmap0_cell0.size();
477 for (std::size_t i = 0; i < dmap0_cell1.size(); ++i)
478 for (int k = 0; k < bs0; ++k)
479 b[bs0 * dmap0_cell1[i] + k] += be[offset_be + bs0 * i + k];
480 }
481}
488template <typename T, int _bs = -1>
489void assemble_cells(
490 const std::function<void(const std::span<T>&,
491 const std::span<const std::uint32_t>&,
492 std::int32_t, int)>& dof_transform,
493 std::span<T> b, const mesh::Geometry& geometry,
494 std::span<const std::int32_t> cells,
495 const graph::AdjacencyList<std::int32_t>& dofmap, int bs,
496 FEkernel<T> auto kernel, std::span<const T> constants,
497 std::span<const T> coeffs, int cstride,
498 std::span<const std::uint32_t> cell_info)
499{
500 assert(_bs < 0 or _bs == bs);
501
502 if (cells.empty())
503 return;
504
505 // Prepare cell geometry
506 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
507 const std::size_t num_dofs_g = geometry.cmap().dim();
508 std::span<const double> x = geometry.x();
509
510 // FIXME: Add proper interface for num_dofs
511 // Create data structures used in assembly
512 const int num_dofs = dofmap.links(0).size();
513 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * num_dofs_g);
514 std::vector<T> be(bs * num_dofs);
515 std::span<T> _be(be);
516
517 // Iterate over active cells
518 for (std::size_t index = 0; index < cells.size(); ++index)
519 {
520 std::int32_t c = cells[index];
521
522 // Get cell coordinates/geometry
523 auto x_dofs = x_dofmap.links(c);
524 for (std::size_t i = 0; i < x_dofs.size(); ++i)
525 {
526 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
527 std::next(coordinate_dofs.begin(), 3 * i));
528 }
529
530 // Tabulate vector for cell
531 std::fill(be.begin(), be.end(), 0);
532 kernel(be.data(), coeffs.data() + index * cstride, constants.data(),
533 coordinate_dofs.data(), nullptr, nullptr);
534 dof_transform(_be, cell_info, c, 1);
535
536 // Scatter cell vector to 'global' vector array
537 auto dofs = dofmap.links(c);
538 if constexpr (_bs > 0)
539 {
540 for (int i = 0; i < num_dofs; ++i)
541 for (int k = 0; k < _bs; ++k)
542 b[_bs * dofs[i] + k] += be[_bs * i + k];
543 }
544 else
545 {
546 for (int i = 0; i < num_dofs; ++i)
547 for (int k = 0; k < bs; ++k)
548 b[bs * dofs[i] + k] += be[bs * i + k];
549 }
550 }
551}
552
559template <typename T, int _bs = -1>
560void assemble_exterior_facets(
561 const std::function<void(const std::span<T>&,
562 const std::span<const std::uint32_t>&,
563 std::int32_t, int)>& dof_transform,
564 std::span<T> b, const mesh::Geometry& geometry,
565 std::span<const std::int32_t> facets,
566 const graph::AdjacencyList<std::int32_t>& dofmap, int bs,
567 FEkernel<T> auto fn, std::span<const T> constants,
568 std::span<const T> coeffs, int cstride,
569 std::span<const std::uint32_t> cell_info)
570{
571 assert(_bs < 0 or _bs == bs);
572
573 if (facets.empty())
574 return;
575
576 // Prepare cell geometry
577 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
578 const std::size_t num_dofs_g = geometry.cmap().dim();
579 std::span<const double> x = geometry.x();
580
581 // FIXME: Add proper interface for num_dofs
582 // Create data structures used in assembly
583 const int num_dofs = dofmap.links(0).size();
584 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * num_dofs_g);
585 std::vector<T> be(bs * num_dofs);
586 std::span<T> _be(be);
587 assert(facets.size() % 2 == 0);
588 for (std::size_t index = 0; index < facets.size(); index += 2)
589 {
590 std::int32_t cell = facets[index];
591 std::int32_t local_facet = facets[index + 1];
592
593 // Get cell coordinates/geometry
594 auto x_dofs = x_dofmap.links(cell);
595 for (std::size_t i = 0; i < x_dofs.size(); ++i)
596 {
597 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
598 std::next(coordinate_dofs.begin(), 3 * i));
599 }
600
601 // Tabulate element vector
602 std::fill(be.begin(), be.end(), 0);
603 fn(be.data(), coeffs.data() + index / 2 * cstride, constants.data(),
604 coordinate_dofs.data(), &local_facet, nullptr);
605
606 dof_transform(_be, cell_info, cell, 1);
607
608 // Add element vector to global vector
609 auto dofs = dofmap.links(cell);
610 if constexpr (_bs > 0)
611 {
612 for (int i = 0; i < num_dofs; ++i)
613 for (int k = 0; k < _bs; ++k)
614 b[_bs * dofs[i] + k] += be[_bs * i + k];
615 }
616 else
617 {
618 for (int i = 0; i < num_dofs; ++i)
619 for (int k = 0; k < bs; ++k)
620 b[bs * dofs[i] + k] += be[bs * i + k];
621 }
622 }
623}
624
631template <typename T, int _bs = -1>
632void assemble_interior_facets(
633 const std::function<void(const std::span<T>&,
634 const std::span<const std::uint32_t>&,
635 std::int32_t, int)>& dof_transform,
636 std::span<T> b, const mesh::Geometry& geometry, int num_cell_facets,
637 std::span<const std::int32_t> facets, const fem::DofMap& dofmap,
638 FEkernel<T> auto fn, std::span<const T> constants,
639 std::span<const T> coeffs, int cstride,
640 std::span<const std::uint32_t> cell_info,
641 const std::function<std::uint8_t(std::size_t)>& get_perm)
642{
643 // Prepare cell geometry
644 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
645 const std::size_t num_dofs_g = geometry.cmap().dim();
646 std::span<const double> x = geometry.x();
647
648 // Create data structures used in assembly
649 using X = scalar_value_type_t<T>;
650 std::vector<X> coordinate_dofs(2 * num_dofs_g * 3);
651 std::span<X> cdofs0(coordinate_dofs.data(), num_dofs_g * 3);
652 std::span<X> cdofs1(coordinate_dofs.data() + num_dofs_g * 3, num_dofs_g * 3);
653 std::vector<T> be;
654
655 const int bs = dofmap.bs();
656 assert(_bs < 0 or _bs == bs);
657 assert(facets.size() % 4 == 0);
658 for (std::size_t index = 0; index < facets.size(); index += 4)
659 {
660 std::array<std::int32_t, 2> cells = {facets[index], facets[index + 2]};
661 std::array<std::int32_t, 2> local_facet
662 = {facets[index + 1], facets[index + 3]};
663
664 // Get cell geometry
665 auto x_dofs0 = x_dofmap.links(cells[0]);
666 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
667 {
668 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
669 std::next(cdofs0.begin(), 3 * i));
670 }
671 auto x_dofs1 = x_dofmap.links(cells[1]);
672 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
673 {
674 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
675 std::next(cdofs1.begin(), 3 * i));
676 }
677
678 // Get dofmaps for cells
679 std::span<const std::int32_t> dmap0 = dofmap.cell_dofs(cells[0]);
680 std::span<const std::int32_t> dmap1 = dofmap.cell_dofs(cells[1]);
681
682 // Tabulate element vector
683 be.resize(bs * (dmap0.size() + dmap1.size()));
684 std::fill(be.begin(), be.end(), 0);
685 const std::array perm{
686 get_perm(cells[0] * num_cell_facets + local_facet[0]),
687 get_perm(cells[1] * num_cell_facets + local_facet[1])};
688 fn(be.data(), coeffs.data() + index / 2 * cstride, constants.data(),
689 coordinate_dofs.data(), local_facet.data(), perm.data());
690
691 std::span<T> _be(be);
692 std::span<T> sub_be = _be.subspan(bs * dmap0.size(), bs * dmap1.size());
693
694 dof_transform(be, cell_info, cells[0], 1);
695 dof_transform(sub_be, cell_info, cells[1], 1);
696
697 // Add element vector to global vector
698 if constexpr (_bs > 0)
699 {
700 for (std::size_t i = 0; i < dmap0.size(); ++i)
701 for (int k = 0; k < _bs; ++k)
702 b[_bs * dmap0[i] + k] += be[_bs * i + k];
703 for (std::size_t i = 0; i < dmap1.size(); ++i)
704 for (int k = 0; k < _bs; ++k)
705 b[_bs * dmap1[i] + k] += be[_bs * (i + dmap0.size()) + k];
706 }
707 else
708 {
709 for (std::size_t i = 0; i < dmap0.size(); ++i)
710 for (int k = 0; k < bs; ++k)
711 b[bs * dmap0[i] + k] += be[bs * i + k];
712 for (std::size_t i = 0; i < dmap1.size(); ++i)
713 for (int k = 0; k < bs; ++k)
714 b[bs * dmap1[i] + k] += be[bs * (i + dmap0.size()) + k];
715 }
716 }
717}
718
733template <typename T>
734void lift_bc(std::span<T> b, const Form<T>& a, std::span<const T> constants,
735 const std::map<std::pair<IntegralType, int>,
736 std::pair<std::span<const T>, int>>& coefficients,
737 std::span<const T> bc_values1,
738 std::span<const std::int8_t> bc_markers1, std::span<const T> x0,
739 double scale)
740{
741 std::shared_ptr<const mesh::Mesh> mesh = a.mesh();
742 assert(mesh);
743
744 // Get dofmap for columns and rows of a
745 assert(a.function_spaces().at(0));
746 assert(a.function_spaces().at(1));
747 const graph::AdjacencyList<std::int32_t>& dofmap0
748 = a.function_spaces()[0]->dofmap()->list();
749 const int bs0 = a.function_spaces()[0]->dofmap()->bs();
750 std::shared_ptr<const fem::FiniteElement> element0
751 = a.function_spaces()[0]->element();
752 const graph::AdjacencyList<std::int32_t>& dofmap1
753 = a.function_spaces()[1]->dofmap()->list();
754 const int bs1 = a.function_spaces()[1]->dofmap()->bs();
755 std::shared_ptr<const fem::FiniteElement> element1
756 = a.function_spaces()[1]->element();
757
758 const bool needs_transformation_data
759 = element0->needs_dof_transformations()
760 or element1->needs_dof_transformations()
761 or a.needs_facet_permutations();
762
763 std::span<const std::uint32_t> cell_info;
764 if (needs_transformation_data)
765 {
766 mesh->topology_mutable().create_entity_permutations();
767 cell_info = std::span(mesh->topology().get_cell_permutation_info());
768 }
769
770 const std::function<void(const std::span<T>&,
771 const std::span<const std::uint32_t>&, std::int32_t,
772 int)>
773 dof_transform = element0->get_dof_transformation_function<T>();
774 const std::function<void(const std::span<T>&,
775 const std::span<const std::uint32_t>&, std::int32_t,
776 int)>
777 dof_transform_to_transpose
778 = element1->get_dof_transformation_to_transpose_function<T>();
779
780 for (int i : a.integral_ids(IntegralType::cell))
781 {
782 const auto& kernel = a.kernel(IntegralType::cell, i);
783 const auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
784 const std::vector<std::int32_t>& cells = a.cell_domains(i);
785 if (bs0 == 1 and bs1 == 1)
786 {
787 _lift_bc_cells<T, 1, 1>(b, mesh->geometry(), kernel, cells, dof_transform,
788 dofmap0, bs0, dof_transform_to_transpose, dofmap1,
789 bs1, constants, coeffs, cstride, cell_info,
790 bc_values1, bc_markers1, x0, scale);
791 }
792 else if (bs0 == 3 and bs1 == 3)
793 {
794 _lift_bc_cells<T, 3, 3>(b, mesh->geometry(), kernel, cells, dof_transform,
795 dofmap0, bs0, dof_transform_to_transpose, dofmap1,
796 bs1, constants, coeffs, cstride, cell_info,
797 bc_values1, bc_markers1, x0, scale);
798 }
799 else
800 {
801 _lift_bc_cells(b, mesh->geometry(), kernel, cells, dof_transform, dofmap0,
802 bs0, dof_transform_to_transpose, dofmap1, bs1, constants,
803 coeffs, cstride, cell_info, bc_values1, bc_markers1, x0,
804 scale);
805 }
806 }
807
808 for (int i : a.integral_ids(IntegralType::exterior_facet))
809 {
810 const auto& kernel = a.kernel(IntegralType::exterior_facet, i);
811 const auto& [coeffs, cstride]
812 = coefficients.at({IntegralType::exterior_facet, i});
813 const std::vector<std::int32_t>& facets = a.exterior_facet_domains(i);
814 _lift_bc_exterior_facets(b, mesh->geometry(), kernel, facets, dof_transform,
815 dofmap0, bs0, dof_transform_to_transpose, dofmap1,
816 bs1, constants, coeffs, cstride, cell_info,
817 bc_values1, bc_markers1, x0, scale);
818 }
819
820 if (a.num_integrals(IntegralType::interior_facet) > 0)
821 {
822
823 std::function<std::uint8_t(std::size_t)> get_perm;
824 if (a.needs_facet_permutations())
825 {
826 mesh->topology_mutable().create_entity_permutations();
827 const std::vector<std::uint8_t>& perms
828 = mesh->topology().get_facet_permutations();
829 get_perm = [&perms](std::size_t i) { return perms[i]; };
830 }
831 else
832 get_perm = [](std::size_t) { return 0; };
833
834 int num_cell_facets = mesh::cell_num_entities(mesh->topology().cell_type(),
835 mesh->topology().dim() - 1);
836 for (int i : a.integral_ids(IntegralType::interior_facet))
837 {
838 const auto& kernel = a.kernel(IntegralType::interior_facet, i);
839 const auto& [coeffs, cstride]
840 = coefficients.at({IntegralType::interior_facet, i});
841 const std::vector<std::int32_t>& facets = a.interior_facet_domains(i);
842 _lift_bc_interior_facets(b, mesh->geometry(), num_cell_facets, kernel,
843 facets, dof_transform, dofmap0, bs0,
844 dof_transform_to_transpose, dofmap1, bs1,
845 constants, coeffs, cstride, cell_info, get_perm,
846 bc_values1, bc_markers1, x0, scale);
847 }
848 }
849}
850
870template <typename T>
871void apply_lifting(
872 std::span<T> b, const std::vector<std::shared_ptr<const Form<T>>> a,
873 const std::vector<std::span<const T>>& constants,
874 const std::vector<std::map<std::pair<IntegralType, int>,
875 std::pair<std::span<const T>, int>>>& coeffs,
876 const std::vector<std::vector<std::shared_ptr<const DirichletBC<T>>>>& bcs1,
877 const std::vector<std::span<const T>>& x0, double scale)
878{
879 // FIXME: make changes to reactivate this check
880 if (!x0.empty() and x0.size() != a.size())
881 {
882 throw std::runtime_error(
883 "Mismatch in size between x0 and bilinear form in assembler.");
884 }
885
886 if (a.size() != bcs1.size())
887 {
888 throw std::runtime_error(
889 "Mismatch in size between a and bcs in assembler.");
890 }
891
892 for (std::size_t j = 0; j < a.size(); ++j)
893 {
894 std::vector<std::int8_t> bc_markers1;
895 std::vector<T> bc_values1;
896 if (a[j] and !bcs1[j].empty())
897 {
898 assert(a[j]->function_spaces().at(0));
899
900 auto V1 = a[j]->function_spaces()[1];
901 assert(V1);
902 auto map1 = V1->dofmap()->index_map;
903 const int bs1 = V1->dofmap()->index_map_bs();
904 assert(map1);
905 const int crange = bs1 * (map1->size_local() + map1->num_ghosts());
906 bc_markers1.assign(crange, false);
907 bc_values1.assign(crange, 0.0);
908 for (const std::shared_ptr<const DirichletBC<T>>& bc : bcs1[j])
909 {
910 bc->mark_dofs(bc_markers1);
911 bc->dof_values(bc_values1);
912 }
913
914 if (!x0.empty())
915 {
916 lift_bc<T>(b, *a[j], constants[j], coeffs[j], bc_values1, bc_markers1,
917 x0[j], scale);
918 }
919 else
920 {
921 lift_bc<T>(b, *a[j], constants[j], coeffs[j], bc_values1, bc_markers1,
922 std::span<const T>(), scale);
923 }
924 }
925 }
926}
927
934template <typename T>
935void assemble_vector(
936 std::span<T> b, const Form<T>& L, std::span<const T> constants,
937 const std::map<std::pair<IntegralType, int>,
938 std::pair<std::span<const T>, int>>& coefficients)
939{
940 std::shared_ptr<const mesh::Mesh> mesh = L.mesh();
941 assert(mesh);
942
943 // Get dofmap data
944 assert(L.function_spaces().at(0));
945 std::shared_ptr<const fem::FiniteElement> element
946 = L.function_spaces().at(0)->element();
947 std::shared_ptr<const fem::DofMap> dofmap
948 = L.function_spaces().at(0)->dofmap();
949 assert(dofmap);
950 const graph::AdjacencyList<std::int32_t>& dofs = dofmap->list();
951 const int bs = dofmap->bs();
952
953 const std::function<void(const std::span<T>&,
954 const std::span<const std::uint32_t>&, std::int32_t,
955 int)>
956 dof_transform = element->get_dof_transformation_function<T>();
957
958 const bool needs_transformation_data
959 = element->needs_dof_transformations() or L.needs_facet_permutations();
960 std::span<const std::uint32_t> cell_info;
961 if (needs_transformation_data)
962 {
963 mesh->topology_mutable().create_entity_permutations();
964 cell_info = std::span(mesh->topology().get_cell_permutation_info());
965 }
966
967 for (int i : L.integral_ids(IntegralType::cell))
968 {
969 const auto& fn = L.kernel(IntegralType::cell, i);
970 const auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
971 const std::vector<std::int32_t>& cells = L.cell_domains(i);
972 if (bs == 1)
973 {
974 impl::assemble_cells<T, 1>(dof_transform, b, mesh->geometry(), cells,
975 dofs, bs, fn, constants, coeffs, cstride,
976 cell_info);
977 }
978 else if (bs == 3)
979 {
980 impl::assemble_cells<T, 3>(dof_transform, b, mesh->geometry(), cells,
981 dofs, bs, fn, constants, coeffs, cstride,
982 cell_info);
983 }
984 else
985 {
986 impl::assemble_cells(dof_transform, b, mesh->geometry(), cells, dofs, bs,
987 fn, constants, coeffs, cstride, cell_info);
988 }
989 }
990
991 for (int i : L.integral_ids(IntegralType::exterior_facet))
992 {
993 const auto& fn = L.kernel(IntegralType::exterior_facet, i);
994 const auto& [coeffs, cstride]
995 = coefficients.at({IntegralType::exterior_facet, i});
996 const std::vector<std::int32_t>& facets = L.exterior_facet_domains(i);
997 if (bs == 1)
998 {
999 impl::assemble_exterior_facets<T, 1>(dof_transform, b, mesh->geometry(),
1000 facets, dofs, bs, fn, constants,
1001 coeffs, cstride, cell_info);
1002 }
1003 else if (bs == 3)
1004 {
1005 impl::assemble_exterior_facets<T, 3>(dof_transform, b, mesh->geometry(),
1006 facets, dofs, bs, fn, constants,
1007 coeffs, cstride, cell_info);
1008 }
1009 else
1010 {
1011 impl::assemble_exterior_facets(dof_transform, b, mesh->geometry(), facets,
1012 dofs, bs, fn, constants, coeffs, cstride,
1013 cell_info);
1014 }
1015 }
1016
1017 if (L.num_integrals(IntegralType::interior_facet) > 0)
1018 {
1019 std::function<std::uint8_t(std::size_t)> get_perm;
1020 if (L.needs_facet_permutations())
1021 {
1022 mesh->topology_mutable().create_entity_permutations();
1023 const std::vector<std::uint8_t>& perms
1024 = mesh->topology().get_facet_permutations();
1025 get_perm = [&perms](std::size_t i) { return perms[i]; };
1026 }
1027 else
1028 get_perm = [](std::size_t) { return 0; };
1029
1030 int num_cell_facets = mesh::cell_num_entities(mesh->topology().cell_type(),
1031 mesh->topology().dim() - 1);
1032 for (int i : L.integral_ids(IntegralType::interior_facet))
1033 {
1034 const auto& fn = L.kernel(IntegralType::interior_facet, i);
1035 const auto& [coeffs, cstride]
1036 = coefficients.at({IntegralType::interior_facet, i});
1037 const std::vector<std::int32_t>& facets = L.interior_facet_domains(i);
1038 if (bs == 1)
1039 {
1040 impl::assemble_interior_facets<T, 1>(
1041 dof_transform, b, mesh->geometry(), num_cell_facets, facets,
1042 *dofmap, fn, constants, coeffs, cstride, cell_info, get_perm);
1043 }
1044 else if (bs == 3)
1045 {
1046 impl::assemble_interior_facets<T, 3>(
1047 dof_transform, b, mesh->geometry(), num_cell_facets, facets,
1048 *dofmap, fn, constants, coeffs, cstride, cell_info, get_perm);
1049 }
1050 else
1051 {
1052 impl::assemble_interior_facets(
1053 dof_transform, b, mesh->geometry(), num_cell_facets, facets,
1054 *dofmap, fn, constants, coeffs, cstride, cell_info, get_perm);
1055 }
1056 }
1057 }
1058}
1059} // namespace dolfinx::fem::impl
Degree-of-freedom map representations and tools.
void cells(la::SparsityPattern &pattern, const mesh::Topology &topology, const std::array< const std::reference_wrapper< const DofMap >, 2 > &dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition: sparsitybuild.cpp:18
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
int cell_num_entities(CellType type, int dim)
Number of entities of dimension dim.
Definition: cell_types.cpp:139