Reference documentation for deal.II version 8.5.1
function_lib.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2016 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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/function.h>
22 #include <deal.II/base/point.h>
23 #include <deal.II/base/table.h>
24 
25 #include <deal.II/base/std_cxx11/array.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
37 namespace Functions
38 {
39 
40 
51  template<int dim>
52  class SquareFunction : public Function<dim>
53  {
54  public:
55  virtual double value (const Point<dim> &p,
56  const unsigned int component = 0) const;
57  virtual void vector_value (const Point<dim> &p,
58  Vector<double> &values) const;
59  virtual void value_list (const std::vector<Point<dim> > &points,
60  std::vector<double> &values,
61  const unsigned int component = 0) const;
62  virtual Tensor<1,dim> gradient (const Point<dim> &p,
63  const unsigned int component = 0) const;
64  virtual void vector_gradient (const Point<dim> &p,
65  std::vector<Tensor<1,dim> > &gradient) const;
66  virtual void gradient_list (const std::vector<Point<dim> > &points,
67  std::vector<Tensor<1,dim> > &gradients,
68  const unsigned int component = 0) const;
69  virtual double laplacian (const Point<dim> &p,
70  const unsigned int component = 0) const;
71  virtual void laplacian_list (const std::vector<Point<dim> > &points,
72  std::vector<double> &values,
73  const unsigned int component = 0) const;
74  };
75 
76 
77 
85  template<int dim>
86  class Q1WedgeFunction : public Function<dim>
87  {
88  public:
89  virtual double value (const Point<dim> &p,
90  const unsigned int component = 0) const;
91 
92  virtual void value_list (const std::vector<Point<dim> > &points,
93  std::vector<double> &values,
94  const unsigned int component = 0) const;
95 
96  virtual void vector_value_list (const std::vector<Point<dim> > &points,
97  std::vector<Vector<double> > &values) const;
98 
99  virtual Tensor<1,dim> gradient (const Point<dim> &p,
100  const unsigned int component = 0) const;
101 
102  virtual void gradient_list (const std::vector<Point<dim> > &points,
103  std::vector<Tensor<1,dim> > &gradients,
104  const unsigned int component = 0) const;
105 
106  virtual void vector_gradient_list (const std::vector<Point<dim> > &,
107  std::vector<std::vector<Tensor<1,dim> > > &) const;
108 
112  virtual double laplacian (const Point<dim> &p,
113  const unsigned int component = 0) const;
114 
118  virtual void laplacian_list (const std::vector<Point<dim> > &points,
119  std::vector<double> &values,
120  const unsigned int component = 0) const;
121  };
122 
123 
124 
140  template<int dim>
141  class PillowFunction : public Function<dim>
142  {
143  public:
148  PillowFunction (const double offset=0.);
149 
153  virtual double value (const Point<dim> &p,
154  const unsigned int component = 0) const;
155 
159  virtual void value_list (const std::vector<Point<dim> > &points,
160  std::vector<double> &values,
161  const unsigned int component = 0) const;
162 
166  virtual Tensor<1,dim> gradient (const Point<dim> &p,
167  const unsigned int component = 0) const;
168 
172  virtual void gradient_list (const std::vector<Point<dim> > &points,
173  std::vector<Tensor<1,dim> > &gradients,
174  const unsigned int component = 0) const;
175 
179  virtual double laplacian (const Point<dim> &p,
180  const unsigned int component = 0) const;
181 
185  virtual void laplacian_list (const std::vector<Point<dim> > &points,
186  std::vector<double> &values,
187  const unsigned int component = 0) const;
188  private:
189  const double offset;
190  };
191 
192 
193 
202  template<int dim>
203  class CosineFunction : public Function<dim>
204  {
205  public:
210  CosineFunction (const unsigned int n_components = 1);
211 
212  virtual double value (const Point<dim> &p,
213  const unsigned int component = 0) const;
214 
215  virtual void value_list (const std::vector<Point<dim> > &points,
216  std::vector<double> &values,
217  const unsigned int component = 0) const;
218 
219  virtual void vector_value_list (const std::vector<Point<dim> > &points,
220  std::vector<Vector<double> > &values) const;
221 
222  virtual Tensor<1,dim> gradient (const Point<dim> &p,
223  const unsigned int component = 0) const;
224 
225  virtual void gradient_list (const std::vector<Point<dim> > &points,
226  std::vector<Tensor<1,dim> > &gradients,
227  const unsigned int component = 0) const;
228 
229  virtual double laplacian (const Point<dim> &p,
230  const unsigned int component = 0) const;
231 
232  virtual void laplacian_list (const std::vector<Point<dim> > &points,
233  std::vector<double> &values,
234  const unsigned int component = 0) const;
235 
239  virtual SymmetricTensor<2,dim> hessian (const Point<dim> &p,
240  const unsigned int component = 0) const;
241 
245  virtual void hessian_list (const std::vector<Point<dim> > &points,
246  std::vector<SymmetricTensor<2,dim> > &hessians,
247  const unsigned int component = 0) const;
248  };
249 
250 
251 
263  template<int dim>
264  class CosineGradFunction : public Function<dim>
265  {
266  public:
271 
272  virtual double value (const Point<dim> &p,
273  const unsigned int component) const;
274  virtual void vector_value (const Point<dim> &p,
275  Vector<double> &values) const;
276  virtual void value_list (const std::vector<Point<dim> > &points,
277  std::vector<double> &values,
278  const unsigned int component) const;
279 
280  virtual void vector_value_list (const std::vector<Point<dim> > &points,
281  std::vector<Vector<double> > &values) const;
282 
283  virtual Tensor<1,dim> gradient (const Point<dim> &p,
284  const unsigned int component) const;
285 
286  virtual void gradient_list (const std::vector<Point<dim> > &points,
287  std::vector<Tensor<1,dim> > &gradients,
288  const unsigned int component) const;
289 
290  virtual void vector_gradient_list (const std::vector<Point<dim> > &points,
291  std::vector<std::vector<Tensor<1,dim> > > &gradients) const;
292 
293  virtual double laplacian (const Point<dim> &p,
294  const unsigned int component) const;
295  };
296 
297 
298 
305  template<int dim>
306  class ExpFunction : public Function<dim>
307  {
308  public:
312  virtual double value (const Point<dim> &p,
313  const unsigned int component = 0) const;
314 
318  virtual void value_list (const std::vector<Point<dim> > &points,
319  std::vector<double> &values,
320  const unsigned int component = 0) const;
321 
325  virtual Tensor<1,dim> gradient (const Point<dim> &p,
326  const unsigned int component = 0) const;
327 
331  virtual void gradient_list (const std::vector<Point<dim> > &points,
332  std::vector<Tensor<1,dim> > &gradients,
333  const unsigned int component = 0) const;
334 
338  virtual double laplacian (const Point<dim> &p,
339  const unsigned int component = 0) const;
340 
344  virtual void laplacian_list (const std::vector<Point<dim> > &points,
345  std::vector<double> &values,
346  const unsigned int component = 0) const;
347  };
348 
349 
350 
358  class LSingularityFunction : public Function<2>
359  {
360  public:
361  virtual double value (const Point<2> &p,
362  const unsigned int component = 0) const;
363 
364  virtual void value_list (const std::vector<Point<2> > &points,
365  std::vector<double> &values,
366  const unsigned int component = 0) const;
367 
368  virtual void vector_value_list (const std::vector<Point<2> > &points,
369  std::vector<Vector<double> > &values) const;
370 
371  virtual Tensor<1,2> gradient (const Point<2> &p,
372  const unsigned int component = 0) const;
373 
374  virtual void gradient_list (const std::vector<Point<2> > &points,
375  std::vector<Tensor<1,2> > &gradients,
376  const unsigned int component = 0) const;
377 
378  virtual void vector_gradient_list (const std::vector<Point<2> > &,
379  std::vector<std::vector<Tensor<1,2> > > &) const;
380 
381  virtual double laplacian (const Point<2> &p,
382  const unsigned int component = 0) const;
383 
384  virtual void laplacian_list (const std::vector<Point<2> > &points,
385  std::vector<double> &values,
386  const unsigned int component = 0) const;
387  };
388 
389 
390 
401  {
402  public:
407  virtual double value (const Point<2> &p,
408  const unsigned int component) const;
409 
410  virtual void value_list (const std::vector<Point<2> > &points,
411  std::vector<double> &values,
412  const unsigned int component) const;
413 
414  virtual void vector_value_list (const std::vector<Point<2> > &points,
415  std::vector<Vector<double> > &values) const;
416 
417  virtual Tensor<1,2> gradient (const Point<2> &p,
418  const unsigned int component) const;
419 
420  virtual void gradient_list (const std::vector<Point<2> > &points,
421  std::vector<Tensor<1,2> > &gradients,
422  const unsigned int component) const;
423 
424  virtual void vector_gradient_list (const std::vector<Point<2> > &,
425  std::vector<std::vector<Tensor<1,2> > > &) const;
426 
427  virtual double laplacian (const Point<2> &p,
428  const unsigned int component) const;
429 
430  virtual void laplacian_list (const std::vector<Point<2> > &points,
431  std::vector<double> &values,
432  const unsigned int component) const;
433  };
434 
435 
436 
443  template <int dim>
444  class SlitSingularityFunction : public Function<dim>
445  {
446  public:
447  virtual double value (const Point<dim> &p,
448  const unsigned int component = 0) const;
449 
450  virtual void value_list (const std::vector<Point<dim> > &points,
451  std::vector<double> &values,
452  const unsigned int component = 0) const;
453 
454  virtual void vector_value_list (const std::vector<Point<dim> > &points,
455  std::vector<Vector<double> > &values) const;
456 
457  virtual Tensor<1,dim> gradient (const Point<dim> &p,
458  const unsigned int component = 0) const;
459 
460  virtual void gradient_list (const std::vector<Point<dim> > &points,
461  std::vector<Tensor<1,dim> > &gradients,
462  const unsigned int component = 0) const;
463 
464  virtual void vector_gradient_list (const std::vector<Point<dim> > &,
465  std::vector<std::vector<Tensor<1,dim> > > &) const;
466 
467  virtual double laplacian (const Point<dim> &p,
468  const unsigned int component = 0) const;
469 
470  virtual void laplacian_list (const std::vector<Point<dim> > &points,
471  std::vector<double> &values,
472  const unsigned int component = 0) const;
473  };
474 
475 
483  {
484  public:
485  virtual double value (const Point<2> &p,
486  const unsigned int component = 0) const;
487 
488  virtual void value_list (const std::vector<Point<2> > &points,
489  std::vector<double> &values,
490  const unsigned int component = 0) const;
491 
492  virtual void vector_value_list (const std::vector<Point<2> > &points,
493  std::vector<Vector<double> > &values) const;
494 
495  virtual Tensor<1,2> gradient (const Point<2> &p,
496  const unsigned int component = 0) const;
497 
498  virtual void gradient_list (const std::vector<Point<2> > &points,
499  std::vector<Tensor<1,2> > &gradients,
500  const unsigned int component = 0) const;
501 
502  virtual void vector_gradient_list (const std::vector<Point<2> > &,
503  std::vector<std::vector<Tensor<1,2> > > &) const;
504 
505  virtual double laplacian (const Point<2> &p,
506  const unsigned int component = 0) const;
507 
508  virtual void laplacian_list (const std::vector<Point<2> > &points,
509  std::vector<double> &values,
510  const unsigned int component = 0) const;
511  };
512 
513 
514 
530  template<int dim>
531  class JumpFunction : public Function<dim>
532  {
533  public:
539  const double steepness);
540 
544  virtual double value (const Point<dim> &p,
545  const unsigned int component = 0) const;
546 
550  virtual void value_list (const std::vector<Point<dim> > &points,
551  std::vector<double> &values,
552  const unsigned int component = 0) const;
553 
557  virtual Tensor<1,dim> gradient (const Point<dim> &p,
558  const unsigned int component = 0) const;
559 
563  virtual void gradient_list (const std::vector<Point<dim> > &points,
564  std::vector<Tensor<1,dim> > &gradients,
565  const unsigned int component = 0) const;
566 
570  virtual double laplacian (const Point<dim> &p,
571  const unsigned int component = 0) const;
572 
576  virtual void laplacian_list (const std::vector<Point<dim> > &points,
577  std::vector<double> &values,
578  const unsigned int component = 0) const;
579 
586  std::size_t memory_consumption () const;
587 
588  protected:
593 
597  const double steepness;
598 
602  double angle;
603 
607  double sine;
608 
612  double cosine;
613  };
614 
615 
616 
629  template <int dim>
630  class FourierCosineFunction : public Function<dim>
631  {
632  public:
638 
645  virtual double value (const Point<dim> &p,
646  const unsigned int component = 0) const;
647 
652  virtual Tensor<1,dim> gradient (const Point<dim> &p,
653  const unsigned int component = 0) const;
654 
658  virtual double laplacian (const Point<dim> &p,
659  const unsigned int component = 0) const;
660  private:
665  };
666 
667 
668 
681  template <int dim>
682  class FourierSineFunction : public Function<dim>
683  {
684  public:
690 
697  virtual double value (const Point<dim> &p,
698  const unsigned int component = 0) const;
699 
704  virtual Tensor<1,dim> gradient (const Point<dim> &p,
705  const unsigned int component = 0) const;
706 
710  virtual double laplacian (const Point<dim> &p,
711  const unsigned int component = 0) const;
712  private:
717  };
718 
719 
729  template <int dim>
730  class FourierSineSum : public Function<dim>
731  {
732  public:
737  FourierSineSum (const std::vector<Point<dim> > &fourier_coefficients,
738  const std::vector<double> &weights);
739 
746  virtual double value (const Point<dim> &p,
747  const unsigned int component = 0) const;
748 
753  virtual Tensor<1,dim> gradient (const Point<dim> &p,
754  const unsigned int component = 0) const;
755 
759  virtual double laplacian (const Point<dim> &p,
760  const unsigned int component = 0) const;
761  private:
765  const std::vector<Point<dim> > fourier_coefficients;
766  const std::vector<double> weights;
767  };
768 
769 
770 
781  template <int dim>
782  class FourierCosineSum : public Function<dim>
783  {
784  public:
789  FourierCosineSum (const std::vector<Point<dim> > &fourier_coefficients,
790  const std::vector<double> &weights);
791 
798  virtual double value (const Point<dim> &p,
799  const unsigned int component = 0) const;
800 
805  virtual Tensor<1,dim> gradient (const Point<dim> &p,
806  const unsigned int component = 0) const;
807 
811  virtual double laplacian (const Point<dim> &p,
812  const unsigned int component = 0) const;
813 
814  private:
818  const std::vector<Point<dim> > fourier_coefficients;
819  const std::vector<double> weights;
820  };
821 
822 
831  template <int dim>
832  class CutOffFunctionBase : public Function<dim>
833  {
834  public:
839  static const unsigned int no_component = numbers::invalid_unsigned_int;
840 
847  CutOffFunctionBase (const double radius = 1.,
848  const Point<dim> = Point<dim>(),
849  const unsigned int n_components = 1,
850  const unsigned int select = CutOffFunctionBase<dim>::no_component);
851 
855  void new_center (const Point<dim> &p);
856 
860  void new_radius (const double r);
861 
862  protected:
867 
871  double radius;
872 
877  const unsigned int selected;
878  };
879 
880 
881 
891  template<int dim>
893  {
894  public:
901  CutOffFunctionLinfty (const double radius = 1.,
902  const Point<dim> = Point<dim>(),
903  const unsigned int n_components = 1,
904  const unsigned int select = CutOffFunctionBase<dim>::no_component);
905 
909  virtual double value (const Point<dim> &p,
910  const unsigned int component = 0) const;
911 
915  virtual void value_list (const std::vector<Point<dim> > &points,
916  std::vector<double> &values,
917  const unsigned int component = 0) const;
918 
922  virtual void vector_value_list (const std::vector<Point<dim> > &points,
923  std::vector<Vector<double> > &values) const;
924  };
925 
926 
936  template<int dim>
938  {
939  public:
947  CutOffFunctionW1 (const double radius = 1.,
948  const Point<dim> = Point<dim>(),
949  const unsigned int n_components = 1,
950  const unsigned int select = CutOffFunctionBase<dim>::no_component);
951 
955  virtual double value (const Point<dim> &p,
956  const unsigned int component = 0) const;
957 
961  virtual void value_list (const std::vector<Point<dim> > &points,
962  std::vector<double> &values,
963  const unsigned int component = 0) const;
964 
968  virtual void vector_value_list (const std::vector<Point<dim> > &points,
969  std::vector<Vector<double> > &values) const;
970  };
971 
972 
983  template<int dim>
985  {
986  public:
994  CutOffFunctionCinfty (const double radius = 1.,
995  const Point<dim> = Point<dim>(),
996  const unsigned int n_components = 1,
997  const unsigned int select = CutOffFunctionBase<dim>::no_component);
998 
1002  virtual double value (const Point<dim> &p,
1003  const unsigned int component = 0) const;
1004 
1008  virtual void value_list (const std::vector<Point<dim> > &points,
1009  std::vector<double> &values,
1010  const unsigned int component = 0) const;
1011 
1015  virtual void vector_value_list (const std::vector<Point<dim> > &points,
1016  std::vector<Vector<double> > &values) const;
1017 
1021  virtual Tensor<1,dim> gradient (const Point<dim> &p,
1022  const unsigned int component = 0) const;
1023  };
1024 
1025 
1026 
1040  template <int dim>
1041  class Monomial : public Function<dim>
1042  {
1043  public:
1051  const unsigned int n_components = 1);
1052 
1056  virtual double value (const Point<dim> &p,
1057  const unsigned int component = 0) const;
1058 
1065  virtual void vector_value (const Point<dim> &p,
1066  Vector<double> &values) const;
1067 
1071  virtual void value_list (const std::vector<Point<dim> > &points,
1072  std::vector<double> &values,
1073  const unsigned int component = 0) const;
1074 
1078  virtual Tensor<1,dim> gradient (const Point<dim> &p,
1079  const unsigned int component = 0) const;
1080 
1081  private:
1086  };
1087 
1088 
1089 
1122  template <int dim>
1124  {
1125  public:
1139  InterpolatedTensorProductGridData (const std_cxx11::array<std::vector<double>,dim> &coordinate_values,
1141 
1152  virtual
1153  double
1154  value (const Point<dim> &p,
1155  const unsigned int component = 0) const;
1156 
1168  virtual
1170  gradient (const Point<dim> &p,
1171  const unsigned int component = 0) const;
1172 
1173  protected:
1178 
1182  const std_cxx11::array<std::vector<double>,dim> coordinate_values;
1183 
1188  };
1189 
1190 
1226  template <int dim>
1228  {
1229  public:
1244  InterpolatedUniformGridData (const std_cxx11::array<std::pair<double,double>,dim> &interval_endpoints,
1245  const std_cxx11::array<unsigned int,dim> &n_subintervals,
1247 
1258  virtual
1259  double
1260  value (const Point<dim> &p,
1261  const unsigned int component = 0) const;
1262 
1263  private:
1267  const std_cxx11::array<std::pair<double,double>,dim> interval_endpoints;
1268 
1272  const std_cxx11::array<unsigned int,dim> n_subintervals;
1273 
1278  };
1279 
1280 
1293  template <int dim>
1294  class Polynomial : public Function<dim>
1295  {
1296  public:
1307  const std::vector<double> &coefficients);
1308 
1312  virtual double value (const Point<dim> &p,
1313  const unsigned int component = 0) const;
1314 
1315 
1319  virtual void value_list (const std::vector<Point<dim> > &points,
1320  std::vector<double> &values,
1321  const unsigned int component = 0) const;
1322 
1326  virtual Tensor<1,dim> gradient (const Point<dim> &p,
1327  const unsigned int component = 0) const;
1328 
1329  private:
1330 
1335 
1339  const std::vector<double> coefficients;
1340  };
1341 
1342 
1343 
1344 }
1345 DEAL_II_NAMESPACE_CLOSE
1346 
1347 #endif
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
const Tensor< 1, dim > fourier_coefficients
Definition: function_lib.h:664
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
static const unsigned int invalid_unsigned_int
Definition: types.h:170
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component) const
virtual double value(const Point< dim > &p, const unsigned int component) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
TableIndices< dim > table_index_of_point(const Point< dim > &p) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const
Monomial(const Tensor< 1, dim > &exponents, const unsigned int n_components=1)
Polynomial(const Table< 2, double > &exponents, const std::vector< double > &coefficients)
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
CutOffFunctionLinfty(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component)
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
InterpolatedUniformGridData(const std_cxx11::array< std::pair< double, double >, dim > &interval_endpoints, const std_cxx11::array< unsigned int, dim > &n_subintervals, const Table< dim, double > &data_values)
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component) const
const Tensor< 1, dim > fourier_coefficients
Definition: function_lib.h:716
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
Definition: function_lib.cc:53
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const
Definition: function_lib.cc:43
const std_cxx11::array< unsigned int, dim > n_subintervals
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const
CosineFunction(const unsigned int n_components=1)
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
Definition: function_lib.cc:34
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
const unsigned int n_components
Definition: function.h:144
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual SymmetricTensor< 2, dim > hessian(const Point< dim > &p, const unsigned int component=0) const
void new_center(const Point< dim > &p)
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
std::size_t memory_consumption() const
const std_cxx11::array< std::vector< double >, dim > coordinate_values
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
const Table< 2, double > exponents
CutOffFunctionBase(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component)
const std::vector< Point< dim > > fourier_coefficients
Definition: function_lib.h:818
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
virtual void hessian_list(const std::vector< Point< dim > > &points, std::vector< SymmetricTensor< 2, dim > > &hessians, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
const std::vector< double > coefficients
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
const unsigned int selected
Definition: function_lib.h:877
virtual double value(const Point< dim > &p, const unsigned int component=0) const
const Point< dim > direction
Definition: function_lib.h:592
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
const Table< dim, double > data_values
const std::vector< Point< dim > > fourier_coefficients
Definition: function_lib.h:765
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
PillowFunction(const double offset=0.)
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual double value(const Point< dim > &p, const unsigned int component=0) const
static const unsigned int no_component
Definition: function_lib.h:839
virtual double value(const Point< dim > &p, const unsigned int component=0) const
FourierCosineFunction(const Tensor< 1, dim > &fourier_coefficients)
virtual double value(const Point< dim > &p, const unsigned int component=0) const
FourierCosineSum(const std::vector< Point< dim > > &fourier_coefficients, const std::vector< double > &weights)
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
FourierSineFunction(const Tensor< 1, dim > &fourier_coefficients)
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
const Tensor< 1, dim > exponents
CutOffFunctionCinfty(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component)
CutOffFunctionW1(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component)
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
Definition: function_lib.cc:94
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
Definition: function_lib.cc:79
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
const std_cxx11::array< std::pair< double, double >, dim > interval_endpoints
virtual double value(const Point< dim > &p, const unsigned int component=0) const
InterpolatedTensorProductGridData(const std_cxx11::array< std::vector< double >, dim > &coordinate_values, const Table< dim, double > &data_values)
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
Definition: function_lib.cc:70
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
JumpFunction(const Point< dim > &direction, const double steepness)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
FourierSineSum(const std::vector< Point< dim > > &fourier_coefficients, const std::vector< double > &weights)