]> https://gitweb.dealii.org/ - dealii.git/commitdiff
vector valued cut-off functions
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Thu, 10 Jan 2002 12:37:53 +0000 (12:37 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Thu, 10 Jan 2002 12:37:53 +0000 (12:37 +0000)
git-svn-id: https://svn.dealii.org/trunk@5376 0785d39b-7218-0410-832d-ea1e28bc413d

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

index def9951f8469a3ad35307355a43049152c70cb2b..30374dbe67a2315c22f6d4eb9c8efc942cbe589a 100644 (file)
@@ -788,9 +788,10 @@ namespace Functions
 /**
  * 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}.
+ * with a specified @p{radius}. If vector valued, it can be restricted
+ * to a single component.
  *
- * @author Guido Kanschat, 2001
+ * @author Guido Kanschat, 2001, 2002
  */
   template<int dim>
   class CutOffFunctionLinfty : public Function<dim>
@@ -800,9 +801,17 @@ namespace Functions
                                        * Constructor. Arguments are the
                                        * center of the ball and its
                                        * radius.
+                                       *
+                                       * If an argument @p{select} is
+                                       * given and not -1, the
+                                       * cut-off function will be
+                                       * non-zero for this component
+                                       * only.
                                        */
       CutOffFunctionLinfty (const double radius = 1.,
-                           Point<dim> = Point<dim>());
+                           Point<dim> = Point<dim>(),
+                           const unsigned int n_components = 1,
+                           const unsigned int select = no_component);
     
                                       /**
                                        * Function value at one point.
@@ -817,6 +826,12 @@ namespace Functions
                               std::vector<double>            &values,
                               const unsigned int              component = 0) const;
 
+                                      /**
+                                       * Function values at multiple points.
+                                       */
+      virtual void vector_value_list (const typename std::vector<Point<dim> > &points,
+                                     std::vector<Vector<double> >           &values) const;
+
     private:
                                       /**
                                        * Center of the integration ball.
@@ -827,15 +842,27 @@ namespace Functions
                                        * Radius of the ball.
                                        */
       const double radius;
+
+                                      /**
+                                       * Component selected. If
+                                       * @p{no_component}, the function is
+                                       * the same in all components.
+                                       */
+      const unsigned int selected;
+                                      /**
+                                       * Value for no selected component.
+                                       */
+      static const unsigned int no_component = static_cast<unsigned int>(-1);
   };
   
   
 /**
  * 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.
+ * maximum value is 1. If vector valued, it can be restricted
+ * to a single component.
  *
- * @author Guido Kanschat, 2001
+ * @author Guido Kanschat, 2001, 2002
  */
   template<int dim>
   class CutOffFunctionW1 : public Function<dim>
@@ -845,9 +872,17 @@ namespace Functions
                                        * Constructor. Arguments are the
                                        * center of the ball and its
                                        * radius.
+                                       * radius.
+                                       *
+                                       * If an argument @p{select} is
+                                       * given, the cut-off function
+                                       * will be non-zero for this
+                                       * component only.
                                        */
       CutOffFunctionW1 (const double radius = 1.,
-                       Point<dim> = Point<dim>());
+                       Point<dim> = Point<dim>(),
+                       const unsigned int n_components = 1,
+                       const unsigned int select = no_component);
     
                                       /**
                                        * Function value at one point.
@@ -862,6 +897,12 @@ namespace Functions
                               std::vector<double>            &values,
                               const unsigned int              component = 0) const;
 
+                                      /**
+                                       * Function values at multiple points.
+                                       */
+      virtual void vector_value_list (const typename std::vector<Point<dim> > &points,
+                                     std::vector<Vector<double> >           &values) const;
+
     private:
                                       /**
                                        * Center of the integration ball.
@@ -872,6 +913,18 @@ namespace Functions
                                        * Radius of the ball.
                                        */
       const double radius;
+
+                                      /**
+                                       * Component selected. If
+                                       * @p{no_component}, the
+                                       * function is the same in all
+                                       * components.
+                                       */
+      const unsigned int selected;
+                                      /**
+                                       * Value for no selected component.
+                                       */
+      static const unsigned int no_component = static_cast<unsigned int>(-1);
   };
   
   
@@ -879,9 +932,10 @@ namespace Functions
  * 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}, $f(r)=exp(1-1/(1-r**2/s**2))$, where $r$ is the
- * distance to the center, and $s$ is the radius of the sphere.
+ * distance to the center, and $s$ is the radius of the sphere. If
+ * vector valued, it can be restricted to a single component.
  *
- * @author Guido Kanschat, 2001
+ * @author Guido Kanschat, 2001, 2002
  */
   template<int dim>
   class CutOffFunctionCinfty : public Function<dim>
@@ -891,9 +945,17 @@ namespace Functions
                                        * Constructor. Arguments are the
                                        * center of the ball and its
                                        * radius.
+                                       * radius.
+                                       *
+                                       * If an argument @p{select} is
+                                       * given, the cut-off function
+                                       * will be non-zero for this
+                                       * component only.
                                        */
       CutOffFunctionCinfty (const double radius = 1.,
-                           Point<dim> = Point<dim>());
+                           Point<dim> = Point<dim>(),
+                           const unsigned int n_components = 1,
+                           const unsigned int select = no_component);
     
                                       /**
                                        * Function value at one point.
@@ -908,6 +970,12 @@ namespace Functions
                               std::vector<double>            &values,
                               const unsigned int              component = 0) const;
 
+                                      /**
+                                       * Function values at multiple points.
+                                       */
+      virtual void vector_value_list (const typename std::vector<Point<dim> > &points,
+                                     std::vector<Vector<double> >           &values) const;
+
                                       /**
                                        * Function gradient at one point.
                                        */
@@ -924,6 +992,18 @@ namespace Functions
                                        * Radius of the ball.
                                        */
       const double radius;
+
+                                      /**
+                                       * Component selected. If
+                                       * @p{no_component}, the function is
+                                       * the same in all components.
+                                       */
+      const unsigned int selected;
+
+                                      /**
+                                       * Value for no selected component.
+                                       */
+      static const unsigned int no_component = static_cast<unsigned int>(-1);      
   };
   
 };
index 1aaaad1712f4de1d974a9151a6258ad5ba80d4e5..776c566318e7752edb55b5ec205496e0b54aea85 100644 (file)
@@ -15,6 +15,7 @@
 #include <base/tensor.h>
 #include <base/point.h>
 #include <base/function_lib.h>
+#include <lac/vector.h>
 
 #include <cmath>
 
@@ -41,20 +42,25 @@ namespace Functions
   
   template<int dim>
   CutOffFunctionLinfty<dim>::CutOffFunctionLinfty (const double r,
-                                                  const Point<dim> p)
+                                                  const Point<dim> p,
+                                                  const unsigned int n_components,
+                                                  unsigned int select)
                  :
-                 Function<dim> (1),
-                  center(p),
-                  radius(r)
+                 Function<dim> (n_components),
+    center(p),
+    radius(r),
+    selected(select)
   {}
 
 
   template<int dim>
   double
   CutOffFunctionLinfty<dim>::value (const Point<dim>   &p,
-                                   const unsigned int) const
+                                   const unsigned int component) const
   {
-    return ((center.distance(p)<radius) ? 1. : 0.);
+    if (selected==no_component || component==selected)
+      return ((center.distance(p)<radius) ? 1. : 0.);
+    return 0.;
   }
 
 
@@ -62,34 +68,69 @@ namespace Functions
   void 
   CutOffFunctionLinfty<dim>::value_list (const typename std::vector<Point<dim> > &points,
                                         std::vector<double>                     &values,
-                                        const unsigned int) const
+                                        const unsigned int component) const
   {
     Assert (values.size() == points.size(),
            ExcDimensionMismatch(values.size(), points.size()));
+    Assert (component < n_components,
+           ExcIndexRange(component,0,n_components));
+    
 
-    for (unsigned int i=0;i<values.size();++i)
-      values[i] = (center.distance(points[i])<radius) ? 1. : 0.;
+    if (selected==no_component || component==selected)
+      for (unsigned int k=0;k<values.size();++k)
+       values[k] = (center.distance(points[k])<radius) ? 1. : 0.;
+    else
+      fill (values.begin(), values.end(), 0.);
   }
 
+  
+  template<int dim>
+  void
+  CutOffFunctionLinfty<dim>::vector_value_list (
+    const typename std::vector<Point<dim> > &points,
+    std::vector<Vector<double> >           &values) const
+  {
+    Assert (values.size() == points.size(),
+           ExcDimensionMismatch(values.size(), points.size()));
 
+    for (unsigned int k=0;k<values.size();++k)
+      {
+       double val = (center.distance(points[k])<radius) ? 1. : 0.;
+       if (selected==no_component)
+         values[k] = val;
+       else
+         {
+           values[k] = 0.;
+           values[k](selected) = val;
+         }
+      }
+  }
+  
 
   template<int dim>
   CutOffFunctionW1<dim>::CutOffFunctionW1 (const double     r,
-                                          const Point<dim> p)
+                                          const Point<dim> p,
+                                          const unsigned int n_components,
+                                          unsigned int select)
                  :
-                 Function<dim> (1),
+                 Function<dim> (n_components),
                   center(p),
-                  radius(r)
+                  radius(r),
+    selected(select)
   {}
 
 
   template<int dim>
   double
   CutOffFunctionW1<dim>::value (const Point<dim>   &p,
-                               const unsigned int) const
+                               const unsigned int component) const
   {
-    const double d = center.distance(p);
-    return ((d<radius) ? (radius-d) : 0.);
+    if (selected==no_component || component==selected)
+      {
+       const double d = center.distance(p);
+       return ((d<radius) ? (radius-d) : 0.);
+      }
+    return 0.;  
   }
 
 
@@ -97,68 +138,138 @@ namespace Functions
   void 
   CutOffFunctionW1<dim>::value_list (const typename std::vector<Point<dim> > &points,
                                     std::vector<double>                     &values,
-                                    const unsigned int) const
+                                    const unsigned int component) 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.);
-      }
+    if (selected==no_component || component==selected)
+      for (unsigned int i=0;i<values.size();++i)
+       {
+         const double d = center.distance(points[i]);
+         values[i] = ((d<radius) ? (radius-d) : 0.);
+       }
+    else
+      fill (values.begin(), values.end(), 0.);
   }
 
 
 
+  template<int dim>
+  void
+  CutOffFunctionW1<dim>::vector_value_list (
+    const typename std::vector<Point<dim> > &points,
+    std::vector<Vector<double> >           &values) const
+  {
+    Assert (values.size() == points.size(),
+           ExcDimensionMismatch(values.size(), points.size()));
+
+    for (unsigned int k=0;k<values.size();++k)
+      {
+       const double d = center.distance(points[k]);
+       double val = (d<radius) ? (radius-d) : 0.;
+       if (selected==no_component)
+         values[k] = val;
+       else
+         {
+           values[k] = 0.;
+           values[k](selected) = val;
+         }
+      }
+  }
+  
+
   template<int dim>
   CutOffFunctionCinfty<dim>::CutOffFunctionCinfty (const double     r,
-                                                  const Point<dim> p)
+                                                  const Point<dim> p,
+                                                  const unsigned int n_components,
+                                                  unsigned int select)
                  :
-                 Function<dim> (1),
+                 Function<dim> (n_components),
                   center(p),
-                  radius(r)
+                  radius(r),
+    selected(select)
   {}
 
 
   template<int dim>
   double
   CutOffFunctionCinfty<dim>::value (const Point<dim>   &p,
-                                   const unsigned int) const
+                                   const unsigned int component) 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));
+    if (selected==no_component || component==selected)
+      {
+       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));
+      }
+    return 0.;
   }
-
+  
 
   template<int dim>
   void 
   CutOffFunctionCinfty<dim>::value_list (const typename std::vector<Point<dim> > &points,
                                         std::vector<double>                     &values,
-                                        const unsigned int) const
+                                        const unsigned int component) const
   {
     Assert (values.size() == points.size(),
            ExcDimensionMismatch(values.size(), points.size()));
     
     const double r = radius;
 
-    for (unsigned int i=0;i<values.size();++i)
+    if (selected==no_component || component==selected)
+      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);
+           }
+       }
+    else
+      fill (values.begin(), values.end(), 0.);
+  }
+
+
+  template<int dim>
+  void
+  CutOffFunctionCinfty<dim>::vector_value_list (
+    const typename std::vector<Point<dim> > &points,
+    std::vector<Vector<double> >           &values) const
+  {
+    Assert (values.size() == points.size(),
+           ExcDimensionMismatch(values.size(), points.size()));
+
+    for (unsigned int k=0;k<values.size();++k)
       {
-       const double d = center.distance(points[i]);
-       if (d>=r)
+       const double d = center.distance(points[k]);
+       const double r = radius;
+       double e = 0.;
+       double val = 0.;
+       if (d<radius)
+         {
+           e = -r*r/(r*r-d*d);
+           if (e>-50)
+             val = M_E * exp(e);
+         }
+
+       if (selected==no_component)
+         values[k] = val;
+       else
          {
-           values[i] = 0.;
-         } else {
-           const double e = -r*r/(r*r-d*d);
-           values[i] = (e<-50) ? 0. : M_E * exp(e);
+           values[k] = 0.;
+           values[k](selected) = val;
          }
       }
   }
-
+  
 
 
   template<int dim>
index b585ed0fb4c02013c66f3c445a7fd9d7077405db..11fd597f4a893ee9eff196090b6802c9adeda49b 100644 (file)
@@ -107,11 +107,22 @@ documentation, etc</a>.
 <h3>base</h3>
 
 <ol>
+  <li> <p>
+       Improved: The cut-off functions <code
+       class="class">Functios::CutOffFunctionLinfty</code>, <code
+       class="class">Functios::CutOffFunctionW1</code>, and <code
+       class="class">Functios::CutOffFunctionCinfty</code> can be
+       vector-valued now and optionally only a single componente can
+       be selected.
+       <br>
+       (GK 2002/01/10)
+       </p>
+
   <li> <p>
        New: the <code
        class="member">deal_II_exceptions::set_additional_assert_output</code>
        function allows to set additional output to be printed upon
-       triggering an <code class="member">Assert()</code> call. This
+       triggering an <code>Assert()</code> call. This
        is helpful for parallel applications where you only see the
        text of the message but do not know from which cluster node it
        stems.

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.