Reference documentation for deal.II version 8.1.0
function_lib.h
1 // ---------------------------------------------------------------------
2 // @f$Id: function_lib.h 30036 2013-07-18 16:55:32Z maier @f$
3 //
4 // Copyright (C) 1999 - 2013 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 #ifndef __deal2__function_lib_h
18 #define __deal2__function_lib_h
19 
20 
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/function.h>
23 #include <deal.II/base/point.h>
24 
26 
36 namespace Functions
37 {
38 
39 
50  template<int dim>
51  class SquareFunction : public Function<dim>
52  {
53  public:
54  virtual double value (const Point<dim> &p,
55  const unsigned int component = 0) const;
56  virtual void vector_value (const Point<dim> &p,
57  Vector<double> &values) const;
58  virtual void value_list (const std::vector<Point<dim> > &points,
59  std::vector<double> &values,
60  const unsigned int component = 0) const;
61  virtual Tensor<1,dim> gradient (const Point<dim> &p,
62  const unsigned int component = 0) const;
63  virtual void vector_gradient (const Point<dim> &p,
64  std::vector<Tensor<1,dim> > &gradient) const;
65  virtual void gradient_list (const std::vector<Point<dim> > &points,
66  std::vector<Tensor<1,dim> > &gradients,
67  const unsigned int component = 0) const;
68  virtual double laplacian (const Point<dim> &p,
69  const unsigned int component = 0) const;
70  virtual void laplacian_list (const std::vector<Point<dim> > &points,
71  std::vector<double> &values,
72  const unsigned int component = 0) const;
73  };
74 
75 
76 
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:
149  PillowFunction (const double offset=0.);
150 
154  virtual double value (const Point<dim> &p,
155  const unsigned int component = 0) const;
156 
160  virtual void value_list (const std::vector<Point<dim> > &points,
161  std::vector<double> &values,
162  const unsigned int component = 0) const;
163 
167  virtual Tensor<1,dim> gradient (const Point<dim> &p,
168  const unsigned int component = 0) const;
169 
173  virtual void gradient_list (const std::vector<Point<dim> > &points,
174  std::vector<Tensor<1,dim> > &gradients,
175  const unsigned int component = 0) const;
176 
180  virtual double laplacian (const Point<dim> &p,
181  const unsigned int component = 0) const;
182 
186  virtual void laplacian_list (const std::vector<Point<dim> > &points,
187  std::vector<double> &values,
188  const unsigned int component = 0) const;
189  private:
190  const double offset;
191  };
192 
193 
194 
203  template<int dim>
204  class CosineFunction : public Function<dim>
205  {
206  public:
214  CosineFunction (const unsigned int n_components = 1);
215 
216  virtual double value (const Point<dim> &p,
217  const unsigned int component = 0) const;
218 
219  virtual void value_list (const std::vector<Point<dim> > &points,
220  std::vector<double> &values,
221  const unsigned int component = 0) const;
222 
223  virtual void vector_value_list (const std::vector<Point<dim> > &points,
224  std::vector<Vector<double> > &values) const;
225 
226  virtual Tensor<1,dim> gradient (const Point<dim> &p,
227  const unsigned int component = 0) const;
228 
229  virtual void gradient_list (const std::vector<Point<dim> > &points,
230  std::vector<Tensor<1,dim> > &gradients,
231  const unsigned int component = 0) const;
232 
233  virtual double laplacian (const Point<dim> &p,
234  const unsigned int component = 0) const;
235 
236  virtual void laplacian_list (const std::vector<Point<dim> > &points,
237  std::vector<double> &values,
238  const unsigned int component = 0) const;
239 
244  virtual Tensor<2,dim> hessian (const Point<dim> &p,
245  const unsigned int component = 0) const;
246 
251  virtual void hessian_list (const std::vector<Point<dim> > &points,
252  std::vector<Tensor<2,dim> > &hessians,
253  const unsigned int component = 0) const;
254  };
255 
256 
257 
270  template<int dim>
271  class CosineGradFunction : public Function<dim>
272  {
273  public:
279 
280  virtual double value (const Point<dim> &p,
281  const unsigned int component) const;
282  virtual void vector_value (const Point<dim> &p,
283  Vector<double> &values) const;
284  virtual void value_list (const std::vector<Point<dim> > &points,
285  std::vector<double> &values,
286  const unsigned int component) const;
287 
288  virtual void vector_value_list (const std::vector<Point<dim> > &points,
289  std::vector<Vector<double> > &values) const;
290 
291  virtual Tensor<1,dim> gradient (const Point<dim> &p,
292  const unsigned int component) const;
293 
294  virtual void gradient_list (const std::vector<Point<dim> > &points,
295  std::vector<Tensor<1,dim> > &gradients,
296  const unsigned int component) const;
297 
298  virtual void vector_gradient_list (const std::vector<Point<dim> > &points,
299  std::vector<std::vector<Tensor<1,dim> > > &gradients) const;
300 
301  virtual double laplacian (const Point<dim> &p,
302  const unsigned int component) const;
303  };
304 
305 
306 
313  template<int dim>
314  class ExpFunction : public Function<dim>
315  {
316  public:
320  virtual double value (const Point<dim> &p,
321  const unsigned int component = 0) const;
322 
326  virtual void value_list (const std::vector<Point<dim> > &points,
327  std::vector<double> &values,
328  const unsigned int component = 0) const;
329 
333  virtual Tensor<1,dim> gradient (const Point<dim> &p,
334  const unsigned int component = 0) const;
335 
339  virtual void gradient_list (const std::vector<Point<dim> > &points,
340  std::vector<Tensor<1,dim> > &gradients,
341  const unsigned int component = 0) const;
342 
346  virtual double laplacian (const Point<dim> &p,
347  const unsigned int component = 0) const;
348 
352  virtual void laplacian_list (const std::vector<Point<dim> > &points,
353  std::vector<double> &values,
354  const unsigned int component = 0) const;
355  };
356 
357 
358 
366  class LSingularityFunction : public Function<2>
367  {
368  public:
369  virtual double value (const Point<2> &p,
370  const unsigned int component = 0) const;
371 
372  virtual void value_list (const std::vector<Point<2> > &points,
373  std::vector<double> &values,
374  const unsigned int component = 0) const;
375 
376  virtual void vector_value_list (const std::vector<Point<2> > &points,
377  std::vector<Vector<double> > &values) const;
378 
379  virtual Tensor<1,2> gradient (const Point<2> &p,
380  const unsigned int component = 0) const;
381 
382  virtual void gradient_list (const std::vector<Point<2> > &points,
383  std::vector<Tensor<1,2> > &gradients,
384  const unsigned int component = 0) const;
385 
386  virtual void vector_gradient_list (const std::vector<Point<2> > &,
387  std::vector<std::vector<Tensor<1,2> > > &) const;
388 
389  virtual double laplacian (const Point<2> &p,
390  const unsigned int component = 0) const;
391 
392  virtual void laplacian_list (const std::vector<Point<2> > &points,
393  std::vector<double> &values,
394  const unsigned int component = 0) const;
395  };
396 
397 
398 
409  {
410  public:
416  virtual double value (const Point<2> &p,
417  const unsigned int component) const;
418 
419  virtual void value_list (const std::vector<Point<2> > &points,
420  std::vector<double> &values,
421  const unsigned int component) const;
422 
423  virtual void vector_value_list (const std::vector<Point<2> > &points,
424  std::vector<Vector<double> > &values) const;
425 
426  virtual Tensor<1,2> gradient (const Point<2> &p,
427  const unsigned int component) const;
428 
429  virtual void gradient_list (const std::vector<Point<2> > &points,
430  std::vector<Tensor<1,2> > &gradients,
431  const unsigned int component) const;
432 
433  virtual void vector_gradient_list (const std::vector<Point<2> > &,
434  std::vector<std::vector<Tensor<1,2> > > &) const;
435 
436  virtual double laplacian (const Point<2> &p,
437  const unsigned int component) const;
438 
439  virtual void laplacian_list (const std::vector<Point<2> > &points,
440  std::vector<double> &values,
441  const unsigned int component) const;
442  };
443 
444 
445 
452  template <int dim>
453  class SlitSingularityFunction : public Function<dim>
454  {
455  public:
456  virtual double value (const Point<dim> &p,
457  const unsigned int component = 0) const;
458 
459  virtual void value_list (const std::vector<Point<dim> > &points,
460  std::vector<double> &values,
461  const unsigned int component = 0) const;
462 
463  virtual void vector_value_list (const std::vector<Point<dim> > &points,
464  std::vector<Vector<double> > &values) const;
465 
466  virtual Tensor<1,dim> gradient (const Point<dim> &p,
467  const unsigned int component = 0) const;
468 
469  virtual void gradient_list (const std::vector<Point<dim> > &points,
470  std::vector<Tensor<1,dim> > &gradients,
471  const unsigned int component = 0) const;
472 
473  virtual void vector_gradient_list (const std::vector<Point<dim> > &,
474  std::vector<std::vector<Tensor<1,dim> > > &) const;
475 
476  virtual double laplacian (const Point<dim> &p,
477  const unsigned int component = 0) const;
478 
479  virtual void laplacian_list (const std::vector<Point<dim> > &points,
480  std::vector<double> &values,
481  const unsigned int component = 0) const;
482  };
483 
484 
492  {
493  public:
494  virtual double value (const Point<2> &p,
495  const unsigned int component = 0) const;
496 
497  virtual void value_list (const std::vector<Point<2> > &points,
498  std::vector<double> &values,
499  const unsigned int component = 0) const;
500 
501  virtual void vector_value_list (const std::vector<Point<2> > &points,
502  std::vector<Vector<double> > &values) const;
503 
504  virtual Tensor<1,2> gradient (const Point<2> &p,
505  const unsigned int component = 0) const;
506 
507  virtual void gradient_list (const std::vector<Point<2> > &points,
508  std::vector<Tensor<1,2> > &gradients,
509  const unsigned int component = 0) const;
510 
511  virtual void vector_gradient_list (const std::vector<Point<2> > &,
512  std::vector<std::vector<Tensor<1,2> > > &) const;
513 
514  virtual double laplacian (const Point<2> &p,
515  const unsigned int component = 0) const;
516 
517  virtual void laplacian_list (const std::vector<Point<2> > &points,
518  std::vector<double> &values,
519  const unsigned int component = 0) const;
520  };
521 
522 
523 
539  template<int dim>
540  class JumpFunction : public Function<dim>
541  {
542  public:
549  const double steepness);
550 
554  virtual double value (const Point<dim> &p,
555  const unsigned int component = 0) const;
556 
561  virtual void value_list (const std::vector<Point<dim> > &points,
562  std::vector<double> &values,
563  const unsigned int component = 0) const;
564 
568  virtual Tensor<1,dim> gradient (const Point<dim> &p,
569  const unsigned int component = 0) const;
570 
574  virtual void gradient_list (const std::vector<Point<dim> > &points,
575  std::vector<Tensor<1,dim> > &gradients,
576  const unsigned int component = 0) const;
577 
581  virtual double laplacian (const Point<dim> &p,
582  const unsigned int component = 0) const;
583 
587  virtual void laplacian_list (const std::vector<Point<dim> > &points,
588  std::vector<double> &values,
589  const unsigned int component = 0) const;
590 
607  std::size_t memory_consumption () const;
608 
609  protected:
614 
619  const double steepness;
620 
624  double angle;
625 
629  double sine;
630 
634  double cosine;
635  };
636 
637 
638 
651  template <int dim>
652  class FourierCosineFunction : public Function<dim>
653  {
654  public:
661 
673  virtual double value (const Point<dim> &p,
674  const unsigned int component = 0) const;
675 
681  virtual Tensor<1,dim> gradient (const Point<dim> &p,
682  const unsigned int component = 0) const;
683 
688  virtual double laplacian (const Point<dim> &p,
689  const unsigned int component = 0) const;
690  private:
695  };
696 
697 
698 
711  template <int dim>
712  class FourierSineFunction : public Function<dim>
713  {
714  public:
721 
733  virtual double value (const Point<dim> &p,
734  const unsigned int component = 0) const;
735 
741  virtual Tensor<1,dim> gradient (const Point<dim> &p,
742  const unsigned int component = 0) const;
743 
748  virtual double laplacian (const Point<dim> &p,
749  const unsigned int component = 0) const;
750  private:
755  };
756 
757 
768  template <int dim>
769  class FourierSineSum : public Function<dim>
770  {
771  public:
777  FourierSineSum (const std::vector<Point<dim> > &fourier_coefficients,
778  const std::vector<double> &weights);
779 
791  virtual double value (const Point<dim> &p,
792  const unsigned int component = 0) const;
793 
799  virtual Tensor<1,dim> gradient (const Point<dim> &p,
800  const unsigned int component = 0) const;
801 
806  virtual double laplacian (const Point<dim> &p,
807  const unsigned int component = 0) const;
808  private:
813  const std::vector<Point<dim> > fourier_coefficients;
814  const std::vector<double> weights;
815  };
816 
817 
818 
829  template <int dim>
830  class FourierCosineSum : public Function<dim>
831  {
832  public:
838  FourierCosineSum (const std::vector<Point<dim> > &fourier_coefficients,
839  const std::vector<double> &weights);
840 
852  virtual double value (const Point<dim> &p,
853  const unsigned int component = 0) const;
854 
860  virtual Tensor<1,dim> gradient (const Point<dim> &p,
861  const unsigned int component = 0) const;
862 
867  virtual double laplacian (const Point<dim> &p,
868  const unsigned int component = 0) const;
869 
870  private:
875  const std::vector<Point<dim> > fourier_coefficients;
876  const std::vector<double> weights;
877  };
878 
879 
889  template <int dim>
890  class CutOffFunctionBase : public Function<dim>
891  {
892  public:
900  static const unsigned int no_component = numbers::invalid_unsigned_int;
901 
913  CutOffFunctionBase (const double radius = 1.,
914  const Point<dim> = Point<dim>(),
915  const unsigned int n_components = 1,
916  const unsigned int select = CutOffFunctionBase<dim>::no_component);
917 
922  void new_center (const Point<dim> &p);
923 
927  void new_radius (const double r);
928 
929  protected:
934 
938  double radius;
939 
945  const unsigned int selected;
946  };
947 
948 
949 
963  template<int dim>
965  {
966  public:
978  CutOffFunctionLinfty (const double radius = 1.,
979  const Point<dim> = Point<dim>(),
980  const unsigned int n_components = 1,
981  const unsigned int select = CutOffFunctionBase<dim>::no_component);
982 
986  virtual double value (const Point<dim> &p,
987  const unsigned int component = 0) const;
988 
992  virtual void value_list (const std::vector<Point<dim> > &points,
993  std::vector<double> &values,
994  const unsigned int component = 0) const;
995 
999  virtual void vector_value_list (const std::vector<Point<dim> > &points,
1000  std::vector<Vector<double> > &values) const;
1001  };
1002 
1003 
1013  template<int dim>
1015  {
1016  public:
1028  CutOffFunctionW1 (const double radius = 1.,
1029  const Point<dim> = Point<dim>(),
1030  const unsigned int n_components = 1,
1031  const unsigned int select = CutOffFunctionBase<dim>::no_component);
1032 
1036  virtual double value (const Point<dim> &p,
1037  const unsigned int component = 0) const;
1038 
1042  virtual void value_list (const std::vector<Point<dim> > &points,
1043  std::vector<double> &values,
1044  const unsigned int component = 0) const;
1045 
1049  virtual void vector_value_list (const std::vector<Point<dim> > &points,
1050  std::vector<Vector<double> > &values) const;
1051  };
1052 
1053 
1064  template<int dim>
1066  {
1067  public:
1079  CutOffFunctionCinfty (const double radius = 1.,
1080  const Point<dim> = Point<dim>(),
1081  const unsigned int n_components = 1,
1082  const unsigned int select = CutOffFunctionBase<dim>::no_component);
1083 
1087  virtual double value (const Point<dim> &p,
1088  const unsigned int component = 0) const;
1089 
1093  virtual void value_list (const std::vector<Point<dim> > &points,
1094  std::vector<double> &values,
1095  const unsigned int component = 0) const;
1096 
1100  virtual void vector_value_list (const std::vector<Point<dim> > &points,
1101  std::vector<Vector<double> > &values) const;
1102 
1106  virtual Tensor<1,dim> gradient (const Point<dim> &p,
1107  const unsigned int component = 0) const;
1108  };
1109 
1110 
1111 
1124  template <int dim>
1125  class Monomial : public Function<dim>
1126  {
1127  public:
1138  const unsigned int n_components = 1);
1139 
1143  virtual double value (const Point<dim> &p,
1144  const unsigned int component = 0) const;
1145 
1155  virtual void vector_value (const Point<dim> &p,
1156  Vector<double> &values) const;
1157 
1161  virtual void value_list (const std::vector<Point<dim> > &points,
1162  std::vector<double> &values,
1163  const unsigned int component = 0) const;
1164 
1168  virtual Tensor<1,dim> gradient (const Point<dim> &p,
1169  const unsigned int component = 0) const;
1170 
1171  private:
1176  };
1177 
1178 }
1179 DEAL_II_NAMESPACE_CLOSE
1180 
1181 #endif
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
const Point< dim > fourier_coefficients
Definition: function_lib.h:694
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
static const unsigned int invalid_unsigned_int
Definition: types.h:191
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, 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 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
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) const
Monomial(const Tensor< 1, dim > &exponents, const unsigned int n_components=1)
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, 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 double laplacian(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 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
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 Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_gradient(const Point< dim > &p, std::vector< Tensor< 1, dim > > &gradient) const
virtual double value(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
std::size_t memory_consumption() const
virtual void vector_gradient_list(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const
virtual double value(const Point< dim > &p, 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
CosineFunction(const unsigned int n_components=1)
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 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 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 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
void new_center(const Point< dim > &p)
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
virtual void vector_gradient_list(const std::vector< Point< dim > > &, std::vector< std::vector< Tensor< 1, dim > > > &) const
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 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 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
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component) 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 void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
CutOffFunctionBase(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
FourierCosineFunction(const Point< dim > &fourier_coefficients)
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
const std::vector< Point< dim > > fourier_coefficients
Definition: function_lib.h:875
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
const unsigned int n_components
Definition: function.h:130
const Point< dim > fourier_coefficients
Definition: function_lib.h:754
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
FourierSineFunction(const Point< 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
virtual double value(const Point< dim > &p, const unsigned int component=0) const
const unsigned int selected
Definition: function_lib.h:945
virtual void vector_gradient_list(const std::vector< Point< dim > > &, std::vector< std::vector< Tensor< 1, dim > > > &) 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
const Point< dim > direction
Definition: function_lib.h:613
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const
const std::vector< Point< dim > > fourier_coefficients
Definition: function_lib.h:813
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, 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
PillowFunction(const double offset=0.)
virtual Tensor< 2, dim > hessian(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 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
static const unsigned int no_component
Definition: function_lib.h:900
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 laplacian(const Point< dim > &p, const unsigned int component) const
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 hessian_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 2, dim > > &hessians, 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
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 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) 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 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
void new_radius(const double r)
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component) const
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 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
JumpFunction(const Point< dim > &direction, const double steepness)
FourierSineSum(const std::vector< Point< dim > > &fourier_coefficients, const std::vector< double > &weights)