]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
FourierSine and FourierCosineFunction
authorwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 18 Jul 2001 09:50:43 +0000 (09:50 +0000)
committerwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 18 Jul 2001 09:50:43 +0000 (09:50 +0000)
git-svn-id: https://svn.dealii.org/trunk@4849 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/base/include/base/function_lib.h
deal.II/base/source/function_lib.cc
deal.II/doc/news/2001/c-3-1.html

index 1172b3a2e82841a02da1ec47ff18488af247f266..10c42f8e52c1fc9012bd13ff414baab104f2c24a 100644 (file)
@@ -525,5 +525,123 @@ class JumpFunction : public Function<dim>
 
 
 
+/**
+ * Given a wavenumber vector generate a cosine function. The
+ * wavenumber coefficient is given as a @p{d}-dimensional point @p{k}
+ * in Fourier space, and the function is then recovered as @p{f(x) =
+ * \prod_i cos(k_i x_i) = Re(\exp(i k.x))}.
+ *
+ * The class has its name from the fact that it resembles one
+ * component of a Fourier cosine decomposition.
+ *
+ * @author Wolfgang Bangerth, 2001
+ */
+template <int dim>
+class FourierCosineFunction : public Function<dim> 
+{
+  public:
+                                    /**
+                                     * Constructor. Take the Fourier
+                                     * coefficients in each space
+                                     * direction as argument.
+                                     */
+    FourierCosineFunction (const Point<dim> &fourier_coefficients);
+    
+                                    /**
+                                     * Return the value of the
+                                     * function at the given
+                                     * point. Unless there is only
+                                     * one component (i.e. the
+                                     * function is scalar), you
+                                     * should state the component you
+                                     * want to have evaluated; it
+                                     * defaults to zero, i.e. the
+                                     * first component.
+                                     */
+    virtual double value (const Point<dim>   &p,
+                         const unsigned int  component = 0) const;
+
+                                    /**
+                                     * Return the gradient of the
+                                     * specified component of the
+                                     * function at the given point.
+                                     */
+    virtual Tensor<1,dim> gradient (const Point<dim>   &p,
+                                   const unsigned int  component = 0) const;
+
+                                    /**
+                                     * Compute the Laplacian of a
+                                     * given component at point @p{p}.
+                                     */
+    virtual double laplacian (const Point<dim>   &p,
+                             const unsigned int  component = 0) const;
+  private:
+                                    /**
+                                     * Stored Fourier coefficients.
+                                     */
+    const Point<dim> fourier_coefficients;
+};
+
+
+
+/**
+ * Given a wavenumber vector generate a sine function. The
+ * wavenumber coefficient is given as a @p{d}-dimensional point @p{k}
+ * in Fourier space, and the function is then recovered as @p{f(x) =
+ * \prod_i sin(k_i x_i) = Im(\exp(i k.x))}.
+ *
+ * The class has its name from the fact that it resembles one
+ * component of a Fourier sine decomposition.
+ *
+ * @author Wolfgang Bangerth, 2001
+ */
+template <int dim>
+class FourierSineFunction : public Function<dim> 
+{
+  public:
+                                    /**
+                                     * Constructor. Take the Fourier
+                                     * coefficients in each space
+                                     * direction as argument.
+                                     */
+    FourierSineFunction (const Point<dim> &fourier_coefficients);
+    
+                                    /**
+                                     * Return the value of the
+                                     * function at the given
+                                     * point. Unless there is only
+                                     * one component (i.e. the
+                                     * function is scalar), you
+                                     * should state the component you
+                                     * want to have evaluated; it
+                                     * defaults to zero, i.e. the
+                                     * first component.
+                                     */
+    virtual double value (const Point<dim>   &p,
+                         const unsigned int  component = 0) const;
+
+                                    /**
+                                     * Return the gradient of the
+                                     * specified component of the
+                                     * function at the given point.
+                                     */
+    virtual Tensor<1,dim> gradient (const Point<dim>   &p,
+                                   const unsigned int  component = 0) const;
+
+                                    /**
+                                     * Compute the Laplacian of a
+                                     * given component at point @p{p}.
+                                     */
+    virtual double laplacian (const Point<dim>   &p,
+                             const unsigned int  component = 0) const;
+  private:
+                                    /**
+                                     * Stored Fourier coefficients.
+                                     */
+    const Point<dim> fourier_coefficients;
+};
+
+
+
 
 #endif
index 36c52c679ae6ef4ba3e5794cc66de1c2b818fcf1..5da39ec8f9dd40dd466b29f9223ebbaf1cfce83e 100644 (file)
@@ -1112,6 +1112,148 @@ JumpFunction<dim>::memory_consumption () const
 };
 
 
+
+
+
+/* ---------------------- FourierSineFunction ----------------------- */
+
+
+template <int dim>
+FourierCosineFunction<dim>::
+FourierCosineFunction (const Point<dim> &fourier_coefficients)
+               :
+               Function<dim> (1),
+                fourier_coefficients (fourier_coefficients)
+{};
+
+
+
+template <int dim>
+double
+FourierCosineFunction<dim>::value (const Point<dim>   &p,
+                                  const unsigned int  component) const
+{
+  Assert (component==0, ExcIndexRange(component,0,1));
+  double val=1;
+  for (unsigned int i=0; i<dim; ++i)
+    val *= std::cos(fourier_coefficients[i] * p[i]);
+  return val;
+};
+
+
+
+template <int dim>
+Tensor<1,dim>
+FourierCosineFunction<dim>::gradient (const Point<dim>   &p,
+                                     const unsigned int  component) const
+{
+  Assert (component==0, ExcIndexRange(component,0,1));
+  Tensor<1,dim> grad;
+  for (unsigned int i=0; i<dim; ++i)
+    grad[i] = 1;
+
+  for (unsigned int i=0; i<dim; ++i)
+    {
+      const double cos_i = std::cos(fourier_coefficients[i] * p[i]);
+      const double sin_i = std::sin(fourier_coefficients[i] * p[i]);
+      
+      for (unsigned int d=0; d<dim; ++d)
+       if (d==i)
+         grad[d] *= - fourier_coefficients[i] * sin_i;
+       else
+         grad[d] *= cos_i;
+    };
+  
+  return grad;
+};
+
+
+
+template <int dim>
+double
+FourierCosineFunction<dim>::laplacian (const Point<dim>   &p,
+                                      const unsigned int  component) const
+{
+  Assert (component==0, ExcIndexRange(component,0,1));
+  double val = -(fourier_coefficients*fourier_coefficients);
+  for (unsigned int i=0; i<dim; ++i)
+    val *= std::cos(fourier_coefficients[i] * p[i]);
+  return val;
+};
+
+
+
+
+/* ---------------------- FourierSineFunction ----------------------- */
+
+
+
+template <int dim>
+FourierSineFunction<dim>::
+FourierSineFunction (const Point<dim> &fourier_coefficients)
+               :
+               Function<dim> (1),
+                fourier_coefficients (fourier_coefficients)
+{};
+
+
+
+template <int dim>
+double
+FourierSineFunction<dim>::value (const Point<dim>   &p,
+                                const unsigned int  component) const
+{
+  Assert (component==0, ExcIndexRange(component,0,1));
+  double val=1;
+  for (unsigned int i=0; i<dim; ++i)
+    val *= std::sin(fourier_coefficients[i] * p[i]);
+  return val;
+};
+
+
+
+template <int dim>
+Tensor<1,dim>
+FourierSineFunction<dim>::gradient (const Point<dim>   &p,
+                                   const unsigned int  component) const
+{
+  Assert (component==0, ExcIndexRange(component,0,1));
+  Tensor<1,dim> grad;
+  for (unsigned int i=0; i<dim; ++i)
+    grad[i] = 1;
+
+  for (unsigned int i=0; i<dim; ++i)
+    {
+      const double cos_i = std::cos(fourier_coefficients[i] * p[i]);
+      const double sin_i = std::sin(fourier_coefficients[i] * p[i]);
+      
+      for (unsigned int d=0; d<dim; ++d)
+       if (d==i)
+         grad[d] *= fourier_coefficients[i] * cos_i;
+       else
+         grad[d] *= sin_i;
+    };
+  
+  return grad;
+};
+
+
+
+template <int dim>
+double
+FourierSineFunction<dim>::laplacian (const Point<dim>   &p,
+                                    const unsigned int  component) const
+{
+  Assert (component==0, ExcIndexRange(component,0,1));
+  double val = -(fourier_coefficients*fourier_coefficients);
+  for (unsigned int i=0; i<dim; ++i)
+    val *= std::sin(fourier_coefficients[i] * p[i]);
+  return val;
+};
+
+
+
+
 template class SquareFunction<1>;
 template class SquareFunction<2>;
 template class SquareFunction<3>;
@@ -1130,3 +1272,10 @@ template class ExpFunction<3>;
 template class JumpFunction<1>;
 template class JumpFunction<2>;
 template class JumpFunction<3>;
+template class FourierCosineFunction<1>;
+template class FourierCosineFunction<2>;
+template class FourierCosineFunction<3>;
+template class FourierSineFunction<1>;
+template class FourierSineFunction<2>;
+template class FourierSineFunction<3>;
+
index 548eecb0f2091a22609928e24642595a126ba0c0..d1595b8c4472656320131142cbaf3cb60548381a 100644 (file)
@@ -158,9 +158,14 @@ documentation, etc</a>.
 <h3>base</h3>
 
 <ol>
-
   <li> <p>
+       New: classes <code class="class">FourierSineFunction</code> and
+       <code class="class">FourierCosineFunction</code>, resembling
+       one mode of a Fourier decomposition.
+       <br>
+       (GK 2001/07/18)
 
+  <li> <p>
        New: class <code class="class">vector2d</code> was introduced
        in analogy to STL class <code class="class">vector</code>. The
        bew class provides a two-dimensional array and replaces the use

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.