Reference documentation for deal.II version 9.5.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// Copyright (C) 1999 - 2023 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_function_lib_h
17#define dealii_function_lib_h
18
19
20#include <deal.II/base/config.h>
21
23#include <deal.II/base/point.h>
25#include <deal.II/base/table.h>
26
27#include <array>
28
30
39namespace Functions
40{
50 template <int dim>
51 class SquareFunction : public Function<dim>
52 {
53 public:
54 virtual double
55 value(const Point<dim> &p, const unsigned int component = 0) const override;
56 virtual void
57 vector_value(const Point<dim> &p, Vector<double> &values) const override;
58 virtual void
59 value_list(const std::vector<Point<dim>> &points,
60 std::vector<double> & values,
61 const unsigned int component = 0) const override;
62 virtual Tensor<1, dim>
63 gradient(const Point<dim> & p,
64 const unsigned int component = 0) const override;
65 virtual void
67 std::vector<Tensor<1, dim>> &gradient) const override;
68 virtual void
69 gradient_list(const std::vector<Point<dim>> &points,
70 std::vector<Tensor<1, dim>> & gradients,
71 const unsigned int component = 0) const override;
72 virtual double
73 laplacian(const Point<dim> & p,
74 const unsigned int component = 0) const override;
75 virtual void
76 laplacian_list(const std::vector<Point<dim>> &points,
77 std::vector<double> & values,
78 const unsigned int component = 0) const override;
79 };
80
81
82
89 template <int dim>
90 class Q1WedgeFunction : public Function<dim>
91 {
92 public:
93 virtual double
94 value(const Point<dim> &p, const unsigned int component = 0) const override;
95
96 virtual void
97 value_list(const std::vector<Point<dim>> &points,
98 std::vector<double> & values,
99 const unsigned int component = 0) const override;
100
101 virtual void
102 vector_value_list(const std::vector<Point<dim>> &points,
103 std::vector<Vector<double>> & values) const override;
104
105 virtual Tensor<1, dim>
106 gradient(const Point<dim> & p,
107 const unsigned int component = 0) const override;
108
109 virtual void
110 gradient_list(const std::vector<Point<dim>> &points,
111 std::vector<Tensor<1, dim>> & gradients,
112 const unsigned int component = 0) const override;
113
114 virtual void
116 const std::vector<Point<dim>> &,
117 std::vector<std::vector<Tensor<1, dim>>> &) const override;
118
122 virtual double
123 laplacian(const Point<dim> & p,
124 const unsigned int component = 0) const override;
125
129 virtual void
130 laplacian_list(const std::vector<Point<dim>> &points,
131 std::vector<double> & values,
132 const unsigned int component = 0) const override;
133 };
134
135
136
151 template <int dim>
152 class PillowFunction : public Function<dim>
153 {
154 public:
159 PillowFunction(const double offset = 0.);
160
164 virtual double
165 value(const Point<dim> &p, const unsigned int component = 0) const override;
166
170 virtual void
171 value_list(const std::vector<Point<dim>> &points,
172 std::vector<double> & values,
173 const unsigned int component = 0) const override;
174
178 virtual Tensor<1, dim>
179 gradient(const Point<dim> & p,
180 const unsigned int component = 0) const override;
181
185 virtual void
186 gradient_list(const std::vector<Point<dim>> &points,
187 std::vector<Tensor<1, dim>> & gradients,
188 const unsigned int component = 0) const override;
189
193 virtual double
194 laplacian(const Point<dim> & p,
195 const unsigned int component = 0) const override;
196
200 virtual void
201 laplacian_list(const std::vector<Point<dim>> &points,
202 std::vector<double> & values,
203 const unsigned int component = 0) const override;
204
205 private:
206 const double offset;
207 };
208
209
210
218 template <int dim>
219 class CosineFunction : public Function<dim>
220 {
221 public:
226 CosineFunction(const unsigned int n_components = 1);
227
228 virtual double
229 value(const Point<dim> &p, const unsigned int component = 0) const override;
230
231 virtual void
232 value_list(const std::vector<Point<dim>> &points,
233 std::vector<double> & values,
234 const unsigned int component = 0) const override;
235
236 virtual void
237 vector_value_list(const std::vector<Point<dim>> &points,
238 std::vector<Vector<double>> & values) const override;
239
240 virtual Tensor<1, dim>
241 gradient(const Point<dim> & p,
242 const unsigned int component = 0) const override;
243
244 virtual void
245 gradient_list(const std::vector<Point<dim>> &points,
246 std::vector<Tensor<1, dim>> & gradients,
247 const unsigned int component = 0) const override;
248
249 virtual double
250 laplacian(const Point<dim> & p,
251 const unsigned int component = 0) const override;
252
253 virtual void
254 laplacian_list(const std::vector<Point<dim>> &points,
255 std::vector<double> & values,
256 const unsigned int component = 0) const override;
257
262 hessian(const Point<dim> & p,
263 const unsigned int component = 0) const override;
264
268 virtual void
269 hessian_list(const std::vector<Point<dim>> & points,
270 std::vector<SymmetricTensor<2, dim>> &hessians,
271 const unsigned int component = 0) const override;
272 };
273
274
275
286 template <int dim>
287 class CosineGradFunction : public Function<dim>
288 {
289 public:
294
295 virtual double
296 value(const Point<dim> &p, const unsigned int component) const override;
297 virtual void
298 vector_value(const Point<dim> &p, Vector<double> &values) const override;
299 virtual void
300 value_list(const std::vector<Point<dim>> &points,
301 std::vector<double> & values,
302 const unsigned int component) const override;
303
304 virtual void
305 vector_value_list(const std::vector<Point<dim>> &points,
306 std::vector<Vector<double>> & values) const override;
307
308 virtual Tensor<1, dim>
309 gradient(const Point<dim> &p, const unsigned int component) const override;
310
311 virtual void
312 gradient_list(const std::vector<Point<dim>> &points,
313 std::vector<Tensor<1, dim>> & gradients,
314 const unsigned int component) const override;
315
316 virtual void
318 const std::vector<Point<dim>> & points,
319 std::vector<std::vector<Tensor<1, dim>>> &gradients) const override;
320
321 virtual double
322 laplacian(const Point<dim> &p, const unsigned int component) const override;
323 };
324
325
326
332 template <int dim>
333 class ExpFunction : public Function<dim>
334 {
335 public:
339 virtual double
340 value(const Point<dim> &p, const unsigned int component = 0) const override;
341
345 virtual void
346 value_list(const std::vector<Point<dim>> &points,
347 std::vector<double> & values,
348 const unsigned int component = 0) const override;
349
353 virtual Tensor<1, dim>
354 gradient(const Point<dim> & p,
355 const unsigned int component = 0) const override;
356
360 virtual void
361 gradient_list(const std::vector<Point<dim>> &points,
362 std::vector<Tensor<1, dim>> & gradients,
363 const unsigned int component = 0) const override;
364
368 virtual double
369 laplacian(const Point<dim> & p,
370 const unsigned int component = 0) const override;
371
375 virtual void
376 laplacian_list(const std::vector<Point<dim>> &points,
377 std::vector<double> & values,
378 const unsigned int component = 0) const override;
379 };
380
381
382
411 {
412 public:
413 virtual double
414 value(const Point<2> &p, const unsigned int component = 0) const override;
415
416 virtual void
417 value_list(const std::vector<Point<2>> &points,
418 std::vector<double> & values,
419 const unsigned int component = 0) const override;
420
421 virtual void
422 vector_value_list(const std::vector<Point<2>> &points,
423 std::vector<Vector<double>> &values) const override;
424
425 virtual Tensor<1, 2>
426 gradient(const Point<2> & p,
427 const unsigned int component = 0) const override;
428
429 virtual void
430 gradient_list(const std::vector<Point<2>> &points,
431 std::vector<Tensor<1, 2>> & gradients,
432 const unsigned int component = 0) const override;
433
434 virtual void
436 const std::vector<Point<2>> &,
437 std::vector<std::vector<Tensor<1, 2>>> &) const override;
438
439 virtual double
440 laplacian(const Point<2> & p,
441 const unsigned int component = 0) const override;
442
443 virtual void
444 laplacian_list(const std::vector<Point<2>> &points,
445 std::vector<double> & values,
446 const unsigned int component = 0) const override;
447 };
448
449
450
460 {
461 public:
466 virtual double
467 value(const Point<2> &p, const unsigned int component) const override;
468
469 virtual void
470 value_list(const std::vector<Point<2>> &points,
471 std::vector<double> & values,
472 const unsigned int component) const override;
473
474 virtual void
475 vector_value_list(const std::vector<Point<2>> &points,
476 std::vector<Vector<double>> &values) const override;
477
478 virtual Tensor<1, 2>
479 gradient(const Point<2> &p, const unsigned int component) const override;
480
481 virtual void
482 gradient_list(const std::vector<Point<2>> &points,
483 std::vector<Tensor<1, 2>> & gradients,
484 const unsigned int component) const override;
485
486 virtual void
488 const std::vector<Point<2>> &,
489 std::vector<std::vector<Tensor<1, 2>>> &) const override;
490
491 virtual double
492 laplacian(const Point<2> &p, const unsigned int component) const override;
493
494 virtual void
495 laplacian_list(const std::vector<Point<2>> &points,
496 std::vector<double> & values,
497 const unsigned int component) const override;
498 };
499
500
501
507 template <int dim>
509 {
510 public:
511 virtual double
512 value(const Point<dim> &p, const unsigned int component = 0) const override;
513
514 virtual void
515 value_list(const std::vector<Point<dim>> &points,
516 std::vector<double> & values,
517 const unsigned int component = 0) const override;
518
519 virtual void
520 vector_value_list(const std::vector<Point<dim>> &points,
521 std::vector<Vector<double>> & values) const override;
522
523 virtual Tensor<1, dim>
524 gradient(const Point<dim> & p,
525 const unsigned int component = 0) const override;
526
527 virtual void
528 gradient_list(const std::vector<Point<dim>> &points,
529 std::vector<Tensor<1, dim>> & gradients,
530 const unsigned int component = 0) const override;
531
532 virtual void
534 const std::vector<Point<dim>> &,
535 std::vector<std::vector<Tensor<1, dim>>> &) const override;
536
537 virtual double
538 laplacian(const Point<dim> & p,
539 const unsigned int component = 0) const override;
540
541 virtual void
542 laplacian_list(const std::vector<Point<dim>> &points,
543 std::vector<double> & values,
544 const unsigned int component = 0) const override;
545 };
546
547
554 {
555 public:
556 virtual double
557 value(const Point<2> &p, const unsigned int component = 0) const override;
558
559 virtual void
560 value_list(const std::vector<Point<2>> &points,
561 std::vector<double> & values,
562 const unsigned int component = 0) const override;
563
564 virtual void
565 vector_value_list(const std::vector<Point<2>> &points,
566 std::vector<Vector<double>> &values) const override;
567
568 virtual Tensor<1, 2>
569 gradient(const Point<2> & p,
570 const unsigned int component = 0) const override;
571
572 virtual void
573 gradient_list(const std::vector<Point<2>> &points,
574 std::vector<Tensor<1, 2>> & gradients,
575 const unsigned int component = 0) const override;
576
577 virtual void
579 const std::vector<Point<2>> &,
580 std::vector<std::vector<Tensor<1, 2>>> &) const override;
581
582 virtual double
583 laplacian(const Point<2> & p,
584 const unsigned int component = 0) const override;
585
586 virtual void
587 laplacian_list(const std::vector<Point<2>> &points,
588 std::vector<double> & values,
589 const unsigned int component = 0) const override;
590 };
591
592
593
608 template <int dim>
609 class JumpFunction : public Function<dim>
610 {
611 public:
616 JumpFunction(const Point<dim> &direction, const double steepness);
617
621 virtual double
622 value(const Point<dim> &p, const unsigned int component = 0) const override;
623
627 virtual void
628 value_list(const std::vector<Point<dim>> &points,
629 std::vector<double> & values,
630 const unsigned int component = 0) const override;
631
635 virtual Tensor<1, dim>
636 gradient(const Point<dim> & p,
637 const unsigned int component = 0) const override;
638
642 virtual void
643 gradient_list(const std::vector<Point<dim>> &points,
644 std::vector<Tensor<1, dim>> & gradients,
645 const unsigned int component = 0) const override;
646
650 virtual double
651 laplacian(const Point<dim> & p,
652 const unsigned int component = 0) const override;
653
657 virtual void
658 laplacian_list(const std::vector<Point<dim>> &points,
659 std::vector<double> & values,
660 const unsigned int component = 0) const override;
661
668 virtual std::size_t
669 memory_consumption() const override;
670
671 protected:
676
680 const double steepness;
681
685 double angle;
686
690 double sine;
691
695 double cosine;
696 };
697
698
699
711 template <int dim>
712 class FourierCosineFunction : public Function<dim>
713 {
714 public:
720
727 virtual double
728 value(const Point<dim> &p, const unsigned int component = 0) const override;
729
734 virtual Tensor<1, dim>
735 gradient(const Point<dim> & p,
736 const unsigned int component = 0) const override;
737
741 virtual double
742 laplacian(const Point<dim> & p,
743 const unsigned int component = 0) const override;
744
745 private:
750 };
751
752
753
765 template <int dim>
766 class FourierSineFunction : public Function<dim>
767 {
768 public:
774
781 virtual double
782 value(const Point<dim> &p, const unsigned int component = 0) const override;
783
788 virtual Tensor<1, dim>
789 gradient(const Point<dim> & p,
790 const unsigned int component = 0) const override;
791
795 virtual double
796 laplacian(const Point<dim> & p,
797 const unsigned int component = 0) const override;
798
799 private:
804 };
805
806
815 template <int dim>
816 class FourierSineSum : public Function<dim>
817 {
818 public:
824 const std::vector<double> & weights);
825
832 virtual double
833 value(const Point<dim> &p, const unsigned int component = 0) const override;
834
839 virtual Tensor<1, dim>
840 gradient(const Point<dim> & p,
841 const unsigned int component = 0) const override;
842
846 virtual double
847 laplacian(const Point<dim> & p,
848 const unsigned int component = 0) const override;
849
850 private:
854 const std::vector<Point<dim>> fourier_coefficients;
855 const std::vector<double> weights;
856 };
857
858
859
869 template <int dim>
870 class FourierCosineSum : public Function<dim>
871 {
872 public:
878 const std::vector<double> & weights);
879
886 virtual double
887 value(const Point<dim> &p, const unsigned int component = 0) const override;
888
893 virtual Tensor<1, dim>
894 gradient(const Point<dim> & p,
895 const unsigned int component = 0) const override;
896
900 virtual double
901 laplacian(const Point<dim> & p,
902 const unsigned int component = 0) const override;
903
904 private:
908 const std::vector<Point<dim>> fourier_coefficients;
909 const std::vector<double> weights;
910 };
911
912
924 template <int dim>
925 class CutOffFunctionBase : public Function<dim>
926 {
927 public:
932 static const unsigned int no_component = numbers::invalid_unsigned_int;
933
950 const double radius = 1.,
951 const Point<dim> center = Point<dim>(),
952 const unsigned int n_components = 1,
953 const unsigned int select = CutOffFunctionBase<dim>::no_component,
954 const bool integrate_to_one = false,
955 const double unitary_integral_value = 1.0);
956
960 virtual ~CutOffFunctionBase() = default;
961
965 virtual void
966 set_center(const Point<dim> &p);
967
971 virtual void
972 set_radius(const double r);
973
977 const Point<dim> &
978 get_center() const;
979
983 double
984 get_radius() const;
985
989 bool
990 integrates_to_one() const;
991
992 protected:
997
1001 double radius;
1002
1007 const unsigned int selected;
1008
1013
1019
1024 };
1025
1026
1036 template <int dim>
1038 {
1039 public:
1050 double radius = 1.0,
1051 const Point<dim> & center = Point<dim>(),
1052 const unsigned int n_components = 1,
1053 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1054 const bool integrate_to_one = false);
1055
1060 template <template <int> class CutOffFunctionBaseType>
1061 void
1063
1067 virtual void
1068 set_center(const Point<dim> &center) override;
1069
1073 virtual void
1074 set_radius(const double radius) override;
1075
1079 virtual double
1080 value(const Point<dim> &p, const unsigned int component = 0) const override;
1081
1085 virtual Tensor<1, dim>
1086 gradient(const Point<dim> & p,
1087 const unsigned int component = 0) const override;
1088
1089 private:
1090 std::array<std::unique_ptr<CutOffFunctionBase<1>>, dim> base;
1091
1093 };
1094
1095
1096
1105 template <int dim>
1107 {
1108 public:
1116 const double radius = 1.,
1117 const Point<dim> = Point<dim>(),
1118 const unsigned int n_components = 1,
1119 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1120 const bool integrate_to_one = false);
1121
1125 virtual double
1126 value(const Point<dim> &p, const unsigned int component = 0) const override;
1127
1131 virtual void
1132 value_list(const std::vector<Point<dim>> &points,
1133 std::vector<double> & values,
1134 const unsigned int component = 0) const override;
1135
1139 virtual void
1140 vector_value_list(const std::vector<Point<dim>> &points,
1141 std::vector<Vector<double>> & values) const override;
1142 };
1143
1144
1153 template <int dim>
1155 {
1156 public:
1164 const double radius = 1.,
1165 const Point<dim> = Point<dim>(),
1166 const unsigned int n_components = 1,
1167 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1168 const bool integrate_to_one = false);
1169
1173 virtual double
1174 value(const Point<dim> &p, const unsigned int component = 0) const override;
1175
1179 virtual void
1180 value_list(const std::vector<Point<dim>> &points,
1181 std::vector<double> & values,
1182 const unsigned int component = 0) const override;
1183
1187 virtual void
1188 vector_value_list(const std::vector<Point<dim>> &points,
1189 std::vector<Vector<double>> & values) const override;
1190 };
1191
1192
1205 template <int dim>
1207 {
1208 public:
1213 const double radius = 1.,
1214 const Point<dim> = Point<dim>(),
1215 const unsigned int n_components = 1,
1216 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1217 bool integrate_to_one = false);
1218
1222 virtual double
1223 value(const Point<dim> &p, const unsigned int component = 0) const override;
1224
1228 virtual void
1229 value_list(const std::vector<Point<dim>> &points,
1230 std::vector<double> & values,
1231 const unsigned int component = 0) const override;
1232
1236 virtual void
1237 vector_value_list(const std::vector<Point<dim>> &points,
1238 std::vector<Vector<double>> & values) const override;
1239
1243 virtual Tensor<1, dim>
1244 gradient(const Point<dim> & p,
1245 const unsigned int component = 0) const override;
1246 };
1247
1248
1258 template <int dim>
1260 {
1261 public:
1269 const double radius = 1.,
1270 const Point<dim> = Point<dim>(),
1271 const unsigned int n_components = 1,
1272 const unsigned int select = CutOffFunctionBase<dim>::no_component,
1273 bool integrate_to_one = false);
1274
1278 virtual double
1279 value(const Point<dim> &p, const unsigned int component = 0) const override;
1280
1284 virtual void
1285 value_list(const std::vector<Point<dim>> &points,
1286 std::vector<double> & values,
1287 const unsigned int component = 0) const override;
1288
1292 virtual void
1293 vector_value_list(const std::vector<Point<dim>> &points,
1294 std::vector<Vector<double>> & values) const override;
1295
1299 virtual Tensor<1, dim>
1300 gradient(const Point<dim> & p,
1301 const unsigned int component = 0) const override;
1302 };
1303
1304
1305
1319 template <int dim, typename Number = double>
1320 class Monomial : public Function<dim, Number>
1321 {
1322 public:
1330 const unsigned int n_components = 1);
1331
1335 virtual Number
1336 value(const Point<dim> &p, const unsigned int component = 0) const override;
1337
1344 virtual void
1345 vector_value(const Point<dim> &p, Vector<Number> &values) const override;
1346
1350 virtual void
1351 value_list(const std::vector<Point<dim>> &points,
1352 std::vector<Number> & values,
1353 const unsigned int component = 0) const override;
1354
1359 gradient(const Point<dim> & p,
1360 const unsigned int component = 0) const override;
1361
1362 private:
1367 };
1368
1369
1370
1452 template <int dim>
1454 {
1455 public:
1474 const std::array<std::vector<double>, dim> &coordinate_values,
1476
1487 std::array<std::vector<double>, dim> &&coordinate_values,
1489
1500 virtual double
1501 value(const Point<dim> &p, const unsigned int component = 0) const override;
1502
1514 virtual Tensor<1, dim>
1515 gradient(const Point<dim> & p,
1516 const unsigned int component = 0) const override;
1517
1521 virtual std::size_t
1522 memory_consumption() const override;
1523
1527 const Table<dim, double> &
1528 get_data() const;
1529
1530 protected:
1535 table_index_of_point(const Point<dim> &p) const;
1536
1540 const std::array<std::vector<double>, dim> coordinate_values;
1541
1546 };
1547
1548
1592 template <int dim>
1594 {
1595 public:
1611 const std::array<std::pair<double, double>, dim> &interval_endpoints,
1612 const std::array<unsigned int, dim> & n_subintervals,
1614
1625 std::array<std::pair<double, double>, dim> &&interval_endpoints,
1626 std::array<unsigned int, dim> && n_subintervals,
1628
1639 virtual double
1640 value(const Point<dim> &p, const unsigned int component = 0) const override;
1641
1653 virtual Tensor<1, dim>
1654 gradient(const Point<dim> & p,
1655 const unsigned int component = 0) const override;
1656
1660 virtual std::size_t
1661 memory_consumption() const override;
1662
1666 const Table<dim, double> &
1667 get_data() const;
1668
1669 private:
1673 const std::array<std::pair<double, double>, dim> interval_endpoints;
1674
1678 const std::array<unsigned int, dim> n_subintervals;
1679
1684 };
1685
1686
1699 template <int dim>
1700 class Polynomial : public Function<dim>
1701 {
1702 public:
1713 const std::vector<double> &coefficients);
1714
1718 virtual double
1719 value(const Point<dim> &p, const unsigned int component = 0) const override;
1720
1721
1725 virtual void
1726 value_list(const std::vector<Point<dim>> &points,
1727 std::vector<double> & values,
1728 const unsigned int component = 0) const override;
1729
1733 virtual Tensor<1, dim>
1734 gradient(const Point<dim> & p,
1735 const unsigned int component = 0) const override;
1736
1740 virtual std::size_t
1741 memory_consumption() const override;
1742
1743 private:
1748
1752 const std::vector<double> coefficients;
1753 };
1754
1755#ifndef DOXYGEN
1756
1757
1758
1759 // Template definitions
1760 template <int dim>
1761 template <template <int> class CutOffFunctionBaseType>
1762 void
1764 {
1765 initialized = true;
1766 static_assert(
1767 std::is_base_of<CutOffFunctionBase<1>, CutOffFunctionBaseType<1>>::value,
1768 "You can only construct a CutOffFunctionTensorProduct function from "
1769 "a class derived from CutOffFunctionBase.");
1770
1771 for (unsigned int i = 0; i < dim; ++i)
1772 base[i].reset(new CutOffFunctionBaseType<1>(this->radius,
1773 Point<1>(this->center[i]),
1774 this->n_components,
1775 this->selected,
1776 this->integrate_to_one));
1777 }
1778
1779
1780
1781#endif
1782
1783} // namespace Functions
1785
1786#endif
const unsigned int n_components
Definition function.h:164
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
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
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 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
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
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
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
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
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
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
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
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
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
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
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
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
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
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
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:112
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
Point< 3 > center
static const unsigned int invalid_unsigned_int
Definition types.h:213