Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/df/df5/Function_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Function.h
1// Copyright (C) 2003-2022 Anders Logg, Garth N. Wells and Massimiliano Leoni
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 "DofMap.h"
10#include "FiniteElement.h"
11#include "FunctionSpace.h"
12#include "interpolate.h"
13#include <concepts>
14#include <dolfinx/common/IndexMap.h>
15#include <dolfinx/common/types.h>
16#include <dolfinx/la/Vector.h>
17#include <dolfinx/mesh/Geometry.h>
18#include <dolfinx/mesh/Mesh.h>
19#include <dolfinx/mesh/Topology.h>
20#include <functional>
21#include <memory>
22#include <numeric>
23#include <span>
24#include <string>
25#include <utility>
26#include <vector>
27
28namespace dolfinx::fem
29{
30template <dolfinx::scalar T, std::floating_point U>
31class Expression;
32
42template <dolfinx::scalar T,
43 std::floating_point U = dolfinx::scalar_value_type_t<T>>
45{
46public:
49 using value_type = T;
51 using geometry_type = U;
52
55 explicit Function(std::shared_ptr<const FunctionSpace<geometry_type>> V)
56 : _function_space(V),
57 _x(std::make_shared<la::Vector<value_type>>(
58 V->dofmap()->index_map, V->dofmap()->index_map_bs()))
59 {
60 if (!V->component().empty())
61 {
62 throw std::runtime_error("Cannot create Function from subspace. Consider "
63 "collapsing the function space");
64 }
65 }
66
75 Function(std::shared_ptr<const FunctionSpace<geometry_type>> V,
76 std::shared_ptr<la::Vector<value_type>> x)
77 : _function_space(V), _x(x)
78 {
79 // We do not check for a subspace since this constructor is used for
80 // creating subfunctions
81
82 // Assertion uses '<=' to deal with sub-functions
83 assert(V->dofmap());
84 assert(V->dofmap()->index_map->size_global() * V->dofmap()->index_map_bs()
85 <= _x->bs() * _x->index_map()->size_global());
86 }
87
88 // Copy constructor
89 Function(const Function& v) = delete;
90
92 Function(Function&& v) = default;
93
95 ~Function() = default;
96
98 Function& operator=(Function&& v) = default;
99
100 // Assignment
101 Function& operator=(const Function& v) = delete;
102
106 Function sub(int i) const
107 {
108 auto sub_space = std::make_shared<FunctionSpace<geometry_type>>(
109 _function_space->sub({i}));
110 assert(sub_space);
111 return Function(sub_space, _x);
112 }
113
118 {
119 // Create new collapsed FunctionSpace
120 auto [V, map] = _function_space->collapse();
121
122 // Create new vector
123 auto x = std::make_shared<la::Vector<value_type>>(
124 V.dofmap()->index_map, V.dofmap()->index_map_bs());
125
126 // Copy values into new vector
127 std::span<const value_type> x_old = _x->array();
128 std::span<value_type> x_new = x->mutable_array();
129 for (std::size_t i = 0; i < map.size(); ++i)
130 {
131 assert((int)i < x_new.size());
132 assert(map[i] < x_old.size());
133 x_new[i] = x_old[map[i]];
134 }
135
136 return Function(
137 std::make_shared<FunctionSpace<geometry_type>>(std::move(V)), x);
138 }
139
142 std::shared_ptr<const FunctionSpace<geometry_type>> function_space() const
143 {
144 return _function_space;
145 }
146
148 std::shared_ptr<const la::Vector<value_type>> x() const { return _x; }
149
151 std::shared_ptr<la::Vector<value_type>> x() { return _x; }
152
163 std::span<const std::int32_t> cells,
164 std::span<const std::int32_t> cell_map,
165 const std::tuple<std::span<const std::int32_t>,
166 std::span<const std::int32_t>,
167 std::span<const geometry_type>,
168 std::span<const std::int32_t>>& nmm_interpolation_data
169 = {})
170 {
171 fem::interpolate(*this, v, cells, cell_map, nmm_interpolation_data);
172 }
173
182 std::span<const std::int32_t> cell_map = std::span<const std::int32_t>()
183 = {},
184 const std::tuple<std::span<const std::int32_t>,
185 std::span<const std::int32_t>,
186 std::span<const geometry_type>,
187 std::span<const std::int32_t>>& nmm_interpolation_data
188 = {})
189 {
190 assert(_function_space);
191 assert(_function_space->mesh());
192 int tdim = _function_space->mesh()->topology()->dim();
193 auto cell_imap = _function_space->mesh()->topology()->index_map(tdim);
194 assert(cell_imap);
195 std::int32_t num_cells = cell_imap->size_local() + cell_imap->num_ghosts();
196 std::vector<std::int32_t> cells(num_cells, 0);
197 std::iota(cells.begin(), cells.end(), 0);
198 interpolate(v, cells, cell_map, nmm_interpolation_data);
199 }
200
205 const std::function<
206 std::pair<std::vector<value_type>, std::vector<std::size_t>>(
207 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
208 const geometry_type,
209 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
210 std::size_t, 3,
211 MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>)>& f,
212 std::span<const std::int32_t> cells)
213 {
214 assert(_function_space);
215 assert(_function_space->element());
216 assert(_function_space->mesh());
217 const std::vector<geometry_type> x
218 = fem::interpolation_coords<geometry_type>(
219 *_function_space->element(), _function_space->mesh()->geometry(),
220 cells);
221 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
222 const geometry_type,
223 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
224 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>
225 _x(x.data(), 3, x.size() / 3);
226
227 const auto [fx, fshape] = f(_x);
228 assert(fshape.size() <= 2);
229 if (int vs = _function_space->value_size(); vs == 1 and fshape.size() == 1)
230 {
231 // Check for scalar-valued functions
232 if (fshape.front() != x.size() / 3)
233 throw std::runtime_error("Data returned by callable has wrong length");
234 }
235 else
236 {
237 // Check for vector/tensor value
238 if (fshape.size() != 2)
239 throw std::runtime_error("Expected 2D array of data");
240
241 if (fshape[0] != vs)
242 {
243 throw std::runtime_error(
244 "Data returned by callable has wrong shape(0) size");
245 }
246
247 if (fshape[1] != x.size() / 3)
248 {
249 throw std::runtime_error(
250 "Data returned by callable has wrong shape(1) size");
251 }
252 }
253
254 std::array<std::size_t, 2> _fshape;
255 if (fshape.size() == 1)
256 _fshape = {1, fshape[0]};
257 else
258 _fshape = {fshape[0], fshape[1]};
259
260 fem::interpolate(*this, std::span<const value_type>(fx.data(), fx.size()),
261 _fshape, cells);
262 }
263
266 void
267 interpolate(const std::function<
268 std::pair<std::vector<value_type>, std::vector<std::size_t>>(
269 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
270 const geometry_type,
271 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
272 std::size_t, 3,
273 MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>)>& f)
274 {
275 assert(_function_space);
276 assert(_function_space->mesh());
277 const int tdim = _function_space->mesh()->topology()->dim();
278 auto cell_map = _function_space->mesh()->topology()->index_map(tdim);
279 assert(cell_map);
280 std::int32_t num_cells = cell_map->size_local() + cell_map->num_ghosts();
281 std::vector<std::int32_t> cells(num_cells, 0);
282 std::iota(cells.begin(), cells.end(), 0);
283 interpolate(f, cells);
284 }
285
295 std::span<const std::int32_t> cells,
296 const dolfinx::mesh::Mesh<geometry_type>& expr_mesh,
297 std::span<const std::int32_t> cell_map
298 = std::span<const std::int32_t>())
299 {
300 // Check that spaces are compatible
301 assert(_function_space);
302 assert(_function_space->element());
303 std::size_t value_size = e.value_size();
304 if (e.argument_function_space())
305 throw std::runtime_error("Cannot interpolate Expression with Argument");
306
307 if (value_size != _function_space->value_size())
308 {
309 throw std::runtime_error(
310 "Function value size not equal to Expression value size");
311 }
312
313 {
314 // Compatibility check
315 auto [X0, shape0] = e.X();
316 auto [X1, shape1] = _function_space->element()->interpolation_points();
317 if (shape0 != shape1)
318 {
319 throw std::runtime_error(
320 "Function element interpolation points has different shape to "
321 "Expression interpolation points");
322 }
323
324 for (std::size_t i = 0; i < X0.size(); ++i)
325 {
326 if (std::abs(X0[i] - X1[i]) > 1.0e-10)
327 {
328 throw std::runtime_error("Function element interpolation points not "
329 "equal to Expression interpolation points");
330 }
331 }
332 }
333
334 // Array to hold evaluated Expression
335 std::size_t num_cells = cells.size();
336 std::size_t num_points = e.X().second[0];
337 std::vector<value_type> fdata(num_cells * num_points * value_size);
338 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
339 const value_type,
340 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
341 f(fdata.data(), num_cells, num_points, value_size);
342
343 // Evaluate Expression at points
344 assert(_function_space->mesh());
345
346 std::vector<std::int32_t> cells_expr;
347 cells_expr.reserve(num_cells);
348 // Get mesh and check if mesh is shared
349 if (auto mesh_v = _function_space->mesh();
350 expr_mesh.topology() == mesh_v->topology())
351 {
352 cells_expr.insert(cells_expr.end(), cells.begin(), cells.end());
353 }
354 // If meshes are different and input mapping is given
355 else if (!cell_map.empty())
356 {
357 std::transform(cells.begin(), cells.end(), std::back_inserter(cells_expr),
358 [&cell_map](std::int32_t c) { return cell_map[c]; });
359 }
360 else
361 std::runtime_error("Meshes are different and no cell map is provided");
362
363 e.eval(expr_mesh, cells_expr, fdata, {num_cells, num_points * value_size});
364
365 // Reshape evaluated data to fit interpolate
366 // Expression returns matrix of shape (num_cells, num_points *
367 // value_size), i.e. xyzxyz ordering of dof values per cell per
368 // point. The interpolation uses xxyyzz input, ordered for all
369 // points of each cell, i.e. (value_size, num_cells*num_points)
370 std::vector<value_type> fdata1(num_cells * num_points * value_size);
371 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
372 value_type, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
373 f1(fdata1.data(), value_size, num_cells, num_points);
374 for (std::size_t i = 0; i < f.extent(0); ++i)
375 for (std::size_t j = 0; j < f.extent(1); ++j)
376 for (std::size_t k = 0; k < f.extent(2); ++k)
377 f1(k, i, j) = f(i, j, k);
378
379 // Interpolate values into appropriate space
380 fem::interpolate(*this,
381 std::span<const value_type>(fdata1.data(), fdata1.size()),
382 {value_size, num_cells * num_points}, cells);
383 }
384
391 const dolfinx::mesh::Mesh<geometry_type>& expr_mesh,
392 std::span<const std::int32_t> cell_map
393 = std::span<const std::int32_t>() = {})
394 {
395 assert(_function_space);
396 assert(_function_space->mesh());
397 const int tdim = _function_space->mesh()->topology()->dim();
398 auto cell_imap = _function_space->mesh()->topology()->index_map(tdim);
399 assert(cell_imap);
400 std::int32_t num_cells = cell_imap->size_local() + cell_imap->num_ghosts();
401 std::vector<std::int32_t> cells(num_cells, 0);
402 std::iota(cells.begin(), cells.end(), 0);
403
404 if (cell_map.size() == 0)
405 interpolate(e, cells, expr_mesh, cell_map);
406 else
407 interpolate(e, cells, expr_mesh, std::span(cells.data(), cells.size()));
408 }
409
422 void eval(std::span<const geometry_type> x, std::array<std::size_t, 2> xshape,
423 std::span<const std::int32_t> cells, std::span<value_type> u,
424 std::array<std::size_t, 2> ushape) const
425 {
426 if (cells.empty())
427 return;
428
429 assert(x.size() == xshape[0] * xshape[1]);
430 assert(u.size() == ushape[0] * ushape[1]);
431
432 // TODO: This could be easily made more efficient by exploiting points
433 // being ordered by the cell to which they belong.
434
435 if (xshape[0] != cells.size())
436 {
437 throw std::runtime_error(
438 "Number of points and number of cells must be equal.");
439 }
440
441 if (xshape[0] != ushape[0])
442 {
443 throw std::runtime_error(
444 "Length of array for Function values must be the "
445 "same as the number of points.");
446 }
447
448 // Get mesh
449 assert(_function_space);
450 auto mesh = _function_space->mesh();
451 assert(mesh);
452 const std::size_t gdim = mesh->geometry().dim();
453 const std::size_t tdim = mesh->topology()->dim();
454 auto map = mesh->topology()->index_map(tdim);
455
456 // Get coordinate map
457 const CoordinateElement<geometry_type>& cmap = mesh->geometry().cmap();
458
459 // Get geometry data
460 auto x_dofmap = mesh->geometry().dofmap();
461 const std::size_t num_dofs_g = cmap.dim();
462 auto x_g = mesh->geometry().x();
463
464 // Get element
465 auto element = _function_space->element();
466 assert(element);
467 const int bs_element = element->block_size();
468 const std::size_t reference_value_size
469 = element->reference_value_size() / bs_element;
470 const std::size_t value_size = _function_space->value_size() / bs_element;
471 const std::size_t space_dimension = element->space_dimension() / bs_element;
472
473 // If the space has sub elements, concatenate the evaluations on the
474 // sub elements
475 const int num_sub_elements = element->num_sub_elements();
476 if (num_sub_elements > 1 and num_sub_elements != bs_element)
477 {
478 throw std::runtime_error("Function::eval is not supported for mixed "
479 "elements. Extract subspaces.");
480 }
481
482 // Create work vector for expansion coefficients
483 std::vector<value_type> coefficients(space_dimension * bs_element);
484
485 // Get dofmap
486 std::shared_ptr<const DofMap> dofmap = _function_space->dofmap();
487 assert(dofmap);
488 const int bs_dof = dofmap->bs();
489
490 std::span<const std::uint32_t> cell_info;
491 if (element->needs_dof_transformations())
492 {
493 mesh->topology_mutable()->create_entity_permutations();
494 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
495 }
496
497 std::vector<geometry_type> coord_dofs_b(num_dofs_g * gdim);
498 impl::mdspan_t<geometry_type, 2> coord_dofs(coord_dofs_b.data(), num_dofs_g,
499 gdim);
500 std::vector<geometry_type> xp_b(1 * gdim);
501 impl::mdspan_t<geometry_type, 2> xp(xp_b.data(), 1, gdim);
502
503 // Loop over points
504 std::fill(u.data(), u.data() + u.size(), 0.0);
505 std::span<const value_type> _v = _x->array();
506
507 // Evaluate geometry basis at point (0, 0, 0) on the reference cell.
508 // Used in affine case.
509 std::array<std::size_t, 4> phi0_shape = cmap.tabulate_shape(1, 1);
510 std::vector<geometry_type> phi0_b(std::reduce(
511 phi0_shape.begin(), phi0_shape.end(), 1, std::multiplies{}));
512 impl::mdspan_t<const geometry_type, 4> phi0(phi0_b.data(), phi0_shape);
513 cmap.tabulate(1, std::vector<geometry_type>(tdim), {1, tdim}, phi0_b);
514 auto dphi0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
515 phi0, std::pair(1, tdim + 1), 0,
516 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
517
518 // Data structure for evaluating geometry basis at specific points.
519 // Used in non-affine case.
520 std::array<std::size_t, 4> phi_shape = cmap.tabulate_shape(1, 1);
521 std::vector<geometry_type> phi_b(
522 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
523 impl::mdspan_t<const geometry_type, 4> phi(phi_b.data(), phi_shape);
524 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
525 phi, std::pair(1, tdim + 1), 0,
526 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
527
528 // Reference coordinates for each point
529 std::vector<geometry_type> Xb(xshape[0] * tdim);
530 impl::mdspan_t<geometry_type, 2> X(Xb.data(), xshape[0], tdim);
531
532 // Geometry data at each point
533 std::vector<geometry_type> J_b(xshape[0] * gdim * tdim);
534 impl::mdspan_t<geometry_type, 3> J(J_b.data(), xshape[0], gdim, tdim);
535 std::vector<geometry_type> K_b(xshape[0] * tdim * gdim);
536 impl::mdspan_t<geometry_type, 3> K(K_b.data(), xshape[0], tdim, gdim);
537 std::vector<geometry_type> detJ(xshape[0]);
538 std::vector<geometry_type> det_scratch(2 * gdim * tdim);
539
540 // Prepare geometry data in each cell
541 for (std::size_t p = 0; p < cells.size(); ++p)
542 {
543 const int cell_index = cells[p];
544
545 // Skip negative cell indices
546 if (cell_index < 0)
547 continue;
548
549 // Get cell geometry (coordinate dofs)
550 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
551 x_dofmap, cell_index, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
552 assert(x_dofs.size() == num_dofs_g);
553 for (std::size_t i = 0; i < num_dofs_g; ++i)
554 {
555 const int pos = 3 * x_dofs[i];
556 for (std::size_t j = 0; j < gdim; ++j)
557 coord_dofs(i, j) = x_g[pos + j];
558 }
559
560 for (std::size_t j = 0; j < gdim; ++j)
561 xp(0, j) = x[p * xshape[1] + j];
562
563 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
564 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
565 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
566 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
567 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
568 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
569
570 std::array<geometry_type, 3> Xpb = {0, 0, 0};
571 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
573 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
574 std::size_t, 1, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>
575 Xp(Xpb.data(), 1, tdim);
576
577 // Compute reference coordinates X, and J, detJ and K
578 if (cmap.is_affine())
579 {
581 _J);
583 std::array<geometry_type, 3> x0 = {0, 0, 0};
584 for (std::size_t i = 0; i < coord_dofs.extent(1); ++i)
585 x0[i] += coord_dofs(0, i);
587 detJ[p]
589 _J, det_scratch);
590 }
591 else
592 {
593 // Pull-back physical point xp to reference coordinate Xp
594 cmap.pull_back_nonaffine(Xp, xp, coord_dofs);
595
596 cmap.tabulate(1, std::span(Xpb.data(), tdim), {1, tdim}, phi_b);
598 _J);
600 detJ[p]
602 _J, det_scratch);
603 }
604
605 for (std::size_t j = 0; j < X.extent(1); ++j)
606 X(p, j) = Xpb[j];
607 }
608
609 // Prepare basis function data structures
610 std::vector<geometry_type> basis_derivatives_reference_values_b(
611 1 * xshape[0] * space_dimension * reference_value_size);
612 impl::mdspan_t<const geometry_type, 4> basis_derivatives_reference_values(
613 basis_derivatives_reference_values_b.data(), 1, xshape[0],
614 space_dimension, reference_value_size);
615 std::vector<geometry_type> basis_values_b(space_dimension * value_size);
616 impl::mdspan_t<geometry_type, 2> basis_values(basis_values_b.data(),
617 space_dimension, value_size);
618
619 // Compute basis on reference element
620 element->tabulate(basis_derivatives_reference_values_b, Xb,
621 {X.extent(0), X.extent(1)}, 0);
622
623 using xu_t = impl::mdspan_t<geometry_type, 2>;
624 using xU_t = impl::mdspan_t<const geometry_type, 2>;
625 using xJ_t = impl::mdspan_t<const geometry_type, 2>;
626 using xK_t = impl::mdspan_t<const geometry_type, 2>;
627 auto push_forward_fn
628 = element->basix_element().template map_fn<xu_t, xU_t, xJ_t, xK_t>();
629
630 // Transformation function for basis function values
631 auto apply_dof_transformation
632 = element->template dof_transformation_fn<geometry_type>(
634
635 // Size of tensor for symmetric elements, unused in non-symmetric case, but placed outside the loop
636 // for pre-computation.
637 int matrix_size;
638 if (element->symmetric())
639 {
640 matrix_size = 0;
641 while (matrix_size * matrix_size < ushape[1])
642 ++matrix_size;
643 }
644 const std::size_t num_basis_values = space_dimension * reference_value_size;
645 for (std::size_t p = 0; p < cells.size(); ++p)
646 {
647 const int cell_index = cells[p];
648 if (cell_index < 0) // Skip negative cell indices
649 continue;
650
651 // Permute the reference basis function values to account for the
652 // cell's orientation
653 apply_dof_transformation(
654 std::span(basis_derivatives_reference_values_b.data()
655 + p * num_basis_values,
656 num_basis_values),
657 cell_info, cell_index, reference_value_size);
658
659 {
660 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
661 basis_derivatives_reference_values, 0, p,
662 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
663 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
664 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
665 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
666 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
667 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
668 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
669 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
670 push_forward_fn(basis_values, _U, _J, detJ[p], _K);
671 }
672
673 // Get degrees of freedom for current cell
674 std::span<const std::int32_t> dofs = dofmap->cell_dofs(cell_index);
675 for (std::size_t i = 0; i < dofs.size(); ++i)
676 for (int k = 0; k < bs_dof; ++k)
677 coefficients[bs_dof * i + k] = _v[bs_dof * dofs[i] + k];
678
679 if (element->symmetric())
680 {
681 int row = 0;
682 int rowstart = 0;
683 // Compute expansion
684 for (int k = 0; k < bs_element; ++k)
685 {
686 if (k - rowstart > row)
687 {
688 row++;
689 rowstart = k;
690 }
691 for (std::size_t i = 0; i < space_dimension; ++i)
692 {
693 for (std::size_t j = 0; j < value_size; ++j)
694 {
695 u[p * ushape[1]
696 + (j * bs_element + row * matrix_size + k - rowstart)]
697 += coefficients[bs_element * i + k] * basis_values(i, j);
698 if (k - rowstart != row)
699 {
700 u[p * ushape[1]
701 + (j * bs_element + row + matrix_size * (k - rowstart))]
702 += coefficients[bs_element * i + k] * basis_values(i, j);
703 }
704 }
705 }
706 }
707 }
708 else
709 {
710 // Compute expansion
711 for (int k = 0; k < bs_element; ++k)
712 {
713 for (std::size_t i = 0; i < space_dimension; ++i)
714 {
715 for (std::size_t j = 0; j < value_size; ++j)
716 {
717 u[p * ushape[1] + (j * bs_element + k)]
718 += coefficients[bs_element * i + k] * basis_values(i, j);
719 }
720 }
721 }
722 }
723 }
724 }
725
727 std::string name = "u";
728
729private:
730 // The function space
731 std::shared_ptr<const FunctionSpace<geometry_type>> _function_space;
732
733 // The vector of expansion coefficients (local)
734 std::shared_ptr<la::Vector<value_type>> _x;
735};
736
737} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Definition CoordinateElement.h:38
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Definition CoordinateElement.h:105
static void pull_back_affine(U &&X, const V &K, const std::array< T, 3 > &x0, const W &x)
Compute reference coordinates X for physical coordinates x for an affine map. For the affine case,...
Definition CoordinateElement.h:188
void tabulate(int nd, std::span< const T > X, std::array< std::size_t, 2 > shape, std::span< T > basis) const
Evaluate basis values and derivatives at set of points.
Definition CoordinateElement.cpp:54
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:114
std::array< std::size_t, 4 > tabulate_shape(std::size_t nd, std::size_t num_points) const
Shape of array to fill when calling tabulate.
Definition CoordinateElement.cpp:47
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:194
void pull_back_nonaffine(mdspan2_t< T > X, mdspan2_t< const T > x, mdspan2_t< const T > cell_geometry, double tol=1.0e-6, int maxit=15) const
Compute reference coordinates X for physical coordinates x for a non-affine map.
Definition CoordinateElement.cpp:71
bool is_affine() const noexcept
Definition CoordinateElement.h:240
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:131
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell.
Definition Expression.h:41
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
Definition Function.h:45
U geometry_type
Geometry type of the Mesh that the Function is defined on.
Definition Function.h:51
std::shared_ptr< const FunctionSpace< geometry_type > > function_space() const
Access the function space.
Definition Function.h:142
void interpolate(const std::function< std::pair< std::vector< value_type >, std::vector< std::size_t > >(MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const geometry_type, MDSPAN_IMPL_STANDARD_NAMESPACE::extents< std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent > >)> &f)
Interpolate an expression function on the whole domain.
Definition Function.h:267
void interpolate(const Expression< value_type, geometry_type > &e, const dolfinx::mesh::Mesh< geometry_type > &expr_mesh, std::span< const std::int32_t > cell_map=std::span< const std::int32_t >()={})
Definition Function.h:390
Function(std::shared_ptr< const FunctionSpace< geometry_type > > V)
Definition Function.h:55
Function collapse() const
Collapse a subfunction (view into a Function) to a stand-alone Function.
Definition Function.h:117
Function sub(int i) const
Extract a sub-function (a view into the Function).
Definition Function.h:106
void interpolate(const Function< value_type, geometry_type > &v, std::span< const std::int32_t > cells, std::span< const std::int32_t > cell_map, const std::tuple< std::span< const std::int32_t >, std::span< const std::int32_t >, std::span< const geometry_type >, std::span< const std::int32_t > > &nmm_interpolation_data={})
Interpolate a provided Function.
Definition Function.h:161
std::string name
Name.
Definition Function.h:727
std::shared_ptr< la::Vector< value_type > > x()
Underlying vector.
Definition Function.h:151
void interpolate(const Expression< value_type, geometry_type > &e, std::span< const std::int32_t > cells, const dolfinx::mesh::Mesh< geometry_type > &expr_mesh, std::span< const std::int32_t > cell_map=std::span< const std::int32_t >())
Interpolate an Expression (based on UFL)
Definition Function.h:294
~Function()=default
Destructor.
void interpolate(const std::function< std::pair< std::vector< value_type >, std::vector< std::size_t > >(MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const geometry_type, MDSPAN_IMPL_STANDARD_NAMESPACE::extents< std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent > >)> &f, std::span< const std::int32_t > cells)
Definition Function.h:204
void eval(std::span< const geometry_type > x, std::array< std::size_t, 2 > xshape, std::span< const std::int32_t > cells, std::span< value_type > u, std::array< std::size_t, 2 > ushape) const
Evaluate the Function at points.
Definition Function.h:422
Function(Function &&v)=default
Move constructor.
std::shared_ptr< const la::Vector< value_type > > x() const
Underlying vector.
Definition Function.h:148
Function & operator=(Function &&v)=default
Move assignment.
Function(std::shared_ptr< const FunctionSpace< geometry_type > > V, std::shared_ptr< la::Vector< value_type > > x)
Create function on given function space with a given vector.
Definition Function.h:75
void interpolate(const Function< value_type, geometry_type > &v, std::span< const std::int32_t > cell_map=std::span< const std::int32_t >()={}, const std::tuple< std::span< const std::int32_t >, std::span< const std::int32_t >, std::span< const geometry_type >, std::span< const std::int32_t > > &nmm_interpolation_data={})
Interpolate a provided Function.
Definition Function.h:180
T value_type
Definition Function.h:49
Definition Vector.h:31
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
std::shared_ptr< Topology > topology()
Get mesh topology.
Definition Mesh.h:64
Definition types.h:20
Finite element method functionality.
Definition assemble_matrix_impl.h:26
void interpolate(Function< T, U > &u, std::span< const T > f, std::array< std::size_t, 2 > fshape, std::span< const std::int32_t > cells)
Interpolate an expression f(x) in a finite element space.
Definition interpolate.h:739