DOLFINx 0.9.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Function.h
1// Copyright (C) 2003-2024 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 <algorithm>
14#include <concepts>
15#include <dolfinx/common/IndexMap.h>
16#include <dolfinx/common/types.h>
17#include <dolfinx/la/Vector.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 <numeric>
24#include <span>
25#include <string>
26#include <utility>
27#include <vector>
28
29namespace dolfinx::fem
30{
31template <dolfinx::scalar T, std::floating_point U>
33
43template <dolfinx::scalar T,
44 std::floating_point U = dolfinx::scalar_value_type_t<T>>
45class Function
46{
47public:
50 using value_type = T;
52 using geometry_type = U;
53
56 explicit Function(std::shared_ptr<const FunctionSpace<geometry_type>> V)
57 : _function_space(V),
58 _x(std::make_shared<la::Vector<value_type>>(
59 V->dofmap()->index_map, V->dofmap()->index_map_bs()))
60 {
61 if (!V->component().empty())
62 {
63 throw std::runtime_error("Cannot create Function from subspace. Consider "
64 "collapsing the function space");
65 }
66 }
67
76 Function(std::shared_ptr<const FunctionSpace<geometry_type>> V,
77 std::shared_ptr<la::Vector<value_type>> x)
78 : _function_space(V), _x(x)
79 {
80 // NOTE: We do not check for a subspace since this constructor is
81 // used for creating subfunctions
82
83 // Assertion uses '<=' to deal with sub-functions
84 assert(V->dofmap());
85 assert(V->dofmap()->index_map->size_global() * V->dofmap()->index_map_bs()
86 <= _x->bs() * _x->index_map()->size_global());
87 }
88
89 // Copy constructor
90 Function(const Function& v) = delete;
91
93 Function(Function&& v) = default;
94
96 ~Function() = default;
97
99 Function& operator=(Function&& v) = default;
100
101 // Assignment
102 Function& operator=(const Function& v) = delete;
103
107 Function sub(int i) const
108 {
109 auto sub_space = std::make_shared<FunctionSpace<geometry_type>>(
110 _function_space->sub({i}));
111 assert(sub_space);
112 return Function(sub_space, _x);
113 }
114
119 {
120 // Create new collapsed FunctionSpace
121 auto [V, map] = _function_space->collapse();
122
123 // Create new vector
124 auto x = std::make_shared<la::Vector<value_type>>(
125 V.dofmap()->index_map, V.dofmap()->index_map_bs());
126
127 // Copy values into new vector
128 std::span<const value_type> x_old = _x->array();
129 std::span<value_type> x_new = x->mutable_array();
130 for (std::size_t i = 0; i < map.size(); ++i)
131 {
132 assert((int)i < x_new.size());
133 assert(map[i] < x_old.size());
134 x_new[i] = x_old[map[i]];
135 }
136
137 return Function(
138 std::make_shared<FunctionSpace<geometry_type>>(std::move(V)), x);
139 }
140
143 std::shared_ptr<const FunctionSpace<geometry_type>> function_space() const
144 {
145 return _function_space;
146 }
147
149 std::shared_ptr<const la::Vector<value_type>> x() const { return _x; }
150
152 std::shared_ptr<la::Vector<value_type>> x() { return _x; }
153
156 void
157 interpolate(const std::function<
158 std::pair<std::vector<value_type>, std::vector<std::size_t>>(
159 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
160 const geometry_type,
161 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
162 std::size_t, 3,
163 MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>)>& f)
164 {
165 assert(_function_space);
166 assert(_function_space->mesh());
167 int tdim = _function_space->mesh()->topology()->dim();
168 auto cmap = _function_space->mesh()->topology()->index_map(tdim);
169 assert(cmap);
170 std::vector<std::int32_t> cells(cmap->size_local() + cmap->num_ghosts(), 0);
171 std::iota(cells.begin(), cells.end(), 0);
172 interpolate(f, cells);
173 }
174
179 const std::function<
180 std::pair<std::vector<value_type>, std::vector<std::size_t>>(
181 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
182 const geometry_type,
183 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
184 std::size_t, 3,
185 MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>)>& f,
186 std::span<const std::int32_t> cells)
187 {
188 assert(_function_space);
189 assert(_function_space->element());
190 assert(_function_space->mesh());
191 const std::vector<geometry_type> x
193 *_function_space->element(), _function_space->mesh()->geometry(),
194 cells);
195 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
196 const geometry_type,
197 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
198 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>
199 _x(x.data(), 3, x.size() / 3);
200
201 const auto [fx, fshape] = f(_x);
202 assert(fshape.size() <= 2);
203 if (int vs = _function_space->value_size(); vs == 1 and fshape.size() == 1)
204 {
205 // Check for scalar-valued functions
206 if (fshape.front() != x.size() / 3)
207 throw std::runtime_error("Data returned by callable has wrong length");
208 }
209 else
210 {
211 // Check for vector/tensor value
212 if (fshape.size() != 2)
213 throw std::runtime_error("Expected 2D array of data");
214
215 if (fshape[0] != vs)
216 {
217 throw std::runtime_error(
218 "Data returned by callable has wrong shape(0) size");
219 }
220
221 if (fshape[1] != x.size() / 3)
222 {
223 throw std::runtime_error(
224 "Data returned by callable has wrong shape(1) size");
225 }
226 }
227
228 std::array<std::size_t, 2> _fshape;
229 if (fshape.size() == 1)
230 _fshape = {1, fshape[0]};
231 else
232 _fshape = {fshape[0], fshape[1]};
233
234 fem::interpolate(*this, std::span<const value_type>(fx.data(), fx.size()),
235 _fshape, cells);
236 }
237
244 {
245 assert(_function_space);
246 assert(_function_space->mesh());
247 int tdim = _function_space->mesh()->topology()->dim();
248 auto cmap = _function_space->mesh()->topology()->index_map(tdim);
249 assert(cmap);
250 std::vector<std::int32_t> cells(cmap->size_local() + cmap->num_ghosts(), 0);
251 std::iota(cells.begin(), cells.end(), 0);
252 interpolate(u, cells, cells);
253 }
254
272 std::span<const std::int32_t> cells0,
273 std::span<const std::int32_t> cells1 = {})
274 {
275 if (cells1.empty())
276 cells1 = cells0;
277 fem::interpolate(*this, cells1, u0, cells0);
278 }
279
286 {
287 assert(_function_space);
288 assert(_function_space->mesh());
289 int tdim = _function_space->mesh()->topology()->dim();
290 auto cmap = _function_space->mesh()->topology()->index_map(tdim);
291 assert(cmap);
292 std::vector<std::int32_t> cells(cmap->size_local() + cmap->num_ghosts(), 0);
293 std::iota(cells.begin(), cells.end(), 0);
294 interpolate(e, cells);
295 }
296
314 std::span<const std::int32_t> cells0,
315 std::span<const std::int32_t> cells1 = {})
316 {
317 // Extract mesh
318 const mesh::Mesh<geometry_type>* mesh0 = nullptr;
319 for (auto& c : e0.coefficients())
320 {
321 assert(c);
322 assert(c->function_space());
323 assert(c->function_space()->mesh());
324 if (const mesh::Mesh<geometry_type>* mesh
325 = c->function_space()->mesh().get();
326 !mesh0)
327 {
328 mesh0 = mesh;
329 }
330 else if (mesh != mesh0)
331 {
332 throw std::runtime_error(
333 "Expression coefficient Functions have different meshes.");
334 }
335 }
336
337 // If Expression has no Function coefficients take mesh from `this`.
338 assert(_function_space);
339 assert(_function_space->mesh());
340 if (!mesh0)
341 mesh0 = _function_space->mesh().get();
342
343 // If cells1 is empty and Function and Expression meshes are the
344 // same, make cells1 the same as cells0. Otherwise check that
345 // lengths of cells0 and cells1 are the same.
346 if (cells1.empty() and mesh0 == _function_space->mesh().get())
347 cells1 = cells0;
348 else if (cells0.size() != cells1.size())
349 throw std::runtime_error("Cells lists have different lengths.");
350
351 // Check that Function and Expression spaces are compatible
352 assert(_function_space->element());
353 std::size_t value_size = e0.value_size();
355 throw std::runtime_error("Cannot interpolate Expression with Argument.");
356 if (value_size != _function_space->value_size())
357 {
358 throw std::runtime_error(
359 "Function value size not equal to Expression value size.");
360 }
361
362 // Compatibility check
363 {
364 auto [X0, shape0] = e0.X();
365 auto [X1, shape1] = _function_space->element()->interpolation_points();
366 if (shape0 != shape1)
367 {
368 throw std::runtime_error(
369 "Function element interpolation points has different shape to "
370 "Expression interpolation points");
371 }
372
373 for (std::size_t i = 0; i < X0.size(); ++i)
374 {
375 if (std::abs(X0[i] - X1[i]) > 1.0e-10)
376 {
377 throw std::runtime_error("Function element interpolation points not "
378 "equal to Expression interpolation points");
379 }
380 }
381 }
382
383 // Array to hold evaluated Expression
384 std::size_t num_cells = cells0.size();
385 std::size_t num_points = e0.X().second[0];
386 std::vector<value_type> fdata(num_cells * num_points * value_size);
387 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
388 const value_type,
389 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
390 f(fdata.data(), num_cells, num_points, value_size);
391
392 // Evaluate Expression at points
393 e0.eval(*mesh0, cells0, fdata, {num_cells, num_points * value_size});
394
395 // Reshape evaluated data to fit interpolate.
396 // Expression returns matrix of shape (num_cells, num_points *
397 // value_size), i.e. xyzxyz ordering of dof values per cell per
398 // point. The interpolation uses xxyyzz input, ordered for all
399 // points of each cell, i.e. (value_size, num_cells*num_points).
400 std::vector<value_type> fdata1(num_cells * num_points * value_size);
401 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
402 value_type, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
403 f1(fdata1.data(), value_size, num_cells, num_points);
404 for (std::size_t i = 0; i < f.extent(0); ++i)
405 for (std::size_t j = 0; j < f.extent(1); ++j)
406 for (std::size_t k = 0; k < f.extent(2); ++k)
407 f1(k, i, j) = f(i, j, k);
408
409 // Interpolate values into appropriate space
410 fem::interpolate(*this,
411 std::span<const value_type>(fdata1.data(), fdata1.size()),
412 {value_size, num_cells * num_points}, cells1);
413 }
414
424 std::span<const std::int32_t> cells,
425 const geometry::PointOwnershipData<U>& interpolation_data)
426 {
427 fem::interpolate(*this, v, cells, interpolation_data);
428 }
429
442 void eval(std::span<const geometry_type> x, std::array<std::size_t, 2> xshape,
443 std::span<const std::int32_t> cells, std::span<value_type> u,
444 std::array<std::size_t, 2> ushape) const
445 {
446 if (cells.empty())
447 return;
448
449 assert(x.size() == xshape[0] * xshape[1]);
450 assert(u.size() == ushape[0] * ushape[1]);
451
452 // TODO: This could be easily made more efficient by exploiting
453 // points being ordered by the cell to which they belong.
454
455 if (xshape[0] != cells.size())
456 {
457 throw std::runtime_error(
458 "Number of points and number of cells must be equal.");
459 }
460
461 if (xshape[0] != ushape[0])
462 {
463 throw std::runtime_error(
464 "Length of array for Function values must be the "
465 "same as the number of points.");
466 }
467
468 // Get mesh
469 assert(_function_space);
470 auto mesh = _function_space->mesh();
471 assert(mesh);
472 const std::size_t gdim = mesh->geometry().dim();
473 const std::size_t tdim = mesh->topology()->dim();
474 auto map = mesh->topology()->index_map(tdim);
475
476 // Get coordinate map
477 const CoordinateElement<geometry_type>& cmap = mesh->geometry().cmap();
478
479 // Get geometry data
480 auto x_dofmap = mesh->geometry().dofmap();
481 const std::size_t num_dofs_g = cmap.dim();
482 auto x_g = mesh->geometry().x();
483
484 // Get element
485 auto element = _function_space->element();
486 assert(element);
487 const int bs_element = element->block_size();
488 const std::size_t reference_value_size
489 = element->reference_value_size() / bs_element;
490 const std::size_t value_size = _function_space->value_size() / bs_element;
491 const std::size_t space_dimension = element->space_dimension() / bs_element;
492
493 // If the space has sub elements, concatenate the evaluations on the
494 // sub elements
495 const int num_sub_elements = element->num_sub_elements();
496 if (num_sub_elements > 1 and num_sub_elements != bs_element)
497 {
498 throw std::runtime_error("Function::eval is not supported for mixed "
499 "elements. Extract subspaces.");
500 }
501
502 // Create work vector for expansion coefficients
503 std::vector<value_type> coefficients(space_dimension * bs_element);
504
505 // Get dofmap
506 std::shared_ptr<const DofMap> dofmap = _function_space->dofmap();
507 assert(dofmap);
508 const int bs_dof = dofmap->bs();
509
510 std::span<const std::uint32_t> cell_info;
511 if (element->needs_dof_transformations())
512 {
513 mesh->topology_mutable()->create_entity_permutations();
514 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
515 }
516
517 std::vector<geometry_type> coord_dofs_b(num_dofs_g * gdim);
518 impl::mdspan_t<geometry_type, 2> coord_dofs(coord_dofs_b.data(), num_dofs_g,
519 gdim);
520 std::vector<geometry_type> xp_b(1 * gdim);
521 impl::mdspan_t<geometry_type, 2> xp(xp_b.data(), 1, gdim);
522
523 // Loop over points
524 std::ranges::fill(u, 0.0);
525 std::span<const value_type> _v = _x->array();
526
527 // Evaluate geometry basis at point (0, 0, 0) on the reference cell.
528 // Used in affine case.
529 std::array<std::size_t, 4> phi0_shape = cmap.tabulate_shape(1, 1);
530 std::vector<geometry_type> phi0_b(std::reduce(
531 phi0_shape.begin(), phi0_shape.end(), 1, std::multiplies{}));
532 impl::mdspan_t<const geometry_type, 4> phi0(phi0_b.data(), phi0_shape);
533 cmap.tabulate(1, std::vector<geometry_type>(tdim), {1, tdim}, phi0_b);
534 auto dphi0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
535 phi0, std::pair(1, tdim + 1), 0,
536 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
537
538 // Data structure for evaluating geometry basis at specific points.
539 // Used in non-affine case.
540 std::array<std::size_t, 4> phi_shape = cmap.tabulate_shape(1, 1);
541 std::vector<geometry_type> phi_b(
542 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
543 impl::mdspan_t<const geometry_type, 4> phi(phi_b.data(), phi_shape);
544 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
545 phi, std::pair(1, tdim + 1), 0,
546 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
547
548 // Reference coordinates for each point
549 std::vector<geometry_type> Xb(xshape[0] * tdim);
550 impl::mdspan_t<geometry_type, 2> X(Xb.data(), xshape[0], tdim);
551
552 // Geometry data at each point
553 std::vector<geometry_type> J_b(xshape[0] * gdim * tdim);
554 impl::mdspan_t<geometry_type, 3> J(J_b.data(), xshape[0], gdim, tdim);
555 std::vector<geometry_type> K_b(xshape[0] * tdim * gdim);
556 impl::mdspan_t<geometry_type, 3> K(K_b.data(), xshape[0], tdim, gdim);
557 std::vector<geometry_type> detJ(xshape[0]);
558 std::vector<geometry_type> det_scratch(2 * gdim * tdim);
559
560 // Prepare geometry data in each cell
561 for (std::size_t p = 0; p < cells.size(); ++p)
562 {
563 const int cell_index = cells[p];
564
565 // Skip negative cell indices
566 if (cell_index < 0)
567 continue;
568
569 // Get cell geometry (coordinate dofs)
570 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
571 x_dofmap, cell_index, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
572 assert(x_dofs.size() == num_dofs_g);
573 for (std::size_t i = 0; i < num_dofs_g; ++i)
574 {
575 const int pos = 3 * x_dofs[i];
576 for (std::size_t j = 0; j < gdim; ++j)
577 coord_dofs(i, j) = x_g[pos + j];
578 }
579
580 for (std::size_t j = 0; j < gdim; ++j)
581 xp(0, j) = x[p * xshape[1] + j];
582
583 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
584 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
585 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
586 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
587 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
588 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
589
590 std::array<geometry_type, 3> Xpb = {0, 0, 0};
591 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
593 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
594 std::size_t, 1, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>
595 Xp(Xpb.data(), 1, tdim);
596
597 // Compute reference coordinates X, and J, detJ and K
598 if (cmap.is_affine())
599 {
601 _J);
603 std::array<geometry_type, 3> x0 = {0, 0, 0};
604 for (std::size_t i = 0; i < coord_dofs.extent(1); ++i)
605 x0[i] += coord_dofs(0, i);
607 detJ[p]
609 _J, det_scratch);
610 }
611 else
612 {
613 // Pull-back physical point xp to reference coordinate Xp
614 cmap.pull_back_nonaffine(Xp, xp, coord_dofs);
615 cmap.tabulate(1, std::span(Xpb.data(), tdim), {1, tdim}, phi_b);
617 _J);
619 detJ[p]
621 _J, det_scratch);
622 }
623
624 for (std::size_t j = 0; j < X.extent(1); ++j)
625 X(p, j) = Xpb[j];
626 }
627
628 // Prepare basis function data structures
629 std::vector<geometry_type> basis_derivatives_reference_values_b(
630 1 * xshape[0] * space_dimension * reference_value_size);
631 impl::mdspan_t<const geometry_type, 4> basis_derivatives_reference_values(
632 basis_derivatives_reference_values_b.data(), 1, xshape[0],
633 space_dimension, reference_value_size);
634 std::vector<geometry_type> basis_values_b(space_dimension * value_size);
635 impl::mdspan_t<geometry_type, 2> basis_values(basis_values_b.data(),
636 space_dimension, value_size);
637
638 // Compute basis on reference element
639 element->tabulate(basis_derivatives_reference_values_b, Xb,
640 {X.extent(0), X.extent(1)}, 0);
641
642 using xu_t = impl::mdspan_t<geometry_type, 2>;
643 using xU_t = impl::mdspan_t<const geometry_type, 2>;
644 using xJ_t = impl::mdspan_t<const geometry_type, 2>;
645 using xK_t = impl::mdspan_t<const geometry_type, 2>;
646 auto push_forward_fn
647 = element->basix_element().template map_fn<xu_t, xU_t, xJ_t, xK_t>();
648
649 // Transformation function for basis function values
650 auto apply_dof_transformation
651 = element->template dof_transformation_fn<geometry_type>(
653
654 // Size of tensor for symmetric elements, unused in non-symmetric case, but
655 // placed outside the loop for pre-computation.
656 int matrix_size;
657 if (element->symmetric())
658 {
659 matrix_size = 0;
660 while (matrix_size * matrix_size < ushape[1])
661 ++matrix_size;
662 }
663
664 const std::size_t num_basis_values = space_dimension * reference_value_size;
665 for (std::size_t p = 0; p < cells.size(); ++p)
666 {
667 const int cell_index = cells[p];
668 if (cell_index < 0) // Skip negative cell indices
669 continue;
670
671 // Permute the reference basis function values to account for the
672 // cell's orientation
673 apply_dof_transformation(
674 std::span(basis_derivatives_reference_values_b.data()
675 + p * num_basis_values,
676 num_basis_values),
677 cell_info, cell_index, reference_value_size);
678
679 {
680 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
681 basis_derivatives_reference_values, 0, p,
682 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
683 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
684 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
685 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
686 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
687 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
688 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
689 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
690 push_forward_fn(basis_values, _U, _J, detJ[p], _K);
691 }
692
693 // Get degrees of freedom for current cell
694 std::span<const std::int32_t> dofs = dofmap->cell_dofs(cell_index);
695 for (std::size_t i = 0; i < dofs.size(); ++i)
696 for (int k = 0; k < bs_dof; ++k)
697 coefficients[bs_dof * i + k] = _v[bs_dof * dofs[i] + k];
698
699 if (element->symmetric())
700 {
701 int row = 0;
702 int rowstart = 0;
703 // Compute expansion
704 for (int k = 0; k < bs_element; ++k)
705 {
706 if (k - rowstart > row)
707 {
708 row++;
709 rowstart = k;
710 }
711 for (std::size_t i = 0; i < space_dimension; ++i)
712 {
713 for (std::size_t j = 0; j < value_size; ++j)
714 {
715 u[p * ushape[1]
716 + (j * bs_element + row * matrix_size + k - rowstart)]
717 += coefficients[bs_element * i + k] * basis_values(i, j);
718 if (k - rowstart != row)
719 {
720 u[p * ushape[1]
721 + (j * bs_element + row + matrix_size * (k - rowstart))]
722 += coefficients[bs_element * i + k] * basis_values(i, j);
723 }
724 }
725 }
726 }
727 }
728 else
729 {
730 // Compute expansion
731 for (int k = 0; k < bs_element; ++k)
732 {
733 for (std::size_t i = 0; i < space_dimension; ++i)
734 {
735 for (std::size_t j = 0; j < value_size; ++j)
736 {
737 u[p * ushape[1] + (j * bs_element + k)]
738 += coefficients[bs_element * i + k] * basis_values(i, j);
739 }
740 }
741 }
742 }
743 }
744 }
745
747 std::string name = "u";
748
749private:
750 // Function space
751 std::shared_ptr<const FunctionSpace<geometry_type>> _function_space;
752
753 // Vector of expansion coefficients (local)
754 std::shared_ptr<la::Vector<value_type>> _x;
755};
756
757} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Definition XDMFFile.h:27
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Definition CoordinateElement.h:126
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:209
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:55
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:135
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:48
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:206
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:83
bool is_affine() const noexcept
Definition CoordinateElement.h:261
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:152
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell.
Definition Function.h:32
std::pair< std::vector< geometry_type >, std::array< std::size_t, 2 > > X() const
Evaluation points on the reference cell.
Definition Expression.h:276
const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > & coefficients() const
Get coefficients.
Definition Expression.h:114
void eval(const mesh::Mesh< geometry_type > &mesh, std::span< const std::int32_t > entities, std::span< scalar_type > values, std::array< std::size_t, 2 > vshape) const
Evaluate Expression on cells or facets.
Definition Expression.h:156
int value_size() const
Get value size.
Definition Expression.h:264
std::shared_ptr< const FunctionSpace< geometry_type > > argument_function_space() const
Get argument function space.
Definition Expression.h:105
This class represents a finite element function space defined by a mesh, a finite element,...
Definition vtk_utils.h:32
Definition XDMFFile.h:29
U geometry_type
Geometry type of the Mesh that the Function is defined on.
Definition Function.h:52
void interpolate(const Function< value_type, geometry_type > &u)
Interpolate a Function over all cells.
Definition Function.h:243
std::shared_ptr< const FunctionSpace< geometry_type > > function_space() const
Access the function space.
Definition Function.h:143
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 f(x) on the whole domain.
Definition Function.h:157
void interpolate(const Function< value_type, geometry_type > &v, std::span< const std::int32_t > cells, const geometry::PointOwnershipData< U > &interpolation_data)
Interpolate a Function defined on a different mesh.
Definition Function.h:423
void interpolate(const Expression< value_type, geometry_type > &e0, std::span< const std::int32_t > cells0, std::span< const std::int32_t > cells1={})
Interpolate an Expression over a subset of cells.
Definition Function.h:313
Function(std::shared_ptr< const FunctionSpace< geometry_type > > V)
Create function on given function space.
Definition Function.h:56
void interpolate(const Function< value_type, geometry_type > &u0, std::span< const std::int32_t > cells0, std::span< const std::int32_t > cells1={})
Interpolate a Function over a subset of cells.
Definition Function.h:271
Function collapse() const
Collapse a subfunction (view into a Function) to a stand-alone Function.
Definition Function.h:118
Function sub(int i) const
Extract a sub-function (a view into the Function).
Definition Function.h:107
std::string name
Name.
Definition Function.h:747
std::shared_ptr< la::Vector< value_type > > x()
Underlying vector.
Definition Function.h:152
void interpolate(const Expression< value_type, geometry_type > &e)
Interpolate an Expression on all cells.
Definition Function.h:285
~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)
Interpolate an expression f(x) over a set of cells.
Definition Function.h:178
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:442
Function(Function &&v)=default
Move constructor.
std::shared_ptr< const la::Vector< value_type > > x() const
Underlying vector (const version).
Definition Function.h:149
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:76
T value_type
Definition Function.h:50
Definition Vector.h:32
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
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 evaluated expression f(x) in a finite element space.
Definition interpolate.h:703
std::vector< T > interpolation_coords(const fem::FiniteElement< T > &element, const mesh::Geometry< T > &geometry, std::span< const std::int32_t > cells)
Compute the evaluation points in the physical space at which an expression should be computed to inte...
Definition interpolate.h:50
Information on the ownership of points distributed across processes.
Definition utils.h:30