]> https://gitweb.dealii.org/ - dealii.git/commitdiff
cutt-off functions
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Thu, 11 Oct 2001 21:34:45 +0000 (21:34 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Thu, 11 Oct 2001 21:34:45 +0000 (21:34 +0000)
git-svn-id: https://svn.dealii.org/trunk@5144 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/base/include/base/function_lib.h
deal.II/base/source/function_lib_cutoff.cc [new file with mode: 0644]

index 3371177c1c9e40ce102cca52c4dca6a3cdf465f6..40c33ed38c193b7108f1eac528101ce17f613b37 100644 (file)
@@ -782,6 +782,142 @@ namespace Functions
       const typename std::vector<Point<dim> > fourier_coefficients;
       const std::vector<double>      weights;
   };
+
+
+/**
+ * Cut-off function in L-infinity for an arbitrary ball.  This
+ * function is the characteristic function of a ball around @p{center}
+ * with a specified @p{radius}.
+ *
+ * @author Guido Kanschat, 2001
+ */
+  template<int dim>
+  class CutOffFunctionLinfty : public Function<dim>
+  {
+    public:
+                                    /**
+                                     * Constructor. Arguments are the
+                                     * center of the ball and its
+                                     * radius.
+                                     */
+    CutOffFunctionLinfty (const double radius = 1.,
+                     Point<dim> = Point<dim>());
+    
+                                      /**
+                                       * Function value at one point.
+                                       */
+      virtual double value (const Point<dim>   &p,
+                           const unsigned int  component = 0) const;
+      
+                                      /**
+                                       * Function values at multiple points.
+                                       */
+      virtual void value_list (const typename std::vector<Point<dim> > &points,
+                              std::vector<double>            &values,
+                              const unsigned int              component = 0) const;
+
+  private:
+                                    /**
+                                     * Center of the integration ball.
+                                     */
+    const Point<dim> center;
+
+                                    /**
+                                     * Radius of the ball.
+                                     */
+    const double radius;
+  };
+  
+  
+/**
+ * Cut-off function for an arbitrary ball. This function is a cone
+ * with support in a ball of certain @p{radius} around @p{center}. The
+ * maximum value is 1.
+ *
+ * @author Guido Kanschat, 2001
+ */
+  template<int dim>
+  class CutOffFunctionW1 : public Function<dim>
+  {
+    public:
+                                    /**
+                                     * Constructor. Arguments are the
+                                     * center of the ball and its
+                                     * radius.
+                                     */
+    CutOffFunctionW1 (const double radius = 1.,
+                     Point<dim> = Point<dim>());
+    
+                                      /**
+                                       * Function value at one point.
+                                       */
+      virtual double value (const Point<dim>   &p,
+                           const unsigned int  component = 0) const;
+      
+                                      /**
+                                       * Function values at multiple points.
+                                       */
+      virtual void value_list (const typename std::vector<Point<dim> > &points,
+                              std::vector<double>            &values,
+                              const unsigned int              component = 0) const;
+
+  private:
+                                    /**
+                                     * Center of the integration ball.
+                                     */
+    const Point<dim> center;
+
+                                    /**
+                                     * Radius of the ball.
+                                     */
+    const double radius;
+  };
+  
+  
+/**
+ * Cut-off function for an arbitrary ball. This is the traditional
+ * cut-off function in C-infinity for a ball of certain @p{radius}
+ * around @p{center}.
+ *
+ * @author Guido Kanschat, 2001
+ */
+  template<int dim>
+  class CutOffFunctionCinfty : public Function<dim>
+  {
+    public:
+                                    /**
+                                     * Constructor. Arguments are the
+                                     * center of the ball and its
+                                     * radius.
+                                     */
+    CutOffFunctionCinfty (const double radius = 1.,
+                         Point<dim> = Point<dim>());
+    
+                                      /**
+                                       * Function value at one point.
+                                       */
+      virtual double value (const Point<dim>   &p,
+                           const unsigned int  component = 0) const;
+      
+                                      /**
+                                       * Function values at multiple points.
+                                       */
+      virtual void value_list (const typename std::vector<Point<dim> > &points,
+                              std::vector<double>            &values,
+                              const unsigned int              component = 0) const;
+
+  private:
+                                    /**
+                                     * Center of the integration ball.
+                                     */
+    const Point<dim> center;
+
+                                    /**
+                                     * Radius of the ball.
+                                     */
+    const double radius;
+  };
+  
   
   
 };
diff --git a/deal.II/base/source/function_lib_cutoff.cc b/deal.II/base/source/function_lib_cutoff.cc
new file mode 100644 (file)
index 0000000..7712347
--- /dev/null
@@ -0,0 +1,172 @@
+//----------------------------  function_lib.cc  ---------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 1998, 1999, 2000, 2001 by the deal authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  function_lib.cc  ---------------------------
+
+
+#include <base/tensor.h>
+#include <base/point.h>
+#include <base/function_lib.h>
+
+#include <cmath>
+
+
+// in strict ANSI C mode, the following constants are not defined by
+// default, so we do it ourselves
+#ifndef M_PI
+#  define      M_PI            3.14159265358979323846
+#endif
+
+#ifndef M_PI_2
+#  define      M_PI_2          1.57079632679489661923
+#endif
+
+#ifndef M_E
+#  define       M_E             2.7182818284590452354
+#endif
+
+namespace Functions
+{
+
+  using namespace std;
+  
+  template<int dim>
+  CutOffFunctionLinfty<dim>::CutOffFunctionLinfty (const double r,
+                                                  const Point<dim> p)
+    : Function<dim> (1),
+    center(p),
+    radius(r)
+  {}
+
+
+  template<int dim>
+  double
+  CutOffFunctionLinfty<dim>::value (const Point<dim>   &p,
+                                   const unsigned int) const
+  {
+    return ((center.distance(p)<radius) ? 1. : 0.);
+  }
+
+
+  template<int dim>
+  void 
+  CutOffFunctionLinfty<dim>::value_list (const typename std::vector<Point<dim> > &points,
+                                        std::vector<double>                     &values,
+                                        const unsigned int) const
+  {
+    Assert (values.size() == points.size(),
+           ExcDimensionMismatch(values.size(), points.size()));
+
+    for (unsigned int i=0;i<values.size();++i)
+      values[i] = (center.distance(points[i])<radius) ? 1. : 0.;
+  }
+
+
+
+  template<int dim>
+  CutOffFunctionW1<dim>::CutOffFunctionW1 (const double     r,
+                                          const Point<dim> p)
+    : Function<dim> (1),
+    center(p),
+    radius(r)
+  {}
+
+
+  template<int dim>
+  double
+  CutOffFunctionW1<dim>::value (const Point<dim>   &p,
+                             const unsigned int) const
+  {
+    const double d = center.distance(p);
+    return ((d<radius) ? (radius-d) : 0.);
+  }
+
+
+  template<int dim>
+  void 
+  CutOffFunctionW1<dim>::value_list (const typename std::vector<Point<dim> > &points,
+                                    std::vector<double>                     &values,
+                                    const unsigned int) const
+  {
+    Assert (values.size() == points.size(),
+           ExcDimensionMismatch(values.size(), points.size()));
+    
+    for (unsigned int i=0;i<values.size();++i)
+      {
+       const double d = center.distance(points[i]);
+       values[i] = ((d<radius) ? (radius-d) : 0.);
+      }
+  }
+
+
+
+  template<int dim>
+  CutOffFunctionCinfty<dim>::CutOffFunctionCinfty (const double     r,
+                                                  const Point<dim> p)
+    : Function<dim> (1),
+    center(p),
+    radius(r)
+  {}
+
+
+  template<int dim>
+  double
+  CutOffFunctionCinfty<dim>::value (const Point<dim>   &p,
+                             const unsigned int) const
+  {
+    const double d = center.distance(p);
+    const double r = radius;
+    if (d>=r)
+      return 0.;
+    const double e = -r*r/(r*r-d*d);
+    return  ((e<-50) ? 0. : M_E * exp(e));
+  }
+
+
+  template<int dim>
+  void 
+  CutOffFunctionCinfty<dim>::value_list (const typename std::vector<Point<dim> > &points,
+                                        std::vector<double>                     &values,
+                                        const unsigned int) const
+  {
+    Assert (values.size() == points.size(),
+           ExcDimensionMismatch(values.size(), points.size()));
+    
+    const double r = radius;
+
+    for (unsigned int i=0;i<values.size();++i)
+      {
+       const double d = center.distance(points[i]);
+       if (d>=r)
+         {
+           values[i] = 0.;
+         } else {
+           const double e = -r*r/(r*r-d*d);
+           values[i] = (e<-50) ? 0. : M_E * exp(e);
+         }
+      }
+  }
+
+
+  template class CutOffFunctionLinfty <1>;
+  template class CutOffFunctionLinfty <2>;
+  template class CutOffFunctionLinfty <3>;
+  
+  template class CutOffFunctionW1 <1>;
+  template class CutOffFunctionW1 <2>;
+  template class CutOffFunctionW1 <3>;
+  
+  template class CutOffFunctionCinfty <1>;
+  template class CutOffFunctionCinfty <2>;
+  template class CutOffFunctionCinfty <3>;
+
+
+}

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.