26 namespace SignedDistance
41 const unsigned int component)
const
54 const unsigned int component)
const
69 const unsigned int component)
const
75 const double distance = center_to_point.
norm();
78 unit_symmetric_tensor<dim>() / distance -
80 Utilities::fixed_power<3>(distance);
89 : point_in_plane(point)
100 const unsigned int component)
const
105 return normal * (point - point_in_plane);
136 const std::array<double, dim> &radii,
137 const double tolerance,
138 const unsigned int max_iter)
141 , tolerance(tolerance)
144 for (
unsigned int d = 0; d < dim; ++d)
153 const unsigned int component)
const
161 return compute_signed_distance_ellipse(point);
173 const unsigned int component)
const
183 const Point<dim> point_in_centered_coordinate_system =
185 grad = compute_analyical_normal_vector_on_ellipse(
186 point_in_centered_coordinate_system);
191 if (grad.
norm() > 1e-12)
192 return grad / grad.
norm();
204 for (
unsigned int d = 0; d < dim; ++d)
205 val += Utilities::fixed_power<2>((point[d] -
center[d]) / radii[d]);
234 const double sign_px = std::copysign(1.0, point[0] -
center[0]);
235 const double sign_py = std::copysign(1.0, point[1] -
center[1]);
237 const double &a = radii[0];
238 const double &b = radii[1];
244 unsigned int iter = 0;
250 (a * a - b * b) * Utilities::fixed_power<3>(
std::cos(t)) / a;
252 (b * b - a * a) * Utilities::fixed_power<3>(
std::sin(t)) / b;
254 const double rx = x - ex;
255 const double ry = y - ey;
257 const double qx = px - ex;
258 const double qy = py - ey;
260 const double r = std::hypot(rx, ry);
262 const double q = std::hypot(qx, qy);
264 const double delta_c = r * std::asin((rx * qy - ry * qx) / (r * q));
266 delta_t = delta_c /
std::sqrt(a * a + b * b - x * x - y * y);
274 while (
std::abs(delta_t) > tolerance && iter < max_iter);
301 const auto &a = radii[0];
302 const auto &b = radii[1];
303 const auto &x = point[0];
304 const auto &y = point[1];
326 return *std::min_element(radii.begin(), radii.end()) * -1.;
328 const Point<2> &closest_point = compute_closest_point_ellipse(point);
330 const double distance =
331 std::hypot(closest_point[0] - point[0], closest_point[1] - point[1]);
333 return evaluate_ellipsoid(point) < 0.0 ? -distance : distance;
341 : bounding_box({bottom_left, top_right})
348 : bounding_box(bounding_box)
356 const unsigned int component)
const
361 return bounding_box.signed_distance(p);
369 const double notch_width,
370 const double notch_height)
377 std::numeric_limits<double>::lowest()) :
381 std::numeric_limits<double>::lowest()),
387 center[1] + notch_height - radius) :
391 center[2] + notch_height - radius))
394 notch_width <= 2 * radius,
396 "The width of the notch must be less than the circle diameter."));
398 notch_height <= 2 * radius,
400 "The height of the notch must be less than the circle diameter."));
408 const unsigned int component)
const
415 return std::max(sphere.value(p), -notch.value(p));
420#include "function_signed_distance.inst"
double value(const Point< dim > &point, const unsigned int component=0) const override
Tensor< 1, dim > gradient(const Point< dim > &, const unsigned int component=0) const override
double compute_signed_distance_ellipse(const Point< dim > &point) const
Point< dim > compute_closest_point_ellipse(const Point< dim > &point) const
Ellipsoid(const Point< dim > ¢er, const std::array< double, dim > &radii, const double tolerance=1e-14, const unsigned int max_iter=10)
double evaluate_ellipsoid(const Point< dim > &point) const
const std::array< double, dim > radii
Tensor< 1, dim, double > compute_analyical_normal_vector_on_ellipse(const Point< dim > &point) const
Plane(const Point< dim > &point, const Tensor< 1, dim > &normal)
double value(const Point< dim > &point, const unsigned int component=0) const override
SymmetricTensor< 2, dim > hessian(const Point< dim > &, const unsigned int component=0) const override
Tensor< 1, dim > gradient(const Point< dim > &, const unsigned int component=0) const override
const Tensor< 1, dim > normal
double value(const Point< dim > &p, const unsigned int component=0) const override
Rectangle(const Point< dim > &bottom_left, const Point< dim > &top_right)
Sphere(const Point< dim > ¢er=Point< dim >(), const double radius=1)
SymmetricTensor< 2, dim > hessian(const Point< dim > &point, const unsigned int component=0) const override
double value(const Point< dim > &point, const unsigned int component=0) const override
Tensor< 1, dim > gradient(const Point< dim > &point, const unsigned int component=0) const override
ZalesakDisk(const Point< dim > ¢er, const double radius, const double notch_width, const double notch_height)
double value(const Point< dim > &p, const unsigned int component=0) const override
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
numbers::NumberTraits< Number >::real_type norm() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIsFinite(number)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
static constexpr double PI_2
static constexpr double PI_4
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
DEAL_II_HOST constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)