Reference documentation for deal.II version 9.1.0-pre
function_lib.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2018 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/table.h>
25 
26 #include <array>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
38 namespace Functions
39 {
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
66  vector_gradient(const Point<dim> & p,
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 
90  template <int dim>
91  class Q1WedgeFunction : public Function<dim>
92  {
93  public:
94  virtual double
95  value(const Point<dim> &p, const unsigned int component = 0) const override;
96 
97  virtual void
98  value_list(const std::vector<Point<dim>> &points,
99  std::vector<double> & values,
100  const unsigned int component = 0) const override;
101 
102  virtual void
103  vector_value_list(const std::vector<Point<dim>> &points,
104  std::vector<Vector<double>> & values) const override;
105 
106  virtual Tensor<1, dim>
107  gradient(const Point<dim> & p,
108  const unsigned int component = 0) const override;
109 
110  virtual void
111  gradient_list(const std::vector<Point<dim>> &points,
112  std::vector<Tensor<1, dim>> & gradients,
113  const unsigned int component = 0) const override;
114 
115  virtual void
117  const std::vector<Point<dim>> &,
118  std::vector<std::vector<Tensor<1, dim>>> &) const override;
119 
123  virtual double
124  laplacian(const Point<dim> & p,
125  const unsigned int component = 0) const override;
126 
130  virtual void
131  laplacian_list(const std::vector<Point<dim>> &points,
132  std::vector<double> & values,
133  const unsigned int component = 0) const override;
134  };
135 
136 
137 
153  template <int dim>
154  class PillowFunction : public Function<dim>
155  {
156  public:
161  PillowFunction(const double offset = 0.);
162 
166  virtual double
167  value(const Point<dim> &p, const unsigned int component = 0) const override;
168 
172  virtual void
173  value_list(const std::vector<Point<dim>> &points,
174  std::vector<double> & values,
175  const unsigned int component = 0) const override;
176 
180  virtual Tensor<1, dim>
181  gradient(const Point<dim> & p,
182  const unsigned int component = 0) const override;
183 
187  virtual void
188  gradient_list(const std::vector<Point<dim>> &points,
189  std::vector<Tensor<1, dim>> & gradients,
190  const unsigned int component = 0) const override;
191 
195  virtual double
196  laplacian(const Point<dim> & p,
197  const unsigned int component = 0) const override;
198 
202  virtual void
203  laplacian_list(const std::vector<Point<dim>> &points,
204  std::vector<double> & values,
205  const unsigned int component = 0) const override;
206 
207  private:
208  const double offset;
209  };
210 
211 
212 
221  template <int dim>
222  class CosineFunction : public Function<dim>
223  {
224  public:
229  CosineFunction(const unsigned int n_components = 1);
230 
231  virtual double
232  value(const Point<dim> &p, const unsigned int component = 0) const override;
233 
234  virtual void
235  value_list(const std::vector<Point<dim>> &points,
236  std::vector<double> & values,
237  const unsigned int component = 0) const override;
238 
239  virtual void
240  vector_value_list(const std::vector<Point<dim>> &points,
241  std::vector<Vector<double>> & values) const override;
242 
243  virtual Tensor<1, dim>
244  gradient(const Point<dim> & p,
245  const unsigned int component = 0) const override;
246 
247  virtual void
248  gradient_list(const std::vector<Point<dim>> &points,
249  std::vector<Tensor<1, dim>> & gradients,
250  const unsigned int component = 0) const override;
251 
252  virtual double
253  laplacian(const Point<dim> & p,
254  const unsigned int component = 0) const override;
255 
256  virtual void
257  laplacian_list(const std::vector<Point<dim>> &points,
258  std::vector<double> & values,
259  const unsigned int component = 0) const override;
260 
265  hessian(const Point<dim> & p,
266  const unsigned int component = 0) const override;
267 
271  virtual void
272  hessian_list(const std::vector<Point<dim>> & points,
273  std::vector<SymmetricTensor<2, dim>> &hessians,
274  const unsigned int component = 0) const override;
275  };
276 
277 
278 
290  template <int dim>
291  class CosineGradFunction : public Function<dim>
292  {
293  public:
298 
299  virtual double
300  value(const Point<dim> &p, const unsigned int component) const override;
301  virtual void
302  vector_value(const Point<dim> &p, Vector<double> &values) const override;
303  virtual void
304  value_list(const std::vector<Point<dim>> &points,
305  std::vector<double> & values,
306  const unsigned int component) const override;
307 
308  virtual void
309  vector_value_list(const std::vector<Point<dim>> &points,
310  std::vector<Vector<double>> & values) const override;
311 
312  virtual Tensor<1, dim>
313  gradient(const Point<dim> &p, const unsigned int component) const override;
314 
315  virtual void
316  gradient_list(const std::vector<Point<dim>> &points,
317  std::vector<Tensor<1, dim>> & gradients,
318  const unsigned int component) const override;
319 
320  virtual void
322  const std::vector<Point<dim>> & points,
323  std::vector<std::vector<Tensor<1, dim>>> &gradients) const override;
324 
325  virtual double
326  laplacian(const Point<dim> &p, const unsigned int component) const override;
327  };
328 
329 
330 
337  template <int dim>
338  class ExpFunction : public Function<dim>
339  {
340  public:
344  virtual double
345  value(const Point<dim> &p, const unsigned int component = 0) const override;
346 
350  virtual void
351  value_list(const std::vector<Point<dim>> &points,
352  std::vector<double> & values,
353  const unsigned int component = 0) const override;
354 
358  virtual Tensor<1, dim>
359  gradient(const Point<dim> & p,
360  const unsigned int component = 0) const override;
361 
365  virtual void
366  gradient_list(const std::vector<Point<dim>> &points,
367  std::vector<Tensor<1, dim>> & gradients,
368  const unsigned int component = 0) const override;
369 
373  virtual double
374  laplacian(const Point<dim> & p,
375  const unsigned int component = 0) const override;
376 
380  virtual void
381  laplacian_list(const std::vector<Point<dim>> &points,
382  std::vector<double> & values,
383  const unsigned int component = 0) const override;
384  };
385 
386 
387 
399  class LSingularityFunction : public Function<2>
400  {
401  public:
402  virtual double
403  value(const Point<2> &p, const unsigned int component = 0) const override;
404 
405  virtual void
406  value_list(const std::vector<Point<2>> &points,
407  std::vector<double> & values,
408  const unsigned int component = 0) const override;
409 
410  virtual void
411  vector_value_list(const std::vector<Point<2>> &points,
412  std::vector<Vector<double>> &values) const override;
413 
414  virtual Tensor<1, 2>
415  gradient(const Point<2> & p,
416  const unsigned int component = 0) const override;
417 
418  virtual void
419  gradient_list(const std::vector<Point<2>> &points,
420  std::vector<Tensor<1, 2>> & gradients,
421  const unsigned int component = 0) const override;
422 
423  virtual void
425  const std::vector<Point<2>> &,
426  std::vector<std::vector<Tensor<1, 2>>> &) const override;
427 
428  virtual double
429  laplacian(const Point<2> & p,
430  const unsigned int component = 0) const override;
431 
432  virtual void
433  laplacian_list(const std::vector<Point<2>> &points,
434  std::vector<double> & values,
435  const unsigned int component = 0) const override;
436  };
437 
438 
439 
450  {
451  public:
456  virtual double
457  value(const Point<2> &p, const unsigned int component) const override;
458 
459  virtual void
460  value_list(const std::vector<Point<2>> &points,
461  std::vector<double> & values,
462  const unsigned int component) const override;
463 
464  virtual void
465  vector_value_list(const std::vector<Point<2>> &points,
466  std::vector<Vector<double>> &values) const override;
467 
468  virtual Tensor<1, 2>
469  gradient(const Point<2> &p, const unsigned int component) const override;
470 
471  virtual void
472  gradient_list(const std::vector<Point<2>> &points,
473  std::vector<Tensor<1, 2>> & gradients,
474  const unsigned int component) const override;
475 
476  virtual void
478  const std::vector<Point<2>> &,
479  std::vector<std::vector<Tensor<1, 2>>> &) const override;
480 
481  virtual double
482  laplacian(const Point<2> &p, const unsigned int component) const override;
483 
484  virtual void
485  laplacian_list(const std::vector<Point<2>> &points,
486  std::vector<double> & values,
487  const unsigned int component) const override;
488  };
489 
490 
491 
498  template <int dim>
499  class SlitSingularityFunction : public Function<dim>
500  {
501  public:
502  virtual double
503  value(const Point<dim> &p, const unsigned int component = 0) const override;
504 
505  virtual void
506  value_list(const std::vector<Point<dim>> &points,
507  std::vector<double> & values,
508  const unsigned int component = 0) const override;
509 
510  virtual void
511  vector_value_list(const std::vector<Point<dim>> &points,
512  std::vector<Vector<double>> & values) const override;
513 
514  virtual Tensor<1, dim>
515  gradient(const Point<dim> & p,
516  const unsigned int component = 0) const override;
517 
518  virtual void
519  gradient_list(const std::vector<Point<dim>> &points,
520  std::vector<Tensor<1, dim>> & gradients,
521  const unsigned int component = 0) const override;
522 
523  virtual void
525  const std::vector<Point<dim>> &,
526  std::vector<std::vector<Tensor<1, dim>>> &) const override;
527 
528  virtual double
529  laplacian(const Point<dim> & p,
530  const unsigned int component = 0) const override;
531 
532  virtual void
533  laplacian_list(const std::vector<Point<dim>> &points,
534  std::vector<double> & values,
535  const unsigned int component = 0) const override;
536  };
537 
538 
546  {
547  public:
548  virtual double
549  value(const Point<2> &p, const unsigned int component = 0) const override;
550 
551  virtual void
552  value_list(const std::vector<Point<2>> &points,
553  std::vector<double> & values,
554  const unsigned int component = 0) const override;
555 
556  virtual void
557  vector_value_list(const std::vector<Point<2>> &points,
558  std::vector<Vector<double>> &values) const override;
559 
560  virtual Tensor<1, 2>
561  gradient(const Point<2> & p,
562  const unsigned int component = 0) const override;
563 
564  virtual void
565  gradient_list(const std::vector<Point<2>> &points,
566  std::vector<Tensor<1, 2>> & gradients,
567  const unsigned int component = 0) const override;
568 
569  virtual void
571  const std::vector<Point<2>> &,
572  std::vector<std::vector<Tensor<1, 2>>> &) const override;
573 
574  virtual double
575  laplacian(const Point<2> & p,
576  const unsigned int component = 0) const override;
577 
578  virtual void
579  laplacian_list(const std::vector<Point<2>> &points,
580  std::vector<double> & values,
581  const unsigned int component = 0) const override;
582  };
583 
584 
585 
601  template <int dim>
602  class JumpFunction : public Function<dim>
603  {
604  public:
609  JumpFunction(const Point<dim> &direction, const double steepness);
610 
614  virtual double
615  value(const Point<dim> &p, const unsigned int component = 0) const override;
616 
620  virtual void
621  value_list(const std::vector<Point<dim>> &points,
622  std::vector<double> & values,
623  const unsigned int component = 0) const override;
624 
628  virtual Tensor<1, dim>
629  gradient(const Point<dim> & p,
630  const unsigned int component = 0) const override;
631 
635  virtual void
636  gradient_list(const std::vector<Point<dim>> &points,
637  std::vector<Tensor<1, dim>> & gradients,
638  const unsigned int component = 0) const override;
639 
643  virtual double
644  laplacian(const Point<dim> & p,
645  const unsigned int component = 0) const override;
646 
650  virtual void
651  laplacian_list(const std::vector<Point<dim>> &points,
652  std::vector<double> & values,
653  const unsigned int component = 0) const override;
654 
661  std::size_t
662  memory_consumption() const;
663 
664  protected:
669 
673  const double steepness;
674 
678  double angle;
679 
683  double sine;
684 
688  double cosine;
689  };
690 
691 
692 
705  template <int dim>
706  class FourierCosineFunction : public Function<dim>
707  {
708  public:
713  FourierCosineFunction(const Tensor<1, dim> &fourier_coefficients);
714 
721  virtual double
722  value(const Point<dim> &p, const unsigned int component = 0) const override;
723 
728  virtual Tensor<1, dim>
729  gradient(const Point<dim> & p,
730  const unsigned int component = 0) const override;
731 
735  virtual double
736  laplacian(const Point<dim> & p,
737  const unsigned int component = 0) const override;
738 
739  private:
744  };
745 
746 
747 
760  template <int dim>
761  class FourierSineFunction : public Function<dim>
762  {
763  public:
768  FourierSineFunction(const Tensor<1, dim> &fourier_coefficients);
769 
776  virtual double
777  value(const Point<dim> &p, const unsigned int component = 0) const override;
778 
783  virtual Tensor<1, dim>
784  gradient(const Point<dim> & p,
785  const unsigned int component = 0) const override;
786 
790  virtual double
791  laplacian(const Point<dim> & p,
792  const unsigned int component = 0) const override;
793 
794  private:
799  };
800 
801 
811  template <int dim>
812  class FourierSineSum : public Function<dim>
813  {
814  public:
819  FourierSineSum(const std::vector<Point<dim>> &fourier_coefficients,
820  const std::vector<double> & weights);
821 
828  virtual double
829  value(const Point<dim> &p, const unsigned int component = 0) const override;
830 
835  virtual Tensor<1, dim>
836  gradient(const Point<dim> & p,
837  const unsigned int component = 0) const override;
838 
842  virtual double
843  laplacian(const Point<dim> & p,
844  const unsigned int component = 0) const override;
845 
846  private:
850  const std::vector<Point<dim>> fourier_coefficients;
851  const std::vector<double> weights;
852  };
853 
854 
855 
866  template <int dim>
867  class FourierCosineSum : public Function<dim>
868  {
869  public:
874  FourierCosineSum(const std::vector<Point<dim>> &fourier_coefficients,
875  const std::vector<double> & weights);
876 
883  virtual double
884  value(const Point<dim> &p, const unsigned int component = 0) const override;
885 
890  virtual Tensor<1, dim>
891  gradient(const Point<dim> & p,
892  const unsigned int component = 0) const override;
893 
897  virtual double
898  laplacian(const Point<dim> & p,
899  const unsigned int component = 0) const override;
900 
901  private:
905  const std::vector<Point<dim>> fourier_coefficients;
906  const std::vector<double> weights;
907  };
908 
909 
918  template <int dim>
919  class CutOffFunctionBase : public Function<dim>
920  {
921  public:
926  static const unsigned int no_component = numbers::invalid_unsigned_int;
927 
935  const double radius = 1.,
936  const Point<dim> = Point<dim>(),
937  const unsigned int n_components = 1,
938  const unsigned int select = CutOffFunctionBase<dim>::no_component);
939 
943  void
944  new_center(const Point<dim> &p);
945 
949  void
950  new_radius(const double r);
951 
952  protected:
957 
961  double radius;
962 
967  const unsigned int selected;
968  };
969 
970 
971 
981  template <int dim>
983  {
984  public:
992  const double radius = 1.,
993  const Point<dim> = Point<dim>(),
994  const unsigned int n_components = 1,
995  const unsigned int select = CutOffFunctionBase<dim>::no_component);
996 
1000  virtual double
1001  value(const Point<dim> &p, const unsigned int component = 0) const override;
1002 
1006  virtual void
1007  value_list(const std::vector<Point<dim>> &points,
1008  std::vector<double> & values,
1009  const unsigned int component = 0) const override;
1010 
1014  virtual void
1015  vector_value_list(const std::vector<Point<dim>> &points,
1016  std::vector<Vector<double>> & values) const override;
1017  };
1018 
1019 
1029  template <int dim>
1031  {
1032  public:
1040  const double radius = 1.,
1041  const Point<dim> = Point<dim>(),
1042  const unsigned int n_components = 1,
1043  const unsigned int select = CutOffFunctionBase<dim>::no_component);
1044 
1048  virtual double
1049  value(const Point<dim> &p, const unsigned int component = 0) const override;
1050 
1054  virtual void
1055  value_list(const std::vector<Point<dim>> &points,
1056  std::vector<double> & values,
1057  const unsigned int component = 0) const override;
1058 
1062  virtual void
1063  vector_value_list(const std::vector<Point<dim>> &points,
1064  std::vector<Vector<double>> & values) const override;
1065  };
1066 
1067 
1078  template <int dim>
1080  {
1081  public:
1089  const double radius = 1.,
1090  const Point<dim> = Point<dim>(),
1091  const unsigned int n_components = 1,
1092  const unsigned int select = CutOffFunctionBase<dim>::no_component);
1093 
1097  virtual double
1098  value(const Point<dim> &p, const unsigned int component = 0) const override;
1099 
1103  virtual void
1104  value_list(const std::vector<Point<dim>> &points,
1105  std::vector<double> & values,
1106  const unsigned int component = 0) const override;
1107 
1111  virtual void
1112  vector_value_list(const std::vector<Point<dim>> &points,
1113  std::vector<Vector<double>> & values) const override;
1114 
1118  virtual Tensor<1, dim>
1119  gradient(const Point<dim> & p,
1120  const unsigned int component = 0) const override;
1121  };
1122 
1123 
1124 
1138  template <int dim>
1139  class Monomial : public Function<dim>
1140  {
1141  public:
1148  Monomial(const Tensor<1, dim> &exponents,
1149  const unsigned int n_components = 1);
1150 
1154  virtual double
1155  value(const Point<dim> &p, const unsigned int component = 0) const override;
1156 
1163  virtual void
1164  vector_value(const Point<dim> &p, Vector<double> &values) const override;
1165 
1169  virtual void
1170  value_list(const std::vector<Point<dim>> &points,
1171  std::vector<double> & values,
1172  const unsigned int component = 0) const override;
1173 
1177  virtual Tensor<1, dim>
1178  gradient(const Point<dim> & p,
1179  const unsigned int component = 0) const override;
1180 
1181  private:
1186  };
1187 
1188 
1189 
1222  template <int dim>
1224  {
1225  public:
1240  const std::array<std::vector<double>, dim> &coordinate_values,
1241  const Table<dim, double> & data_values);
1242 
1253  virtual double
1254  value(const Point<dim> &p, const unsigned int component = 0) const override;
1255 
1267  virtual Tensor<1, dim>
1268  gradient(const Point<dim> & p,
1269  const unsigned int component = 0) const override;
1270 
1271  protected:
1276  table_index_of_point(const Point<dim> &p) const;
1277 
1281  const std::array<std::vector<double>, dim> coordinate_values;
1282 
1287  };
1288 
1289 
1325  template <int dim>
1327  {
1328  public:
1344  const std::array<std::pair<double, double>, dim> &interval_endpoints,
1345  const std::array<unsigned int, dim> & n_subintervals,
1346  const Table<dim, double> & data_values);
1347 
1358  virtual double
1359  value(const Point<dim> &p, const unsigned int component = 0) const override;
1360 
1361  private:
1365  const std::array<std::pair<double, double>, dim> interval_endpoints;
1366 
1370  const std::array<unsigned int, dim> n_subintervals;
1371 
1376  };
1377 
1378 
1391  template <int dim>
1392  class Polynomial : public Function<dim>
1393  {
1394  public:
1404  Polynomial(const Table<2, double> & exponents,
1405  const std::vector<double> &coefficients);
1406 
1410  virtual double
1411  value(const Point<dim> &p, const unsigned int component = 0) const override;
1412 
1413 
1417  virtual void
1418  value_list(const std::vector<Point<dim>> &points,
1419  std::vector<double> & values,
1420  const unsigned int component = 0) const override;
1421 
1425  virtual Tensor<1, dim>
1426  gradient(const Point<dim> & p,
1427  const unsigned int component = 0) const override;
1428 
1429  private:
1434 
1438  const std::vector<double> coefficients;
1439  };
1440 
1441 
1442 
1443 } // namespace Functions
1444 DEAL_II_NAMESPACE_CLOSE
1445 
1446 #endif
const Tensor< 1, dim > exponents
static const unsigned int invalid_unsigned_int
Definition: types.h:173
const unsigned int n_components
Definition: function.h:160
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_gradient_list(const std::vector< Point< dim >> &points, std::vector< std::vector< Tensor< 1, dim, double >>> &gradients) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
Definition: function_lib.cc:34
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
Definition: function_lib.cc:92
const Tensor< 1, dim > fourier_coefficients
Definition: function_lib.h:743
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const override
Definition: function_lib.cc:42
virtual SymmetricTensor< 2, dim, double > hessian(const Point< dim > &p, const unsigned int component=0) const
std::size_t memory_consumption() const
const Table< dim, double > data_values
const std::vector< Point< dim > > fourier_coefficients
Definition: function_lib.h:905
virtual void hessian_list(const std::vector< Point< dim >> &points, std::vector< SymmetricTensor< 2, dim, double >> &values, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
Definition: function_lib.cc:69
const Table< 2, double > exponents
const std::vector< double > coefficients
const unsigned int selected
Definition: function_lib.h:967
const Point< dim > direction
Definition: function_lib.h:668
const std::vector< Point< dim > > fourier_coefficients
Definition: function_lib.h:850
const std::array< std::vector< double >, dim > coordinate_values
const std::array< unsigned int, dim > n_subintervals
const Tensor< 1, dim > fourier_coefficients
Definition: function_lib.h:798
virtual void vector_value_list(const std::vector< Point< dim >> &points, std::vector< Vector< double >> &values) const
const std::array< std::pair< double, double >, dim > interval_endpoints
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