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