Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
function_lib.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2019 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 
22 #include <deal.II/base/function.h>
23 #include <deal.II/base/point.h>
24 #include <deal.II/base/smartpointer.h>
25 #include <deal.II/base/table.h>
26 
27 #include <array>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
39 namespace Functions
40 {
51  template <int dim>
52  class SquareFunction : public Function<dim>
53  {
54  public:
55  virtual double
56  value(const Point<dim> &p, const unsigned int component = 0) const override;
57  virtual void
58  vector_value(const Point<dim> &p, Vector<double> &values) const override;
59  virtual void
60  value_list(const std::vector<Point<dim>> &points,
61  std::vector<double> & values,
62  const unsigned int component = 0) const override;
63  virtual Tensor<1, dim>
64  gradient(const Point<dim> & p,
65  const unsigned int component = 0) const override;
66  virtual void
67  vector_gradient(const Point<dim> & p,
68  std::vector<Tensor<1, dim>> &gradient) const override;
69  virtual void
70  gradient_list(const std::vector<Point<dim>> &points,
71  std::vector<Tensor<1, dim>> & gradients,
72  const unsigned int component = 0) const override;
73  virtual double
74  laplacian(const Point<dim> & p,
75  const unsigned int component = 0) const override;
76  virtual void
77  laplacian_list(const std::vector<Point<dim>> &points,
78  std::vector<double> & values,
79  const unsigned int component = 0) const override;
80  };
81 
82 
83 
91  template <int dim>
92  class Q1WedgeFunction : public Function<dim>
93  {
94  public:
95  virtual double
96  value(const Point<dim> &p, const unsigned int component = 0) const override;
97 
98  virtual void
99  value_list(const std::vector<Point<dim>> &points,
100  std::vector<double> & values,
101  const unsigned int component = 0) const override;
102 
103  virtual void
104  vector_value_list(const std::vector<Point<dim>> &points,
105  std::vector<Vector<double>> & values) const override;
106 
107  virtual Tensor<1, dim>
108  gradient(const Point<dim> & p,
109  const unsigned int component = 0) const override;
110 
111  virtual void
112  gradient_list(const std::vector<Point<dim>> &points,
113  std::vector<Tensor<1, dim>> & gradients,
114  const unsigned int component = 0) const override;
115 
116  virtual void
117  vector_gradient_list(
118  const std::vector<Point<dim>> &,
119  std::vector<std::vector<Tensor<1, dim>>> &) const override;
120 
124  virtual double
125  laplacian(const Point<dim> & p,
126  const unsigned int component = 0) const override;
127 
131  virtual void
132  laplacian_list(const std::vector<Point<dim>> &points,
133  std::vector<double> & values,
134  const unsigned int component = 0) const override;
135  };
136 
137 
138 
154  template <int dim>
155  class PillowFunction : public Function<dim>
156  {
157  public:
162  PillowFunction(const double offset = 0.);
163 
167  virtual double
168  value(const Point<dim> &p, const unsigned int component = 0) const override;
169 
173  virtual void
174  value_list(const std::vector<Point<dim>> &points,
175  std::vector<double> & values,
176  const unsigned int component = 0) const override;
177 
181  virtual Tensor<1, dim>
182  gradient(const Point<dim> & p,
183  const unsigned int component = 0) const override;
184 
188  virtual void
189  gradient_list(const std::vector<Point<dim>> &points,
190  std::vector<Tensor<1, dim>> & gradients,
191  const unsigned int component = 0) const override;
192 
196  virtual double
197  laplacian(const Point<dim> & p,
198  const unsigned int component = 0) const override;
199 
203  virtual void
204  laplacian_list(const std::vector<Point<dim>> &points,
205  std::vector<double> & values,
206  const unsigned int component = 0) const override;
207 
208  private:
209  const double offset;
210  };
211 
212 
213 
222  template <int dim>
223  class CosineFunction : public Function<dim>
224  {
225  public:
230  CosineFunction(const unsigned int n_components = 1);
231 
232  virtual double
233  value(const Point<dim> &p, const unsigned int component = 0) const override;
234 
235  virtual void
236  value_list(const std::vector<Point<dim>> &points,
237  std::vector<double> & values,
238  const unsigned int component = 0) const override;
239 
240  virtual void
241  vector_value_list(const std::vector<Point<dim>> &points,
242  std::vector<Vector<double>> & values) const override;
243 
244  virtual Tensor<1, dim>
245  gradient(const Point<dim> & p,
246  const unsigned int component = 0) const override;
247 
248  virtual void
249  gradient_list(const std::vector<Point<dim>> &points,
250  std::vector<Tensor<1, dim>> & gradients,
251  const unsigned int component = 0) const override;
252 
253  virtual double
254  laplacian(const Point<dim> & p,
255  const unsigned int component = 0) const override;
256 
257  virtual void
258  laplacian_list(const std::vector<Point<dim>> &points,
259  std::vector<double> & values,
260  const unsigned int component = 0) const override;
261 
266  hessian(const Point<dim> & p,
267  const unsigned int component = 0) const override;
268 
272  virtual void
273  hessian_list(const std::vector<Point<dim>> & points,
274  std::vector<SymmetricTensor<2, dim>> &hessians,
275  const unsigned int component = 0) const override;
276  };
277 
278 
279 
291  template <int dim>
292  class CosineGradFunction : public Function<dim>
293  {
294  public:
299 
300  virtual double
301  value(const Point<dim> &p, const unsigned int component) const override;
302  virtual void
303  vector_value(const Point<dim> &p, Vector<double> &values) const override;
304  virtual void
305  value_list(const std::vector<Point<dim>> &points,
306  std::vector<double> & values,
307  const unsigned int component) const override;
308 
309  virtual void
310  vector_value_list(const std::vector<Point<dim>> &points,
311  std::vector<Vector<double>> & values) const override;
312 
313  virtual Tensor<1, dim>
314  gradient(const Point<dim> &p, const unsigned int component) const override;
315 
316  virtual void
317  gradient_list(const std::vector<Point<dim>> &points,
318  std::vector<Tensor<1, dim>> & gradients,
319  const unsigned int component) const override;
320 
321  virtual void
322  vector_gradient_list(
323  const std::vector<Point<dim>> & points,
324  std::vector<std::vector<Tensor<1, dim>>> &gradients) const override;
325 
326  virtual double
327  laplacian(const Point<dim> &p, const unsigned int component) const override;
328  };
329 
330 
331 
338  template <int dim>
339  class ExpFunction : public Function<dim>
340  {
341  public:
345  virtual double
346  value(const Point<dim> &p, const unsigned int component = 0) const override;
347 
351  virtual void
352  value_list(const std::vector<Point<dim>> &points,
353  std::vector<double> & values,
354  const unsigned int component = 0) const override;
355 
359  virtual Tensor<1, dim>
360  gradient(const Point<dim> & p,
361  const unsigned int component = 0) const override;
362 
366  virtual void
367  gradient_list(const std::vector<Point<dim>> &points,
368  std::vector<Tensor<1, dim>> & gradients,
369  const unsigned int component = 0) const override;
370 
374  virtual double
375  laplacian(const Point<dim> & p,
376  const unsigned int component = 0) const override;
377 
381  virtual void
382  laplacian_list(const std::vector<Point<dim>> &points,
383  std::vector<double> & values,
384  const unsigned int component = 0) const override;
385  };
386 
387 
388 
400  class LSingularityFunction : public Function<2>
401  {
402  public:
403  virtual double
404  value(const Point<2> &p, const unsigned int component = 0) const override;
405 
406  virtual void
407  value_list(const std::vector<Point<2>> &points,
408  std::vector<double> & values,
409  const unsigned int component = 0) const override;
410 
411  virtual void
412  vector_value_list(const std::vector<Point<2>> &points,
413  std::vector<Vector<double>> &values) const override;
414 
415  virtual Tensor<1, 2>
416  gradient(const Point<2> & p,
417  const unsigned int component = 0) const override;
418 
419  virtual void
420  gradient_list(const std::vector<Point<2>> &points,
421  std::vector<Tensor<1, 2>> & gradients,
422  const unsigned int component = 0) const override;
423 
424  virtual void
425  vector_gradient_list(
426  const std::vector<Point<2>> &,
427  std::vector<std::vector<Tensor<1, 2>>> &) const override;
428 
429  virtual double
430  laplacian(const Point<2> & p,
431  const unsigned int component = 0) const override;
432 
433  virtual void
434  laplacian_list(const std::vector<Point<2>> &points,
435  std::vector<double> & values,
436  const unsigned int component = 0) const override;
437  };
438 
439 
440 
451  {
452  public:
457  virtual double
458  value(const Point<2> &p, const unsigned int component) const override;
459 
460  virtual void
461  value_list(const std::vector<Point<2>> &points,
462  std::vector<double> & values,
463  const unsigned int component) const override;
464 
465  virtual void
466  vector_value_list(const std::vector<Point<2>> &points,
467  std::vector<Vector<double>> &values) const override;
468 
469  virtual Tensor<1, 2>
470  gradient(const Point<2> &p, const unsigned int component) const override;
471 
472  virtual void
473  gradient_list(const std::vector<Point<2>> &points,
474  std::vector<Tensor<1, 2>> & gradients,
475  const unsigned int component) const override;
476 
477  virtual void
478  vector_gradient_list(
479  const std::vector<Point<2>> &,
480  std::vector<std::vector<Tensor<1, 2>>> &) const override;
481 
482  virtual double
483  laplacian(const Point<2> &p, const unsigned int component) const override;
484 
485  virtual void
486  laplacian_list(const std::vector<Point<2>> &points,
487  std::vector<double> & values,
488  const unsigned int component) const override;
489  };
490 
491 
492 
499  template <int dim>
500  class SlitSingularityFunction : public Function<dim>
501  {
502  public:
503  virtual double
504  value(const Point<dim> &p, const unsigned int component = 0) const override;
505 
506  virtual void
507  value_list(const std::vector<Point<dim>> &points,
508  std::vector<double> & values,
509  const unsigned int component = 0) const override;
510 
511  virtual void
512  vector_value_list(const std::vector<Point<dim>> &points,
513  std::vector<Vector<double>> & values) const override;
514 
515  virtual Tensor<1, dim>
516  gradient(const Point<dim> & p,
517  const unsigned int component = 0) const override;
518 
519  virtual void
520  gradient_list(const std::vector<Point<dim>> &points,
521  std::vector<Tensor<1, dim>> & gradients,
522  const unsigned int component = 0) const override;
523 
524  virtual void
525  vector_gradient_list(
526  const std::vector<Point<dim>> &,
527  std::vector<std::vector<Tensor<1, dim>>> &) const override;
528 
529  virtual double
530  laplacian(const Point<dim> & p,
531  const unsigned int component = 0) const override;
532 
533  virtual void
534  laplacian_list(const std::vector<Point<dim>> &points,
535  std::vector<double> & values,
536  const unsigned int component = 0) const override;
537  };
538 
539 
547  {
548  public:
549  virtual double
550  value(const Point<2> &p, const unsigned int component = 0) const override;
551 
552  virtual void
553  value_list(const std::vector<Point<2>> &points,
554  std::vector<double> & values,
555  const unsigned int component = 0) const override;
556 
557  virtual void
558  vector_value_list(const std::vector<Point<2>> &points,
559  std::vector<Vector<double>> &values) const override;
560 
561  virtual Tensor<1, 2>
562  gradient(const Point<2> & p,
563  const unsigned int component = 0) const override;
564 
565  virtual void
566  gradient_list(const std::vector<Point<2>> &points,
567  std::vector<Tensor<1, 2>> & gradients,
568  const unsigned int component = 0) const override;
569 
570  virtual void
571  vector_gradient_list(
572  const std::vector<Point<2>> &,
573  std::vector<std::vector<Tensor<1, 2>>> &) const override;
574 
575  virtual double
576  laplacian(const Point<2> & p,
577  const unsigned int component = 0) const override;
578 
579  virtual void
580  laplacian_list(const std::vector<Point<2>> &points,
581  std::vector<double> & values,
582  const unsigned int component = 0) const override;
583  };
584 
585 
586 
602  template <int dim>
603  class JumpFunction : public Function<dim>
604  {
605  public:
610  JumpFunction(const Point<dim> &direction, const double steepness);
611 
615  virtual double
616  value(const Point<dim> &p, const unsigned int component = 0) const override;
617 
621  virtual void
622  value_list(const std::vector<Point<dim>> &points,
623  std::vector<double> & values,
624  const unsigned int component = 0) const override;
625 
629  virtual Tensor<1, dim>
630  gradient(const Point<dim> & p,
631  const unsigned int component = 0) const override;
632 
636  virtual void
637  gradient_list(const std::vector<Point<dim>> &points,
638  std::vector<Tensor<1, dim>> & gradients,
639  const unsigned int component = 0) const override;
640 
644  virtual double
645  laplacian(const Point<dim> & p,
646  const unsigned int component = 0) const override;
647 
651  virtual void
652  laplacian_list(const std::vector<Point<dim>> &points,
653  std::vector<double> & values,
654  const unsigned int component = 0) const override;
655 
662  std::size_t
663  memory_consumption() const;
664 
665  protected:
670 
674  const double steepness;
675 
679  double angle;
680 
684  double sine;
685 
689  double cosine;
690  };
691 
692 
693 
706  template <int dim>
707  class FourierCosineFunction : public Function<dim>
708  {
709  public:
715 
722  virtual double
723  value(const Point<dim> &p, const unsigned int component = 0) const override;
724 
729  virtual Tensor<1, dim>
730  gradient(const Point<dim> & p,
731  const unsigned int component = 0) const override;
732 
736  virtual double
737  laplacian(const Point<dim> & p,
738  const unsigned int component = 0) const override;
739 
740  private:
745  };
746 
747 
748 
761  template <int dim>
762  class FourierSineFunction : public Function<dim>
763  {
764  public:
770 
777  virtual double
778  value(const Point<dim> &p, const unsigned int component = 0) const override;
779 
784  virtual Tensor<1, dim>
785  gradient(const Point<dim> & p,
786  const unsigned int component = 0) const override;
787 
791  virtual double
792  laplacian(const Point<dim> & p,
793  const unsigned int component = 0) const override;
794 
795  private:
800  };
801 
802 
812  template <int dim>
813  class FourierSineSum : public Function<dim>
814  {
815  public:
820  FourierSineSum(const std::vector<Point<dim>> &fourier_coefficients,
821  const std::vector<double> & weights);
822 
829  virtual double
830  value(const Point<dim> &p, const unsigned int component = 0) const override;
831 
836  virtual Tensor<1, dim>
837  gradient(const Point<dim> & p,
838  const unsigned int component = 0) const override;
839 
843  virtual double
844  laplacian(const Point<dim> & p,
845  const unsigned int component = 0) const override;
846 
847  private:
851  const std::vector<Point<dim>> fourier_coefficients;
852  const std::vector<double> weights;
853  };
854 
855 
856 
867  template <int dim>
868  class FourierCosineSum : public Function<dim>
869  {
870  public:
876  const std::vector<double> & weights);
877 
884  virtual double
885  value(const Point<dim> &p, const unsigned int component = 0) const override;
886 
891  virtual Tensor<1, dim>
892  gradient(const Point<dim> & p,
893  const unsigned int component = 0) const override;
894 
898  virtual double
899  laplacian(const Point<dim> & p,
900  const unsigned int component = 0) const override;
901 
902  private:
906  const std::vector<Point<dim>> fourier_coefficients;
907  const std::vector<double> weights;
908  };
909 
910 
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 
966  DEAL_II_DEPRECATED void
967  new_center(const Point<dim> &p);
968 
974  DEAL_II_DEPRECATED void
975  new_radius(const double r);
976 
980  virtual void
981  set_center(const Point<dim> &p);
982 
986  virtual void
987  set_radius(const double r);
988 
992  const Point<dim> &
993  get_center() const;
994 
998  double
999  get_radius() const;
1000 
1004  bool
1005  integrates_to_one() const;
1006 
1007  protected:
1012 
1016  double radius;
1017 
1022  const unsigned int selected;
1023 
1028 
1034 
1038  double rescaling;
1039  };
1040 
1041 
1052  template <int dim>
1054  {
1055  public:
1066  double radius = 1.0,
1067  const Point<dim> & center = Point<dim>(),
1068  const unsigned int n_components = 1,
1069  const unsigned int select = CutOffFunctionBase<dim>::no_component,
1070  const bool integrate_to_one = false);
1071 
1076  template <template <int> class CutOffFunctionBaseType>
1077  void
1078  set_base();
1079 
1083  virtual void
1084  set_center(const Point<dim> &center) override;
1085 
1089  virtual void
1090  set_radius(const double radius) override;
1091 
1095  virtual double
1096  value(const Point<dim> &p, const unsigned int component = 0) const override;
1097 
1101  virtual Tensor<1, dim>
1102  gradient(const Point<dim> & p,
1103  const unsigned int component = 0) const override;
1104 
1105  private:
1106  std::array<std::unique_ptr<CutOffFunctionBase<1>>, dim> base;
1107 
1108  bool initialized;
1109  };
1110 
1111 
1112 
1122  template <int dim>
1124  {
1125  public:
1133  const double radius = 1.,
1134  const Point<dim> = Point<dim>(),
1135  const unsigned int n_components = 1,
1136  const unsigned int select = CutOffFunctionBase<dim>::no_component,
1137  const bool integrate_to_one = false);
1138 
1142  virtual double
1143  value(const Point<dim> &p, const unsigned int component = 0) const override;
1144 
1148  virtual void
1149  value_list(const std::vector<Point<dim>> &points,
1150  std::vector<double> & values,
1151  const unsigned int component = 0) const override;
1152 
1156  virtual void
1157  vector_value_list(const std::vector<Point<dim>> &points,
1158  std::vector<Vector<double>> & values) const override;
1159  };
1160 
1161 
1171  template <int dim>
1173  {
1174  public:
1182  const double radius = 1.,
1183  const Point<dim> = Point<dim>(),
1184  const unsigned int n_components = 1,
1185  const unsigned int select = CutOffFunctionBase<dim>::no_component,
1186  const bool integrate_to_one = false);
1187 
1191  virtual double
1192  value(const Point<dim> &p, const unsigned int component = 0) const override;
1193 
1197  virtual void
1198  value_list(const std::vector<Point<dim>> &points,
1199  std::vector<double> & values,
1200  const unsigned int component = 0) const override;
1201 
1205  virtual void
1206  vector_value_list(const std::vector<Point<dim>> &points,
1207  std::vector<Vector<double>> & values) const override;
1208  };
1209 
1210 
1224  template <int dim>
1226  {
1227  public:
1232  const double radius = 1.,
1233  const Point<dim> = Point<dim>(),
1234  const unsigned int n_components = 1,
1235  const unsigned int select = CutOffFunctionBase<dim>::no_component,
1236  bool integrate_to_one = false);
1237 
1241  virtual double
1242  value(const Point<dim> &p, const unsigned int component = 0) const override;
1243 
1247  virtual void
1248  value_list(const std::vector<Point<dim>> &points,
1249  std::vector<double> & values,
1250  const unsigned int component = 0) const override;
1251 
1255  virtual void
1256  vector_value_list(const std::vector<Point<dim>> &points,
1257  std::vector<Vector<double>> & values) const override;
1258 
1262  virtual Tensor<1, dim>
1263  gradient(const Point<dim> & p,
1264  const unsigned int component = 0) const override;
1265  };
1266 
1267 
1278  template <int dim>
1280  {
1281  public:
1289  const double radius = 1.,
1290  const Point<dim> = Point<dim>(),
1291  const unsigned int n_components = 1,
1292  const unsigned int select = CutOffFunctionBase<dim>::no_component,
1293  bool integrate_to_one = false);
1294 
1298  virtual double
1299  value(const Point<dim> &p, const unsigned int component = 0) const override;
1300 
1304  virtual void
1305  value_list(const std::vector<Point<dim>> &points,
1306  std::vector<double> & values,
1307  const unsigned int component = 0) const override;
1308 
1312  virtual void
1313  vector_value_list(const std::vector<Point<dim>> &points,
1314  std::vector<Vector<double>> & values) const override;
1315 
1319  virtual Tensor<1, dim>
1320  gradient(const Point<dim> & p,
1321  const unsigned int component = 0) const override;
1322  };
1323 
1324 
1325 
1340  template <int dim>
1341  class Monomial : public Function<dim>
1342  {
1343  public:
1351  const unsigned int n_components = 1);
1352 
1356  virtual double
1357  value(const Point<dim> &p, const unsigned int component = 0) const override;
1358 
1365  virtual void
1366  vector_value(const Point<dim> &p, Vector<double> &values) const override;
1367 
1371  virtual void
1372  value_list(const std::vector<Point<dim>> &points,
1373  std::vector<double> & values,
1374  const unsigned int component = 0) const override;
1375 
1379  virtual Tensor<1, dim>
1380  gradient(const Point<dim> & p,
1381  const unsigned int component = 0) const override;
1382 
1383  private:
1388  };
1389 
1390 
1391 
1425  template <int dim>
1427  {
1428  public:
1443  const std::array<std::vector<double>, dim> &coordinate_values,
1445 
1456  virtual double
1457  value(const Point<dim> &p, const unsigned int component = 0) const override;
1458 
1470  virtual Tensor<1, dim>
1471  gradient(const Point<dim> & p,
1472  const unsigned int component = 0) const override;
1473 
1474  protected:
1479  table_index_of_point(const Point<dim> &p) const;
1480 
1484  const std::array<std::vector<double>, dim> coordinate_values;
1485 
1490  };
1491 
1492 
1529  template <int dim>
1531  {
1532  public:
1548  const std::array<std::pair<double, double>, dim> &interval_endpoints,
1549  const std::array<unsigned int, dim> & n_subintervals,
1551 
1562  virtual double
1563  value(const Point<dim> &p, const unsigned int component = 0) const override;
1564 
1565  private:
1569  const std::array<std::pair<double, double>, dim> interval_endpoints;
1570 
1574  const std::array<unsigned int, dim> n_subintervals;
1575 
1580  };
1581 
1582 
1596  template <int dim>
1597  class Polynomial : public Function<dim>
1598  {
1599  public:
1610  const std::vector<double> &coefficients);
1611 
1615  virtual double
1616  value(const Point<dim> &p, const unsigned int component = 0) const override;
1617 
1618 
1622  virtual void
1623  value_list(const std::vector<Point<dim>> &points,
1624  std::vector<double> & values,
1625  const unsigned int component = 0) const override;
1626 
1630  virtual Tensor<1, dim>
1631  gradient(const Point<dim> & p,
1632  const unsigned int component = 0) const override;
1633 
1634  private:
1639 
1643  const std::vector<double> coefficients;
1644  };
1645 
1646 #ifndef DOXYGEN
1647 
1648 
1649 
1650  // Template definitions
1651  template <int dim>
1652  template <template <int> class CutOffFunctionBaseType>
1653  void
1655  {
1656  initialized = true;
1657  static_assert(
1658  std::is_base_of<CutOffFunctionBase<1>, CutOffFunctionBaseType<1>>::value,
1659  "You can only construct a CutOffFunctionTensorProduct function from "
1660  "a class derived from CutOffFunctionBase.");
1661 
1662  for (unsigned int i = 0; i < dim; ++i)
1663  base[i].reset(new CutOffFunctionBaseType<1>(this->radius,
1664  Point<1>(this->center[i]),
1665  this->n_components,
1666  this->selected,
1667  this->integrate_to_one));
1668  }
1669 
1670 
1671 
1672 #endif
1673 
1674 } // namespace Functions
1675 DEAL_II_NAMESPACE_CLOSE
1676 
1677 #endif
const Tensor< 1, dim > exponents
virtual Tensor< 1, dim > gradient(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) const override
virtual SymmetricTensor< 2, dim > hessian(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
static const unsigned int invalid_unsigned_int
Definition: types.h:173
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
const unsigned int n_components
Definition: function.h:162
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
TableIndices< dim > table_index_of_point(const Point< dim > &p) const
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 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) const override
Monomial(const Tensor< 1, dim > &exponents, const unsigned int n_components=1)
Polynomial(const Table< 2, double > &exponents, const std::vector< double > &coefficients)
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)
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 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 value_list(const std::vector< Point< dim >> &points, std::vector< double > &values, const unsigned int component=0) const override
Definition: function_lib.cc:52
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const override
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_center(const Point< dim > &center) override
FourierCosineSum(const std::vector< Point< dim >> &fourier_coefficients, const std::vector< double > &weights)
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
Definition: function_lib.cc:34
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 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
CosineFunction(const unsigned int n_components=1)
virtual double laplacian(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 Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void set_center(const Point< dim > &p)
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
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)
const Point< dim > & get_center() const
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 Tensor< 1, dim > gradient(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
Definition: function_lib.cc:92
virtual void set_radius(const double radius) override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
void new_center(const Point< dim > &p)
const Tensor< 1, dim > fourier_coefficients
Definition: function_lib.h:744
std::size_t memory_consumption() const
virtual void value_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 Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const override
Definition: function_lib.cc:42
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
virtual ~CutOffFunctionBase()=default
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) 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)
const Table< dim, double > data_values
virtual Tensor< 1, dim > gradient(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(const Point< dim > &p, Vector< double > &values) 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
const std::vector< Point< dim > > fourier_coefficients
Definition: function_lib.h:906
virtual double value(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 Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component) const override
virtual void value_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
Definition: function_lib.cc:69
const Table< 2, double > exponents
virtual void laplacian_list(const std::vector< Point< dim >> &points, std::vector< double > &values, 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
const std::vector< double > coefficients
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
const unsigned int selected
virtual void vector_value_list(const std::vector< Point< dim >> &points, std::vector< Vector< double >> &values) const override
virtual void vector_value_list(const std::vector< Point< dim >> &points, std::vector< Vector< double >> &values) const override
const Point< dim > direction
Definition: function_lib.h:669
virtual void value_list(const std::vector< Point< dim >> &points, std::vector< double > &values, const unsigned int component=0) const override
const std::vector< Point< dim > > fourier_coefficients
Definition: function_lib.h:851
virtual Tensor< 1, dim > gradient(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 double value(const Point< dim > &p, const unsigned int component=0) const override
PillowFunction(const double offset=0.)
const std::array< std::vector< double >, dim > coordinate_values
virtual Tensor< 1, dim > gradient(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 std::array< unsigned int, dim > n_subintervals
static const unsigned int no_component
Definition: function_lib.h:931
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)
FourierCosineFunction(const Tensor< 1, dim > &fourier_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 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
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
FourierSineFunction(const Tensor< 1, dim > &fourier_coefficients)
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 void vector_value_list(const std::vector< Point< dim >> &points, std::vector< Vector< double >> &values) 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
const Tensor< 1, dim > fourier_coefficients
Definition: function_lib.h:799
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 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< dim >> &points, std::vector< Tensor< 1, dim >> &gradients, 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
const std::array< std::pair< double, double >, dim > interval_endpoints
virtual void value_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) 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 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 double laplacian(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
Definition: function_lib.cc:77
virtual void set_radius(const double r)
virtual void value_list(const std::vector< Point< dim >> &points, std::vector< double > &values, const unsigned int component=0) const override
JumpFunction(const Point< dim > &direction, const double steepness)