From c15bcf897feeaa200dd3221d1535c1e50609c874 Mon Sep 17 00:00:00 2001 From: Guido Kanschat Date: Thu, 11 Oct 2001 21:34:45 +0000 Subject: [PATCH] cutt-off functions git-svn-id: https://svn.dealii.org/trunk@5144 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/base/include/base/function_lib.h | 136 ++++++++++++++++ deal.II/base/source/function_lib_cutoff.cc | 172 +++++++++++++++++++++ 2 files changed, 308 insertions(+) create mode 100644 deal.II/base/source/function_lib_cutoff.cc diff --git a/deal.II/base/include/base/function_lib.h b/deal.II/base/include/base/function_lib.h index 3371177c1c..40c33ed38c 100644 --- a/deal.II/base/include/base/function_lib.h +++ b/deal.II/base/include/base/function_lib.h @@ -782,6 +782,142 @@ namespace Functions const typename std::vector > fourier_coefficients; const std::vector 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 + class CutOffFunctionLinfty : public Function + { + public: + /** + * Constructor. Arguments are the + * center of the ball and its + * radius. + */ + CutOffFunctionLinfty (const double radius = 1., + Point = Point()); + + /** + * Function value at one point. + */ + virtual double value (const Point &p, + const unsigned int component = 0) const; + + /** + * Function values at multiple points. + */ + virtual void value_list (const typename std::vector > &points, + std::vector &values, + const unsigned int component = 0) const; + + private: + /** + * Center of the integration ball. + */ + const Point 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 + class CutOffFunctionW1 : public Function + { + public: + /** + * Constructor. Arguments are the + * center of the ball and its + * radius. + */ + CutOffFunctionW1 (const double radius = 1., + Point = Point()); + + /** + * Function value at one point. + */ + virtual double value (const Point &p, + const unsigned int component = 0) const; + + /** + * Function values at multiple points. + */ + virtual void value_list (const typename std::vector > &points, + std::vector &values, + const unsigned int component = 0) const; + + private: + /** + * Center of the integration ball. + */ + const Point 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 + class CutOffFunctionCinfty : public Function + { + public: + /** + * Constructor. Arguments are the + * center of the ball and its + * radius. + */ + CutOffFunctionCinfty (const double radius = 1., + Point = Point()); + + /** + * Function value at one point. + */ + virtual double value (const Point &p, + const unsigned int component = 0) const; + + /** + * Function values at multiple points. + */ + virtual void value_list (const typename std::vector > &points, + std::vector &values, + const unsigned int component = 0) const; + + private: + /** + * Center of the integration ball. + */ + const Point 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 index 0000000000..7712347050 --- /dev/null +++ b/deal.II/base/source/function_lib_cutoff.cc @@ -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 +#include +#include + +#include + + +// 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 + CutOffFunctionLinfty::CutOffFunctionLinfty (const double r, + const Point p) + : Function (1), + center(p), + radius(r) + {} + + + template + double + CutOffFunctionLinfty::value (const Point &p, + const unsigned int) const + { + return ((center.distance(p) + void + CutOffFunctionLinfty::value_list (const typename std::vector > &points, + std::vector &values, + const unsigned int) const + { + Assert (values.size() == points.size(), + ExcDimensionMismatch(values.size(), points.size())); + + for (unsigned int i=0;i + CutOffFunctionW1::CutOffFunctionW1 (const double r, + const Point p) + : Function (1), + center(p), + radius(r) + {} + + + template + double + CutOffFunctionW1::value (const Point &p, + const unsigned int) const + { + const double d = center.distance(p); + return ((d + void + CutOffFunctionW1::value_list (const typename std::vector > &points, + std::vector &values, + const unsigned int) const + { + Assert (values.size() == points.size(), + ExcDimensionMismatch(values.size(), points.size())); + + for (unsigned int i=0;i + CutOffFunctionCinfty::CutOffFunctionCinfty (const double r, + const Point p) + : Function (1), + center(p), + radius(r) + {} + + + template + double + CutOffFunctionCinfty::value (const Point &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 + void + CutOffFunctionCinfty::value_list (const typename std::vector > &points, + std::vector &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=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>; + + +} -- 2.39.5