File: Oracle_base.h

package info (click to toggle)
cgal 6.1.1-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 144,952 kB
  • sloc: cpp: 811,597; ansic: 208,576; sh: 493; python: 411; makefile: 286; javascript: 174
file content (368 lines) | stat: -rw-r--r-- 11,273 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
// Copyright (c) 2019-2022 Google LLC (USA).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
//
// $URL: https://github.com/CGAL/cgal/blob/v6.1.1/Alpha_wrap_3/include/CGAL/Alpha_wrap_3/internal/Oracle_base.h $
// $Id: include/CGAL/Alpha_wrap_3/internal/Oracle_base.h 08b27d3db14 $
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
//
// Author(s)     : Mael Rouxel-Labbé
//
#ifndef CGAL_ALPHA_WRAP_3_INTERNAL_ORACLE_BASE_H
#define CGAL_ALPHA_WRAP_3_INTERNAL_ORACLE_BASE_H

#include <CGAL/license/Alpha_wrap_3.h>

#include <CGAL/Alpha_wrap_3/internal/offset_intersection.h>

#include <CGAL/AABB_tree/internal/AABB_traversal_traits.h>
#include <CGAL/Default.h>

#include <algorithm>
#include <memory>

namespace CGAL {
namespace Alpha_wraps_3 {
namespace internal {

template <typename AABBTraits>
struct Default_traversal_traits
{
  using Projection_traits = CGAL::internal::AABB_tree::Projection_traits<AABBTraits>;

  template <typename Query>
  using Do_intersect_traits = CGAL::internal::AABB_tree::Do_intersect_traits<AABBTraits, Query>;

  template <typename Query>
  using First_intersection_traits = CGAL::internal::AABB_tree::First_intersection_traits<AABBTraits, Query>;
};

// Factorize the implementation of the functions calling the AABB tree
template <typename AABBTree,
          typename AABBTraversalTraits>
struct AABB_tree_oracle_helper
{
  using Self = AABB_tree_oracle_helper<AABBTree, AABBTraversalTraits>;

  using AABB_traits = typename AABBTree::AABB_traits;
  using GT = typename AABB_traits::Geom_traits;

  using FT = typename AABB_traits::FT;
  using Point_3 = typename AABB_traits::Point;

  template <typename Query>
  static bool do_intersect(const Query& query,
                           const AABBTree& tree)
  {
    CGAL_precondition(!tree.empty());

    using Do_intersect_traits = typename AABBTraversalTraits::template Do_intersect_traits<Query>;

    Do_intersect_traits traversal_traits(tree.traits());
    tree.traversal(query, traversal_traits);
    return traversal_traits.is_intersection_found();
  }

  static Point_3 closest_point(const Point_3& p,
                               const AABBTree& tree)
  {
    CGAL_precondition(!tree.empty());

    using Projection_traits = typename AABBTraversalTraits::Projection_traits;

    const auto& hint = tree.best_hint(p);

    Projection_traits projection_traits(hint.first, hint.second, tree.traits());
    tree.traversal(p, projection_traits);
    return projection_traits.closest_point();
  }

  static FT squared_distance(const Point_3& p,
                             const AABBTree& tree)
  {
    CGAL_precondition(!tree.empty());

    const Point_3 closest = Self::closest_point(p, tree);
    return tree.traits().squared_distance_object()(p, closest);
  }

  static bool first_intersection(const Point_3& p, const Point_3& q, Point_3& o,
                                 const FT offset_size,
                                 const FT intersection_precision,
                                 const AABBTree& tree)
  {
    CGAL_precondition(!tree.empty());

    using AABB_distance_oracle = internal::AABB_distance_oracle<AABBTree, AABBTraversalTraits>;
    using Offset_intersection = internal::Offset_intersection<GT, AABB_distance_oracle>;

    AABB_distance_oracle dist_oracle(tree);
    Offset_intersection offset_intersection(dist_oracle, offset_size, intersection_precision, 1 /*lip*/);
    return offset_intersection.first_intersection(p, q, o);
  }
};

template <typename GT,
          typename AABBTree,
          typename AABBTraversalTraits = CGAL::Default,
          typename BaseOracle = int> // base oracle
class AABB_tree_oracle
  : public BaseOracle
{
protected:
  using Geom_traits = GT;
  using FT = typename GT::FT;
  using Point_3 = typename GT::Point_3;

  // When building oracle stacks, there are copies of (empty) trees, which isn't possible, thus pointer
  using AABB_tree = AABBTree;
  using AABB_tree_ptr = std::shared_ptr<AABB_tree>;
  using AABB_traits = typename AABB_tree::AABB_traits;
  using AABB_traversal_traits = typename Default::Get<AABBTraversalTraits,
                                                      Default_traversal_traits<AABB_traits> >::type;

  using AABB_helper = AABB_tree_oracle_helper<AABB_tree, AABB_traversal_traits>;

public:
  AABB_tree_oracle(const BaseOracle& base,
                   const GT& gt)
    : BaseOracle(base),
      m_gt(gt),
      m_tree_ptr(std::make_shared<AABB_tree>())
  { }

  AABB_tree_oracle(const AABB_tree_oracle&) = default;

public:
  const Geom_traits& geom_traits() const { return m_gt; }

  AABB_tree& tree() { return *m_tree_ptr; }
  const AABB_tree& tree() const { return *m_tree_ptr; }
  BaseOracle& base() { return static_cast<BaseOracle&>(*this); }
  const BaseOracle& base() const { return static_cast<const BaseOracle&>(*this); }

  bool empty() const { return m_tree_ptr->empty(); }
  bool do_call() const { return (!empty() || base().do_call()); }

  void clear() { m_tree_ptr->clear() && base().clear(); }

public:
  typename AABB_tree::Bounding_box bbox() const
  {
    CGAL_precondition(do_call());

    typename AABB_tree::Bounding_box bb;

    if(base().do_call())
      bb += base().bbox();

    if(!empty())
      bb += tree().bbox();

    return bb;
  }

  template <typename T>
  bool do_intersect(const T& t) const
  {
    CGAL_precondition(do_call());

    if(base().do_call() && base().do_intersect(t))
      return true;

    if(!empty())
      return AABB_helper::do_intersect(t, tree());

    return false;
  }

  FT squared_distance(const Point_3& p) const
  {
    CGAL_precondition(do_call());

    if(base().do_call())
    {
      if(!empty()) // both non empty
      {
        const FT base_sqd = base().squared_distance(p);
        // @speed could do a smarter traversal, no need to search deeper than the current best
        const FT this_sqd = AABB_helper::squared_distance(p, tree());
        return (std::min)(base_sqd, this_sqd);
      }
      else // this level is empty
      {
        return base().squared_distance(p);
      }
    }
    else // empty base
    {
      return AABB_helper::squared_distance(p, tree());
    }
  }

  Point_3 closest_point(const Point_3& p) const
  {
    CGAL_precondition(do_call());

    if(base().do_call())
    {
      if(!empty()) // both non empty
      {
        const Point_3 base_c = base().closest_point(p);
        // @speed could do a smarter traversal, no need to search deeper than the current best
        const Point_3 this_c = AABB_helper::closest_point(p, tree());
        return (compare_distance_to_point(p, base_c, this_c) == CGAL::SMALLER) ? base_c : this_c;
      }
      else // this level is empty
      {
        return base().closest_point(p);
      }
    }
    else // empty base
    {
      return AABB_helper::closest_point(p, tree());
    }
  }

  bool first_intersection(const Point_3& p, const Point_3& q,
                          Point_3& o,
                          const FT offset_size,
                          const FT intersection_precision) const
  {
    CGAL_precondition(do_call());

    if(base().do_call())
    {
      if(!empty()) // both non empty
      {
        Point_3 base_o;
        bool base_b = base().first_intersection(p, q, base_o, offset_size, intersection_precision);

        if(base_b) // intersection found in base
        {
          // @speed could do a smarter traversal, no need to search deeper than the current best
          Point_3 this_o;
          bool this_b = AABB_helper::first_intersection(p, q, this_o, offset_size, intersection_precision, tree());
          if(this_b)
            o = (compare_distance_to_point(p, base_o, this_o) == SMALLER) ? base_o : this_o;
          else
            o = base_o;

          return true;
        }
        else // no intersection found in non-empty base
        {
          return AABB_helper::first_intersection(p, q, o, offset_size, intersection_precision, tree());
        }
      }
      else // this level is empty
      {
        return base().first_intersection(p, q, o, offset_size, intersection_precision);
      }
    }
    else // empty base
    {
      return AABB_helper::first_intersection(p, q, o, offset_size, intersection_precision, tree());
    }
  }

  bool first_intersection(const Point_3& p, const Point_3& q,
                          Point_3& o,
                          const FT offset_size) const
  {
    return first_intersection(p, q, o, offset_size, 1e-2 * offset_size);
  }

protected:
  Geom_traits m_gt;
  AABB_tree_ptr m_tree_ptr;
};

// partial specialization, when there is no further oracle beneath in the stack.
//
// `int` is used to denote the absence of base rather than `void`,
// as to use the same constructor for both versions (thus requires a default construction)
template <typename GT,
          typename AABBTree,
          typename AABBTraversalTraits>
class AABB_tree_oracle<GT, AABBTree, AABBTraversalTraits, int>
{
protected:
  using Geom_traits = GT;
  using FT = typename GT::FT;
  using Point_3 = typename GT::Point_3;

  using AABB_tree = AABBTree;
  using AABB_tree_ptr = std::shared_ptr<AABB_tree>;
  using AABB_traits = typename AABB_tree::AABB_traits;
  using AABB_traversal_traits = typename Default::Get<AABBTraversalTraits,
                                                      Default_traversal_traits<AABB_traits> >::type;

  using AABB_helper = AABB_tree_oracle_helper<AABB_tree, AABB_traversal_traits>;

public:
  AABB_tree_oracle(const int, // to have a common constructor API between both versions
                   const GT& gt)
    : m_gt(gt), m_tree_ptr(std::make_shared<AABB_tree>())
  { }

public:
  const Geom_traits& geom_traits() const { return m_gt; }
  AABB_tree& tree() { return *m_tree_ptr; }
  const AABB_tree& tree() const { return *m_tree_ptr; }

  bool empty() const { return m_tree_ptr->empty(); }
  bool do_call() const { return !empty(); }

  void clear() { m_tree_ptr->clear(); }

public:
  typename AABB_tree::Bounding_box bbox() const
  {
    CGAL_precondition(!empty());
    return tree().bbox();
  }

  template <typename T>
  bool do_intersect(const T& t) const
  {
    CGAL_precondition(!empty());
    return AABB_helper::do_intersect(t, tree());
  }

  FT squared_distance(const Point_3& p) const
  {
    CGAL_precondition(!empty());
    return AABB_helper::squared_distance(p, tree());
  }

  Point_3 closest_point(const Point_3& p) const
  {
    CGAL_precondition(!empty());
    return AABB_helper::closest_point(p, tree());
  }

  bool first_intersection(const Point_3& p, const Point_3& q, Point_3& o,
                          const FT offset_size, const FT intersection_precision) const
  {
    CGAL_precondition(!empty());
    return AABB_helper::first_intersection(p, q, o, offset_size, intersection_precision, tree());
  }

  bool first_intersection(const Point_3& p, const Point_3& q, Point_3& o, const FT offset_size) const
  {
    CGAL_precondition(!empty());
    return AABB_helper::first_intersection(p, q, o, offset_size, 1e-2 * offset_size, tree());
  }

private:
  Geom_traits m_gt;
  AABB_tree_ptr m_tree_ptr;
};

} // namespace internal
} // namespace Alpha_wraps_3
} // namespace CGAL

#endif // CGAL_ALPHA_WRAP_3_INTERNAL_ORACLE_BASE_H