+/**
+ * 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
};
+
+
+
+/* ---------------------- 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>;
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>;
+
<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