Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
function_lib.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1999 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_function_lib_h
16#define dealii_function_lib_h
17
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/point.h>
24#include <deal.II/base/table.h>
25
26#include <array>
27
29
38namespace Functions
39{
49 template <int dim>
50 class SquareFunction : public Function<dim>
51 {
52 public:
53 virtual double
54 value(const Point<dim> &p, const unsigned int component = 0) const override;
55 virtual void
56 vector_value(const Point<dim> &p, Vector<double> &values) const override;
57 virtual void
58 value_list(const std::vector<Point<dim>> &points,
59 std::vector<double> &values,
60 const unsigned int component = 0) const override;
61 virtual Tensor<1, dim>
62 gradient(const Point<dim> &p,
63 const unsigned int component = 0) const override;
64 virtual void
66 std::vector<Tensor<1, dim>> &gradient) const override;
67 virtual void
68 gradient_list(const std::vector<Point<dim>> &points,
69 std::vector<Tensor<1, dim>> &gradients,
70 const unsigned int component = 0) const override;
71 virtual double
72 laplacian(const Point<dim> &p,
73 const unsigned int component = 0) const override;
74 virtual void
75 laplacian_list(const std::vector<Point<dim>> &points,
76 std::vector<double> &values,
77 const unsigned int component = 0) const override;
78 };
79
80
81
88 template <int dim>
89 class Q1WedgeFunction : public Function<dim>
90 {
91 public:
92 virtual double
93 value(const Point<dim> &p, const unsigned int component = 0) const override;
94
95 virtual void
96 value_list(const std::vector<Point<dim>> &points,
97 std::vector<double> &values,
98 const unsigned int component = 0) const override;
99
100 virtual void
101 vector_value_list(const std::vector<Point<dim>> &points,
102 std::vector<Vector<double>> &values) const override;
103
104 virtual Tensor<1, dim>
105 gradient(const Point<dim> &p,
106 const unsigned int component = 0) const override;
107
108 virtual void
109 gradient_list(const std::vector<Point<dim>> &points,
110 std::vector<Tensor<1, dim>> &gradients,
111 const unsigned int component = 0) const override;
112
113 virtual void
115 const std::vector<Point<dim>> &,
116 std::vector<std::vector<Tensor<1, dim>>> &) const override;
117
121 virtual double
122 laplacian(const Point<dim> &p,
123 const unsigned int component = 0) const override;
124
128 virtual void
129 laplacian_list(const std::vector<Point<dim>> &points,
130 std::vector<double> &values,
131 const unsigned int component = 0) const override;
132 };
133
134
135
150 template <int dim>
151 class PillowFunction : public Function<dim>
152 {
153 public:
158 PillowFunction(const double offset = 0.);
159
163 virtual double
164 value(const Point<dim> &p, const unsigned int component = 0) const override;
165
169 virtual void
170 value_list(const std::vector<Point<dim>> &points,
171 std::vector<double> &values,
172 const unsigned int component = 0) const override;
173
177 virtual Tensor<1, dim>
178 gradient(const Point<dim> &p,
179 const unsigned int component = 0) const override;
180
184 virtual void
185 gradient_list(const std::vector<Point<dim>> &points,
186 std::vector<Tensor<1, dim>> &gradients,
187 const unsigned int component = 0) const override;
188
192 virtual double
193 laplacian(const Point<dim> &p,
194 const unsigned int component = 0) const override;
195
199 virtual void
200 laplacian_list(const std::vector<Point<dim>> &points,
201 std::vector<double> &values,
202 const unsigned int component = 0) const override;
203
204 private:
205 const double offset;
206 };
207
208
209
217 template <int dim>
218 class CosineFunction : public Function<dim>
219 {
220 public:
225 CosineFunction(const unsigned int n_components = 1);
226
227 virtual double
228 value(const Point<dim> &p, const unsigned int component = 0) const override;
229
230 virtual void
231 value_list(const std::vector<Point<dim>> &points,
232 std::vector<double> &values,
233 const unsigned int component = 0) const override;
234
235 virtual void
236 vector_value_list(const std::vector<Point<dim>> &points,
237 std::vector<Vector<double>> &values) const override;
238
239 virtual Tensor<1, dim>
240 gradient(const Point<dim> &p,
241 const unsigned int component = 0) const override;
242
243 virtual void
244 gradient_list(const std::vector<Point<dim>> &points,
245 std::vector<Tensor<1, dim>> &gradients,
246 const unsigned int component = 0) const override;
247
248 virtual double
249 laplacian(const Point<dim> &p,
250 const unsigned int component = 0) const override;
251
252 virtual void
253 laplacian_list(const std::vector<Point<dim>> &points,
254 std::vector<double> &values,
255 const unsigned int component = 0) const override;
256
261 hessian(const Point<dim> &p,
262 const unsigned int component = 0) const override;
263
267 virtual void
268 hessian_list(const std::vector<Point<dim>> &points,
269 std::vector<SymmetricTensor<2, dim>> &hessians,
270 const unsigned int component = 0) const override;
271 };
272
273
274
285 template <int dim>
286 class CosineGradFunction : public Function<dim>
287 {
288 public:
293
294 virtual double
295 value(const Point<dim> &p, const unsigned int component) const override;
296 virtual void
297 vector_value(const Point<dim> &p, Vector<double> &values) const override;
298 virtual void
299 value_list(const std::vector<Point<dim>> &points,
300 std::vector<double> &values,
301 const unsigned int component) const override;
302
303 virtual void
304 vector_value_list(const std::vector<Point<dim>> &points,
305 std::vector<Vector<double>> &values) const override;
306
307 virtual Tensor<1, dim>
308 gradient(const Point<dim> &p, const unsigned int component) const override;
309
310 virtual void
311 gradient_list(const std::vector<Point<dim>> &points,
312 std::vector<Tensor<1, dim>> &gradients,
313 const unsigned int component) const override;
314
315 virtual void
317 const std::vector<Point<dim>> &points,
318 std::vector<std::vector<Tensor<1, dim>>> &gradients) const override;
319
320 virtual double
321 laplacian(const Point<dim> &p, const unsigned int component) const override;
322 };
323
324
325
331 template <int dim>
332 class ExpFunction : public Function<dim>
333 {
334 public:
338 virtual double
339 value(const Point<dim> &p, const unsigned int component = 0) const override;
340
344 virtual void
345 value_list(const std::vector<Point<dim>> &points,
346 std::vector<double> &values,
347 const unsigned int component = 0) const override;
348
352 virtual Tensor<1, dim>
353 gradient(const Point<dim> &p,
354 const unsigned int component = 0) const override;
355
359 virtual void
360 gradient_list(const std::vector<Point<dim>> &points,
361 std::vector<Tensor<1, dim>> &gradients,
362 const unsigned int component = 0) const override;
363
367 virtual double
368 laplacian(const Point<dim> &p,
369 const unsigned int component = 0) const override;
370
374 virtual void
375 laplacian_list(const std::vector<Point<dim>> &points,
376 std::vector<double> &values,
377 const unsigned int component = 0) const override;
378 };
379
380
381
410 {
411 public:
412 virtual double
413 value(const Point<2> &p, const unsigned int component = 0) const override;
414
415 virtual void
416 value_list(const std::vector<Point<2>> &points,
417 std::vector<double> &values,
418 const unsigned int component = 0) const override;
419
420 virtual void
421 vector_value_list(const std::vector<Point<2>> &points,
422 std::vector<Vector<double>> &values) const override;
423
424 virtual Tensor<1, 2>
425 gradient(const Point<2> &p,
426 const unsigned int component = 0) const override;
427
428 virtual void
429 gradient_list(const std::vector<Point<2>> &points,
430 std::vector<Tensor<1, 2>> &gradients,
431 const unsigned int component = 0) const override;
432
433 virtual void
435 const std::vector<Point<2>> &,
436 std::vector<std::vector<Tensor<1, 2>>> &) const override;
437
438 virtual double
439 laplacian(const Point<2> &p,
440 const unsigned int component = 0) const override;
441
442 virtual void
443 laplacian_list(const std::vector<Point<2>> &points,
444 std::vector<double> &values,
445 const unsigned int component = 0) const override;
446 };
447
448
449
459 {
460 public:
465 virtual double
466 value(const Point<2> &p, const unsigned int component) const override;
467
468 virtual void
469 value_list(const std::vector<Point<2>> &points,
470 std::vector<double> &values,
471 const unsigned int component) const override;
472
473 virtual void
474 vector_value_list(const std::vector<Point<2>> &points,
475 std::vector<Vector<double>> &values) const override;
476
477 virtual Tensor<1, 2>
478 gradient(const Point<2> &p, const unsigned int component) const override;
479
480 virtual void
481 gradient_list(const std::vector<Point<2>> &points,
482 std::vector<Tensor<1, 2>> &gradients,
483 const unsigned int component) const override;
484
485 virtual void
487 const std::vector<Point<2>> &,
488 std::vector<std::vector<Tensor<1, 2>>> &) const override;
489
490 virtual double
491 laplacian(const Point<2> &p, const unsigned int component) const override;
492
493 virtual void
494 laplacian_list(const std::vector<Point<2>> &points,
495 std::vector<double> &values,
496 const unsigned int component) const override;
497 };
498
499
500
506 template <int dim>
508 {
509 public:
510 virtual double
511 value(const Point<dim> &p, const unsigned int component = 0) const override;
512
513 virtual void
514 value_list(const std::vector<Point<dim>> &points,
515 std::vector<double> &values,
516 const unsigned int component = 0) const override;
517
518 virtual void
519 vector_value_list(const std::vector<Point<dim>> &points,
520 std::vector<Vector<double>> &values) const override;
521
522 virtual Tensor<1, dim>
523 gradient(const Point<dim> &p,
524 const unsigned int component = 0) const override;
525
526 virtual void
527 gradient_list(const std::vector<Point<dim>> &points,
528 std::vector<Tensor<1, dim>> &gradients,
529 const unsigned int component = 0) const override;
530
531 virtual void
533 const std::vector<Point<dim>> &,
534 std::vector<std::vector<Tensor<1, dim>>> &) const override;
535
536 virtual double
537 laplacian(const Point<dim> &p,
538 const unsigned int component = 0) const override;
539
540 virtual void
541 laplacian_list(const std::vector<Point<dim>> &points,
542 std::vector<double> &values,
543 const unsigned int component = 0) const override;
544 };
545
546
553 {
554 public:
555 virtual double
556 value(const Point<2> &p, const unsigned int component = 0) const override;
557
558 virtual void
559 value_list(const std::vector<Point<2>> &points,
560 std::vector<double> &values,
561 const unsigned int component = 0) const override;
562
563 virtual void
564 vector_value_list(const std::vector<Point<2>> &points,
565 std::vector<Vector<double>> &values) const override;
566
567 virtual Tensor<1, 2>
568 gradient(const Point<2> &p,
569 const unsigned int component = 0) const override;
570
571 virtual void
572 gradient_list(const std::vector<Point<2>> &points,
573 std::vector<Tensor<1, 2>> &gradients,
574 const unsigned int component = 0) const override;
575
576 virtual void
578 const std::vector<Point<2>> &,
579 std::vector<std::vector<Tensor<1, 2>>> &) const override;
580
581 virtual double
582 laplacian(const Point<2> &p,
583 const unsigned int component = 0) const override;
584
585 virtual void
586 laplacian_list(const std::vector<Point<2>> &points,
587 std::vector<double> &values,
588 const unsigned int component = 0) const override;
589 };
590
591
592
607 template <int dim>
608 class JumpFunction : public Function<dim>
609 {
610 public:
615 JumpFunction(const Point<dim> &direction, const double steepness);
616
620 virtual double
621 value(const Point<dim> &p, const unsigned int component = 0) const override;
622
626 virtual void
627 value_list(const std::vector<Point<dim>> &points,
628 std::vector<double> &values,
629 const unsigned int component = 0) const override;
630
634 virtual Tensor<1, dim>
635 gradient(const Point<dim> &p,
636 const unsigned int component = 0) const override;
637
641 virtual void
642 gradient_list(const std::vector<Point<dim>> &points,
643 std::vector<Tensor<1, dim>> &gradients,
644 const unsigned int component = 0) const override;
645
649 virtual double
650 laplacian(const Point<dim> &p,
651 const unsigned int component = 0) const override;
652
656 virtual void
657 laplacian_list(const std::vector<Point<dim>> &points,
658 std::vector<double> &values,
659 const unsigned int component = 0) const override;
660
667 virtual std::size_t
668 memory_consumption() const override;
669
670 protected:
675
679 const double steepness;
680
684 double angle;
685
689 double sine;
690
694 double cosine;
695 };
696
697
698
710 template <int dim>
711 class FourierCosineFunction : public Function<dim>
712 {
713 public:
719
726 virtual double
727 value(const Point<dim> &p, const unsigned int component = 0) const override;
728
733 virtual Tensor<1, dim>
734 gradient(const Point<dim> &p,
735 const unsigned int component = 0) const override;
736
740 virtual double
741 laplacian(const Point<dim> &p,
742 const unsigned int component = 0) const override;
743
744 private:
749 };
750
751
752
764 template <int dim>
765 class FourierSineFunction : public Function<dim>
766 {
767 public:
773
780 virtual double
781 value(const Point<dim> &p, const unsigned int component = 0) const override;
782
787 virtual Tensor<1, dim>
788 gradient(const Point<dim> &p,
789 const unsigned int component = 0) const override;
790
794 virtual double
795 laplacian(const Point<dim> &p,
796 const unsigned int component = 0) const override;
797
798 private:
803 };
804
805
814 template <int dim>
815 class FourierSineSum : public Function<dim>
816 {
817 public:
823 const std::vector<double> &weights);
824
831 virtual double
832 value(const Point<dim> &p, const unsigned int component = 0) const override;
833
838 virtual Tensor<1, dim>
839 gradient(const Point<dim> &p,
840 const unsigned int component = 0) const override;
841
845 virtual double
846 laplacian(const Point<dim> &p,
847 const unsigned int component = 0) const override;
848
849 private:
853 const std::vector<Point<dim>> fourier_coefficients;
854 const std::vector<double> weights;
855 };
856
857
858
868 template <int dim>
869 class FourierCosineSum : public Function<dim>
870 {
871 public:
877 const std::vector<double> &weights);
878
885 virtual double
886 value(const Point<dim> &p, const unsigned int component = 0) const override;
887
892 virtual Tensor<1, dim>
893 gradient(const Point<dim> &p,
894 const unsigned int component = 0) const override;
895
899 virtual double
900 laplacian(const Point<dim> &p,
901 const unsigned int component = 0) const override;
902
903 private:
907 const std::vector<Point<dim>> fourier_coefficients;
908 const std::vector<double> weights;
909 };
910
911
923 template <int dim>
924 class CutOffFunctionBase : public Function<dim>
925 {
926 public:
931 static const unsigned int no_component = numbers::invalid_unsigned_int;
932
949 const double radius = 1.,
950 const Point<dim> center = Point<dim>(),
951 const unsigned int n_components = 1,
952 const unsigned int select = CutOffFunctionBase<dim>::no_component,
953 const bool integrate_to_one = false,
954 const double unitary_integral_value = 1.0);
955
959 virtual ~CutOffFunctionBase() = default;
960
964 virtual void
965 set_center(const Point<dim> &p);
966
970 virtual void
971 set_radius(const double r);
972
976 const Point<dim> &
977 get_center() const;
978
982 double
983 get_radius() const;
984
988 bool
989 integrates_to_one() const;
990
991 protected:
996
1000 double radius;
1001
1006 const unsigned int selected;
1007
1012
1018
1023 };
1024
1025
1035 template <int dim>
1037 {
1038 public:
1049 double radius = 1.0,
1050 const Point<dim> &center = Point<dim>(),
1051 const unsigned int n_components = 1,
1052 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1053 const bool integrate_to_one = false);
1054
1059 template <template <int> class CutOffFunctionBaseType>
1060 void
1062
1066 virtual void
1067 set_center(const Point<dim> &center) override;
1068
1072 virtual void
1073 set_radius(const double radius) override;
1074
1078 virtual double
1079 value(const Point<dim> &p, const unsigned int component = 0) const override;
1080
1084 virtual Tensor<1, dim>
1085 gradient(const Point<dim> &p,
1086 const unsigned int component = 0) const override;
1087
1088 private:
1089 std::array<std::unique_ptr<CutOffFunctionBase<1>>, dim> base;
1090
1092 };
1093
1094
1095
1104 template <int dim>
1106 {
1107 public:
1115 const double radius = 1.,
1116 const Point<dim> = Point<dim>(),
1117 const unsigned int n_components = 1,
1118 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1119 const bool integrate_to_one = false);
1120
1124 virtual double
1125 value(const Point<dim> &p, const unsigned int component = 0) const override;
1126
1130 virtual void
1131 value_list(const std::vector<Point<dim>> &points,
1132 std::vector<double> &values,
1133 const unsigned int component = 0) const override;
1134
1138 virtual void
1139 vector_value_list(const std::vector<Point<dim>> &points,
1140 std::vector<Vector<double>> &values) const override;
1141 };
1142
1143
1152 template <int dim>
1154 {
1155 public:
1163 const double radius = 1.,
1164 const Point<dim> = Point<dim>(),
1165 const unsigned int n_components = 1,
1166 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1167 const bool integrate_to_one = false);
1168
1172 virtual double
1173 value(const Point<dim> &p, const unsigned int component = 0) const override;
1174
1178 virtual void
1179 value_list(const std::vector<Point<dim>> &points,
1180 std::vector<double> &values,
1181 const unsigned int component = 0) const override;
1182
1186 virtual void
1187 vector_value_list(const std::vector<Point<dim>> &points,
1188 std::vector<Vector<double>> &values) const override;
1189 };
1190
1191
1204 template <int dim>
1206 {
1207 public:
1212 const double radius = 1.,
1213 const Point<dim> = Point<dim>(),
1214 const unsigned int n_components = 1,
1215 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1216 bool integrate_to_one = false);
1217
1221 virtual double
1222 value(const Point<dim> &p, const unsigned int component = 0) const override;
1223
1227 virtual void
1228 value_list(const std::vector<Point<dim>> &points,
1229 std::vector<double> &values,
1230 const unsigned int component = 0) const override;
1231
1235 virtual void
1236 vector_value_list(const std::vector<Point<dim>> &points,
1237 std::vector<Vector<double>> &values) const override;
1238
1242 virtual Tensor<1, dim>
1243 gradient(const Point<dim> &p,
1244 const unsigned int component = 0) const override;
1245 };
1246
1247
1257 template <int dim>
1259 {
1260 public:
1268 const double radius = 1.,
1269 const Point<dim> = Point<dim>(),
1270 const unsigned int n_components = 1,
1271 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1272 bool integrate_to_one = false);
1273
1277 virtual double
1278 value(const Point<dim> &p, const unsigned int component = 0) const override;
1279
1283 virtual void
1284 value_list(const std::vector<Point<dim>> &points,
1285 std::vector<double> &values,
1286 const unsigned int component = 0) const override;
1287
1291 virtual void
1292 vector_value_list(const std::vector<Point<dim>> &points,
1293 std::vector<Vector<double>> &values) const override;
1294
1298 virtual Tensor<1, dim>
1299 gradient(const Point<dim> &p,
1300 const unsigned int component = 0) const override;
1301 };
1302
1303
1304
1318 template <int dim, typename Number = double>
1319 class Monomial : public Function<dim, Number>
1320 {
1321 public:
1329 const unsigned int n_components = 1);
1330
1334 virtual Number
1335 value(const Point<dim> &p, const unsigned int component = 0) const override;
1336
1343 virtual void
1344 vector_value(const Point<dim> &p, Vector<Number> &values) const override;
1345
1349 virtual void
1350 value_list(const std::vector<Point<dim>> &points,
1351 std::vector<Number> &values,
1352 const unsigned int component = 0) const override;
1353
1358 gradient(const Point<dim> &p,
1359 const unsigned int component = 0) const override;
1360
1361 private:
1366 };
1367
1368
1369
1451 template <int dim>
1453 {
1454 public:
1473 const std::array<std::vector<double>, dim> &coordinate_values,
1475
1493 std::array<std::vector<double>, dim> &&coordinate_values,
1495
1506 virtual double
1507 value(const Point<dim> &p, const unsigned int component = 0) const override;
1508
1520 virtual Tensor<1, dim>
1521 gradient(const Point<dim> &p,
1522 const unsigned int component = 0) const override;
1523
1527 virtual std::size_t
1528 memory_consumption() const override;
1529
1533 const Table<dim, double> &
1534 get_data() const;
1535
1536 protected:
1541 table_index_of_point(const Point<dim> &p) const;
1542
1546 const std::array<std::vector<double>, dim> coordinate_values;
1547
1552 };
1553
1554
1598 template <int dim>
1600 {
1601 public:
1617 const std::array<std::pair<double, double>, dim> &interval_endpoints,
1618 const std::array<unsigned int, dim> &n_subintervals,
1620
1638 std::array<std::pair<double, double>, dim> &&interval_endpoints,
1639 std::array<unsigned int, dim> &&n_subintervals,
1641
1652 virtual double
1653 value(const Point<dim> &p, const unsigned int component = 0) const override;
1654
1666 virtual Tensor<1, dim>
1667 gradient(const Point<dim> &p,
1668 const unsigned int component = 0) const override;
1669
1673 virtual std::size_t
1674 memory_consumption() const override;
1675
1679 const Table<dim, double> &
1680 get_data() const;
1681
1682 private:
1686 const std::array<std::pair<double, double>, dim> interval_endpoints;
1687
1691 const std::array<unsigned int, dim> n_subintervals;
1692
1697 };
1698
1699
1712 template <int dim>
1713 class Polynomial : public Function<dim>
1714 {
1715 public:
1726 const std::vector<double> &coefficients);
1727
1731 virtual double
1732 value(const Point<dim> &p, const unsigned int component = 0) const override;
1733
1734
1738 virtual void
1739 value_list(const std::vector<Point<dim>> &points,
1740 std::vector<double> &values,
1741 const unsigned int component = 0) const override;
1742
1746 virtual Tensor<1, dim>
1747 gradient(const Point<dim> &p,
1748 const unsigned int component = 0) const override;
1749
1753 virtual std::size_t
1754 memory_consumption() const override;
1755
1756 private:
1761
1765 const std::vector<double> coefficients;
1766 };
1767
1768
1811 template <int dim>
1812 class RayleighKotheVortex : public Function<dim>
1813 {
1814 public:
1818 RayleighKotheVortex(const double T = 1.0);
1819
1823 virtual void
1824 vector_value(const Point<dim> &point,
1825 Vector<double> &values) const override;
1826
1827 private:
1831 const double T;
1832 };
1833
1834#ifndef DOXYGEN
1835
1836
1837
1838 // Template definitions
1839 template <int dim>
1840 template <template <int> class CutOffFunctionBaseType>
1841 void
1843 {
1844 initialized = true;
1845 static_assert(
1846 std::is_base_of_v<CutOffFunctionBase<1>, CutOffFunctionBaseType<1>>,
1847 "You can only construct a CutOffFunctionTensorProduct function from "
1848 "a class derived from CutOffFunctionBase.");
1849
1850 for (unsigned int i = 0; i < dim; ++i)
1851 base[i].reset(new CutOffFunctionBaseType<1>(this->radius,
1852 Point<1>(this->center[i]),
1853 this->n_components,
1854 this->selected,
1855 this->integrate_to_one));
1856 }
1857
1858
1859
1860#endif
1861
1862} // namespace Functions
1864
1865#endif
const unsigned int n_components
Definition function.h:163
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
CosineFunction(const unsigned int n_components=1)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual SymmetricTensor< 2, dim > hessian(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual void hessian_list(const std::vector< Point< dim > > &points, std::vector< SymmetricTensor< 2, dim > > &hessians, const unsigned int component=0) const override
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component) const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component) const override
virtual double value(const Point< dim > &p, const unsigned int component) const override
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component) const override
virtual void vector_gradient_list(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const override
virtual void set_radius(const double r)
virtual ~CutOffFunctionBase()=default
virtual void set_center(const Point< dim > &p)
const Point< dim > & get_center() const
static const unsigned int no_component
CutOffFunctionBase(const double radius=1., const Point< dim > center=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false, const double unitary_integral_value=1.0)
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
CutOffFunctionC1(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, bool integrate_to_one=false)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
CutOffFunctionCinfty(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, bool integrate_to_one=false)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
CutOffFunctionLinfty(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false)
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
std::array< std::unique_ptr< CutOffFunctionBase< 1 > >, dim > base
CutOffFunctionTensorProduct(double radius=1.0, const Point< dim > &center=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void set_radius(const double radius) override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void set_center(const Point< dim > &center) override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
CutOffFunctionW1(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false)
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
FourierCosineFunction(const Tensor< 1, dim > &fourier_coefficients)
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
const Tensor< 1, dim > fourier_coefficients
const std::vector< Point< dim > > fourier_coefficients
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
FourierCosineSum(const std::vector< Point< dim > > &fourier_coefficients, const std::vector< double > &weights)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
const std::vector< double > weights
FourierSineFunction(const Tensor< 1, dim > &fourier_coefficients)
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
const Tensor< 1, dim > fourier_coefficients
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
const std::vector< Point< dim > > fourier_coefficients
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
FourierSineSum(const std::vector< Point< dim > > &fourier_coefficients, const std::vector< double > &weights)
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
const std::vector< double > weights
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
InterpolatedTensorProductGridData(const std::array< std::vector< double >, dim > &coordinate_values, const Table< dim, double > &data_values)
const Table< dim, double > & get_data() const
TableIndices< dim > table_index_of_point(const Point< dim > &p) const
virtual std::size_t memory_consumption() const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
const std::array< std::vector< double >, dim > coordinate_values
const Table< dim, double > & get_data() const
const std::array< unsigned int, dim > n_subintervals
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
InterpolatedUniformGridData(const std::array< std::pair< double, double >, dim > &interval_endpoints, const std::array< unsigned int, dim > &n_subintervals, const Table< dim, double > &data_values)
const Table< dim, double > data_values
const std::array< std::pair< double, double >, dim > interval_endpoints
virtual std::size_t memory_consumption() const override
const Point< dim > direction
virtual std::size_t memory_consumption() const override
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
JumpFunction(const Point< dim > &direction, const double steepness)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void laplacian_list(const std::vector< Point< 2 > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< 2 > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual Tensor< 1, 2 > gradient(const Point< 2 > &p, const unsigned int component=0) const override
virtual double laplacian(const Point< 2 > &p, const unsigned int component=0) const override
virtual double value(const Point< 2 > &p, const unsigned int component=0) const override
virtual void vector_gradient_list(const std::vector< Point< 2 > > &, std::vector< std::vector< Tensor< 1, 2 > > > &) const override
virtual void gradient_list(const std::vector< Point< 2 > > &points, std::vector< Tensor< 1, 2 > > &gradients, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< 2 > > &points, std::vector< Vector< double > > &values) const override
virtual void vector_gradient_list(const std::vector< Point< 2 > > &, std::vector< std::vector< Tensor< 1, 2 > > > &) const override
virtual double laplacian(const Point< 2 > &p, const unsigned int component) const override
virtual void value_list(const std::vector< Point< 2 > > &points, std::vector< double > &values, const unsigned int component) const override
virtual void vector_value_list(const std::vector< Point< 2 > > &points, std::vector< Vector< double > > &values) const override
virtual void laplacian_list(const std::vector< Point< 2 > > &points, std::vector< double > &values, const unsigned int component) const override
virtual void gradient_list(const std::vector< Point< 2 > > &points, std::vector< Tensor< 1, 2 > > &gradients, const unsigned int component) const override
virtual Tensor< 1, 2 > gradient(const Point< 2 > &p, const unsigned int component) const override
virtual double value(const Point< 2 > &p, const unsigned int component) const override
virtual Number value(const Point< dim > &p, const unsigned int component=0) const override
virtual Tensor< 1, dim, Number > gradient(const Point< dim > &p, const unsigned int component=0) const override
const Tensor< 1, dim, Number > exponents
Monomial(const Tensor< 1, dim, Number > &exponents, const unsigned int n_components=1)
virtual void vector_value(const Point< dim > &p, Vector< Number > &values) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< Number > &values, const unsigned int component=0) const override
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
PillowFunction(const double offset=0.)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
const Table< 2, double > exponents
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual std::size_t memory_consumption() const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
const std::vector< double > coefficients
Polynomial(const Table< 2, double > &exponents, const std::vector< double > &coefficients)
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_gradient_list(const std::vector< Point< dim > > &, std::vector< std::vector< Tensor< 1, dim > > > &) const override
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
RayleighKotheVortex(const double T=1.0)
virtual void vector_value(const Point< dim > &point, Vector< double > &values) const override
virtual void gradient_list(const std::vector< Point< 2 > > &points, std::vector< Tensor< 1, 2 > > &gradients, const unsigned int component=0) const override
virtual void laplacian_list(const std::vector< Point< 2 > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< 2 > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void vector_gradient_list(const std::vector< Point< 2 > > &, std::vector< std::vector< Tensor< 1, 2 > > > &) const override
virtual double value(const Point< 2 > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< 2 > > &points, std::vector< Vector< double > > &values) const override
virtual double laplacian(const Point< 2 > &p, const unsigned int component=0) const override
virtual Tensor< 1, 2 > gradient(const Point< 2 > &p, const unsigned int component=0) const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void vector_gradient_list(const std::vector< Point< dim > > &, std::vector< std::vector< Tensor< 1, dim > > > &) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const override
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void vector_gradient(const Point< dim > &p, std::vector< Tensor< 1, dim > > &gradient) const override
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
Point< 3 > center
static const unsigned int invalid_unsigned_int
Definition types.h:220