33 const unsigned int n_components,
34 const unsigned int select,
35 const bool integrate_to_one,
36 const double unitary_integral_value)
41 , integrate_to_one(integrate_to_one)
42 , unitary_integral_value(unitary_integral_value)
43 , rescaling(integrate_to_one ? 1. / (unitary_integral_value *
78 1. / (unitary_integral_value * Utilities::fixed_power<dim>(radius));
98 return integrate_to_one;
107 const unsigned int n_components,
108 const unsigned int select,
109 const bool integrate_to_one)
125 for (
unsigned int i = 0; i < dim; ++i)
126 base[i]->set_center(
Point<1>(p[i]));
137 for (
unsigned int i = 0; i < dim; ++i)
138 base[i]->set_radius(r);
147 const unsigned int component)
const
151 for (
unsigned int i = 0; i < dim; ++i)
152 ret *= base[i]->value(
Point<1>(p[i]), component);
161 const unsigned int component)
const
165 for (
unsigned int d = 0; d < dim; ++d)
167 ret[d] = base[d]->gradient(
Point<1>(p[d]), component)[0];
168 for (
unsigned int i = 0; i < dim; ++i)
170 ret[d] *= base[i]->value(
Point<1>(p[i]), component);
182 const double integral_Linfty[] = {2.0,
183 3.14159265358979323846264338328,
184 4.18879020478639098461685784437};
188 const double integral_W1[] = {1.0,
189 1.04719755119659774615421446109,
190 1.04719755119659774615421446109};
194 const double integral_Cinfty[] = {1.20690032243787617533623799633,
195 1.26811216112759608094632335664,
196 1.1990039070192139033798473858};
200 const double integral_C1[] = {1.0,
201 0.93417655442731527615578663815,
202 0.821155557658032806157358815206};
210 const unsigned int n_components,
211 const unsigned int select,
212 const bool integrate_to_one)
218 integral_Linfty[dim - 1])
225 const unsigned int component)
const
228 component == this->selected)
229 return ((this->
center.
distance(p) < this->radius) ? this->rescaling : 0.);
237 std::vector<double> & values,
238 const unsigned int component)
const
240 Assert(values.size() == points.size(),
246 component == this->selected)
247 for (
unsigned int k = 0; k < values.size(); ++k)
252 std::fill(values.begin(), values.end(), 0.);
262 Assert(values.size() == points.size(),
265 for (
unsigned int k = 0; k < values.size(); ++k)
267 const double val = (this->
center.
distance(points[k]) < this->radius) ?
275 values[k](this->selected) = val;
283 const unsigned int n_components,
284 const unsigned int select,
285 const bool integrate_to_one)
291 integral_W1[dim - 1])
298 const unsigned int component)
const
301 component == this->selected)
304 return ((d < this->radius) ?
305 (this->radius - d) / this->radius * this->rescaling :
315 std::vector<double> & values,
316 const unsigned int component)
const
318 Assert(values.size() == points.size(),
322 component == this->selected)
323 for (
unsigned int i = 0; i < values.size(); ++i)
326 values[i] = ((d < this->radius) ?
327 (this->radius - d) / this->radius * this->rescaling :
331 std::fill(values.begin(), values.end(), 0.);
342 Assert(values.size() == points.size(),
345 for (
unsigned int k = 0; k < values.size(); ++k)
350 (this->radius - d) / this->radius * this->rescaling :
357 values[k](this->selected) = val;
367 const unsigned int n_components,
368 const unsigned int select,
369 bool integrate_to_one)
375 integral_Cinfty[dim - 1])
382 const unsigned int component)
const
385 component == this->selected)
388 const double r = this->radius;
391 const double e = -r * r / (r * r - d * d);
401 std::vector<double> & values,
402 const unsigned int component)
const
404 Assert(values.size() == points.size(),
407 const double r = this->radius;
410 component == this->selected)
411 for (
unsigned int i = 0; i < values.size(); ++i)
420 const double e = -r * r / (r * r - d * d);
426 std::fill(values.begin(), values.end(), 0.);
436 Assert(values.size() == points.size(),
439 for (
unsigned int k = 0; k < values.size(); ++k)
442 const double r = this->radius;
444 if (d < this->radius)
446 const double e = -r * r / (r * r - d * d);
456 values[k](this->selected) = val;
466 const unsigned int)
const
469 const double r = this->radius;
472 const double e = -d * d / (r - d) / (r + d);
476 (-2.0 * r * r / Utilities::fixed_power<2>(-r * r + d * d) * d *
486 const unsigned int n_components,
487 const unsigned int select,
488 bool integrate_to_one)
494 integral_C1[dim - 1])
501 const unsigned int component)
const
504 component == this->selected)
507 const double r = this->radius;
519 std::vector<double> & values,
520 const unsigned int component)
const
522 Assert(values.size() == points.size(),
525 const double r = this->radius;
528 component == this->selected)
529 for (
unsigned int i = 0; i < values.size(); ++i)
543 std::fill(values.begin(), values.end(), 0.);
553 Assert(values.size() == points.size(),
556 for (
unsigned int k = 0; k < values.size(); ++k)
559 const double r = this->radius;
561 if (d < this->radius)
571 values[k](this->selected) = val;
583 const double r = this->radius;
587 (p - this->
center) / d * this->rescaling;
virtual void set_radius(const double r)
bool integrates_to_one() const
double get_radius() const
virtual void set_center(const Point< dim > &p)
const Point< dim > & get_center() const
CutOffFunctionBase(const double radius=1., const Point< dim > center=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false, const double unitary_integral_value=1.0)
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
CutOffFunctionC1(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, bool integrate_to_one=false)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
CutOffFunctionCinfty(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, bool integrate_to_one=false)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
CutOffFunctionLinfty(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false)
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
CutOffFunctionTensorProduct(double radius=1.0, const Point< dim > ¢er=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void set_radius(const double radius) override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
virtual void set_center(const Point< dim > ¢er) override
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
CutOffFunctionW1(const double radius=1., const Point< dim >=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false)
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
static constexpr double E
static constexpr double PI
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)