32 const unsigned int n_components,
33 const unsigned int select,
34 const bool integrate_to_one,
35 const double unitary_integral_value)
40 , integrate_to_one(integrate_to_one)
41 , unitary_integral_value(unitary_integral_value)
42 , rescaling(integrate_to_one ? 1. / (unitary_integral_value *
77 1. / (unitary_integral_value * Utilities::fixed_power<dim>(radius));
97 return integrate_to_one;
106 const unsigned int n_components,
107 const unsigned int select,
108 const bool integrate_to_one)
124 for (
unsigned int i = 0; i < dim; ++i)
125 base[i]->set_center(
Point<1>(p[i]));
136 for (
unsigned int i = 0; i < dim; ++i)
137 base[i]->set_radius(r);
146 const unsigned int component)
const
150 for (
unsigned int i = 0; i < dim; ++i)
151 ret *= base[i]->value(
Point<1>(p[i]), component);
160 const unsigned int component)
const
164 for (
unsigned int d = 0; d < dim; ++d)
166 ret[d] = base[d]->gradient(
Point<1>(p[d]), component)[0];
167 for (
unsigned int i = 0; i < dim; ++i)
169 ret[d] *= base[i]->value(
Point<1>(p[i]), component);
181 const double integral_Linfty[] = {2.0,
182 3.14159265358979323846264338328,
183 4.18879020478639098461685784437};
187 const double integral_W1[] = {1.0,
188 1.04719755119659774615421446109,
189 1.04719755119659774615421446109};
193 const double integral_Cinfty[] = {1.20690032243787617533623799633,
194 1.26811216112759608094632335664,
195 1.1990039070192139033798473858};
199 const double integral_C1[] = {1.0,
200 0.93417655442731527615578663815,
201 0.821155557658032806157358815206};
209 const unsigned int n_components,
210 const unsigned int select,
211 const bool integrate_to_one)
217 integral_Linfty[dim - 1])
224 const unsigned int component)
const
227 component == this->selected)
228 return ((this->center.distance(p) < this->radius) ? this->rescaling : 0.);
236 std::vector<double> &values,
237 const unsigned int component)
const
239 Assert(values.size() == points.size(),
245 component == this->selected)
246 for (
unsigned int k = 0; k < values.size(); ++k)
247 values[k] = (this->center.distance(points[k]) < this->radius) ?
251 std::fill(values.begin(), values.end(), 0.);
261 Assert(values.size() == points.size(),
264 for (
unsigned int k = 0; k < values.size(); ++k)
266 const double val = (this->center.distance(points[k]) < this->radius) ?
274 values[k](this->selected) = val;
282 const unsigned int n_components,
283 const unsigned int select,
284 const bool integrate_to_one)
290 integral_W1[dim - 1])
297 const unsigned int component)
const
300 component == this->selected)
302 const double d = this->center.distance(p);
303 return ((d < this->radius) ?
304 (this->radius - d) / this->radius * this->rescaling :
314 std::vector<double> &values,
315 const unsigned int component)
const
317 Assert(values.size() == points.size(),
321 component == this->selected)
322 for (
unsigned int i = 0; i < values.size(); ++i)
324 const double d = this->center.distance(points[i]);
325 values[i] = ((d < this->radius) ?
326 (this->radius - d) / this->radius * this->rescaling :
330 std::fill(values.begin(), values.end(), 0.);
341 Assert(values.size() == points.size(),
344 for (
unsigned int k = 0; k < values.size(); ++k)
346 const double d = this->center.distance(points[k]);
349 (this->radius - d) / this->radius * this->rescaling :
356 values[k](this->selected) = val;
366 const unsigned int n_components,
367 const unsigned int select,
368 bool integrate_to_one)
374 integral_Cinfty[dim - 1])
381 const unsigned int component)
const
384 component == this->selected)
386 const double d = this->center.distance(p);
387 const double r = this->radius;
390 const double e = -r * r / (r * r - d * d);
400 std::vector<double> &values,
401 const unsigned int component)
const
403 Assert(values.size() == points.size(),
406 const double r = this->radius;
409 component == this->selected)
410 for (
unsigned int i = 0; i < values.size(); ++i)
412 const double d = this->center.distance(points[i]);
419 const double e = -r * r / (r * r - d * d);
425 std::fill(values.begin(), values.end(), 0.);
435 Assert(values.size() == points.size(),
438 for (
unsigned int k = 0; k < values.size(); ++k)
440 const double d = this->center.distance(points[k]);
441 const double r = this->radius;
443 if (d < this->radius)
445 const double e = -r * r / (r * r - d * d);
455 values[k](this->selected) = val;
465 const unsigned int)
const
467 const double d = this->center.distance(p);
468 const double r = this->radius;
471 const double e = -d * d / (r - d) / (r + d);
474 (p - this->center) / d *
475 (-2.0 * r * r / Utilities::fixed_power<2>(-r * r + d * d) * d *
485 const unsigned int n_components,
486 const unsigned int select,
487 bool integrate_to_one)
493 integral_C1[dim - 1])
500 const unsigned int component)
const
503 component == this->selected)
505 const double d = this->center.distance(p);
506 const double r = this->radius;
518 std::vector<double> &values,
519 const unsigned int component)
const
521 Assert(values.size() == points.size(),
524 const double r = this->radius;
527 component == this->selected)
528 for (
unsigned int i = 0; i < values.size(); ++i)
530 const double d = this->center.distance(points[i]);
542 std::fill(values.begin(), values.end(), 0.);
552 Assert(values.size() == points.size(),
555 for (
unsigned int k = 0; k < values.size(); ++k)
557 const double d = this->center.distance(points[k]);
558 const double r = this->radius;
560 if (d < this->radius)
570 values[k](this->selected) = val;
581 const double d = this->center.distance(p);
582 const double r = this->radius;
586 (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
#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 > &)