TTK
Loading...
Searching...
No Matches
PeriodicImplicitTriangulation.h
Go to the documentation of this file.
1
12
13#pragma once
14
15// base code includes
17
18#include <array>
19#include <numeric>
20
21namespace ttk {
22
24
25 public:
28
30 = default;
34 = default;
36 = default;
37
38 int getCellEdgeInternal(const SimplexId &cellId,
39 const int &id,
40 SimplexId &edgeId) const override;
41
42 SimplexId getCellEdgeNumberInternal(const SimplexId &cellId) const override;
43
44 const std::vector<std::vector<SimplexId>> *getCellEdgesInternal() override;
45
47 const SimplexId &cellId,
48 const int &localNeighborId,
49 SimplexId &neighborId) const override;
50
52 const SimplexId &cellId) const override;
53
54 const std::vector<std::vector<SimplexId>> *
56
57 int getCellTriangleInternal(const SimplexId &cellId,
58 const int &id,
59 SimplexId &triangleId) const override;
60
62 const SimplexId &ttkNotUsed(cellId)) const override {
63 // NOTE: the output is always 4 here. let's keep the function in there
64 // in case of further generalization to CW-complexes
65 return 4;
66 }
67
68 const std::vector<std::vector<SimplexId>> *
69 getCellTrianglesInternal() override;
70
72 const SimplexId &cellId,
73 const int &localVertexId,
74 SimplexId &vertexId) const override;
75
77 const SimplexId &cellId) const override;
78
80 return dimensionality_;
81 }
82
84 const SimplexId &edgeId) const override;
85
86 const std::vector<std::vector<SimplexId>> *
88
89 const std::vector<std::vector<SimplexId>> *
91
92 const std::vector<std::vector<SimplexId>> *
93 getEdgeTrianglesInternal() override;
94
95 const std::vector<std::array<SimplexId, 2>> *
97
99 return cellNumber_;
100 }
101
103 return edgeNumber_;
104 }
105
107 return triangleNumber_;
108 }
109
111 return vertexNumber_;
112 }
113
114 virtual int getTetrahedronEdge(const SimplexId &tetId,
115 const int &id,
116 SimplexId &edgeId) const = 0;
117
118 int getTetrahedronEdges(std::vector<std::vector<SimplexId>> &edges) const;
119
120 virtual int getTetrahedronTriangle(const SimplexId &tetId,
121 const int &id,
122 SimplexId &triangleId) const = 0;
123
125 std::vector<std::vector<SimplexId>> &triangles) const;
126
127 virtual int getTetrahedronNeighbor(const SimplexId &tetId,
128 const int &localNeighborId,
129 SimplexId &neighborId) const = 0;
130
132
133 int getTetrahedronNeighbors(std::vector<std::vector<SimplexId>> &neighbors);
134
135 virtual int getTetrahedronVertex(const SimplexId &tetId,
136 const int &localVertexId,
137 SimplexId &vertexId) const = 0;
138
140 const SimplexId &ttkNotUsed(triangleId)) const override {
141 // NOTE: the output is always 3 here. let's keep the function in there
142 // in case of further generalization to CW-complexes
143 return 3;
144 }
145
146 const std::vector<std::vector<SimplexId>> *
147 getTriangleEdgesInternal() override;
148
150 std::vector<std::vector<SimplexId>> &edges) const;
151
153 const SimplexId &triangleId) const override;
154
155 const std::vector<std::vector<SimplexId>> *
157
158 virtual int getTriangleNeighbor(const SimplexId &triangleId,
159 const int &localNeighborId,
160 SimplexId &neighborId) const = 0;
161
162 SimplexId getTriangleNeighborNumber(const SimplexId &triangleId) const;
163
164 int getTriangleNeighbors(std::vector<std::vector<SimplexId>> &neighbors);
165
166 const std::vector<std::vector<SimplexId>> *
168
169 const std::vector<std::array<SimplexId, 3>> *
171
173 getVertexEdgeNumberInternal(const SimplexId &vertexId) const override;
174
175 const std::vector<std::vector<SimplexId>> *
176 getVertexEdgesInternal() override;
177
179 const SimplexId &vertexId) const override;
180
181 const std::vector<std::vector<SimplexId>> *
183
185 const SimplexId &vertexId) const override;
186
187 const std::vector<std::vector<SimplexId>> *
189
191 const SimplexId &vertexId) const override;
192
193 const std::vector<std::vector<SimplexId>> *
195
197 getVertexTriangleNumberInternal(const SimplexId &vertexId) const override;
198
199 const std::vector<std::vector<SimplexId>> *
201
203 const SimplexId &edgeId) const override;
204
205 inline bool isEmpty() const override {
206 return !vertexNumber_;
207 }
208
210 const SimplexId &triangleId) const override;
211
213 const SimplexId &vertexId) const override;
214
215 int setInputGrid(const float &xOrigin,
216 const float &yOrigin,
217 const float &zOrigin,
218 const float &xSpacing,
219 const float &ySpacing,
220 const float &zSpacing,
221 const SimplexId &xDim,
222 const SimplexId &yDim,
223 const SimplexId &zDim) override;
224
225 inline const std::array<ttk::SimplexId, 3> &
226 getGridDimensions() const override {
227 return this->dimensions_;
228 }
229
231 int preconditionEdgesInternal() override = 0;
234
236 if(dimensionality_ == 3) {
238 } else if(dimensionality_ == 2 && !hasPreconditionedTriangles_) {
240 return this->preconditionTrianglesInternal();
241 }
242 return 0;
243 }
244
250 }
251 return 0;
252 }
253
254 inline int getCellVTKIDInternal(const int &ttkId,
255 int &vtkId) const override {
256#ifndef TTK_ENABLE_KAMIKAZE
257 if(ttkId < 0) {
258 return -1;
259 }
260#endif // TTK_ENABLE_KAMIKAZE
261 const SimplexId nSimplexPerCell{this->getDimensionality() == 3 ? 6 : 2};
262 vtkId = ttkId / nSimplexPerCell;
263 return 0;
264 }
265
266 protected:
268 float origin_[3]; //
269 float spacing_[3]; //
270 SimplexId nbvoxels_[3]; // nombre de voxels par axe
272
273 // Vertex helper //
274 SimplexId vshift_[2]; // VertexShift
275
276 // Edge helper //
277 SimplexId esetdims_[7]; // EdgeSetDimensions
278 SimplexId esetshift_[7]; // EdgeSetShift
279 SimplexId eshift_[14]; // EdgeShift
280
281 // Triangle helper //
282 SimplexId tsetdims_[6]; // TriangleSetDimensions
283 SimplexId tsetshift_[6]; // TriangleSetShift
284 SimplexId tshift_[12]; // TriangleShift
285
286 // Tetrahedron helper //
287 SimplexId tetshift_[2]; // TetrahedronShift
288
289 SimplexId cellNumber_; // number of cells
290 SimplexId vertexNumber_; // number of vertices
291 SimplexId edgeNumber_; // number of edges
292 SimplexId triangleNumber_; // number of triangles
293 SimplexId tetrahedronNumber_; // number of tetrahedra
294
295 // 2d helpers
298
299 // acceleration variables
303
304 enum class EdgePosition : char {
305 // e--------f
306 // /| /|
307 // / | / |
308 // a--------b |
309 // | g-----|--h
310 // | / | /
311 // |/ |/
312 // c--------d
313
314 // length (ab)
315 L_3D,
316 // height (ac)
317 H_3D,
318 // depth (ae)
319 P_3D,
320 // diagonal1 (bc)
321 D1_3D,
322 // diagonal2 (ag)
323 D2_3D,
324 // diagonal3 (be)
325 D3_3D,
326 // diagonal4 (bg)
327 D4_3D,
328
329 // length (ab)
330 L_2D,
331 // height (ac)
332 H_2D,
333 // diagonal1 (bc)
334 D1_2D,
335
338 CENTER_1D,
339 };
340
341 enum class TrianglePosition : char {
342 // e--------f
343 // /| /|
344 // / | / |
345 // a--------b |
346 // | g-----|--h
347 // | / | /
348 // |/ |/
349 // c--------d
350
351 F_3D, // face (abc, bcd)
352 C_3D, // side (abe, bef)
353 H_3D, // top (acg, aeg)
354 D1_3D, // diagonal1 (bdg, beg)
355 D2_3D, // diagonal2 (abg, bgh)
356 D3_3D, // diagonal3 (bcg, bfg)
357
358 TOP_2D, // abc
359 BOTTOM_2D, // bcd
360 };
361
363
364 // acceleration functions
365 int checkAcceleration();
366 bool isPowerOfTwo(unsigned long long int v, unsigned long long int &r);
367
368 //\cond
369 // 2D //
370 void vertexToPosition2d(const SimplexId vertex,
371 SimplexId p[2]) const override;
372 void
373 edgeToPosition2d(const SimplexId edge, const int k, SimplexId p[2]) const;
374 void triangleToPosition2d(const SimplexId triangle,
375 SimplexId p[2]) const override;
376
377 SimplexId getVertexNeighbor2d(const SimplexId p[2],
378 const SimplexId v,
379 const int id) const;
380 SimplexId getVertexEdge2d(const SimplexId p[2], const int id) const;
381 SimplexId getVertexStar2d(const SimplexId p[2], const int id) const;
382 SimplexId getVertexLink2d(const SimplexId p[2], const int id) const;
383
384 SimplexId getEdgeTriangle2dL(const SimplexId p[3], const int id) const;
385 SimplexId getEdgeTriangle2dH(const SimplexId p[3], const int id) const;
386 SimplexId getEdgeTriangle2dD1(const SimplexId p[3], const int id) const;
387
388 SimplexId getEdgeLink2dL(const SimplexId p[2], const int id) const;
389 SimplexId getEdgeLink2dH(const SimplexId p[2], const int id) const;
390 SimplexId getEdgeLink2dD1(const SimplexId p[2], const int id) const;
391
392 SimplexId getEdgeStar2dL(const SimplexId p[2], const int id) const;
393 SimplexId getEdgeStar2dH(const SimplexId p[2], const int id) const;
394
395 // 3D //
396 void vertexToPosition(const SimplexId vertex,
397 SimplexId p[3]) const override;
398 void
399 edgeToPosition(const SimplexId edge, const int k, SimplexId p[3]) const;
400 void triangleToPosition(const SimplexId triangle,
401 const int k,
402 SimplexId p[3]) const override;
403 void tetrahedronToPosition(const SimplexId tetrahedron,
404 SimplexId p[3]) const override;
405
406 SimplexId getVertexNeighbor3d(const SimplexId p[3],
407 const SimplexId v,
408 const int id) const;
409 SimplexId getVertexEdge3d(const SimplexId p[3], const int id) const;
410 SimplexId getVertexTriangle3d(const SimplexId p[3], const int id) const;
411 SimplexId getVertexLink3d(const SimplexId p[3], const int id) const;
412 SimplexId getVertexStar3d(const SimplexId p[3], const int id) const;
413
414 SimplexId getEdgeTriangle3dL(const SimplexId p[3], const int id) const;
415 SimplexId getEdgeTriangle3dH(const SimplexId p[3], const int id) const;
416 SimplexId getEdgeTriangle3dP(const SimplexId p[3], const int id) const;
417 SimplexId getEdgeTriangle3dD1(const SimplexId p[3], const int id) const;
418 SimplexId getEdgeTriangle3dD2(const SimplexId p[3], const int id) const;
419 SimplexId getEdgeTriangle3dD3(const SimplexId p[3], const int id) const;
420 SimplexId getEdgeTriangle3dD4(const SimplexId p[3], const int id) const;
421
422 SimplexId getEdgeLinkL(const SimplexId p[3], const int id) const;
423 SimplexId getEdgeLinkH(const SimplexId p[3], const int id) const;
424 SimplexId getEdgeLinkP(const SimplexId p[3], const int id) const;
425 SimplexId getEdgeLinkD1(const SimplexId p[3], const int id) const;
426 SimplexId getEdgeLinkD2(const SimplexId p[3], const int id) const;
427 SimplexId getEdgeLinkD3(const SimplexId p[3], const int id) const;
428 SimplexId getEdgeLinkD4(const SimplexId p[3], const int id) const;
429
430 SimplexId getEdgeStarL(const SimplexId p[3], const int id) const;
431 SimplexId getEdgeStarH(const SimplexId p[3], const int id) const;
432 SimplexId getEdgeStarP(const SimplexId p[3], const int id) const;
433 SimplexId getEdgeStarD1(const SimplexId p[3], const int id) const;
434 SimplexId getEdgeStarD2(const SimplexId p[3], const int id) const;
435 SimplexId getEdgeStarD3(const SimplexId p[3], const int id) const;
436
437 SimplexId getTriangleVertexF(const SimplexId p[3], const int id) const;
438 SimplexId getTriangleVertexH(const SimplexId p[3], const int id) const;
439 SimplexId getTriangleVertexC(const SimplexId p[3], const int id) const;
440 SimplexId getTriangleVertexD1(const SimplexId p[3], const int id) const;
441 SimplexId getTriangleVertexD2(const SimplexId p[3], const int id) const;
442 SimplexId getTriangleVertexD3(const SimplexId p[3], const int id) const;
443
444 SimplexId getTriangleEdgeF_0(const SimplexId p[3], const int id) const;
445 SimplexId getTriangleEdgeF_1(const SimplexId p[3], const int id) const;
446 SimplexId getTriangleEdgeH_0(const SimplexId p[3], const int id) const;
447 SimplexId getTriangleEdgeH_1(const SimplexId p[3], const int id) const;
448 SimplexId getTriangleEdgeC_0(const SimplexId p[3], const int id) const;
449 SimplexId getTriangleEdgeC_1(const SimplexId p[3], const int id) const;
450 SimplexId getTriangleEdgeD1_0(const SimplexId p[3], const int id) const;
451 SimplexId getTriangleEdgeD1_1(const SimplexId p[3], const int id) const;
452 SimplexId getTriangleEdgeD2_0(const SimplexId p[3], const int id) const;
453 SimplexId getTriangleEdgeD2_1(const SimplexId p[3], const int id) const;
454 SimplexId getTriangleEdgeD3_0(const SimplexId p[3], const int id) const;
455 SimplexId getTriangleEdgeD3_1(const SimplexId p[3], const int id) const;
456
457 SimplexId getTriangleLinkF(const SimplexId p[3], const int id) const;
458 SimplexId getTriangleLinkH(const SimplexId p[3], const int id) const;
459 SimplexId getTriangleLinkC(const SimplexId p[3], const int id) const;
460 SimplexId getTriangleLinkD1(const SimplexId p[3], const int id) const;
461 SimplexId getTriangleLinkD2(const SimplexId p[3], const int id) const;
462 SimplexId getTriangleLinkD3(const SimplexId p[3], const int id) const;
463
464 SimplexId getTriangleStarF(const SimplexId p[3], const int id) const;
465 SimplexId getTriangleStarH(const SimplexId p[3], const int id) const;
466 SimplexId getTriangleStarC(const SimplexId p[3], const int id) const;
467 SimplexId getTriangleStarD1(const SimplexId p[3], const int id) const;
468 SimplexId getTriangleStarD2(const SimplexId p[3], const int id) const;
469 SimplexId getTriangleStarD3(const SimplexId p[3], const int id) const;
470
471 SimplexId getTetrahedronVertexABCG(const SimplexId p[3],
472 const int id) const;
473 SimplexId getTetrahedronVertexBCDG(const SimplexId p[3],
474 const int id) const;
475 SimplexId getTetrahedronVertexABEG(const SimplexId p[3],
476 const int id) const;
477 SimplexId getTetrahedronVertexBEFG(const SimplexId p[3],
478 const int id) const;
479 SimplexId getTetrahedronVertexBFGH(const SimplexId p[3],
480 const int id) const;
481 SimplexId getTetrahedronVertexBDGH(const SimplexId p[3],
482 const int id) const;
483
484 SimplexId getTetrahedronEdgeABCG(const SimplexId p[3], const int id) const;
485 SimplexId getTetrahedronEdgeBCDG(const SimplexId p[3], const int id) const;
486 SimplexId getTetrahedronEdgeABEG(const SimplexId p[3], const int id) const;
487 SimplexId getTetrahedronEdgeBEFG(const SimplexId p[3], const int id) const;
488 SimplexId getTetrahedronEdgeBFGH(const SimplexId p[3], const int id) const;
489 SimplexId getTetrahedronEdgeBDGH(const SimplexId p[3], const int id) const;
490
491 SimplexId getTetrahedronTriangleABCG(const SimplexId p[3],
492 const int id) const;
493 SimplexId getTetrahedronTriangleBCDG(const SimplexId p[3],
494 const int id) const;
495 SimplexId getTetrahedronTriangleABEG(const SimplexId p[3],
496 const int id) const;
497 SimplexId getTetrahedronTriangleBEFG(const SimplexId p[3],
498 const int id) const;
499 SimplexId getTetrahedronTriangleBFGH(const SimplexId p[3],
500 const int id) const;
501 SimplexId getTetrahedronTriangleBDGH(const SimplexId p[3],
502 const int id) const;
503
504 SimplexId getTetrahedronNeighborABCG(const SimplexId t,
505 const SimplexId p[3],
506 const int id) const;
507 SimplexId getTetrahedronNeighborBCDG(const SimplexId t,
508 const SimplexId p[3],
509 const int id) const;
510 SimplexId getTetrahedronNeighborABEG(const SimplexId t,
511 const SimplexId p[3],
512 const int id) const;
513 SimplexId getTetrahedronNeighborBEFG(const SimplexId t,
514 const SimplexId p[3],
515 const int id) const;
516 SimplexId getTetrahedronNeighborBFGH(const SimplexId t,
517 const SimplexId p[3],
518 const int id) const;
519 SimplexId getTetrahedronNeighborBDGH(const SimplexId t,
520 const SimplexId p[3],
521 const int id) const;
522 //\endcond
523
524#ifdef TTK_ENABLE_MPI
525
526 protected:
527 int preconditionDistributedCells() override;
528 // std::shared_ptr<PeriodicImplicitTriangulation> metaGrid_;
529 std::array<unsigned char, 6> isBoundaryPeriodic{};
530
531 public:
532 void createMetaGrid(const double *const bounds) override;
533 void setIsBoundaryPeriodic(std::array<unsigned char, 6> boundary);
534 int getCellRankInternal(const SimplexId lcid) const override;
535
536 protected:
537 std::array<SimplexId, 3>
538 getVertGlobalCoords(const SimplexId lvid) const override;
539 std::array<SimplexId, 3>
540 getVertLocalCoords(const SimplexId gvid) const override;
541
542#endif // TTK_ENABLE_MPI
543 };
544
545 template <typename Derived>
548 inline Derived &underlying() {
549 return static_cast<Derived &>(*this);
550 }
551 inline Derived const &underlying() const {
552 return static_cast<Derived const &>(*this);
553 }
554
555 public:
557 const SimplexId &vertexId,
558 const int &localNeighborId,
559 SimplexId &neighborId) const override;
560
561 int getVertexEdgeInternal(const SimplexId &vertexId,
562 const int &id,
563 SimplexId &edgeId) const override;
564
565 int getVertexTriangleInternal(const SimplexId &vertexId,
566 const int &id,
567 SimplexId &triangleId) const override;
568
570 const SimplexId &vertexId,
571 const int &localLinkId,
572 SimplexId &linkId) const override;
573
575 const SimplexId &vertexId,
576 const int &localStarId,
577 SimplexId &starId) const override;
578
580 float &x,
581 float &y,
582 float &z) const override;
583
584 int getEdgeVertexInternal(const SimplexId &edgeId,
585 const int &localVertexId,
586 SimplexId &vertexId) const override;
587
589 getEdgeTriangleNumberInternal(const SimplexId &edgeId) const override;
590
591 int getEdgeTriangleInternal(const SimplexId &edgeId,
592 const int &id,
593 SimplexId &triangleId) const override;
594
595 int
597 const int &localLinkId,
598 SimplexId &linkId) const override;
599
601 const SimplexId &edgeId) const override;
602
603 int
605 const int &localStarId,
606 SimplexId &starId) const override;
607
608 int getTriangleVertexInternal(const SimplexId &triangleId,
609 const int &localVertexId,
610 SimplexId &vertexId) const override;
611
612 int getTriangleEdgeInternal(const SimplexId &triangleId,
613 const int &id,
614 SimplexId &edgeId) const override;
615
617 const SimplexId &triangleId,
618 const int &localLinkId,
619 SimplexId &linkId) const override;
620
622 const SimplexId &triangleId,
623 const int &localStarId,
624 SimplexId &starId) const override;
625
627 const SimplexId &triangleId) const override;
628
629 int getTriangleNeighbor(const SimplexId &triangleId,
630 const int &localNeighborId,
631 SimplexId &neighborId) const override;
632
633 int getTetrahedronVertex(const SimplexId &tetId,
634 const int &localVertexId,
635 SimplexId &vertexId) const override;
636
637 int getTetrahedronEdge(const SimplexId &tetId,
638 const int &id,
639 SimplexId &edgeId) const override;
640
641 int getTetrahedronTriangle(const SimplexId &tetId,
642 const int &id,
643 SimplexId &triangleId) const override;
644
645 int getTetrahedronNeighbor(const SimplexId &tetId,
646 const int &localNeighborId,
647 SimplexId &neighborId) const override;
648
652 virtual int getEdgeIncenter(SimplexId edgeId, float incenter[3]) const {
653 SimplexId v0{}, v1{};
654 getEdgeVertexInternal(edgeId, 0, v0);
655 getEdgeVertexInternal(edgeId, 1, v1);
656
657 std::array<float, 3> p0{}, p1{};
658 getVertexPointInternal(v0, p0[0], p0[1], p0[2]);
659 getVertexPointInternal(v1, p1[0], p1[1], p1[2]);
660
661 const auto &ind0 = this->underlying().getVertexCoords(v0);
662 const auto &ind1 = this->underlying().getVertexCoords(v1);
663
664 for(int i = 0; i < dimensionality_; ++i) {
665 if(ind1[i] == nbvoxels_[i]) {
666 p0[i] += (ind0[i] == 0) * dimensions_[i] * spacing_[i];
667 } else if(ind0[i] == nbvoxels_[i]) {
668 p1[i] += (ind1[i] == 0) * dimensions_[i] * spacing_[i];
669 }
670 }
671
672 for(int i = 0; i < 3; ++i) {
673 incenter[i] = 0.5f * (p0[i] + p1[i]);
674 }
675
676 return 0;
677 }
678
683 virtual int getTriangleIncenter(SimplexId triangleId,
684 float incenter[3]) const {
685
686 SimplexId v0{}, v1{}, v2{};
687 getTriangleVertexInternal(triangleId, 0, v0);
688 getTriangleVertexInternal(triangleId, 1, v1);
689 getTriangleVertexInternal(triangleId, 2, v2);
690
691 std::array<float, 3> p0{}, p1{}, p2{};
692 getVertexPointInternal(v0, p0[0], p0[1], p0[2]);
693 getVertexPointInternal(v1, p1[0], p1[1], p1[2]);
694 getVertexPointInternal(v2, p2[0], p2[1], p2[2]);
695
696 const auto &ind0 = this->underlying().getVertexCoords(v0);
697 const auto &ind1 = this->underlying().getVertexCoords(v1);
698 const auto &ind2 = this->underlying().getVertexCoords(v2);
699
700 for(int i = 0; i < dimensionality_; ++i) {
701 if(ind0[i] == nbvoxels_[i]) {
702 p1[i] += (ind1[i] == 0) * dimensions_[i] * spacing_[i];
703 p2[i] += (ind2[i] == 0) * dimensions_[i] * spacing_[i];
704 } else if(ind1[i] == nbvoxels_[i]) {
705 p0[i] += (ind0[i] == 0) * dimensions_[i] * spacing_[i];
706 p2[i] += (ind2[i] == 0) * dimensions_[i] * spacing_[i];
707 } else if(ind2[i] == nbvoxels_[i]) {
708 p0[i] += (ind0[i] == 0) * dimensions_[i] * spacing_[i];
709 p1[i] += (ind1[i] == 0) * dimensions_[i] * spacing_[i];
710 }
711 }
712
713 std::array<float, 3> d{Geometry::distance(p1.data(), p2.data()),
714 Geometry::distance(p2.data(), p0.data()),
715 Geometry::distance(p0.data(), p1.data())};
716 const float sum = d[0] + d[1] + d[2];
717 for(int i = 0; i < 3; ++i) {
718 incenter[i] = (d[0] * p0[i] + d[1] * p1[i] + d[2] * p2[i]) / sum;
719 }
720
721 return 0;
722 }
723
728 virtual int getTetraIncenter(SimplexId tetraId, float incenter[3]) const {
729
730 SimplexId v0{}, v1{}, v2{}, v3{};
731 getCellVertexInternal(tetraId, 0, v0);
732 getCellVertexInternal(tetraId, 1, v1);
733 getCellVertexInternal(tetraId, 2, v2);
734 getCellVertexInternal(tetraId, 3, v3);
735
736 std::array<float, 3> p0{}, p1{}, p2{}, p3{};
737 getVertexPointInternal(v0, p0[0], p0[1], p0[2]);
738 getVertexPointInternal(v1, p1[0], p1[1], p1[2]);
739 getVertexPointInternal(v2, p2[0], p2[1], p2[2]);
740 getVertexPointInternal(v3, p3[0], p3[1], p3[2]);
741
742 const auto &ind0 = this->underlying().getVertexCoords(v0);
743 const auto &ind1 = this->underlying().getVertexCoords(v1);
744 const auto &ind2 = this->underlying().getVertexCoords(v2);
745 const auto &ind3 = this->underlying().getVertexCoords(v3);
746
747 for(int i = 0; i < dimensionality_; ++i) {
748 if(ind0[i] == nbvoxels_[i]) {
749 p1[i] += (ind1[i] == 0) * dimensions_[i] * spacing_[i];
750 p2[i] += (ind2[i] == 0) * dimensions_[i] * spacing_[i];
751 p3[i] += (ind3[i] == 0) * dimensions_[i] * spacing_[i];
752 } else if(ind1[i] == nbvoxels_[i]) {
753 p0[i] += (ind0[i] == 0) * dimensions_[i] * spacing_[i];
754 p2[i] += (ind2[i] == 0) * dimensions_[i] * spacing_[i];
755 p3[i] += (ind3[i] == 0) * dimensions_[i] * spacing_[i];
756 } else if(ind2[i] == nbvoxels_[i]) {
757 p0[i] += (ind0[i] == 0) * dimensions_[i] * spacing_[i];
758 p1[i] += (ind1[i] == 0) * dimensions_[i] * spacing_[i];
759 p3[i] += (ind3[i] == 0) * dimensions_[i] * spacing_[i];
760 } else if(ind3[i] == nbvoxels_[i]) {
761 p0[i] += (ind0[i] == 0) * dimensions_[i] * spacing_[i];
762 p1[i] += (ind1[i] == 0) * dimensions_[i] * spacing_[i];
763 p2[i] += (ind2[i] == 0) * dimensions_[i] * spacing_[i];
764 }
765 }
766
767 for(int i = 0; i < 3; ++i) {
768 incenter[i] = 0.25f * (p0[i] + p1[i] + p2[i] + p3[i]);
769 }
770 return 0;
771 }
772 };
773
774} // namespace ttk
775
777
778inline void
780 SimplexId p[2]) const {
781 if(isAccelerated_) {
782 p[0] = vertex & mod_[0];
783 p[1] = vertex >> div_[0];
784 } else {
785 p[0] = vertex % vshift_[0];
786 p[1] = vertex / vshift_[0];
787 }
788}
789
790inline void ttk::PeriodicImplicitTriangulation::edgeToPosition2d(
791 const SimplexId edge, const int k, SimplexId p[2]) const {
792 const int e = (k) ? edge - esetshift_[k - 1] : edge;
793 p[0] = e % eshift_[2 * k];
794 p[1] = e / eshift_[2 * k];
795}
796
798 const SimplexId triangle, SimplexId p[2]) const {
799 p[0] = triangle % tshift_[0];
800 p[1] = triangle / tshift_[0];
801}
802
803inline ttk::SimplexId ttk::PeriodicImplicitTriangulation::getVertexNeighbor2d(
804 const SimplexId p[2], const SimplexId v, const int id) const {
805 SimplexId wrapXLeft = 0;
806 SimplexId wrapXRight = 0;
807 SimplexId wrapYTop = 0;
808 SimplexId wrapYBottom = 0;
809 if(p[0] == 0)
810 wrapXLeft = wrap_[0];
811 if(p[0] == nbvoxels_[Di_])
812 wrapXRight = -wrap_[0];
813 if(p[1] == 0)
814 wrapYTop = wrap_[1];
815 if(p[1] == nbvoxels_[Dj_])
816 wrapYBottom = -wrap_[1];
817 switch(id) {
818 case 0:
819 return v - 1 + wrapXLeft;
820 case 1:
821 return v - vshift_[0] + wrapYTop;
822 case 2:
823 return v - vshift_[0] + 1 + wrapXRight + wrapYTop;
824 case 3:
825 return v + 1 + wrapXRight;
826 case 4:
827 return v + vshift_[0] + wrapYBottom;
828 case 5:
829 return v + vshift_[0] - 1 + wrapXLeft + wrapYBottom;
830 }
831 return -1;
832}
833
834inline ttk::SimplexId
835 ttk::PeriodicImplicitTriangulation::getVertexEdge2d(const SimplexId p[2],
836 const int id) const {
837 SimplexId wrapXLeft = 0;
838 SimplexId wrapYTop = 0;
839 if(p[0] == 0)
840 wrapXLeft = wrap_[0];
841 if(p[1] == 0)
842 wrapYTop = wrap_[1];
843 switch(id) {
844 case 0:
845 return esetshift_[0] + p[0] + (p[1] - 1) * eshift_[2] + wrapYTop;
846 case 1:
847 return p[0] + p[1] * eshift_[0] - 1 + wrapXLeft;
848 case 2:
849 return esetshift_[1] + p[0] + (p[1] - 1) * eshift_[4] + wrapYTop;
850 case 3:
851 return p[0] + p[1] * eshift_[0];
852 case 4:
853 return esetshift_[0] + p[0] + p[1] * eshift_[2];
854 case 5:
855 return esetshift_[1] + p[0] + p[1] * eshift_[4] - 1 + wrapXLeft;
856 }
857 return -1;
858}
859
860inline ttk::SimplexId
861 ttk::PeriodicImplicitTriangulation::getVertexStar2d(const SimplexId p[2],
862 const int id) const {
863 SimplexId wrapXLeft = 0;
864 SimplexId wrapYTop = 0;
865 if(p[0] == 0)
866 wrapXLeft = 2 * wrap_[0];
867 if(p[1] == 0)
868 wrapYTop = 2 * wrap_[1];
869 switch(id) {
870 case 0:
871 return (p[0] - 1) * 2 + p[1] * tshift_[0] + wrapXLeft;
872 case 1:
873 return (p[0] - 1) * 2 + p[1] * tshift_[0] + 1 + wrapXLeft;
874 case 2:
875 return p[0] * 2 + p[1] * tshift_[0];
876 case 3:
877 return p[0] * 2 + (p[1] - 1) * tshift_[0] + wrapYTop;
878 case 4:
879 return p[0] * 2 + (p[1] - 1) * tshift_[0] + 1 + wrapYTop;
880 case 5:
881 return (p[0] - 1) * 2 + (p[1] - 1) * tshift_[0] + 1 + wrapXLeft
882 + wrapYTop;
883 }
884 return -1;
885}
886
887inline ttk::SimplexId
888 ttk::PeriodicImplicitTriangulation::getVertexLink2d(const SimplexId p[2],
889 const int id) const {
890 SimplexId wrapXLeft = 0;
891 SimplexId wrapXRight = 0;
892 SimplexId wrapYTop = 0;
893 SimplexId wrapYBottom = 0;
894 if(p[0] == 0)
895 wrapXLeft = wrap_[0];
896 if(p[0] == nbvoxels_[Di_])
897 wrapXRight = -wrap_[0];
898 if(p[1] == 0)
899 wrapYTop = wrap_[1];
900 if(p[1] == nbvoxels_[Dj_])
901 wrapYBottom = -wrap_[1];
902 switch(id) {
903 case 0:
904 return esetshift_[0] + p[0] + p[1] * eshift_[2] - 1 + wrapXLeft;
905 case 1:
906 return p[0] + (p[1] + 1) * eshift_[0] - 1 + wrapXLeft + wrapYBottom;
907 case 2:
908 return esetshift_[1] + p[0] + p[1] * eshift_[4];
909 case 3:
910 return esetshift_[0] + p[0] + (p[1] - 1) * eshift_[2] + 1 + wrapXRight
911 + wrapYTop;
912 case 4:
913 return p[0] + (p[1] - 1) * eshift_[0] + wrapYTop;
914 case 5:
915 return esetshift_[1] + p[0] + (p[1] - 1) * eshift_[4] - 1 + wrapXLeft
916 + wrapYTop;
917 }
918 return -1;
919}
920
921inline ttk::SimplexId
922 ttk::PeriodicImplicitTriangulation::getEdgeTriangle2dL(const SimplexId p[3],
923 const int id) const {
924 SimplexId wrapYTop = 0;
925 if(p[1] == 0)
926 wrapYTop = 2 * wrap_[1];
927 switch(id) {
928 case 0:
929 return p[Di_] * 2 + p[Dj_] * tshift_[0];
930 case 1:
931 return p[Di_] * 2 + (p[Dj_] - 1) * tshift_[0] + 1 + wrapYTop;
932 }
933 return -1;
934}
935
936inline ttk::SimplexId
937 ttk::PeriodicImplicitTriangulation::getEdgeTriangle2dH(const SimplexId p[3],
938 const int id) const {
939 SimplexId wrapXLeft = 0;
940 if(p[0] == 0)
941 wrapXLeft = 2 * wrap_[0];
942 switch(id) {
943 case 0:
944 return p[Di_] * 2 + p[Dj_] * tshift_[0];
945 case 1:
946 return (p[Di_] - 1) * 2 + p[Dj_] * tshift_[0] + 1 + wrapXLeft;
947 }
948 return -1;
949}
950
951inline ttk::SimplexId
952 ttk::PeriodicImplicitTriangulation::getEdgeTriangle2dD1(const SimplexId p[3],
953 const int id) const {
954 switch(id) {
955 case 0:
956 return p[Di_] * 2 + p[Dj_] * tshift_[0];
957 case 1:
958 return p[Di_] * 2 + p[Dj_] * tshift_[0] + 1;
959 }
960 return -1;
961}
962
963inline ttk::SimplexId
964 ttk::PeriodicImplicitTriangulation::getEdgeLink2dL(const SimplexId p[2],
965 const int id) const {
966 if(p[1] > 0 and p[1] < nbvoxels_[Dj_]) {
967 switch(id) {
968 case 0:
969 return p[0] + (p[1] + 1) * vshift_[0];
970 case 1:
971 return ((p[0] < nbvoxels_[Di_])
972 ? (p[0] + (p[1] - 1) * vshift_[0] + 1)
973 : (p[0] + (p[1] - 1) * vshift_[0] + 1 - wrap_[0]));
974 }
975 return -1;
976 } else if(p[1] == 0) {
977 switch(id) {
978 case 0:
979 return p[0] + (p[1] + 1) * vshift_[0];
980 case 1:
981 return ((p[0] < nbvoxels_[Di_])
982 ? (p[0] + (p[1] - 1) * vshift_[0] + 1)
983 : (p[0] + (p[1] - 1) * vshift_[0] + 1 - wrap_[0]))
984 + wrap_[1];
985 }
986 return -1;
987 } else {
988 switch(id) {
989 case 0:
990 return p[0] + (p[1] + 1) * vshift_[0] - wrap_[1];
991 case 1:
992 return ((p[0] < nbvoxels_[Di_])
993 ? (p[0] + (p[1] - 1) * vshift_[0] + 1)
994 : (p[0] + (p[1] - 1) * vshift_[0] + 1 - wrap_[0]));
995 }
996 return -1;
997 }
998}
999
1000inline ttk::SimplexId
1001 ttk::PeriodicImplicitTriangulation::getEdgeLink2dH(const SimplexId p[2],
1002 const int id) const {
1003 if(p[0] > 0 and p[0] < nbvoxels_[Di_]) {
1004 switch(id) {
1005 case 0:
1006 return p[0] + p[1] * vshift_[0] + 1;
1007 case 1:
1008 return ((p[1] < nbvoxels_[Dj_])
1009 ? (p[0] + (p[1] + 1) * vshift_[0] - 1)
1010 : (p[0] + (p[1] + 1) * vshift_[0] - 1 - wrap_[1]));
1011 }
1012 return -1;
1013 } else if(p[0] == 0) {
1014 switch(id) {
1015 case 0:
1016 return p[0] + p[1] * vshift_[0] + 1 + wrap_[0];
1017 case 1:
1018 return ((p[1] < nbvoxels_[Dj_])
1019 ? (p[0] + (p[1] + 1) * vshift_[0] - 1)
1020 : (p[0] + (p[1] + 1) * vshift_[0] - 1 - wrap_[1]));
1021 }
1022 return -1;
1023 } else {
1024 switch(id) {
1025 case 0:
1026 return p[0] + p[1] * vshift_[0] + 1;
1027 case 1:
1028 return ((p[1] < nbvoxels_[Dj_])
1029 ? (p[0] + (p[1] + 1) * vshift_[0] - 1)
1030 : (p[0] + (p[1] + 1) * vshift_[0] - 1 - wrap_[1]))
1031 - wrap_[0];
1032 }
1033 return -1;
1034 }
1035}
1036
1037inline ttk::SimplexId
1038 ttk::PeriodicImplicitTriangulation::getEdgeLink2dD1(const SimplexId p[2],
1039 const int id) const {
1040 SimplexId wrapX = (p[0] < nbvoxels_[Di_]) ? 0 : wrap_[0];
1041 SimplexId wrapY = (p[1] < nbvoxels_[Dj_]) ? 0 : wrap_[1];
1042 switch(id) {
1043 case 0:
1044 return p[0] + p[1] * vshift_[0];
1045 case 1:
1046 return p[0] + (p[1] + 1) * vshift_[0] + 1 - wrapX - wrapY;
1047 }
1048 return -1;
1049}
1050
1051inline ttk::SimplexId
1052 ttk::PeriodicImplicitTriangulation::getEdgeStar2dL(const SimplexId p[2],
1053 const int id) const {
1054 if(p[1] == 0) {
1055 switch(id) {
1056 case 0:
1057 return p[0] * 2 + p[1] * tshift_[0];
1058 case 1:
1059 return p[0] * 2 + (p[1] - 1) * tshift_[0] + 1 + 2 * wrap_[1];
1060 }
1061 return -1;
1062 } else {
1063 switch(id) {
1064 case 0:
1065 return p[0] * 2 + p[1] * tshift_[0];
1066 case 1:
1067 return p[0] * 2 + (p[1] - 1) * tshift_[0] + 1;
1068 }
1069 return -1;
1070 }
1071}
1072
1073inline ttk::SimplexId
1074 ttk::PeriodicImplicitTriangulation::getEdgeStar2dH(const SimplexId p[2],
1075 const int id) const {
1076 if(p[0] == 0) {
1077 switch(id) {
1078 case 0:
1079 return p[0] * 2 + p[1] * tshift_[0];
1080 case 1:
1081 return (p[0] - 1) * 2 + p[1] * tshift_[0] + 1 + 2 * wrap_[0];
1082 }
1083 return -1;
1084 } else {
1085 switch(id) {
1086 case 0:
1087 return p[0] * 2 + p[1] * tshift_[0];
1088 case 1:
1089 return (p[0] - 1) * 2 + p[1] * tshift_[0] + 1;
1090 }
1091 return -1;
1092 }
1093}
1094
1095inline void
1097 SimplexId p[3]) const {
1098 if(isAccelerated_) {
1099 p[0] = vertex & mod_[0];
1100 p[1] = (vertex & mod_[1]) >> div_[0];
1101 p[2] = vertex >> div_[1];
1102 } else {
1103 p[0] = vertex % vshift_[0];
1104 p[1] = (vertex % vshift_[1]) / vshift_[0];
1105 p[2] = vertex / vshift_[1];
1106 }
1107}
1108
1109inline void ttk::PeriodicImplicitTriangulation::edgeToPosition(
1110 const SimplexId edge, const int k, SimplexId p[3]) const {
1111 const int e = (k) ? edge - esetshift_[k - 1] : edge;
1112 p[0] = e % eshift_[2 * k];
1113 p[1] = (e % eshift_[2 * k + 1]) / eshift_[2 * k];
1114 p[2] = e / eshift_[2 * k + 1];
1115}
1116
1118 const SimplexId triangle, const int k, SimplexId p[3]) const {
1119 const SimplexId t = (k) ? triangle - tsetshift_[k - 1] : triangle;
1120 p[0] = t % tshift_[2 * k];
1121 p[1] = (t % tshift_[2 * k + 1]) / tshift_[2 * k];
1122 p[2] = t / tshift_[2 * k + 1];
1123}
1124
1126 const SimplexId tetrahedron, SimplexId p[3]) const {
1127 p[0] = (tetrahedron % tetshift_[0]) / 6;
1128 p[1] = (tetrahedron % tetshift_[1]) / tetshift_[0];
1129 p[2] = tetrahedron / tetshift_[1];
1130}
1131
1132inline ttk::SimplexId ttk::PeriodicImplicitTriangulation::getVertexNeighbor3d(
1133 const SimplexId p[3], const SimplexId v, const int id) const {
1134 SimplexId wrapXLeft = 0;
1135 SimplexId wrapXRight = 0;
1136 SimplexId wrapYTop = 0;
1137 SimplexId wrapYBottom = 0;
1138 SimplexId wrapZBack = 0;
1139 SimplexId wrapZFront = 0;
1140 if(p[0] == 0)
1141 wrapXLeft = wrap_[0];
1142 if(p[0] == nbvoxels_[0])
1143 wrapXRight = -wrap_[0];
1144 if(p[1] == 0)
1145 wrapYTop = wrap_[1];
1146 if(p[1] == nbvoxels_[1])
1147 wrapYBottom = -wrap_[1];
1148 if(p[2] == 0)
1149 wrapZBack = wrap_[2];
1150 if(p[2] == nbvoxels_[2])
1151 wrapZFront = -wrap_[2];
1152 switch(id) {
1153 case 0:
1154 return v - vshift_[0] - vshift_[1] + wrapYTop + wrapZBack;
1155 case 1:
1156 return v + 1 - vshift_[0] - vshift_[1] + wrapXRight + wrapYTop
1157 + wrapZBack;
1158 case 2:
1159 return v - vshift_[1] + wrapZBack;
1160 case 3:
1161 return v + 1 - vshift_[1] + wrapXRight + wrapZBack;
1162 case 4:
1163 return v - vshift_[0] + wrapYTop;
1164 case 5:
1165 return v + 1 - vshift_[0] + wrapXRight + wrapYTop;
1166 case 6:
1167 return v + 1 + wrapXRight;
1168 case 7:
1169 return v - 1 + vshift_[1] + wrapXLeft + wrapZFront;
1170 case 8:
1171 return v + vshift_[1] + wrapZFront;
1172 case 9:
1173 return v - 1 + wrapXLeft;
1174 case 10:
1175 return v - 1 + vshift_[0] + wrapXLeft + wrapYBottom;
1176 case 11:
1177 return v + vshift_[0] + wrapYBottom;
1178 case 12:
1179 return v - 1 + vshift_[0] + vshift_[1] + wrapXLeft + wrapYBottom
1180 + wrapZFront;
1181 case 13:
1182 return v + vshift_[0] + vshift_[1] + wrapYBottom + wrapZFront;
1183 }
1184 return -1;
1185}
1186
1187inline ttk::SimplexId
1188 ttk::PeriodicImplicitTriangulation::getVertexEdge3d(const SimplexId p[3],
1189 const int id) const {
1190 SimplexId wrapXLeft = 0;
1191 SimplexId wrapYTop = 0;
1192 SimplexId wrapZBack = 0;
1193 if(p[0] == 0)
1194 wrapXLeft = wrap_[0];
1195 if(p[1] == 0)
1196 wrapYTop = wrap_[1];
1197 if(p[2] == 0)
1198 wrapZBack = wrap_[2];
1199 switch(id) {
1200 case 0:
1201 return esetshift_[3] + p[0] + (p[1] - 1) * eshift_[8]
1202 + (p[2] - 1) * eshift_[9] + wrapYTop + wrapZBack;
1203 case 1:
1204 return esetshift_[5] + p[0] + (p[1] - 1) * eshift_[12]
1205 + (p[2] - 1) * eshift_[13] + wrapYTop + wrapZBack;
1206 case 2:
1207 return esetshift_[1] + p[0] + p[1] * eshift_[4] + (p[2] - 1) * eshift_[5]
1208 + wrapZBack;
1209 case 3:
1210 return esetshift_[4] + p[0] + p[1] * eshift_[10]
1211 + (p[2] - 1) * eshift_[11] + wrapZBack;
1212 case 4:
1213 return esetshift_[0] + p[0] + (p[1] - 1) * eshift_[2] + p[2] * eshift_[3]
1214 + wrapYTop;
1215 case 5:
1216 return esetshift_[2] + p[0] + (p[1] - 1) * eshift_[6] + p[2] * eshift_[7]
1217 + wrapYTop;
1218 case 6:
1219 return p[0] + p[1] * eshift_[0] + p[2] * eshift_[1];
1220 case 7:
1221 return esetshift_[4] + p[0] + p[1] * eshift_[10] + p[2] * eshift_[11] - 1
1222 + wrapXLeft;
1223 case 8:
1224 return esetshift_[1] + p[0] + p[1] * eshift_[4] + p[2] * eshift_[5];
1225 case 9:
1226 return p[0] + p[1] * eshift_[0] + p[2] * eshift_[1] - 1 + wrapXLeft;
1227 case 10:
1228 return esetshift_[2] + p[0] + p[1] * eshift_[6] + p[2] * eshift_[7] - 1
1229 + wrapXLeft;
1230 case 11:
1231 return esetshift_[0] + p[0] + p[1] * eshift_[2] + p[2] * eshift_[3];
1232 case 12:
1233 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13] - 1
1234 + wrapXLeft;
1235 case 13:
1236 return esetshift_[3] + p[0] + p[1] * eshift_[8] + p[2] * eshift_[9];
1237 }
1238 return -1;
1239}
1240
1241inline ttk::SimplexId
1242 ttk::PeriodicImplicitTriangulation::getVertexLink3d(const SimplexId p[3],
1243 const int id) const {
1244 SimplexId wrapXLeft = 0;
1245 SimplexId wrapXRight = 0;
1246 SimplexId wrapYTop = 0;
1247 SimplexId wrapYBottom = 0;
1248 SimplexId wrapZBack = 0;
1249 SimplexId wrapZFront = 0;
1250 if(p[0] == 0)
1251 wrapXLeft = 2 * wrap_[0];
1252 if(p[0] == nbvoxels_[0])
1253 wrapXRight = -2 * wrap_[0];
1254 if(p[1] == 0)
1255 wrapYTop = 2 * wrap_[1];
1256 if(p[1] == nbvoxels_[1])
1257 wrapYBottom = -2 * wrap_[1];
1258 if(p[2] == 0)
1259 wrapZBack = 2 * wrap_[2];
1260 if(p[2] == nbvoxels_[2])
1261 wrapZFront = -2 * wrap_[2];
1262 switch(id) {
1263 case 0:
1264 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10] + p[2] * tshift_[11];
1265 case 1:
1266 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6] + p[2] * tshift_[7]
1267 + 1;
1268 case 2:
1269 return tsetshift_[1] + (p[0] - 1) * 2 + p[1] * tshift_[4]
1270 + p[2] * tshift_[5] + wrapXLeft;
1271 case 3:
1272 return tsetshift_[1] + (p[0] - 1) * 2 + p[1] * tshift_[4]
1273 + p[2] * tshift_[5] + 1 + wrapXLeft;
1274 case 4:
1275 return tsetshift_[0] + (p[0] - 1) * 2 + (p[1] + 1) * tshift_[2]
1276 + p[2] * tshift_[3] + wrapXLeft + wrapYBottom;
1277 case 5:
1278 return tsetshift_[0] + (p[0] - 1) * 2 + (p[1] + 1) * tshift_[2]
1279 + p[2] * tshift_[3] + 1 + wrapXLeft + wrapYBottom;
1280 case 6:
1281 return (p[0] - 1) * 2 + p[1] * tshift_[0] + (p[2] + 1) * tshift_[1]
1282 + wrapXLeft + wrapZFront;
1283 case 7:
1284 return (p[0] - 1) * 2 + p[1] * tshift_[0] + (p[2] + 1) * tshift_[1] + 1
1285 + wrapXLeft + wrapZFront;
1286 case 8:
1287 return tsetshift_[3] + p[0] * 2 + (p[1] - 1) * tshift_[8]
1288 + p[2] * tshift_[9] + wrapYTop;
1289 case 9:
1290 return tsetshift_[2] + p[0] * 2 + (p[1] - 1) * tshift_[6]
1291 + p[2] * tshift_[7] + wrapYTop;
1292 case 10:
1293 return tsetshift_[4] + (p[0] - 1) * 2 + (p[1] - 1) * tshift_[10]
1294 + p[2] * tshift_[11] + wrapXLeft + wrapYTop;
1295 case 11:
1296 return tsetshift_[3] + (p[0] - 1) * 2 + (p[1] - 1) * tshift_[8]
1297 + p[2] * tshift_[9] + 1 + wrapXLeft + wrapYTop;
1298 case 12:
1299 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8]
1300 + (p[2] - 1) * tshift_[9] + wrapZBack;
1301 case 13:
1302 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10]
1303 + (p[2] - 1) * tshift_[11] + 1 + wrapZBack;
1304 case 14:
1305 return tsetshift_[2] + (p[0] - 1) * 2 + p[1] * tshift_[6]
1306 + (p[2] - 1) * tshift_[7] + 1 + wrapXLeft + wrapZBack;
1307 case 15:
1308 return tsetshift_[3] + (p[0] - 1) * 2 + p[1] * tshift_[8]
1309 + (p[2] - 1) * tshift_[9] + 1 + wrapXLeft + wrapZBack;
1310 case 16:
1311 return tsetshift_[1] + (p[0] + 1) * 2 + (p[1] - 1) * tshift_[4]
1312 + (p[2] - 1) * tshift_[5] + wrapXRight + wrapYTop + wrapZBack;
1313 case 17:
1314 return tsetshift_[1] + (p[0] + 1) * 2 + (p[1] - 1) * tshift_[4]
1315 + (p[2] - 1) * tshift_[5] + 1 + wrapXRight + wrapYTop + wrapZBack;
1316 case 18:
1317 return tsetshift_[0] + p[0] * 2 + (p[1] - 1) * tshift_[2]
1318 + (p[2] - 1) * tshift_[3] + wrapYTop + wrapZBack;
1319 case 19:
1320 return tsetshift_[0] + p[0] * 2 + (p[1] - 1) * tshift_[2]
1321 + (p[2] - 1) * tshift_[3] + 1 + wrapYTop + wrapZBack;
1322 case 20:
1323 return p[0] * 2 + (p[1] - 1) * tshift_[0] + (p[2] - 1) * tshift_[1]
1324 + wrapYTop + wrapZBack;
1325 case 21:
1326 return p[0] * 2 + (p[1] - 1) * tshift_[0] + (p[2] - 1) * tshift_[1] + 1
1327 + wrapYTop + wrapZBack;
1328 case 22:
1329 return tsetshift_[4] + (p[0] - 1) * 2 + (p[1] - 1) * tshift_[10]
1330 + (p[2] - 1) * tshift_[11] + 1 + wrapXLeft + wrapYTop + wrapZBack;
1331 case 23:
1332 return tsetshift_[2] + (p[0] - 1) * 2 + (p[1] - 1) * tshift_[6]
1333 + (p[2] - 1) * tshift_[7] + wrapXLeft + wrapYTop + wrapZBack;
1334 }
1335 return -1;
1336}
1337
1338inline ttk::SimplexId
1339 ttk::PeriodicImplicitTriangulation::getVertexTriangle3d(const SimplexId p[3],
1340 const int id) const {
1341 SimplexId wrapXLeft = 0;
1342 SimplexId wrapYTop = 0;
1343 SimplexId wrapZBack = 0;
1344 if(p[0] == 0)
1345 wrapXLeft = 2 * wrap_[0];
1346 if(p[1] == 0)
1347 wrapYTop = 2 * wrap_[1];
1348 if(p[2] == 0)
1349 wrapZBack = 2 * wrap_[2];
1350 switch(id) {
1351 case 0:
1352 return (p[0] - 1) * 2 + p[1] * tshift_[0] + p[2] * tshift_[1] + 1
1353 + wrapXLeft;
1354 case 1:
1355 return tsetshift_[4] + (p[0] - 1) * 2 + p[1] * tshift_[10]
1356 + p[2] * tshift_[11] + wrapXLeft;
1357 case 2:
1358 return tsetshift_[2] + (p[0] - 1) * 2 + p[1] * tshift_[6]
1359 + p[2] * tshift_[7] + wrapXLeft;
1360 case 3:
1361 return tsetshift_[3] + (p[0] - 1) * 2 + p[1] * tshift_[8]
1362 + p[2] * tshift_[9] + 1 + wrapXLeft;
1363 case 4:
1364 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4] + p[2] * tshift_[5];
1365 case 5:
1366 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4] + p[2] * tshift_[5]
1367 + 1;
1368 case 6:
1369 return tsetshift_[4] + (p[0] - 1) * 2 + p[1] * tshift_[10]
1370 + p[2] * tshift_[11] + 1 + wrapXLeft;
1371 case 7:
1372 return tsetshift_[0] + (p[0] - 1) * 2 + p[1] * tshift_[2]
1373 + p[2] * tshift_[3] + 1 + wrapXLeft;
1374 case 8:
1375 return tsetshift_[2] + (p[0] - 1) * 2 + p[1] * tshift_[6]
1376 + p[2] * tshift_[7] + 1 + wrapXLeft;
1377 case 9:
1378 return tsetshift_[3] + (p[0] - 1) * 2 + p[1] * tshift_[8]
1379 + p[2] * tshift_[9] + wrapXLeft;
1380 case 10:
1381 return tsetshift_[0] + (p[0] - 1) * 2 + p[1] * tshift_[2]
1382 + p[2] * tshift_[3] + wrapXLeft;
1383 case 11:
1384 return (p[0] - 1) * 2 + p[1] * tshift_[0] + p[2] * tshift_[1] + wrapXLeft;
1385 case 12:
1386 return p[0] * 2 + tsetshift_[2] + (p[1] - 1) * tshift_[6]
1387 + (p[2] - 1) * tshift_[7] + wrapYTop + wrapZBack;
1388 case 13:
1389 return p[0] * 2 + tsetshift_[2] + (p[1] - 1) * tshift_[6]
1390 + (p[2] - 1) * tshift_[7] + 1 + wrapYTop + wrapZBack;
1391 case 14:
1392 return p[0] * 2 + tsetshift_[1] + (p[1] - 1) * tshift_[4]
1393 + (p[2] - 1) * tshift_[5] + wrapYTop + wrapZBack;
1394 case 15:
1395 return p[0] * 2 + tsetshift_[1] + (p[1] - 1) * tshift_[4]
1396 + (p[2] - 1) * tshift_[5] + 1 + wrapYTop + wrapZBack;
1397 case 16:
1398 return p[0] * 2 + tsetshift_[3] + (p[1] - 1) * tshift_[8]
1399 + (p[2] - 1) * tshift_[9] + wrapYTop + wrapZBack;
1400 case 17:
1401 return p[0] * 2 + tsetshift_[3] + (p[1] - 1) * tshift_[8]
1402 + (p[2] - 1) * tshift_[9] + 1 + wrapYTop + wrapZBack;
1403 case 18:
1404 return p[0] * 2 + tsetshift_[4] + (p[1] - 1) * tshift_[10]
1405 + (p[2] - 1) * tshift_[11] + wrapYTop + wrapZBack;
1406 case 19:
1407 return p[0] * 2 + tsetshift_[4] + (p[1] - 1) * tshift_[10]
1408 + (p[2] - 1) * tshift_[11] + 1 + wrapYTop + wrapZBack;
1409 case 20:
1410 return p[0] * 2 + (p[1] - 1) * tshift_[0] + p[2] * tshift_[1] + wrapYTop;
1411 case 21:
1412 return p[0] * 2 + (p[1] - 1) * tshift_[0] + p[2] * tshift_[1] + 1
1413 + wrapYTop;
1414 case 22:
1415 return p[0] * 2 + tsetshift_[0] + p[1] * tshift_[2]
1416 + (p[2] - 1) * tshift_[3] + wrapZBack;
1417 case 23:
1418 return p[0] * 2 + tsetshift_[0] + p[1] * tshift_[2]
1419 + (p[2] - 1) * tshift_[3] + 1 + wrapZBack;
1420 case 24:
1421 return (p[0] - 1) * 2 + (p[1] - 1) * tshift_[0] + p[2] * tshift_[1] + 1
1422 + wrapXLeft + wrapYTop;
1423 case 25:
1424 return tsetshift_[2] + (p[0] - 1) * 2 + (p[1] - 1) * tshift_[6]
1425 + p[2] * tshift_[7] + wrapXLeft + wrapYTop;
1426 case 26:
1427 return tsetshift_[1] + p[0] * 2 + (p[1] - 1) * tshift_[4]
1428 + p[2] * tshift_[5] + wrapYTop;
1429 case 27:
1430 return p[0] * 2 + tsetshift_[2] + p[1] * tshift_[6]
1431 + (p[2] - 1) * tshift_[7] + 1 + wrapZBack;
1432 case 28:
1433 return p[0] * 2 + tsetshift_[1] + p[1] * tshift_[4]
1434 + (p[2] - 1) * tshift_[5] + 1 + wrapZBack;
1435 case 29:
1436 return p[0] * 2 + p[1] * tshift_[0] + p[2] * tshift_[1];
1437 case 30:
1438 return tsetshift_[3] + (p[0] - 1) * 2 + (p[1] - 1) * tshift_[8]
1439 + (p[2] - 1) * tshift_[9] + 1 + wrapXLeft + wrapYTop + wrapZBack;
1440 case 31:
1441 return tsetshift_[0] + (p[0] - 1) * 2 + p[1] * tshift_[2]
1442 + (p[2] - 1) * tshift_[3] + 1 + wrapXLeft + wrapZBack;
1443 case 32:
1444 return tsetshift_[0] + p[0] * 2 + p[1] * tshift_[2] + p[2] * tshift_[3];
1445 case 33:
1446 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9];
1447 case 34:
1448 return tsetshift_[4] + (p[0] - 1) * 2 + p[1] * tshift_[10]
1449 + (p[2] - 1) * tshift_[11] + 1 + wrapXLeft + wrapZBack;
1450 case 35:
1451 return p[0] * 2 + tsetshift_[4] + (p[1] - 1) * tshift_[10]
1452 + p[2] * tshift_[11] + wrapYTop;
1453 }
1454 return -1;
1455}
1456
1457inline ttk::SimplexId
1458 ttk::PeriodicImplicitTriangulation::getVertexStar3d(const SimplexId p[3],
1459 const int id) const {
1460 SimplexId wrapXLeft = 0;
1461 SimplexId wrapYTop = 0;
1462 SimplexId wrapZBack = 0;
1463 if(p[0] == 0)
1464 wrapXLeft = 6 * wrap_[0];
1465 if(p[1] == 0)
1466 wrapYTop = 6 * wrap_[1];
1467 if(p[2] == 0)
1468 wrapZBack = 6 * wrap_[2];
1469 switch(id) {
1470 case 0:
1471 return (p[0] - 1) * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1]
1472 + wrapXLeft;
1473 case 1:
1474 return (p[0] - 1) * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 1
1475 + wrapXLeft;
1476 case 2:
1477 return (p[0] - 1) * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 2
1478 + wrapXLeft;
1479 case 3:
1480 return (p[0] - 1) * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 3
1481 + wrapXLeft;
1482 case 4:
1483 return (p[0] - 1) * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 4
1484 + wrapXLeft;
1485 case 5:
1486 return (p[0] - 1) * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 5
1487 + wrapXLeft;
1488 case 6:
1489 return p[0] * 6 + (p[1] - 1) * tetshift_[0] + (p[2] - 1) * tetshift_[1]
1490 + wrapYTop + wrapZBack;
1491 case 7:
1492 return p[0] * 6 + (p[1] - 1) * tetshift_[0] + (p[2] - 1) * tetshift_[1]
1493 + 1 + wrapYTop + wrapZBack;
1494 case 8:
1495 return p[0] * 6 + (p[1] - 1) * tetshift_[0] + (p[2] - 1) * tetshift_[1]
1496 + 2 + wrapYTop + wrapZBack;
1497 case 9:
1498 return p[0] * 6 + (p[1] - 1) * tetshift_[0] + (p[2] - 1) * tetshift_[1]
1499 + 3 + wrapYTop + wrapZBack;
1500 case 10:
1501 return p[0] * 6 + (p[1] - 1) * tetshift_[0] + (p[2] - 1) * tetshift_[1]
1502 + 4 + wrapYTop + wrapZBack;
1503 case 11:
1504 return p[0] * 6 + (p[1] - 1) * tetshift_[0] + (p[2] - 1) * tetshift_[1]
1505 + 5 + wrapYTop + wrapZBack;
1506 case 12:
1507 return p[0] * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1];
1508 case 13:
1509 return p[0] * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 2;
1510 case 14:
1511 return p[0] * 6 + (p[1] - 1) * tetshift_[0] + p[2] * tetshift_[1]
1512 + wrapYTop;
1513 case 15:
1514 return p[0] * 6 + (p[1] - 1) * tetshift_[0] + p[2] * tetshift_[1] + 1
1515 + wrapYTop;
1516 case 16:
1517 return (p[0] - 1) * 6 + (p[1] - 1) * tetshift_[0] + p[2] * tetshift_[1]
1518 + 1 + wrapXLeft + wrapYTop;
1519 case 17:
1520 return (p[0] - 1) * 6 + (p[1] - 1) * tetshift_[0] + p[2] * tetshift_[1]
1521 + 5 + wrapXLeft + wrapYTop;
1522 case 18:
1523 return p[0] * 6 + p[1] * tetshift_[0] + (p[2] - 1) * tetshift_[1] + 2
1524 + wrapZBack;
1525 case 19:
1526 return p[0] * 6 + p[1] * tetshift_[0] + (p[2] - 1) * tetshift_[1] + 3
1527 + wrapZBack;
1528 case 20:
1529 return (p[0] - 1) * 6 + p[1] * tetshift_[0] + (p[2] - 1) * tetshift_[1]
1530 + 3 + wrapXLeft + wrapZBack;
1531 case 21:
1532 return (p[0] - 1) * 6 + p[1] * tetshift_[0] + (p[2] - 1) * tetshift_[1]
1533 + 4 + wrapXLeft + wrapZBack;
1534 case 22:
1535 return (p[0] - 1) * 6 + (p[1] - 1) * tetshift_[0]
1536 + (p[2] - 1) * tetshift_[1] + 4 + wrapXLeft + wrapYTop + wrapZBack;
1537 case 23:
1538 return (p[0] - 1) * 6 + (p[1] - 1) * tetshift_[0]
1539 + (p[2] - 1) * tetshift_[1] + 5 + wrapXLeft + wrapYTop + wrapZBack;
1540 }
1541 return -1;
1542}
1543
1544inline ttk::SimplexId
1545 ttk::PeriodicImplicitTriangulation::getEdgeTriangle3dL(const SimplexId p[3],
1546 const int id) const {
1547 SimplexId wrapYTop = 0;
1548 SimplexId wrapZBack = 0;
1549 if(p[1] == 0)
1550 wrapYTop = 2 * wrap_[1];
1551 if(p[2] == 0)
1552 wrapZBack = 2 * wrap_[2];
1553 switch(id) {
1554 case 0:
1555 return p[0] * 2 + (p[1] - 1) * tshift_[0] + p[2] * tshift_[1] + 1
1556 + wrapYTop;
1557 case 1:
1558 return tsetshift_[0] + p[0] * 2 + p[1] * tshift_[2] + p[2] * tshift_[3];
1559 case 2:
1560 return tsetshift_[0] + p[0] * 2 + p[1] * tshift_[2]
1561 + (p[2] - 1) * tshift_[3] + 1 + wrapZBack;
1562 case 3:
1563 return tsetshift_[3] + p[0] * 2 + (p[1] - 1) * tshift_[8]
1564 + (p[2] - 1) * tshift_[9] + 1 + wrapYTop + wrapZBack;
1565 case 4:
1566 return p[0] * 2 + p[1] * tshift_[0] + p[2] * tshift_[1];
1567 case 5:
1568 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9];
1569 }
1570 return -1;
1571}
1572
1573inline ttk::SimplexId
1574 ttk::PeriodicImplicitTriangulation::getEdgeTriangle3dH(const SimplexId p[3],
1575 const int id) const {
1576 SimplexId wrapXLeft = 0;
1577 SimplexId wrapZBack = 0;
1578 if(p[0] == 0)
1579 wrapXLeft = 2 * wrap_[0];
1580 if(p[2] == 0)
1581 wrapZBack = 2 * wrap_[2];
1582 switch(id) {
1583 case 0:
1584 return (p[0] - 1) * 2 + p[1] * tshift_[0] + p[2] * tshift_[1] + 1
1585 + wrapXLeft;
1586 case 1:
1587 return tsetshift_[2] + (p[0] - 1) * 2 + p[1] * tshift_[6]
1588 + p[2] * tshift_[7] + wrapXLeft;
1589 case 2:
1590 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4] + p[2] * tshift_[5];
1591 case 3:
1592 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4]
1593 + (p[2] - 1) * tshift_[5] + 1 + wrapZBack;
1594 case 4:
1595 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6]
1596 + (p[2] - 1) * tshift_[7] + 1 + wrapZBack;
1597 case 5:
1598 return p[0] * 2 + p[1] * tshift_[0] + p[2] * tshift_[1];
1599 }
1600 return -1;
1601}
1602
1603inline ttk::SimplexId
1604 ttk::PeriodicImplicitTriangulation::getEdgeTriangle3dP(const SimplexId p[3],
1605 const int id) const {
1606 SimplexId wrapXLeft = 0;
1607 SimplexId wrapYTop = 0;
1608 if(p[0] == 0)
1609 wrapXLeft = 2 * wrap_[0];
1610 if(p[1] == 0)
1611 wrapYTop = 2 * wrap_[1];
1612 switch(id) {
1613 case 0:
1614 return tsetshift_[0] + (p[0] - 1) * 2 + p[1] * tshift_[2]
1615 + p[2] * tshift_[3] + 1 + wrapXLeft;
1616 case 1:
1617 return tsetshift_[4] + (p[0] - 1) * 2 + p[1] * tshift_[10]
1618 + p[2] * tshift_[11] + 1 + wrapXLeft;
1619 case 2:
1620 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4] + p[2] * tshift_[5]
1621 + 1;
1622 case 3:
1623 return tsetshift_[0] + p[0] * 2 + p[1] * tshift_[2] + p[2] * tshift_[3];
1624 case 4:
1625 return tsetshift_[4] + p[0] * 2 + (p[1] - 1) * tshift_[10]
1626 + p[2] * tshift_[11] + wrapYTop;
1627 case 5:
1628 return tsetshift_[1] + p[0] * 2 + (p[1] - 1) * tshift_[4]
1629 + p[2] * tshift_[5] + wrapYTop;
1630 }
1631 return -1;
1632}
1633
1634inline ttk::SimplexId
1635 ttk::PeriodicImplicitTriangulation::getEdgeTriangle3dD1(const SimplexId p[3],
1636 const int id) const {
1637 SimplexId wrapZBack = 0;
1638 if(p[2] == 0)
1639 wrapZBack = 2 * wrap_[2];
1640 switch(id) {
1641 case 0:
1642 return p[0] * 2 + p[1] * tshift_[0] + p[2] * tshift_[1];
1643 case 1:
1644 return p[0] * 2 + p[1] * tshift_[0] + p[2] * tshift_[1] + 1;
1645 case 2:
1646 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10] + p[2] * tshift_[11];
1647 case 3:
1648 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10]
1649 + (p[2] - 1) * tshift_[11] + 1 + wrapZBack;
1650 }
1651 return -1;
1652}
1653
1654inline ttk::SimplexId
1655 ttk::PeriodicImplicitTriangulation::getEdgeTriangle3dD2(const SimplexId p[3],
1656 const int id) const {
1657 SimplexId wrapXLeft = 0;
1658 if(p[0] == 0)
1659 wrapXLeft = 2 * wrap_[0];
1660 switch(id) {
1661 case 0:
1662 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4] + p[2] * tshift_[5];
1663 case 1:
1664 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4] + p[2] * tshift_[5]
1665 + 1;
1666 case 2:
1667 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9];
1668 case 3:
1669 return tsetshift_[3] + (p[0] - 1) * 2 + p[1] * tshift_[8]
1670 + p[2] * tshift_[9] + 1 + wrapXLeft;
1671 }
1672 return -1;
1673}
1674
1675inline ttk::SimplexId
1676 ttk::PeriodicImplicitTriangulation::getEdgeTriangle3dD3(const SimplexId p[3],
1677 const int id) const {
1678 SimplexId wrapYTop = 0;
1679 if(p[1] == 0)
1680 wrapYTop = 2 * wrap_[1];
1681 switch(id) {
1682 case 0:
1683 return tsetshift_[0] + p[0] * 2 + p[1] * tshift_[2] + p[2] * tshift_[3];
1684 case 1:
1685 return tsetshift_[0] + p[0] * 2 + p[1] * tshift_[2] + p[2] * tshift_[3]
1686 + 1;
1687 case 2:
1688 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6] + p[2] * tshift_[7]
1689 + 1;
1690 case 3:
1691 return tsetshift_[2] + p[0] * 2 + (p[1] - 1) * tshift_[6]
1692 + p[2] * tshift_[7] + wrapYTop;
1693 }
1694 return -1;
1695}
1696
1697inline ttk::SimplexId
1698 ttk::PeriodicImplicitTriangulation::getEdgeTriangle3dD4(const SimplexId p[3],
1699 const int id) const {
1700 switch(id) {
1701 case 0:
1702 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6] + p[2] * tshift_[7];
1703 case 1:
1704 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6] + p[2] * tshift_[7]
1705 + 1;
1706 case 2:
1707 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9];
1708 case 3:
1709 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9]
1710 + 1;
1711 case 4:
1712 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10] + p[2] * tshift_[11];
1713 case 5:
1714 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10] + p[2] * tshift_[11]
1715 + 1;
1716 }
1717 return -1;
1718}
1719
1720inline ttk::SimplexId
1721 ttk::PeriodicImplicitTriangulation::getEdgeLinkL(const SimplexId p[3],
1722 const int id) const {
1723 SimplexId wrapXRight = 0;
1724 SimplexId wrapYTop = 0;
1725 SimplexId wrapYBottom = 0;
1726 SimplexId wrapZBack = 0;
1727 SimplexId wrapZFront = 0;
1728 if(p[0] == nbvoxels_[0])
1729 wrapXRight = -wrap_[0];
1730 if(p[1] == 0)
1731 wrapYTop = wrap_[1];
1732 if(p[1] == nbvoxels_[1])
1733 wrapYBottom = -wrap_[1];
1734 if(p[2] == 0)
1735 wrapZBack = wrap_[2];
1736 if(p[2] == nbvoxels_[2])
1737 wrapZFront = -wrap_[2];
1738 switch(id) {
1739 case 0:
1740 return esetshift_[1] + p[0] + (p[1] + 1) * eshift_[4] + p[2] * eshift_[5]
1741 + wrapYBottom;
1742 case 1:
1743 return esetshift_[0] + p[0] + p[1] * eshift_[2] + (p[2] + 1) * eshift_[3]
1744 + wrapZFront;
1745 case 2:
1746 return esetshift_[5] + p[0] + (p[1] - 1) * eshift_[12]
1747 + p[2] * eshift_[13] + wrapYTop;
1748 case 3:
1749 return esetshift_[1] + p[0] + 1 + (p[1] - 1) * eshift_[4]
1750 + (p[2] - 1) * eshift_[5] + wrapXRight + wrapYTop + wrapZBack;
1751 case 4:
1752 return esetshift_[0] + p[0] + 1 + (p[1] - 1) * eshift_[2]
1753 + (p[2] - 1) * eshift_[3] + wrapXRight + wrapYTop + wrapZBack;
1754 case 5:
1755 return esetshift_[5] + p[0] + p[1] * eshift_[12]
1756 + (p[2] - 1) * eshift_[13] + wrapZBack;
1757 }
1758 return -1;
1759}
1760
1761inline ttk::SimplexId
1762 ttk::PeriodicImplicitTriangulation::getEdgeLinkH(const SimplexId p[3],
1763 const int id) const {
1764 SimplexId wrapXLeft = 0;
1765 SimplexId wrapXRight = 0;
1766 SimplexId wrapYBottom = 0;
1767 SimplexId wrapZBack = 0;
1768 SimplexId wrapZFront = 0;
1769 if(p[0] == 0)
1770 wrapXLeft = wrap_[0];
1771 if(p[0] == nbvoxels_[0])
1772 wrapXRight = -wrap_[0];
1773 if(p[1] == nbvoxels_[1])
1774 wrapYBottom = -wrap_[1];
1775 if(p[2] == 0)
1776 wrapZBack = wrap_[2];
1777 if(p[2] == nbvoxels_[2])
1778 wrapZFront = -wrap_[2];
1779 switch(id) {
1780 case 0:
1781 return p[0] + p[1] * eshift_[0] + (p[2] - 1) * eshift_[1] + wrapZBack;
1782 case 1:
1783 return p[0] + (p[1] + 1) * eshift_[0] + (p[2] + 1) * eshift_[1] - 1
1784 + wrapXLeft + wrapYBottom + wrapZFront;
1785 case 2:
1786 return esetshift_[1] + p[0] - 1 + (p[1] + 1) * eshift_[4]
1787 + p[2] * eshift_[5] + wrapXLeft + wrapYBottom;
1788 case 3:
1789 return esetshift_[1] + p[0] + 1 + p[1] * eshift_[4]
1790 + (p[2] - 1) * eshift_[5] + wrapXRight + wrapZBack;
1791 case 4:
1792 return esetshift_[5] + (p[0] - 1) + p[1] * eshift_[12]
1793 + (p[2] - 1) * eshift_[13] + wrapXLeft + wrapZBack;
1794 case 5:
1795 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13];
1796 }
1797 return -1;
1798}
1799
1800inline ttk::SimplexId
1801 ttk::PeriodicImplicitTriangulation::getEdgeLinkP(const SimplexId p[3],
1802 const int id) const {
1803 SimplexId wrapXLeft = 0;
1804 SimplexId wrapXRight = 0;
1805 SimplexId wrapYTop = 0;
1806 SimplexId wrapYBottom = 0;
1807 SimplexId wrapZFront = 0;
1808 if(p[0] == 0)
1809 wrapXLeft = wrap_[0];
1810 if(p[0] == nbvoxels_[0])
1811 wrapXRight = -wrap_[0];
1812 if(p[1] == 0)
1813 wrapYTop = wrap_[1];
1814 if(p[1] == nbvoxels_[1])
1815 wrapYBottom = -wrap_[1];
1816 if(p[2] == nbvoxels_[2])
1817 wrapZFront = -wrap_[2];
1818 switch(id) {
1819 case 0:
1820 return p[0] + (p[1] - 1) * eshift_[0] + p[2] * eshift_[1] + wrapYTop;
1821 case 1:
1822 return p[0] + (p[1] + 1) * eshift_[0] + (p[2] + 1) * eshift_[1] - 1
1823 + wrapXLeft + wrapYBottom + wrapZFront;
1824 case 2:
1825 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13];
1826 case 3:
1827 return esetshift_[5] + p[0] + (p[1] - 1) * eshift_[12]
1828 + p[2] * eshift_[13] - 1 + wrapXLeft + wrapYTop;
1829 case 4:
1830 return esetshift_[0] + p[0] + p[1] * eshift_[2] + (p[2] + 1) * eshift_[3]
1831 - 1 + wrapXLeft + wrapZFront;
1832 case 5:
1833 return esetshift_[0] + p[0] + (p[1] - 1) * eshift_[2] + p[2] * eshift_[3]
1834 + 1 + wrapXRight + wrapYTop;
1835 }
1836 return -1;
1837}
1838
1839inline ttk::SimplexId
1840 ttk::PeriodicImplicitTriangulation::getEdgeLinkD1(const SimplexId p[3],
1841 const int id) const {
1842 SimplexId wrapXRight = 0;
1843 SimplexId wrapYBottom = 0;
1844 SimplexId wrapZBack = 0;
1845 if(p[0] == nbvoxels_[0])
1846 wrapXRight = -wrap_[0];
1847 if(p[1] == nbvoxels_[1])
1848 wrapYBottom = -wrap_[1];
1849 if(p[2] == 0)
1850 wrapZBack = wrap_[2];
1851 switch(id) {
1852 case 0:
1853 return esetshift_[4] + p[0] + p[1] * eshift_[10]
1854 + (p[2] - 1) * eshift_[11] + wrapZBack;
1855 case 1:
1856 return esetshift_[4] + p[0] + (p[1] + 1) * eshift_[10]
1857 + p[2] * eshift_[11] + wrapYBottom;
1858 case 2:
1859 return esetshift_[3] + p[0] + p[1] * eshift_[8] + p[2] * eshift_[9];
1860 case 3:
1861 return esetshift_[3] + p[0] + p[1] * eshift_[8] + (p[2] - 1) * eshift_[9]
1862 + 1 + wrapXRight + wrapZBack;
1863 }
1864 return -1;
1865}
1866
1867inline ttk::SimplexId
1868 ttk::PeriodicImplicitTriangulation::getEdgeLinkD2(const SimplexId p[3],
1869 const int id) const {
1870 SimplexId wrapXLeft = 0;
1871 SimplexId wrapYBottom = 0;
1872 SimplexId wrapZFront = 0;
1873 if(p[0] == 0)
1874 wrapXLeft = wrap_[0];
1875 if(p[1] == nbvoxels_[1])
1876 wrapYBottom = -wrap_[1];
1877 if(p[2] == nbvoxels_[2])
1878 wrapZFront = -wrap_[2];
1879 switch(id) {
1880 case 0:
1881 return esetshift_[4] + p[0] + p[1] * eshift_[10] + p[2] * eshift_[11];
1882 case 1:
1883 return esetshift_[4] + p[0] + (p[1] + 1) * eshift_[10]
1884 + p[2] * eshift_[11] - 1 + wrapXLeft + wrapYBottom;
1885 case 2:
1886 return esetshift_[2] + p[0] + p[1] * eshift_[6] + p[2] * eshift_[7];
1887 case 3:
1888 return esetshift_[2] + p[0] + p[1] * eshift_[6] + (p[2] + 1) * eshift_[7]
1889 - 1 + wrapXLeft + wrapZFront;
1890 }
1891 return -1;
1892}
1893
1894inline ttk::SimplexId
1895 ttk::PeriodicImplicitTriangulation::getEdgeLinkD3(const SimplexId p[3],
1896 const int id) const {
1897 SimplexId wrapXRight = 0;
1898 SimplexId wrapYTop = 0;
1899 SimplexId wrapZFront = 0;
1900 if(p[0] == nbvoxels_[0])
1901 wrapXRight = -wrap_[0];
1902 if(p[1] == 0)
1903 wrapYTop = wrap_[1];
1904 if(p[2] == nbvoxels_[2])
1905 wrapZFront = -wrap_[2];
1906 switch(id) {
1907 case 0:
1908 return esetshift_[2] + p[0] + (p[1] - 1) * eshift_[6] + p[2] * eshift_[7]
1909 + wrapYTop;
1910 case 1:
1911 return esetshift_[2] + p[0] + p[1] * eshift_[6] + (p[2] + 1) * eshift_[7]
1912 + wrapZFront;
1913 case 2:
1914 return esetshift_[3] + p[0] + p[1] * eshift_[8] + p[2] * eshift_[9];
1915 case 3:
1916 return esetshift_[3] + p[0] + (p[1] - 1) * eshift_[8] + p[2] * eshift_[9]
1917 + 1 + wrapXRight + wrapYTop;
1918 }
1919 return -1;
1920}
1921
1922inline ttk::SimplexId
1923 ttk::PeriodicImplicitTriangulation::getEdgeLinkD4(const SimplexId p[3],
1924 const int id) const {
1925 SimplexId wrapXRight = 0;
1926 SimplexId wrapYBottom = 0;
1927 SimplexId wrapZFront = 0;
1928 if(p[0] == nbvoxels_[0])
1929 wrapXRight = -wrap_[0];
1930 if(p[1] == nbvoxels_[1])
1931 wrapYBottom = -wrap_[1];
1932 if(p[2] == nbvoxels_[2])
1933 wrapZFront = -wrap_[2];
1934 switch(id) {
1935 case 0:
1936 return p[0] + (p[1] + 1) * eshift_[0] + p[2] * eshift_[1] + wrapYBottom;
1937 case 1:
1938 return p[0] + p[1] * eshift_[0] + (p[2] + 1) * eshift_[1] + wrapZFront;
1939 case 2:
1940 return esetshift_[1] + p[0] + p[1] * eshift_[4] + p[2] * eshift_[5];
1941 case 3:
1942 return esetshift_[1] + p[0] + 1 + (p[1] + 1) * eshift_[4]
1943 + p[2] * eshift_[5] + wrapXRight + wrapYBottom;
1944 case 4:
1945 return esetshift_[0] + p[0] + p[1] * eshift_[2] + p[2] * eshift_[3];
1946 case 5:
1947 return esetshift_[0] + p[0] + 1 + p[1] * eshift_[2]
1948 + (p[2] + 1) * eshift_[3] + wrapXRight + wrapZFront;
1949 }
1950 return -1;
1951}
1952
1953inline ttk::SimplexId
1954 ttk::PeriodicImplicitTriangulation::getEdgeStarL(const SimplexId p[3],
1955 const int id) const {
1956 SimplexId wrapYTop = 0;
1957 SimplexId wrapZBack = 0;
1958 if(p[1] == 0)
1959 wrapYTop = 6 * wrap_[1];
1960 if(p[2] == 0)
1961 wrapZBack = 6 * wrap_[2];
1962 switch(id) {
1963 case 0:
1964 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6;
1965 case 1:
1966 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 2;
1967 case 2:
1968 return p[2] * tetshift_[1] + (p[1] - 1) * tetshift_[0] + p[0] * 6 + 1
1969 + wrapYTop;
1970 case 3:
1971 return (p[2] - 1) * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 3
1972 + wrapZBack;
1973 case 4:
1974 return (p[2] - 1) * tetshift_[1] + (p[1] - 1) * tetshift_[0] + p[0] * 6
1975 + 4 + wrapYTop + wrapZBack;
1976 case 5:
1977 return (p[2] - 1) * tetshift_[1] + (p[1] - 1) * tetshift_[0] + p[0] * 6
1978 + 5 + wrapYTop + wrapZBack;
1979 }
1980 return -1;
1981}
1982
1983inline ttk::SimplexId
1984 ttk::PeriodicImplicitTriangulation::getEdgeStarH(const SimplexId p[3],
1985 const int id) const {
1986 SimplexId wrapXLeft = 0;
1987 SimplexId wrapZBack = 0;
1988 if(p[0] == 0)
1989 wrapXLeft = 6 * wrap_[0];
1990 if(p[2] == 0)
1991 wrapZBack = 6 * wrap_[2];
1992 switch(id) {
1993 case 0:
1994 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6;
1995 case 1:
1996 return (p[2] - 1) * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 2
1997 + wrapZBack;
1998 case 2:
1999 return (p[2] - 1) * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 3
2000 + wrapZBack;
2001 case 3:
2002 return (p[2] - 1) * tetshift_[1] + p[1] * tetshift_[0] + (p[0] - 1) * 6
2003 + 4 + wrapXLeft + wrapZBack;
2004 case 4:
2005 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + (p[0] - 1) * 6 + 1
2006 + wrapXLeft;
2007 case 5:
2008 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + (p[0] - 1) * 6 + 5
2009 + wrapXLeft;
2010 }
2011 return -1;
2012}
2013
2014inline ttk::SimplexId
2015 ttk::PeriodicImplicitTriangulation::getEdgeStarP(const SimplexId p[3],
2016 const int id) const {
2017 SimplexId wrapXLeft = 0;
2018 SimplexId wrapYTop = 0;
2019 if(p[0] == 0)
2020 wrapXLeft = 6 * wrap_[0];
2021 if(p[1] == 0)
2022 wrapYTop = 6 * wrap_[1];
2023 switch(id) {
2024 case 0:
2025 return p[2] * tetshift_[1] + (p[1] - 1) * tetshift_[0] + (p[0] - 1) * 6
2026 + 5 + wrapXLeft + wrapYTop;
2027 case 1:
2028 return p[2] * tetshift_[1] + (p[1] - 1) * tetshift_[0] + p[0] * 6
2029 + wrapYTop;
2030 case 2:
2031 return p[2] * tetshift_[1] + (p[1] - 1) * tetshift_[0] + p[0] * 6 + 1
2032 + wrapYTop;
2033 case 3:
2034 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 2;
2035 case 4:
2036 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + (p[0] - 1) * 6 + 3
2037 + wrapXLeft;
2038 case 5:
2039 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + (p[0] - 1) * 6 + 4
2040 + wrapXLeft;
2041 }
2042 return -1;
2043}
2044
2045inline ttk::SimplexId
2046 ttk::PeriodicImplicitTriangulation::getEdgeStarD1(const SimplexId p[3],
2047 const int id) const {
2048 SimplexId wrapZBack = 0;
2049 if(p[2] == 0)
2050 wrapZBack = 6 * wrap_[2];
2051 switch(id) {
2052 case 0:
2053 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6;
2054 case 1:
2055 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 1;
2056 case 2:
2057 return (p[2] - 1) * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 3
2058 + wrapZBack;
2059 case 3:
2060 return (p[2] - 1) * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 4
2061 + wrapZBack;
2062 }
2063 return -1;
2064}
2065
2066inline ttk::SimplexId
2067 ttk::PeriodicImplicitTriangulation::getEdgeStarD2(const SimplexId p[3],
2068 const int id) const {
2069 SimplexId wrapXLeft = 0;
2070 if(p[0] == 0)
2071 wrapXLeft = 6 * wrap_[0];
2072 switch(id) {
2073 case 0:
2074 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6;
2075 case 1:
2076 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 2;
2077 case 2:
2078 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + (p[0] - 1) * 6 + 5
2079 + wrapXLeft;
2080 case 3:
2081 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + (p[0] - 1) * 6 + 4
2082 + wrapXLeft;
2083 }
2084 return -1;
2085}
2086
2087inline ttk::SimplexId
2088 ttk::PeriodicImplicitTriangulation::getEdgeStarD3(const SimplexId p[3],
2089 const int id) const {
2090 SimplexId wrapYTop = 0;
2091 if(p[1] == 0)
2092 wrapYTop = 6 * wrap_[1];
2093 switch(id) {
2094 case 0:
2095 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 2;
2096 case 1:
2097 return p[2] * tetshift_[1] + p[1] * tetshift_[0] + p[0] * 6 + 3;
2098 case 2:
2099 return p[2] * tetshift_[1] + (p[1] - 1) * tetshift_[0] + p[0] * 6 + 1
2100 + wrapYTop;
2101 case 3:
2102 return p[2] * tetshift_[1] + (p[1] - 1) * tetshift_[0] + p[0] * 6 + 5
2103 + wrapYTop;
2104 }
2105 return -1;
2106}
2107
2108inline ttk::SimplexId
2109 ttk::PeriodicImplicitTriangulation::getTriangleVertexF(const SimplexId p[3],
2110 const int id) const {
2111 SimplexId wrapXRight = 0;
2112 SimplexId wrapYBottom = 0;
2113 if(p[0] / 2 == nbvoxels_[0])
2114 wrapXRight = -wrap_[0];
2115 if(p[1] == nbvoxels_[1])
2116 wrapYBottom = -wrap_[1];
2117 if(p[0] % 2) {
2118 switch(id) {
2119 case 0:
2120 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2121 + wrapXRight;
2122 case 1:
2123 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2124 + wrapYBottom;
2125 case 2:
2126 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0] + 1
2127 + wrapXRight + wrapYBottom;
2128 }
2129 } else {
2130 switch(id) {
2131 case 0:
2132 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1];
2133 case 1:
2134 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2135 + wrapXRight;
2136 case 2:
2137 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2138 + wrapYBottom;
2139 }
2140 }
2141 return -1;
2142}
2143
2144inline ttk::SimplexId
2145 ttk::PeriodicImplicitTriangulation::getTriangleVertexH(const SimplexId p[3],
2146 const int id) const {
2147 SimplexId wrapXRight = 0;
2148 SimplexId wrapZFront = 0;
2149 if(p[0] / 2 == nbvoxels_[0])
2150 wrapXRight = -wrap_[0];
2151 if(p[2] == nbvoxels_[2])
2152 wrapZFront = -wrap_[2];
2153 if(p[0] % 2) {
2154 switch(id) {
2155 case 0:
2156 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2157 + wrapXRight;
2158 case 1:
2159 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1] + 1
2160 + wrapXRight + wrapZFront;
2161 case 2:
2162 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2163 + wrapZFront;
2164 }
2165 } else {
2166 switch(id) {
2167 case 0:
2168 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1];
2169 case 1:
2170 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2171 + wrapXRight;
2172 case 2:
2173 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2174 + wrapZFront;
2175 }
2176 }
2177 return -1;
2178}
2179
2180inline ttk::SimplexId
2181 ttk::PeriodicImplicitTriangulation::getTriangleVertexC(const SimplexId p[3],
2182 const int id) const {
2183 SimplexId wrapYBottom = 0;
2184 SimplexId wrapZFront = 0;
2185 if(p[1] == nbvoxels_[1])
2186 wrapYBottom = -wrap_[1];
2187 if(p[2] == nbvoxels_[2])
2188 wrapZFront = -wrap_[2];
2189 if(p[0] % 2) {
2190 switch(id) {
2191 case 0:
2192 return (p[0] / 2) + p[1] * vshift_[0] + p[2] * vshift_[1];
2193 case 1:
2194 return (p[0] / 2) + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2195 + wrapZFront;
2196 case 2:
2197 return (p[0] / 2) + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2198 + vshift_[0] + wrapYBottom + wrapZFront;
2199 }
2200 } else {
2201 switch(id) {
2202 case 0:
2203 return (p[0] / 2) + p[1] * vshift_[0] + p[2] * vshift_[1];
2204 case 1:
2205 return (p[0] / 2) + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2206 + wrapYBottom;
2207 case 2:
2208 return (p[0] / 2) + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2209 + vshift_[0] + wrapYBottom + wrapZFront;
2210 }
2211 }
2212 return -1;
2213}
2214
2215inline ttk::SimplexId
2216 ttk::PeriodicImplicitTriangulation::getTriangleVertexD1(const SimplexId p[3],
2217 const int id) const {
2218 SimplexId wrapXRight = 0;
2219 SimplexId wrapYBottom = 0;
2220 SimplexId wrapZFront = 0;
2221 if(p[0] / 2 == nbvoxels_[0])
2222 wrapXRight = -wrap_[0];
2223 if(p[1] == nbvoxels_[1])
2224 wrapYBottom = -wrap_[1];
2225 if(p[2] == nbvoxels_[2])
2226 wrapZFront = -wrap_[2];
2227 if(p[0] % 2) {
2228 switch(id) {
2229 case 0:
2230 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2231 + wrapXRight;
2232 case 1:
2233 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2234 + wrapZFront;
2235 case 2:
2236 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2237 + vshift_[0] + wrapYBottom + wrapZFront;
2238 }
2239 } else {
2240 switch(id) {
2241 case 0:
2242 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2243 + wrapXRight;
2244 case 1:
2245 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0] + 1
2246 + wrapXRight + wrapYBottom;
2247 case 2:
2248 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2249 + vshift_[0] + wrapYBottom + wrapZFront;
2250 }
2251 }
2252 return -1;
2253}
2254
2255inline ttk::SimplexId
2256 ttk::PeriodicImplicitTriangulation::getTriangleVertexD2(const SimplexId p[3],
2257 const int id) const {
2258 SimplexId wrapXRight = 0;
2259 SimplexId wrapYBottom = 0;
2260 SimplexId wrapZFront = 0;
2261 if(p[0] / 2 == nbvoxels_[0])
2262 wrapXRight = -wrap_[0];
2263 if(p[1] == nbvoxels_[1])
2264 wrapYBottom = -wrap_[1];
2265 if(p[2] == nbvoxels_[2])
2266 wrapZFront = -wrap_[2];
2267 if(p[0] % 2) {
2268 switch(id) {
2269 case 0:
2270 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2271 + wrapXRight;
2272 case 1:
2273 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2274 + vshift_[1] + 1 + wrapXRight + wrapYBottom + wrapZFront;
2275 case 2:
2276 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2277 + vshift_[1] + wrapYBottom + wrapZFront;
2278 }
2279 } else {
2280 switch(id) {
2281 case 0:
2282 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1];
2283 case 1:
2284 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2285 + wrapXRight;
2286 case 2:
2287 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2288 + vshift_[1] + wrapYBottom + wrapZFront;
2289 }
2290 }
2291 return -1;
2292}
2293
2294inline ttk::SimplexId
2295 ttk::PeriodicImplicitTriangulation::getTriangleVertexD3(const SimplexId p[3],
2296 const int id) const {
2297 SimplexId wrapXRight = 0;
2298 SimplexId wrapYBottom = 0;
2299 SimplexId wrapZFront = 0;
2300 if(p[0] / 2 == nbvoxels_[0])
2301 wrapXRight = -wrap_[0];
2302 if(p[1] == nbvoxels_[1])
2303 wrapYBottom = -wrap_[1];
2304 if(p[2] == nbvoxels_[2])
2305 wrapZFront = -wrap_[2];
2306 if(p[0] % 2) {
2307 switch(id) {
2308 case 0:
2309 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2310 + wrapXRight;
2311 case 1:
2312 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1] + 1
2313 + wrapXRight + wrapZFront;
2314 case 2:
2315 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2316 + vshift_[1] + wrapYBottom + wrapZFront;
2317 }
2318 } else {
2319 switch(id) {
2320 case 0:
2321 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1
2322 + wrapXRight;
2323 case 1:
2324 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2325 + wrapYBottom;
2326 case 2:
2327 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2328 + vshift_[1] + wrapYBottom + wrapZFront;
2329 }
2330 }
2331 return -1;
2332}
2333
2334inline ttk::SimplexId
2335 ttk::PeriodicImplicitTriangulation::getTriangleEdgeF_0(const SimplexId p[3],
2336 const int id) const {
2337 switch(id) {
2338 case 0:
2339 return p[0] / 2 + p[1] * eshift_[0] + p[2] * eshift_[1];
2340 case 1:
2341 return esetshift_[0] + p[0] / 2 + p[1] * eshift_[2] + p[2] * eshift_[3];
2342 case 2:
2343 return esetshift_[2] + p[0] / 2 + p[1] * eshift_[6] + p[2] * eshift_[7];
2344 }
2345 return -1;
2346}
2347
2348inline ttk::SimplexId
2349 ttk::PeriodicImplicitTriangulation::getTriangleEdgeF_1(const SimplexId p[3],
2350 const int id) const {
2351 SimplexId wrapXRight = 0;
2352 if(p[0] / 2 == nbvoxels_[0])
2353 wrapXRight = -wrap_[0];
2354 SimplexId wrapYBottom = 0;
2355 if(p[1] == nbvoxels_[1])
2356 wrapYBottom = -wrap_[1];
2357 switch(id) {
2358 case 0:
2359 return p[0] / 2 + (p[1] + 1) * eshift_[0] + p[2] * eshift_[1]
2360 + wrapYBottom;
2361 case 1:
2362 return esetshift_[0] + p[0] / 2 + p[1] * eshift_[2] + p[2] * eshift_[3]
2363 + 1 + wrapXRight;
2364 case 2:
2365 return esetshift_[2] + p[0] / 2 + p[1] * eshift_[6] + p[2] * eshift_[7];
2366 }
2367 return -1;
2368}
2369
2370inline ttk::SimplexId
2371 ttk::PeriodicImplicitTriangulation::getTriangleEdgeH_0(const SimplexId p[3],
2372 const int id) const {
2373 switch(id) {
2374 case 0:
2375 return p[0] / 2 + p[1] * eshift_[0] + p[2] * eshift_[1];
2376 case 1:
2377 return esetshift_[1] + p[0] / 2 + p[1] * eshift_[4] + p[2] * eshift_[5];
2378 case 2:
2379 return esetshift_[4] + p[0] / 2 + p[1] * eshift_[10] + p[2] * eshift_[11];
2380 }
2381 return -1;
2382}
2383
2384inline ttk::SimplexId
2385 ttk::PeriodicImplicitTriangulation::getTriangleEdgeH_1(const SimplexId p[3],
2386 const int id) const {
2387 SimplexId wrapXRight = 0;
2388 if(p[0] / 2 == nbvoxels_[0])
2389 wrapXRight = -wrap_[0];
2390 SimplexId wrapZFront = 0;
2391 if(p[2] == nbvoxels_[2])
2392 wrapZFront = -wrap_[2];
2393 switch(id) {
2394 case 0:
2395 return p[0] / 2 + p[1] * eshift_[0] + (p[2] + 1) * eshift_[1]
2396 + wrapZFront;
2397 case 1:
2398 return esetshift_[1] + p[0] / 2 + p[1] * eshift_[4] + p[2] * eshift_[5]
2399 + 1 + wrapXRight;
2400 case 2:
2401 return esetshift_[4] + p[0] / 2 + p[1] * eshift_[10] + p[2] * eshift_[11];
2402 }
2403 return -1;
2404}
2405
2406inline ttk::SimplexId
2407 ttk::PeriodicImplicitTriangulation::getTriangleEdgeC_0(const SimplexId p[3],
2408 const int id) const {
2409 SimplexId wrapYBottom = 0;
2410 if(p[1] == nbvoxels_[1])
2411 wrapYBottom = -wrap_[1];
2412 switch(id) {
2413 case 0:
2414 return esetshift_[0] + p[0] / 2 + p[1] * eshift_[2] + p[2] * eshift_[3];
2415 case 1:
2416 return esetshift_[1] + p[0] / 2 + (p[1] + 1) * eshift_[4]
2417 + p[2] * eshift_[5] + wrapYBottom;
2418 case 2:
2419 return esetshift_[3] + p[0] / 2 + p[1] * eshift_[8] + p[2] * eshift_[9];
2420 }
2421 return -1;
2422}
2423
2424inline ttk::SimplexId
2425 ttk::PeriodicImplicitTriangulation::getTriangleEdgeC_1(const SimplexId p[3],
2426 const int id) const {
2427 SimplexId wrapZFront = 0;
2428 if(p[2] == nbvoxels_[2])
2429 wrapZFront = -wrap_[2];
2430 switch(id) {
2431 case 0:
2432 return esetshift_[0] + p[0] / 2 + p[1] * eshift_[2]
2433 + (p[2] + 1) * eshift_[3] + wrapZFront;
2434 case 1:
2435 return esetshift_[1] + p[0] / 2 + p[1] * eshift_[4] + p[2] * eshift_[5];
2436 case 2:
2437 return esetshift_[3] + p[0] / 2 + p[1] * eshift_[8] + p[2] * eshift_[9];
2438 }
2439 return -1;
2440}
2441
2442inline ttk::SimplexId
2443 ttk::PeriodicImplicitTriangulation::getTriangleEdgeD1_0(const SimplexId p[3],
2444 const int id) const {
2445 SimplexId wrapXRight = 0;
2446 if(p[0] / 2 == nbvoxels_[0])
2447 wrapXRight = -wrap_[0];
2448 SimplexId wrapYBottom = 0;
2449 if(p[1] == nbvoxels_[1])
2450 wrapYBottom = -wrap_[1];
2451 switch(id) {
2452 case 0:
2453 return esetshift_[0] + p[0] / 2 + p[1] * eshift_[2] + p[2] * eshift_[3]
2454 + 1 + wrapXRight;
2455 case 1:
2456 return esetshift_[4] + p[0] / 2 + (p[1] + 1) * eshift_[10]
2457 + p[2] * eshift_[11] + wrapYBottom;
2458 case 2:
2459 return esetshift_[5] + p[0] / 2 + p[1] * eshift_[12] + p[2] * eshift_[13];
2460 }
2461 return -1;
2462}
2463
2464inline ttk::SimplexId
2465 ttk::PeriodicImplicitTriangulation::getTriangleEdgeD1_1(const SimplexId p[3],
2466 const int id) const {
2467 SimplexId wrapZFront = 0;
2468 if(p[2] == nbvoxels_[2])
2469 wrapZFront = -wrap_[2];
2470 switch(id) {
2471 case 0:
2472 return esetshift_[0] + p[0] / 2 + p[1] * eshift_[2]
2473 + (p[2] + 1) * eshift_[3] + wrapZFront;
2474 case 1:
2475 return esetshift_[4] + p[0] / 2 + p[1] * eshift_[10] + p[2] * eshift_[11];
2476 case 2:
2477 return esetshift_[5] + p[0] / 2 + p[1] * eshift_[12] + p[2] * eshift_[13];
2478 }
2479 return -1;
2480}
2481
2482inline ttk::SimplexId
2483 ttk::PeriodicImplicitTriangulation::getTriangleEdgeD2_0(const SimplexId p[3],
2484 const int id) const {
2485 switch(id) {
2486 case 0:
2487 return p[0] / 2 + p[1] * eshift_[0] + p[2] * eshift_[1];
2488 case 1:
2489 return esetshift_[3] + p[0] / 2 + p[1] * eshift_[8] + p[2] * eshift_[9];
2490 case 2:
2491 return esetshift_[5] + p[0] / 2 + p[1] * eshift_[12] + p[2] * eshift_[13];
2492 }
2493 return -1;
2494}
2495
2496inline ttk::SimplexId
2497 ttk::PeriodicImplicitTriangulation::getTriangleEdgeD2_1(const SimplexId p[3],
2498 const int id) const {
2499 SimplexId wrapXRight = 0;
2500 SimplexId wrapYBottom = 0;
2501 SimplexId wrapZFront = 0;
2502 if(p[0] / 2 == nbvoxels_[0])
2503 wrapXRight = -wrap_[0];
2504 if(p[1] == nbvoxels_[1])
2505 wrapYBottom = -wrap_[1];
2506 if(p[2] == nbvoxels_[2])
2507 wrapZFront = -wrap_[2];
2508 switch(id) {
2509 case 0:
2510 return p[0] / 2 + (p[1] + 1) * eshift_[0] + (p[2] + 1) * eshift_[1]
2511 + wrapYBottom + wrapZFront;
2512 case 1:
2513 return esetshift_[3] + p[0] / 2 + p[1] * eshift_[8] + p[2] * eshift_[9]
2514 + 1 + wrapXRight;
2515 case 2:
2516 return esetshift_[5] + p[0] / 2 + p[1] * eshift_[12] + p[2] * eshift_[13];
2517 }
2518 return -1;
2519}
2520
2521inline ttk::SimplexId
2522 ttk::PeriodicImplicitTriangulation::getTriangleEdgeD3_0(const SimplexId p[3],
2523 const int id) const {
2524 SimplexId wrapYBottom = 0;
2525 if(p[1] == nbvoxels_[1])
2526 wrapYBottom = -wrap_[1];
2527 switch(id) {
2528 case 0:
2529 return esetshift_[1] + p[0] / 2 + (p[1] + 1) * eshift_[4]
2530 + p[2] * eshift_[5] + wrapYBottom;
2531 case 1:
2532 return esetshift_[2] + p[0] / 2 + p[1] * eshift_[6] + p[2] * eshift_[7];
2533 case 2:
2534 return esetshift_[5] + p[0] / 2 + p[1] * eshift_[12] + p[2] * eshift_[13];
2535 }
2536 return -1;
2537}
2538
2539inline ttk::SimplexId
2540 ttk::PeriodicImplicitTriangulation::getTriangleEdgeD3_1(const SimplexId p[3],
2541 const int id) const {
2542 SimplexId wrapXRight = 0;
2543 SimplexId wrapZFront = 0;
2544 if(p[0] / 2 == nbvoxels_[0])
2545 wrapXRight = -wrap_[0];
2546 if(p[2] == nbvoxels_[2])
2547 wrapZFront = -wrap_[2];
2548 switch(id) {
2549 case 0:
2550 return esetshift_[1] + p[0] / 2 + p[1] * eshift_[4] + p[2] * eshift_[5]
2551 + 1 + wrapXRight;
2552 case 1:
2553 return esetshift_[2] + p[0] / 2 + p[1] * eshift_[6]
2554 + (p[2] + 1) * eshift_[7] + wrapZFront;
2555 case 2:
2556 return esetshift_[5] + p[0] / 2 + p[1] * eshift_[12] + p[2] * eshift_[13];
2557 }
2558 return -1;
2559}
2560
2561inline ttk::SimplexId
2562 ttk::PeriodicImplicitTriangulation::getTriangleLinkF(const SimplexId p[3],
2563 const int id) const {
2564 SimplexId wrapXRight = 0;
2565 SimplexId wrapYBottom = 0;
2566 SimplexId wrapZBack = 0;
2567 SimplexId wrapZFront = 0;
2568 if(p[0] / 2 == nbvoxels_[0])
2569 wrapXRight = -wrap_[0];
2570 if(p[1] == nbvoxels_[1])
2571 wrapYBottom = -wrap_[1];
2572 if(p[2] == 0)
2573 wrapZBack = wrap_[2];
2574 if(p[2] == nbvoxels_[2])
2575 wrapZFront = -wrap_[2];
2576 switch(id) {
2577 case 0:
2578 return p[0] / 2 + (p[1] + 1) * vshift_[0] + (p[2] + 1) * vshift_[1]
2579 + wrapYBottom + wrapZFront;
2580 case 1:
2581 return p[0] / 2 + p[1] * vshift_[0] + (p[2] - 1) * vshift_[1] + 1
2582 + wrapXRight + wrapZBack;
2583 }
2584 return -1;
2585}
2586
2587inline ttk::SimplexId
2588 ttk::PeriodicImplicitTriangulation::getTriangleLinkH(const SimplexId p[3],
2589 const int id) const {
2590 SimplexId wrapXRight = 0;
2591 SimplexId wrapYTop = 0;
2592 SimplexId wrapYBottom = 0;
2593 SimplexId wrapZFront = 0;
2594 if(p[0] / 2 == nbvoxels_[0])
2595 wrapXRight = -wrap_[0];
2596 if(p[1] == 0)
2597 wrapYTop = wrap_[1];
2598 if(p[1] == nbvoxels_[1])
2599 wrapYBottom = -wrap_[1];
2600 if(p[2] == nbvoxels_[2])
2601 wrapZFront = -wrap_[2];
2602 switch(id) {
2603 case 0:
2604 return p[0] / 2 + (p[1] + 1) * vshift_[0] + (p[2] + 1) * vshift_[1]
2605 + wrapYBottom + wrapZFront;
2606 case 1:
2607 return p[0] / 2 + (p[1] - 1) * vshift_[0] + p[2] * vshift_[1] + 1
2608 + wrapXRight + wrapYTop;
2609 }
2610 return -1;
2611}
2612
2613inline ttk::SimplexId
2614 ttk::PeriodicImplicitTriangulation::getTriangleLinkC(const SimplexId p[3],
2615 const int id) const {
2616 SimplexId wrapXLeft = 0;
2617 SimplexId wrapXRight = 0;
2618 SimplexId wrapYBottom = 0;
2619 SimplexId wrapZFront = 0;
2620 if(p[0] / 2 == 0)
2621 wrapXLeft = wrap_[0];
2622 if(p[0] / 2 == nbvoxels_[0])
2623 wrapXRight = -wrap_[0];
2624 if(p[1] == nbvoxels_[1])
2625 wrapYBottom = -wrap_[1];
2626 if(p[2] == nbvoxels_[2])
2627 wrapZFront = -wrap_[2];
2628 switch(id) {
2629 case 0:
2630 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1] + 1 + wrapXRight;
2631 case 1:
2632 return p[0] / 2 + (p[1] + 1) * vshift_[0] + (p[2] + 1) * vshift_[1] - 1
2633 + wrapXLeft + wrapYBottom + wrapZFront;
2634 }
2635 return -1;
2636}
2637
2638inline ttk::SimplexId
2639 ttk::PeriodicImplicitTriangulation::getTriangleLinkD1(const SimplexId p[3],
2640 const int id) const {
2641 SimplexId wrapXRight = 0;
2642 SimplexId wrapYBottom = 0;
2643 SimplexId wrapZFront = 0;
2644 if(p[0] / 2 == nbvoxels_[0])
2645 wrapXRight = -wrap_[0];
2646 if(p[1] == nbvoxels_[1])
2647 wrapYBottom = -wrap_[1];
2648 if(p[2] == nbvoxels_[2])
2649 wrapZFront = -wrap_[2];
2650 if(p[0] % 2) {
2651 switch(id) {
2652 case 0:
2653 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1];
2654 case 1:
2655 return p[0] / 2 + p[1] * vshift_[0] + (p[2] + 1) * vshift_[1] + 1
2656 + wrapXRight + wrapZFront;
2657 }
2658 } else {
2659 switch(id) {
2660 case 0:
2661 return p[0] / 2 + (p[1] + 1) * vshift_[0] + p[2] * vshift_[1]
2662 + wrapYBottom;
2663 case 1:
2664 return p[0] / 2 + (p[1] + 1) * vshift_[0] + (p[2] + 1) * vshift_[1] + 1
2665 + wrapXRight + wrapYBottom + wrapZFront;
2666 }
2667 }
2668 return -1;
2669}
2670
2671inline ttk::SimplexId
2672 ttk::PeriodicImplicitTriangulation::getTriangleLinkD2(const SimplexId p[3],
2673 const int id) const {
2674 SimplexId wrapXRight = 0;
2675 SimplexId wrapYBottom = 0;
2676 SimplexId wrapZFront = 0;
2677 if(p[0] / 2 == nbvoxels_[0])
2678 wrapXRight = -wrap_[0];
2679 if(p[1] == nbvoxels_[1])
2680 wrapYBottom = -wrap_[1];
2681 if(p[2] == nbvoxels_[2])
2682 wrapZFront = -wrap_[2];
2683 if(p[0] % 2) {
2684 switch(id) {
2685 case 0:
2686 return p[0] / 2 + (p[1] + 1) * vshift_[0] + p[2] * vshift_[1] + 1
2687 + wrapXRight + wrapYBottom;
2688 case 1:
2689 return p[0] / 2 + p[1] * vshift_[0] + (p[2] + 1) * vshift_[1] + 1
2690 + wrapXRight + wrapZFront;
2691 }
2692 } else {
2693 switch(id) {
2694 case 0:
2695 return p[0] / 2 + (p[1] + 1) * vshift_[0] + p[2] * vshift_[1]
2696 + wrapYBottom;
2697 case 1:
2698 return p[0] / 2 + p[1] * vshift_[0] + (p[2] + 1) * vshift_[1]
2699 + wrapZFront;
2700 }
2701 }
2702 return -1;
2703}
2704
2705inline ttk::SimplexId
2706 ttk::PeriodicImplicitTriangulation::getTriangleLinkD3(const SimplexId p[3],
2707 const int id) const {
2708 SimplexId wrapXRight = 0;
2709 SimplexId wrapYBottom = 0;
2710 SimplexId wrapZFront = 0;
2711 if(p[0] / 2 == nbvoxels_[0])
2712 wrapXRight = -wrap_[0];
2713 if(p[1] == nbvoxels_[1])
2714 wrapYBottom = -wrap_[1];
2715 if(p[2] == nbvoxels_[2])
2716 wrapZFront = -wrap_[2];
2717 if(p[0] % 2) {
2718 switch(id) {
2719 case 0:
2720 return p[0] / 2 + p[1] * vshift_[0] + (p[2] + 1) * vshift_[1]
2721 + wrapZFront;
2722 case 1:
2723 return p[0] / 2 + (p[1] + 1) * vshift_[0] + (p[2] + 1) * vshift_[1] + 1
2724 + wrapXRight + wrapYBottom + wrapZFront;
2725 }
2726 } else {
2727 switch(id) {
2728 case 0:
2729 return p[0] / 2 + p[1] * vshift_[0] + p[2] * vshift_[1];
2730 case 1:
2731 return p[0] / 2 + (p[1] + 1) * vshift_[0] + p[2] * vshift_[1] + 1
2732 + wrapXRight + wrapYBottom;
2733 }
2734 }
2735 return -1;
2736}
2737
2738inline ttk::SimplexId
2739 ttk::PeriodicImplicitTriangulation::getTriangleStarF(const SimplexId p[3],
2740 const int id) const {
2741 SimplexId wrapZBack = 0;
2742 if(p[2] == 0)
2743 wrapZBack = 6 * wrap_[2];
2744 if(p[0] % 2) {
2745 switch(id) {
2746 case 0:
2747 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 1;
2748 case 1:
2749 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + (p[2] - 1) * tetshift_[1]
2750 + 4 + wrapZBack;
2751 }
2752 } else {
2753 switch(id) {
2754 case 0:
2755 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1];
2756 case 1:
2757 return p[0] * 3 + p[1] * tetshift_[0] + (p[2] - 1) * tetshift_[1] + 3
2758 + wrapZBack;
2759 }
2760 }
2761 return -1;
2762}
2763
2764inline ttk::SimplexId
2765 ttk::PeriodicImplicitTriangulation::getTriangleStarH(const SimplexId p[3],
2766 const int id) const {
2767 SimplexId wrapYTop = 0;
2768 if(p[1] == 0)
2769 wrapYTop = 6 * wrap_[1];
2770 if(p[0] % 2) {
2771 switch(id) {
2772 case 0:
2773 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 3;
2774 case 1:
2775 return (p[0] - 1) * 3 + (p[1] - 1) * tetshift_[0] + p[2] * tetshift_[1]
2776 + 5 + wrapYTop;
2777 }
2778 } else {
2779 switch(id) {
2780 case 0:
2781 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 2;
2782 case 1:
2783 return p[0] * 3 + (p[1] - 1) * tetshift_[0] + p[2] * tetshift_[1] + 1
2784 + wrapYTop;
2785 }
2786 }
2787 return -1;
2788}
2789
2790inline ttk::SimplexId
2791 ttk::PeriodicImplicitTriangulation::getTriangleStarC(const SimplexId p[3],
2792 const int id) const {
2793 SimplexId wrapXLeft = 0;
2794 if(p[0] < 2)
2795 wrapXLeft = 6 * wrap_[0];
2796 if(p[0] % 2) {
2797 switch(id) {
2798 case 0:
2799 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 2;
2800 case 1:
2801 return (p[0] / 2 - 1) * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1]
2802 + 4 + wrapXLeft;
2803 }
2804 } else {
2805 switch(id) {
2806 case 0:
2807 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1];
2808 case 1:
2809 return (p[0] / 2 - 1) * 6 + p[1] * tetshift_[0] + p[2] * tetshift_[1]
2810 + 5 + wrapXLeft;
2811 }
2812 }
2813 return -1;
2814}
2815
2816inline ttk::SimplexId
2817 ttk::PeriodicImplicitTriangulation::getTriangleStarD1(const SimplexId p[3],
2818 const int id) const {
2819 if(p[0] % 2) {
2820 switch(id) {
2821 case 0:
2822 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 2;
2823 case 1:
2824 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 3;
2825 }
2826 } else {
2827 switch(id) {
2828 case 0:
2829 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 1;
2830 case 1:
2831 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 5;
2832 }
2833 }
2834 return -1;
2835}
2836
2837inline ttk::SimplexId
2838 ttk::PeriodicImplicitTriangulation::getTriangleStarD2(const SimplexId p[3],
2839 const int id) const {
2840 if(p[0] % 2) {
2841 switch(id) {
2842 case 0:
2843 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 5;
2844 case 1:
2845 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 4;
2846 }
2847 } else {
2848 switch(id) {
2849 case 0:
2850 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1];
2851 case 1:
2852 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 2;
2853 }
2854 }
2855 return -1;
2856}
2857
2858inline ttk::SimplexId
2859 ttk::PeriodicImplicitTriangulation::getTriangleStarD3(const SimplexId p[3],
2860 const int id) const {
2861 if(p[0] % 2) {
2862 switch(id) {
2863 case 0:
2864 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 3;
2865 case 1:
2866 return (p[0] - 1) * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 4;
2867 }
2868 } else {
2869 switch(id) {
2870 case 0:
2871 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1];
2872 case 1:
2873 return p[0] * 3 + p[1] * tetshift_[0] + p[2] * tetshift_[1] + 1;
2874 }
2875 }
2876 return -1;
2877}
2878
2879inline ttk::SimplexId
2880 ttk::PeriodicImplicitTriangulation::getTetrahedronVertexABCG(
2881 const SimplexId p[3], const int id) const {
2882 SimplexId wrapXRight = 0;
2883 SimplexId wrapYBottom = 0;
2884 SimplexId wrapZFront = 0;
2885 if(p[0] == nbvoxels_[0])
2886 wrapXRight = -wrap_[0];
2887 if(p[1] == nbvoxels_[1])
2888 wrapYBottom = -wrap_[1];
2889 if(p[2] == nbvoxels_[2])
2890 wrapZFront = -wrap_[2];
2891 switch(id) {
2892 case 0:
2893 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1];
2894 case 1:
2895 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + 1 + wrapXRight;
2896 case 2:
2897 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2898 + wrapYBottom;
2899 case 3:
2900 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2901 + vshift_[1] + wrapYBottom + wrapZFront;
2902 }
2903 return -1;
2904}
2905
2906inline ttk::SimplexId
2907 ttk::PeriodicImplicitTriangulation::getTetrahedronVertexBCDG(
2908 const SimplexId p[3], const int id) const {
2909 SimplexId wrapXRight = 0;
2910 SimplexId wrapYBottom = 0;
2911 SimplexId wrapZFront = 0;
2912 if(p[0] == nbvoxels_[0])
2913 wrapXRight = -wrap_[0];
2914 if(p[1] == nbvoxels_[1])
2915 wrapYBottom = -wrap_[1];
2916 if(p[2] == nbvoxels_[2])
2917 wrapZFront = -wrap_[2];
2918 switch(id) {
2919 case 0:
2920 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + 1 + wrapXRight;
2921 case 1:
2922 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2923 + wrapYBottom;
2924 case 2:
2925 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0] + 1
2926 + wrapXRight + wrapYBottom;
2927 case 3:
2928 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2929 + vshift_[1] + wrapYBottom + wrapZFront;
2930 }
2931 return -1;
2932}
2933
2934inline ttk::SimplexId
2935 ttk::PeriodicImplicitTriangulation::getTetrahedronVertexABEG(
2936 const SimplexId p[3], const int id) const {
2937 SimplexId wrapXRight = 0;
2938 SimplexId wrapYBottom = 0;
2939 SimplexId wrapZFront = 0;
2940 if(p[0] == nbvoxels_[0])
2941 wrapXRight = -wrap_[0];
2942 if(p[1] == nbvoxels_[1])
2943 wrapYBottom = -wrap_[1];
2944 if(p[2] == nbvoxels_[2])
2945 wrapZFront = -wrap_[2];
2946 switch(id) {
2947 case 0:
2948 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1];
2949 case 1:
2950 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + 1 + wrapXRight;
2951 case 2:
2952 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2953 + wrapZFront;
2954 case 3:
2955 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2956 + vshift_[1] + wrapYBottom + wrapZFront;
2957 }
2958 return -1;
2959}
2960
2961inline ttk::SimplexId
2962 ttk::PeriodicImplicitTriangulation::getTetrahedronVertexBEFG(
2963 const SimplexId p[3], const int id) const {
2964 SimplexId wrapXRight = 0;
2965 SimplexId wrapYBottom = 0;
2966 SimplexId wrapZFront = 0;
2967 if(p[0] == nbvoxels_[0])
2968 wrapXRight = -wrap_[0];
2969 if(p[1] == nbvoxels_[1])
2970 wrapYBottom = -wrap_[1];
2971 if(p[2] == nbvoxels_[2])
2972 wrapZFront = -wrap_[2];
2973 switch(id) {
2974 case 0:
2975 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + 1 + wrapXRight;
2976 case 1:
2977 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1]
2978 + wrapZFront;
2979 case 2:
2980 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1] + 1
2981 + wrapXRight + wrapZFront;
2982 case 3:
2983 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
2984 + vshift_[1] + wrapYBottom + wrapZFront;
2985 }
2986 return -1;
2987}
2988
2989inline ttk::SimplexId
2990 ttk::PeriodicImplicitTriangulation::getTetrahedronVertexBFGH(
2991 const SimplexId p[3], const int id) const {
2992 SimplexId wrapXRight = 0;
2993 SimplexId wrapYBottom = 0;
2994 SimplexId wrapZFront = 0;
2995 if(p[0] == nbvoxels_[0])
2996 wrapXRight = -wrap_[0];
2997 if(p[1] == nbvoxels_[1])
2998 wrapYBottom = -wrap_[1];
2999 if(p[2] == nbvoxels_[2])
3000 wrapZFront = -wrap_[2];
3001 switch(id) {
3002 case 0:
3003 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + 1 + wrapXRight;
3004 case 1:
3005 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[1] + 1
3006 + wrapXRight + wrapZFront;
3007 case 2:
3008 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
3009 + vshift_[1] + wrapYBottom + wrapZFront;
3010 case 3:
3011 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
3012 + vshift_[1] + 1 + wrapXRight + wrapYBottom + wrapZFront;
3013 }
3014 return -1;
3015}
3016
3017inline ttk::SimplexId
3018 ttk::PeriodicImplicitTriangulation::getTetrahedronVertexBDGH(
3019 const SimplexId p[3], const int id) const {
3020 SimplexId wrapXRight = 0;
3021 SimplexId wrapYBottom = 0;
3022 SimplexId wrapZFront = 0;
3023 if(p[0] == nbvoxels_[0])
3024 wrapXRight = -wrap_[0];
3025 if(p[1] == nbvoxels_[1])
3026 wrapYBottom = -wrap_[1];
3027 if(p[2] == nbvoxels_[2])
3028 wrapZFront = -wrap_[2];
3029 switch(id) {
3030 case 0:
3031 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + 1 + wrapXRight;
3032 case 1:
3033 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0] + 1
3034 + wrapXRight + wrapYBottom;
3035 case 2:
3036 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
3037 + vshift_[1] + wrapYBottom + wrapZFront;
3038 case 3:
3039 return p[0] + p[1] * vshift_[0] + p[2] * vshift_[1] + vshift_[0]
3040 + vshift_[1] + 1 + wrapXRight + wrapYBottom + wrapZFront;
3041 }
3042 return -1;
3043}
3044
3045inline ttk::SimplexId
3046 ttk::PeriodicImplicitTriangulation::getTetrahedronEdgeABCG(
3047 const SimplexId p[3], const int id) const {
3048 SimplexId wrapYBottom = 0;
3049 if(p[1] == nbvoxels_[1])
3050 wrapYBottom = -wrap_[1];
3051 switch(id) {
3052 case 0:
3053 return p[0] + p[1] * eshift_[0] + p[2] * eshift_[1];
3054 case 1:
3055 return esetshift_[0] + p[0] + p[1] * eshift_[2] + p[2] * eshift_[3];
3056 case 2:
3057 return esetshift_[1] + p[0] + (p[1] + 1) * eshift_[4] + p[2] * eshift_[5]
3058 + wrapYBottom;
3059 case 3:
3060 return esetshift_[2] + p[0] + p[1] * eshift_[6] + p[2] * eshift_[7];
3061 case 4:
3062 return esetshift_[3] + p[0] + p[1] * eshift_[8] + p[2] * eshift_[9];
3063 case 5:
3064 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13];
3065 }
3066 return -1;
3067}
3068
3069inline ttk::SimplexId
3070 ttk::PeriodicImplicitTriangulation::getTetrahedronEdgeBCDG(
3071 const SimplexId p[3], const int id) const {
3072 SimplexId wrapXRight = 0;
3073 SimplexId wrapYBottom = 0;
3074 if(p[0] == nbvoxels_[0])
3075 wrapXRight = -wrap_[0];
3076 if(p[1] == nbvoxels_[1])
3077 wrapYBottom = -wrap_[1];
3078 switch(id) {
3079 case 0:
3080 return p[0] + (p[1] + 1) * eshift_[0] + p[2] * eshift_[1] + wrapYBottom;
3081 case 1:
3082 return esetshift_[0] + (p[0] + 1) + p[1] * eshift_[2] + p[2] * eshift_[3]
3083 + wrapXRight;
3084 case 2:
3085 return esetshift_[1] + p[0] + (p[1] + 1) * eshift_[4] + p[2] * eshift_[5]
3086 + wrapYBottom;
3087 case 3:
3088 return esetshift_[2] + p[0] + p[1] * eshift_[6] + p[2] * eshift_[7];
3089 case 4:
3090 return esetshift_[4] + p[0] + (p[1] + 1) * eshift_[10]
3091 + p[2] * eshift_[11] + wrapYBottom;
3092 case 5:
3093 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13];
3094 }
3095 return -1;
3096}
3097
3098inline ttk::SimplexId
3099 ttk::PeriodicImplicitTriangulation::getTetrahedronEdgeABEG(
3100 const SimplexId p[3], const int id) const {
3101 SimplexId wrapZFront = 0;
3102 if(p[2] == nbvoxels_[2])
3103 wrapZFront = -wrap_[2];
3104 switch(id) {
3105 case 0:
3106 return p[0] + p[1] * eshift_[0] + p[2] * eshift_[1];
3107 case 1:
3108 return esetshift_[0] + p[0] + p[1] * eshift_[2] + (p[2] + 1) * eshift_[3]
3109 + wrapZFront;
3110 case 2:
3111 return esetshift_[1] + p[0] + p[1] * eshift_[4] + p[2] * eshift_[5];
3112 case 3:
3113 return esetshift_[3] + p[0] + p[1] * eshift_[8] + p[2] * eshift_[9];
3114 case 4:
3115 return esetshift_[4] + p[0] + p[1] * eshift_[10] + p[2] * eshift_[11];
3116 case 5:
3117 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13];
3118 }
3119 return -1;
3120}
3121
3122inline ttk::SimplexId
3123 ttk::PeriodicImplicitTriangulation::getTetrahedronEdgeBEFG(
3124 const SimplexId p[3], const int id) const {
3125 SimplexId wrapXRight = 0;
3126 SimplexId wrapZFront = 0;
3127 if(p[0] == nbvoxels_[0])
3128 wrapXRight = -wrap_[0];
3129 if(p[2] == nbvoxels_[2])
3130 wrapZFront = -wrap_[2];
3131 switch(id) {
3132 case 0:
3133 return p[0] + p[1] * eshift_[0] + (p[2] + 1) * eshift_[1] + wrapZFront;
3134 case 1:
3135 return esetshift_[0] + p[0] + p[1] * eshift_[2] + (p[2] + 1) * eshift_[3]
3136 + wrapZFront;
3137 case 2:
3138 return esetshift_[1] + (p[0] + 1) + p[1] * eshift_[4] + p[2] * eshift_[5]
3139 + wrapXRight;
3140 case 3:
3141 return esetshift_[2] + p[0] + p[1] * eshift_[6] + (p[2] + 1) * eshift_[7]
3142 + wrapZFront;
3143 case 4:
3144 return esetshift_[4] + p[0] + p[1] * eshift_[10] + p[2] * eshift_[11];
3145 case 5:
3146 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13];
3147 }
3148 return -1;
3149}
3150
3151inline ttk::SimplexId
3152 ttk::PeriodicImplicitTriangulation::getTetrahedronEdgeBFGH(
3153 const SimplexId p[3], const int id) const {
3154 SimplexId wrapXRight = 0;
3155 SimplexId wrapYBottom = 0;
3156 SimplexId wrapZFront = 0;
3157 if(p[0] == nbvoxels_[0])
3158 wrapXRight = -wrap_[0];
3159 if(p[1] == nbvoxels_[1])
3160 wrapYBottom = -wrap_[1];
3161 if(p[2] == nbvoxels_[2])
3162 wrapZFront = -wrap_[2];
3163 switch(id) {
3164 case 0:
3165 return p[0] + (p[1] + 1) * eshift_[0] + (p[2] + 1) * eshift_[1]
3166 + wrapYBottom + wrapZFront;
3167 case 1:
3168 return esetshift_[0] + (p[0] + 1) + p[1] * eshift_[2]
3169 + (p[2] + 1) * eshift_[3] + wrapXRight + wrapZFront;
3170 case 2:
3171 return esetshift_[1] + (p[0] + 1) + p[1] * eshift_[4] + p[2] * eshift_[5]
3172 + wrapXRight;
3173 case 3:
3174 return esetshift_[2] + p[0] + p[1] * eshift_[6] + (p[2] + 1) * eshift_[7]
3175 + wrapZFront;
3176 case 4:
3177 return esetshift_[3] + (p[0] + 1) + p[1] * eshift_[8] + p[2] * eshift_[9]
3178 + wrapXRight;
3179 case 5:
3180 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13];
3181 }
3182 return -1;
3183}
3184
3185inline ttk::SimplexId
3186 ttk::PeriodicImplicitTriangulation::getTetrahedronEdgeBDGH(
3187 const SimplexId p[3], const int id) const {
3188 SimplexId wrapXRight = 0;
3189 SimplexId wrapYBottom = 0;
3190 SimplexId wrapZFront = 0;
3191 if(p[0] == nbvoxels_[0])
3192 wrapXRight = -wrap_[0];
3193 if(p[1] == nbvoxels_[1])
3194 wrapYBottom = -wrap_[1];
3195 if(p[2] == nbvoxels_[2])
3196 wrapZFront = -wrap_[2];
3197 switch(id) {
3198 case 0:
3199 return p[0] + (p[1] + 1) * eshift_[0] + (p[2] + 1) * eshift_[1]
3200 + wrapYBottom + wrapZFront;
3201 case 1:
3202 return esetshift_[0] + (p[0] + 1) + p[1] * eshift_[2] + p[2] * eshift_[3]
3203 + wrapXRight;
3204 case 2:
3205 return esetshift_[1] + (p[0] + 1) + (p[1] + 1) * eshift_[4]
3206 + p[2] * eshift_[5] + wrapXRight + wrapYBottom;
3207 case 3:
3208 return esetshift_[3] + (p[0] + 1) + p[1] * eshift_[8] + p[2] * eshift_[9]
3209 + wrapXRight;
3210 case 4:
3211 return esetshift_[4] + p[0] + (p[1] + 1) * eshift_[10]
3212 + p[2] * eshift_[11] + wrapYBottom;
3213 case 5:
3214 return esetshift_[5] + p[0] + p[1] * eshift_[12] + p[2] * eshift_[13];
3215 }
3216 return -1;
3217}
3218
3219inline ttk::SimplexId
3220 ttk::PeriodicImplicitTriangulation::getTetrahedronTriangleABCG(
3221 const SimplexId p[3], const int id) const {
3222 switch(id) {
3223 case 0:
3224 return p[0] * 2 + p[1] * tshift_[0] + p[2] * tshift_[1];
3225 case 1:
3226 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4] + p[2] * tshift_[5];
3227 case 2:
3228 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9];
3229 case 3:
3230 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10] + p[2] * tshift_[11];
3231 }
3232 return -1;
3233}
3234
3235inline ttk::SimplexId
3236 ttk::PeriodicImplicitTriangulation::getTetrahedronTriangleBCDG(
3237 const SimplexId p[3], const int id) const {
3238 switch(id) {
3239 case 0:
3240 return p[0] * 2 + p[1] * tshift_[0] + p[2] * tshift_[1] + 1;
3241 case 1:
3242 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10] + p[2] * tshift_[11];
3243 case 2:
3244 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6] + p[2] * tshift_[7];
3245 case 3:
3246 return (p[1] < nbvoxels_[1])
3247 ? tsetshift_[0] + p[0] * 2 + (p[1] + 1) * tshift_[2]
3248 + p[2] * tshift_[3]
3249 : tsetshift_[0] + p[0] * 2 + (p[1] + 1) * tshift_[2]
3250 + p[2] * tshift_[3] - wrap_[1] * 2;
3251 }
3252 return -1;
3253}
3254
3255inline ttk::SimplexId
3256 ttk::PeriodicImplicitTriangulation::getTetrahedronTriangleABEG(
3257 const SimplexId p[3], const int id) const {
3258 switch(id) {
3259 case 0:
3260 return tsetshift_[0] + p[0] * 2 + p[1] * tshift_[2] + p[2] * tshift_[3];
3261 case 1:
3262 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9];
3263 case 2:
3264 return tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4] + p[2] * tshift_[5]
3265 + 1;
3266 case 3:
3267 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6] + p[2] * tshift_[7]
3268 + 1;
3269 }
3270 return -1;
3271}
3272
3273inline ttk::SimplexId
3274 ttk::PeriodicImplicitTriangulation::getTetrahedronTriangleBEFG(
3275 const SimplexId p[3], const int id) const {
3276 switch(id) {
3277 case 0:
3278 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6] + p[2] * tshift_[7]
3279 + 1;
3280 case 1:
3281 return tsetshift_[0] + p[0] * 2 + p[1] * tshift_[2] + p[2] * tshift_[3]
3282 + 1;
3283 case 2:
3284 return (p[2] < nbvoxels_[2])
3285 ? p[0] * 2 + p[1] * tshift_[0] + (p[2] + 1) * tshift_[1]
3286 : p[0] * 2 + p[1] * tshift_[0] + (p[2] + 1) * tshift_[1]
3287 - wrap_[2] * 2;
3288 case 3:
3289 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10] + p[2] * tshift_[11]
3290 + 1;
3291 }
3292 return -1;
3293}
3294
3295inline ttk::SimplexId
3296 ttk::PeriodicImplicitTriangulation::getTetrahedronTriangleBFGH(
3297 const SimplexId p[3], const int id) const {
3298 switch(id) {
3299 case 0:
3300 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9]
3301 + 1;
3302 case 1:
3303 return tsetshift_[4] + p[0] * 2 + p[1] * tshift_[10] + p[2] * tshift_[11]
3304 + 1;
3305 case 2:
3306 return (p[0] < nbvoxels_[0])
3307 ? tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4]
3308 + p[2] * tshift_[5] + 3
3309 : tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4]
3310 + p[2] * tshift_[5] + 3 - wrap_[0] * 2;
3311 case 3:
3312 return (p[2] < nbvoxels_[2])
3313 ? p[0] * 2 + p[1] * tshift_[0] + (p[2] + 1) * tshift_[1] + 1
3314 : p[0] * 2 + p[1] * tshift_[0] + (p[2] + 1) * tshift_[1] + 1
3315 - wrap_[2] * 2;
3316 }
3317 return -1;
3318}
3319
3320inline ttk::SimplexId
3321 ttk::PeriodicImplicitTriangulation::getTetrahedronTriangleBDGH(
3322 const SimplexId p[3], const int id) const {
3323 switch(id) {
3324 case 0:
3325 return (p[0] < nbvoxels_[0])
3326 ? tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4]
3327 + p[2] * tshift_[5] + 2
3328 : tsetshift_[1] + p[0] * 2 + p[1] * tshift_[4]
3329 + p[2] * tshift_[5] + 2 - wrap_[0] * 2;
3330 case 1:
3331 return tsetshift_[2] + p[0] * 2 + p[1] * tshift_[6] + p[2] * tshift_[7];
3332 case 2:
3333 return tsetshift_[3] + p[0] * 2 + p[1] * tshift_[8] + p[2] * tshift_[9]
3334 + 1;
3335 case 3:
3336 return (p[1] < nbvoxels_[1])
3337 ? tsetshift_[0] + p[0] * 2 + (p[1] + 1) * tshift_[2]
3338 + p[2] * tshift_[3] + 1
3339 : tsetshift_[0] + p[0] * 2 + (p[1] + 1) * tshift_[2]
3340 + p[2] * tshift_[3] + 1 - wrap_[1] * 2;
3341 }
3342 return -1;
3343}
3344
3345inline ttk::SimplexId
3346 ttk::PeriodicImplicitTriangulation::getTetrahedronNeighborABCG(
3347 const SimplexId t, const SimplexId p[3], const int id) const {
3348 switch(id) {
3349 case 0:
3350 return t + 1;
3351 case 1:
3352 return t + 2;
3353 case 2:
3354 return p[0] > 0 ? t - 1 : t - 1 + wrap_[0] * 6;
3355 case 3:
3356 return p[2] > 0 ? t - tetshift_[1] + 3
3357 : t - tetshift_[1] + 3 + wrap_[2] * 6;
3358 }
3359 return -1;
3360}
3361
3362inline ttk::SimplexId
3363 ttk::PeriodicImplicitTriangulation::getTetrahedronNeighborBCDG(
3364 const SimplexId t, const SimplexId p[3], const int id) const {
3365 switch(id) {
3366 case 0:
3367 return t - 1;
3368 case 1:
3369 return t + 4;
3370 case 2:
3371 return p[2] > 0 ? t - tetshift_[1] + 3
3372 : t - tetshift_[1] + 3 + wrap_[2] * 6;
3373 case 3:
3374 return p[1] < nbvoxels_[1] ? t + tetshift_[0] + 1
3375 : t + tetshift_[0] + 1 - wrap_[1] * 6;
3376 }
3377 return -1;
3378}
3379
3380inline ttk::SimplexId
3381 ttk::PeriodicImplicitTriangulation::getTetrahedronNeighborABEG(
3382 const SimplexId t, const SimplexId p[3], const int id) const {
3383 switch(id) {
3384 case 0:
3385 return t - 2;
3386 case 1:
3387 return t + 1;
3388 case 2:
3389 return p[0] > 0 ? t - 4 : t - 4 + wrap_[0] * 6;
3390 case 3:
3391 return p[1] > 0 ? t - tetshift_[0] - 1
3392 : t - tetshift_[0] - 1 + wrap_[1] * 6;
3393 }
3394 return -1;
3395}
3396
3397inline ttk::SimplexId
3398 ttk::PeriodicImplicitTriangulation::getTetrahedronNeighborBEFG(
3399 const SimplexId t, const SimplexId p[3], const int id) const {
3400 switch(id) {
3401 case 0:
3402 return t - 1;
3403 case 1:
3404 return t + 1;
3405 case 2:
3406 return p[1] > 0 ? t - tetshift_[0] + 2
3407 : t - tetshift_[0] + 2 + wrap_[1] * 6;
3408 case 3:
3409 return p[2] < nbvoxels_[2] ? t + tetshift_[1] - 3
3410 : t + tetshift_[1] - 3 - wrap_[2] * 6;
3411 }
3412 return -1;
3413}
3414
3415inline ttk::SimplexId
3416 ttk::PeriodicImplicitTriangulation::getTetrahedronNeighborBFGH(
3417 const SimplexId t, const SimplexId p[3], const int id) const {
3418 switch(id) {
3419 case 0:
3420 return t - 1;
3421 case 1:
3422 return t + 1;
3423 case 2:
3424 return p[0] < nbvoxels_[0] ? t + 4 : t + 4 - wrap_[0] * 6;
3425 case 3:
3426 return p[2] < nbvoxels_[2] ? t + tetshift_[1] - 3
3427 : t + tetshift_[1] - 3 - wrap_[2] * 6;
3428 }
3429 return -1;
3430}
3431
3432inline ttk::SimplexId
3433 ttk::PeriodicImplicitTriangulation::getTetrahedronNeighborBDGH(
3434 const SimplexId t, const SimplexId p[3], const int id) const {
3435 switch(id) {
3436 case 0:
3437 return t - 1;
3438 case 1:
3439 return t - 4;
3440 case 2:
3441 return p[0] < nbvoxels_[0] ? t + 1 : t + 1 - wrap_[0] * 6;
3442 case 3:
3443 return p[1] < nbvoxels_[1] ? t + tetshift_[0] - 2
3444 : t + tetshift_[0] - 2 - wrap_[1] * 6;
3445 }
3446 return -1;
3447}
3448
3449#include <PeriodicPreconditions.h>
3450
#define TTK_TRIANGULATION_INTERNAL(NAME)
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition: BaseClass.h:47
virtual int getVertexPointInternal(const SimplexId &ttkNotUsed(vertexId), float &ttkNotUsed(x), float &ttkNotUsed(y), float &ttkNotUsed(z)) const
virtual int getCellVertexInternal(const SimplexId &ttkNotUsed(cellId), const int &ttkNotUsed(localVertexId), SimplexId &ttkNotUsed(vertexId)) const
virtual int getEdgeIncenter(SimplexId edgeId, float incenter[3]) const
SimplexId TTK_TRIANGULATION_INTERNAL() getEdgeStarNumber(const SimplexId &edgeId) const override
int TTK_TRIANGULATION_INTERNAL() getVertexStar(const SimplexId &vertexId, const int &localStarId, SimplexId &starId) const override
int getTetrahedronVertex(const SimplexId &tetId, const int &localVertexId, SimplexId &vertexId) const override
int getTetrahedronTriangle(const SimplexId &tetId, const int &id, SimplexId &triangleId) const override
int getTriangleVertexInternal(const SimplexId &triangleId, const int &localVertexId, SimplexId &vertexId) const override
int getEdgeVertexInternal(const SimplexId &edgeId, const int &localVertexId, SimplexId &vertexId) const override
virtual int getTetraIncenter(SimplexId tetraId, float incenter[3]) const
int TTK_TRIANGULATION_INTERNAL() getEdgeStar(const SimplexId &edgeId, const int &localStarId, SimplexId &starId) const override
int TTK_TRIANGULATION_INTERNAL() getVertexLink(const SimplexId &vertexId, const int &localLinkId, SimplexId &linkId) const override
int TTK_TRIANGULATION_INTERNAL() getVertexNeighbor(const SimplexId &vertexId, const int &localNeighborId, SimplexId &neighborId) const override
virtual int getTriangleIncenter(SimplexId triangleId, float incenter[3]) const
int getVertexTriangleInternal(const SimplexId &vertexId, const int &id, SimplexId &triangleId) const override
int getTriangleEdgeInternal(const SimplexId &triangleId, const int &id, SimplexId &edgeId) const override
int getTetrahedronNeighbor(const SimplexId &tetId, const int &localNeighborId, SimplexId &neighborId) const override
int getVertexEdgeInternal(const SimplexId &vertexId, const int &id, SimplexId &edgeId) const override
int TTK_TRIANGULATION_INTERNAL() getTriangleLink(const SimplexId &triangleId, const int &localLinkId, SimplexId &linkId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getTriangleStarNumber(const SimplexId &triangleId) const override
int TTK_TRIANGULATION_INTERNAL() getTriangleStar(const SimplexId &triangleId, const int &localStarId, SimplexId &starId) const override
int TTK_TRIANGULATION_INTERNAL() getVertexPoint(const SimplexId &vertexId, float &x, float &y, float &z) const override
SimplexId getEdgeTriangleNumberInternal(const SimplexId &edgeId) const override
int getEdgeTriangleInternal(const SimplexId &edgeId, const int &id, SimplexId &triangleId) const override
int TTK_TRIANGULATION_INTERNAL() getEdgeLink(const SimplexId &edgeId, const int &localLinkId, SimplexId &linkId) const override
int getTetrahedronEdge(const SimplexId &tetId, const int &id, SimplexId &edgeId) const override
int getTriangleNeighbor(const SimplexId &triangleId, const int &localNeighborId, SimplexId &neighborId) const override
TTK triangulation class for grids with periodic boundary conditions implemented in all directions.
virtual int preconditionTetrahedronsInternal()=0
SimplexId TTK_TRIANGULATION_INTERNAL() getTriangleLinkNumber(const SimplexId &triangleId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getEdgeLinks() override
int getTriangleNeighbors(std::vector< std::vector< SimplexId > > &neighbors)
virtual int getTetrahedronNeighbor(const SimplexId &tetId, const int &localNeighborId, SimplexId &neighborId) const =0
int TTK_TRIANGULATION_INTERNAL() getCellVertex(const SimplexId &cellId, const int &localVertexId, SimplexId &vertexId) const override
int getCellTriangleInternal(const SimplexId &cellId, const int &id, SimplexId &triangleId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getTriangleLinks() override
int setInputGrid(const float &xOrigin, const float &yOrigin, const float &zOrigin, const float &xSpacing, const float &ySpacing, const float &zSpacing, const SimplexId &xDim, const SimplexId &yDim, const SimplexId &zDim) override
const std::vector< std::vector< SimplexId > > * getTriangleEdgesInternal() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getTriangleStars() override
SimplexId TTK_TRIANGULATION_INTERNAL() getEdgeLinkNumber(const SimplexId &edgeId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexLinks() override
int getTetrahedronNeighbors(std::vector< std::vector< SimplexId > > &neighbors)
SimplexId getVertexEdgeNumberInternal(const SimplexId &vertexId) const override
SimplexId getCellEdgeNumberInternal(const SimplexId &cellId) const override
const std::vector< std::array< SimplexId, 2 > > *TTK_TRIANGULATION_INTERNAL() getEdges() override
bool TTK_TRIANGULATION_INTERNAL() isTriangleOnBoundary(const SimplexId &triangleId) const override
int TTK_TRIANGULATION_INTERNAL() getDimensionality() const override
int getTetrahedronEdges(std::vector< std::vector< SimplexId > > &edges) const
bool TTK_TRIANGULATION_INTERNAL() isEdgeOnBoundary(const SimplexId &edgeId) const override
SimplexId getTriangleNeighborNumber(const SimplexId &triangleId) const
const std::array< ttk::SimplexId, 3 > & getGridDimensions() const override
virtual int getTetrahedronVertex(const SimplexId &tetId, const int &localVertexId, SimplexId &vertexId) const =0
const std::vector< std::vector< SimplexId > > * getEdgeTrianglesInternal() override
SimplexId TTK_TRIANGULATION_INTERNAL() getNumberOfVertices() const override
int preconditionTrianglesInternal() override=0
SimplexId TTK_TRIANGULATION_INTERNAL() getCellVertexNumber(const SimplexId &cellId) const override
int preconditionEdgesInternal() override=0
virtual int getTriangleNeighbor(const SimplexId &triangleId, const int &localNeighborId, SimplexId &neighborId) const =0
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexStarNumber(const SimplexId &vertexId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getCellNeighborNumber(const SimplexId &cellId) const override
bool isPowerOfTwo(unsigned long long int v, unsigned long long int &r)
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexStars() override
PeriodicImplicitTriangulation(const PeriodicImplicitTriangulation &)=default
virtual int getTetrahedronEdge(const SimplexId &tetId, const int &id, SimplexId &edgeId) const =0
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getEdgeStars() override
virtual int getTetrahedronTriangle(const SimplexId &tetId, const int &id, SimplexId &triangleId) const =0
PeriodicImplicitTriangulation & operator=(const PeriodicImplicitTriangulation &)=default
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getCellNeighbors() override
const std::vector< std::vector< SimplexId > > * getVertexEdgesInternal() override
SimplexId getTriangleEdgeNumberInternal(const SimplexId &ttkNotUsed(triangleId)) const override
PeriodicImplicitTriangulation & operator=(PeriodicImplicitTriangulation &&)=default
SimplexId getCellTriangleNumberInternal(const SimplexId &ttkNotUsed(cellId)) const override
const std::vector< std::array< SimplexId, 3 > > *TTK_TRIANGULATION_INTERNAL() getTriangles() override
SimplexId getVertexTriangleNumberInternal(const SimplexId &vertexId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexLinkNumber(const SimplexId &vertexId) const override
int getTetrahedronTriangles(std::vector< std::vector< SimplexId > > &triangles) const
const std::vector< std::vector< SimplexId > > * getCellEdgesInternal() override
virtual int preconditionVerticesInternal()=0
int getCellEdgeInternal(const SimplexId &cellId, const int &id, SimplexId &edgeId) const override
SimplexId getNumberOfTrianglesInternal() const override
bool TTK_TRIANGULATION_INTERNAL() isVertexOnBoundary(const SimplexId &vertexId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexNeighbors() override
SimplexId getTetrahedronNeighborNumber(const SimplexId &tetId) const
SimplexId TTK_TRIANGULATION_INTERNAL() getNumberOfCells() const override
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexNeighborNumber(const SimplexId &vertexId) const override
int TTK_TRIANGULATION_INTERNAL() getCellNeighbor(const SimplexId &cellId, const int &localNeighborId, SimplexId &neighborId) const override
PeriodicImplicitTriangulation(PeriodicImplicitTriangulation &&)=default
const std::vector< std::vector< SimplexId > > * getCellTrianglesInternal() override
const std::vector< std::vector< SimplexId > > * getVertexTrianglesInternal() override
int getCellVTKIDInternal(const int &ttkId, int &vtkId) const override
RegularGridTriangulation is an abstract subclass of ttk::AbstractTriangulation that exposes a common ...
virtual void tetrahedronToPosition(const SimplexId tetrahedron, SimplexId p[3]) const =0
virtual void vertexToPosition(const SimplexId vertex, SimplexId p[3]) const =0
virtual void triangleToPosition2d(const SimplexId triangle, SimplexId p[2]) const =0
std::array< SimplexId, 3 > dimensions_
virtual void vertexToPosition2d(const SimplexId vertex, SimplexId p[2]) const =0
virtual void triangleToPosition(const SimplexId triangle, const int k, SimplexId p[3]) const =0
T distance(const T *p0, const T *p1, const int &dimension=3)
Definition: Geometry.cpp:344
The Topology ToolKit.
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22