1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736
|
// This file is part of ff3d - http://www.freefem.org/ff3d
// Copyright (C) 2001, 2002, 2003 Pascal Have
// This program is free software; you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation; either version 2, or (at your option)
// any later version.
// This program is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU General Public License for more details.
// You should have received a copy of the GNU General Public License
// along with this program; if not, write to the Free Software Foundation,
// Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
// $Id: triangulation.cpp,v 1.6 2007/05/20 23:02:46 delpinux Exp $
/*
T O D O T O D O T O D O T O D O T O D O T O D O T O D O T O D O T O D O
Ds que Find_P_in_elt sera en quadtree virer tous les 'last_find' d'acclration relative
Drcursivation de insertPoint et surtout scanNeigh
Optimiser la dtection d'une frontire dans le coloriage
Mettre les lignes de contraintes NULL quand les Quadtrees seront fonctionnels
T O D O T O D O T O D O T O D O T O D O T O D O T O D O T O D O T O D O
*/
#include <stack>
#include <fstream>
#include "triangulation.hpp"
#define ALTERNATE_INCIRCLE
#define DEBUG_FIND
//#define STEP_DETAILS
#define STEP_WRITE
#define STANDALONE
#ifdef STANDALONE
template<class FLT>
inline FLT sqr(const FLT &x) { return x * x; }
#else /* STANDALONE */
//#include "timer.hpp"
//#include "communicator.hpp"
//#include "local_math.hpp"
//#include "utils.hpp"
#endif /* STANDALONE */
#include <cstdlib>
int myrand(const int width) {
return (int)(static_cast<double>(width)*rand()/(RAND_MAX+1.0));
}
#include <map>
#include <set>
/************************************************************/
/******************** Fonctions locales *********************/
/************************************************************/
//! Produit vectoriel spcialis 2D
inline double SDet(const Vertex & V1, const Vertex & V2, const Vertex & V3) {
return (V2[0]-V1[0]) * (V3[1]-V1[1]) - (V2[1]-V1[1]) * (V3[0]-V1[0]);
}
// inline double SDet(const Vertex & V1, const Vertex & V2, const Vertex & V3) {
// Triangle tr(&V1,&V2,&V3);
// return Vertex((tr[1]-tr[0]) ^ (tr[2]-tr[0]))[2];
// }
inline bool checkInter(const Vertex &P1, const Vertex &P2, const Vertex &Q1, const Vertex &Q2) {
// retourne true si P1P2 inter Q1Q2 != vide
ASSERT(!(P1 == P2 || P1 == Q1 || P1 == Q2 || P2 == Q1 || P2 == Q2 || Q1 == Q2));
return (SDet(P1,P2,Q1)*SDet(P1,P2,Q2) < 0) && (SDet(Q1,Q2,P1)*SDet(Q1,Q2,P2) < 0); // version optimise
}
inline double Angle(const Vertex &V1, const Vertex &V2, const Vertex &V3) { // cos de l'angle * ||V2-V1|| (V1 point pivot)
const double distance = sqrt(sqr(V3[0]-V1[0])+sqr(V3[1]-V1[1]));
return (((V2[0]-V1[0])*(V3[0]-V1[0]))+((V2[1]-V1[1])*(V3[1]-V1[1])))/distance;
}
/************************************************************/
/********************** Triangulation ***********************/
/************************************************************/
void Triangulation::setBox(const PointNumType xmin, const PointNumType ymin,
const PointNumType xmax, const PointNumType ymax) {
__triangles.clear(); // pour des bases saines
__corners[0] = Vertex(xmin,ymax,0);
__corners[1] = Vertex(xmin,ymin,0);
__corners[2] = Vertex(xmax,ymin,0);
__corners[3] = Vertex(xmax,ymax,0);
__triangles.push_back(ConnectedTriangle(&__corners[0],&__corners[2],&__corners[3]));
__triangles.push_back(ConnectedTriangle(&__corners[0],&__corners[1],&__corners[2]));
__triangles.front().neigh(2) = (&__triangles.back());
__triangles.back().neigh(1) = (&__triangles.front());
#ifdef STEP_WRITE
{ ASSERT(__checkAll()); std::ofstream o("init.mesh"); export_mesh(o); }
#endif /* STEP_WRITE */
}
#ifndef ALTERNATE_INCIRCLE
bool Triangulation::isInCircle(const Vertex & position, const Triangle & tri) const {
double x1,x2,x3,y1,y2,y3;
double a,b,c,d,e,f;
double centre_cercle_x,centre_cercle_y;
double rayon_cercle,distance;
double tolerance = 1.e-13;
x1=tri[0][0];
x2=tri[1][0];
x3=tri[2][0];
y1=tri[0][1];
y2=tri[1][1];
y3=tri[2][1];
a=x2-x3;
b=y2-y3;
c=((y2+y3)*(y2-y3)+(x2+x3)*(x2-x3))/2.;
d=x2-x1;
e=y2-y1;
f=((y2+y1)*(y2-y1)+(x2+x1)*(x2-x1))/2.;
centre_cercle_y = (d*c-a*f)/(b*d-a*e);
centre_cercle_x = (c*e-f*b)/(a*e-b*d);
//printf("centre du cercle :%g %g\n",(d*c-a*f)/(b*d-a*e),(c*e-f*b)/(a*e-b*d));
distance = sqr(centre_cercle_x-position[0]) + sqr(centre_cercle_y-position[1]);
rayon_cercle = sqr(centre_cercle_x-x1) + sqr(centre_cercle_y-y1);
return (distance<=(rayon_cercle+tolerance));
}
#else /* ALTERNATE_INCIRCLE */
bool Triangulation::__isInCircle(const Vertex & position, const Triangle & tri) const {
const double eps = 1e-13;
// On peut viter quelques indirections par des variables locales
/* 9 affectations : 3+3+3*3 (+) et 0+0+3*2 (*) */
const double
BxmAx=tri(1)[0]-tri(0)[0],
CxmAx=tri(2)[0]-tri(0)[0],
DxmAx=position[0]-tri(0)[0];
const double
BymAy=tri(1)[1]-tri(0)[1],
CymAy=tri(2)[1]-tri(0)[1],
DymAy=position[1]-tri(0)[1];
const double
BA2=BxmAx*(tri(1)[0]+tri(0)[0])+BymAy*(tri(1)[1]+tri(0)[1]),
CA2=CxmAx*(tri(2)[0]+tri(0)[0])+CymAy*(tri(2)[1]+tri(0)[1]),
DA2=DxmAx*(position[0]+tri(0)[0])+DymAy*(position[1]+tri(0)[1]);
const double r= // 5 (+) , 3*3 (*)
(BxmAx*CymAy-CxmAx*BymAy)*DA2+(CxmAx*DymAy-DxmAx*CymAy)*BA2+(DxmAx*BymAy-BxmAx*DymAy)*CA2;
// Total 20 (+) et 12 (*)
return (r<eps); // ou <=0
}
#endif /* ALTERNATE_INCIRCLE */
void Triangulation::__permutation(TriangleIndex num_tri1, const unsigned k1) {
// realise la permutation d'une arte commune a deux triangles.
// TriangleIndex = pointeur du triangle.
// k1 = numero du triangle voisin.
// S = numero du sommet
// V = numero des triangles voisins.
// S1 S1
// * *
// V2 * * V1 V2 *** V1
// * * * * *
// * T1 * * * *
// S2*********S4 S2 * T1*T2 *S4
// * T2 * * * *
// * * * * *
// V3 * * V4 V3 *** V4
// * *
// S3 S3
// => modification de V1 qui deivient T2 et V3 qui devient T1
const TriangleIndex num_tri2 = num_tri1->neigh(k1);
ASSERT(num_tri2 != NULL);
const unsigned k2 = num_tri2->localizeNeigh(num_tri1);
const PointIndex
S1=num_tri1->base()[k1],
S2=num_tri1->base()[(k1+1)%3],
S4=num_tri1->base()[(k1+2)%3];
const unsigned tag01 = num_tri1->tag;
const unsigned tag02 = num_tri2->tag;
// say() << *num_tri1 << " " << k1 << " " << tag01 << " " << (tag01 & HSide[k1]) << "\n";
// say() << *num_tri2 << " " << k2 << " " << tag02 << " " << (tag02 & HSide[k2]) << "\n";
ASSERT((tag01 & HSide[k1]) == 0); // Ne flip pas des artes marques
ASSERT((tag02 & HSide[k2]) == 0);
num_tri1->tag = ((tag01&HSide[(k1+2)%3])?HSide[2]:0)|((tag02&HSide[(k2+1)%3])?HSide[0]:0);
num_tri2->tag = ((tag01&HSide[(k1+1)%3])?HSide[1]:0)|((tag02&HSide[(k2+2)%3])?HSide[0]:0);
const TriangleIndex
V1=num_tri1->neigh((k1+1)%3),
V2=num_tri1->neigh((k1+2)%3);
PointIndex S3;
TriangleIndex V3,V4;
// recherche du sommet dans le triangle num_tri2 oppos au cot commun (ie oppos k1 dans num_tri1)
for(unsigned k=0;k<3;++k)
if(num_tri2->neigh(k) == num_tri1) {
S3=num_tri2->base()[k];
V3=num_tri2->neigh((k+1)%3); // Ici, un bug est mort dans d'attroces souffrances
V4=num_tri2->neigh((k+2)%3);
break;
}
for(unsigned k=0;k<3;++k) {
if(V1 != NULL)
if(V1->neigh(k) == num_tri1)
V1->neigh(k) = num_tri2;
if(V3 != NULL)
if(V3->neigh(k) == num_tri2)
V3->neigh(k) = num_tri1;
}
num_tri1->setVertices(S1,S2,S3);
num_tri2->setVertices(S1,S3,S4);
num_tri1->setNeighs(V3,num_tri2,V2);
num_tri2->setNeighs(V4,V1,num_tri1);
}
void Triangulation::__elimination_tri(TriangleIndex num_tri) {
// elimination fu triangle pointe par neigh_index
for(unsigned k=0;k<3;++k) {
TriangleIndex neigh = num_tri->neigh(k);
if(neigh != NULL) {
for(unsigned k2=0;k2<3;++k2)
if(neigh->neigh(k2) == num_tri) {
neigh->neigh(k2) = NULL;
break;
}
}
num_tri->neigh(k) = NULL;
}
// num_tri->tag = -1;
__deleteTriangle(num_tri);
}
// recherche du triangle comportant l'element P
TriangleIndex Triangulation::__find_P_in_elt(const Vertex & P) {
#ifdef DEBUG_FIND
// Version non optimise de debug (en cas de trous de maillages suppos convexe)
for(Triangles::iterator i = __triangles.begin();true;++i) {
ASSERT(i != __triangles.end());
const Vertex & P1 = (*i)(0);
const Vertex & P2 = (*i)(1);
const Vertex & P3 = (*i)(2);
ASSERT(SDet(P1,P2,P3)>=0);
// ffout(4)<<P1<<"\n";
// ffout(4)<<P2<<"\n";
// ffout(4)<<P3<<"\n";
// ffout(4)<<"dans triangles?\n";
// ffout(4)<<SDet(P1,P2,P)<<" "<<SDet(P2,P3,P)<<" "<<SDet(P3,P1,P)<<"\n";
if ( (SDet(P1,P2,P)>=-1e-15) && (SDet(P2,P3,P)>=-1e-15) && (SDet(P3,P1,P)>=-1e-15) ) {
return &*i;
}
}
#else /* DEBUG_FIND */
// static TriangleIndex lastFind = NULL;
// if (lastFind == NULL) lastFind = &*__triangles.begin();
for(TriangleIndex i = &*__triangles.begin();;) {
ASSERT(i != NULL);
const Vertex & P2 = (*i)(1);
const Vertex & P3 = (*i)(2);
if(SDet(P2,P3,P)<0) {
i=i->neigh(0);
} else {
const Vertex & P1 = (*i)(0); // on en a pas besoin avant
if(SDet(P3,P1,P)<0) {
i=i->neigh(1);
} else {
if(SDet(P1,P2,P)<0) {
i=i->neigh(2);
} else {
// lastFind = i; // Li Find_P_in_elt
return i;
}
}
}
}
#endif /* DEBUG_FIND */
}
void Triangulation::__insertPoint(Vertex * newPoint, const TriangleIndex & curTriangle) {
// triangle forcant cette position; vite les problmes dus l'ffacement de triangles
__newTriangle(ConnectedTriangle(&__corners[0],&__corners[0],&__corners[0]));
const Triangles::iterator first = --__triangles.end();
Vertex * last = curTriangle->base()[1];
for(unsigned i=0;i<3;++i) {
if (__scanNeigh(curTriangle, i, newPoint, last) ) {
// au bord de la bulle
TriangleIndex tmpNeigh = curTriangle->neigh(i);
__newTriangle(ConnectedTriangle(last,curTriangle->base()[(i+2)%3],newPoint,NULL,NULL,tmpNeigh));
if (tmpNeigh != NULL)
tmpNeigh->neigh(tmpNeigh->localizeNeigh(curTriangle)) = &__triangles.back();
last = curTriangle->base()[(i+2)%3];
}
}
__deleteTriangle(curTriangle);
Triangles::iterator
LTi = first,
LTf = __triangles.end(),
LTaux;
++LTi;
LTi->neigh(1) = &*(--LTf); // Cas particulier
LTf->neigh(0) = &*LTi;
LTaux = LTi;
++LTaux;
while(LTi != LTf) { // ici LTf = --triangle.end()
LTi->neigh(0) = &*LTaux;
LTaux->neigh(1) = &*LTi;
++LTi; ++LTaux;
}
__triangles.erase(first); // plus utile
}
bool Triangulation::__scanNeigh(const TriangleIndex & start,
const unsigned in,
Vertex * newPoint,
Vertex * &last) {
// Scan partir du triangle Start par le voisin in, la recherche de triangles voisins contenant dans
// leur cercle circonscrit le point newPoint; si OK met jour la liste L entre LP1 et LP2
const TriangleIndex curNeigh = start->neigh(in);
if (curNeigh != NULL && __isInCircle(*newPoint,*curNeigh)) {
const unsigned k = curNeigh->localizeNeigh(&*start); // sommet oppos au cot d'entre
if (__scanNeigh(curNeigh, (k+1)%3, newPoint, last) ) {
// au bord de la bulle
TriangleIndex tmpNeigh = curNeigh->neigh((k+1)%3);
__newTriangle(ConnectedTriangle(last,curNeigh->base()[k],newPoint,NULL,NULL,tmpNeigh));
if (tmpNeigh != NULL)
tmpNeigh->neigh(tmpNeigh->localizeNeigh(curNeigh)) = &__triangles.back();
last = curNeigh->base()[k];
}
// Idem que prcdemment, mais pour pas de boucle, code dupliqu
if (__scanNeigh(curNeigh, (k+2)%3, newPoint, last)) {
// au bord de la bulle
TriangleIndex tmpNeigh = curNeigh->neigh((k+2)%3);
__newTriangle(ConnectedTriangle(last,curNeigh->base()[(k+1)%3],newPoint,NULL,NULL,tmpNeigh));
if (tmpNeigh != NULL)
tmpNeigh->neigh(tmpNeigh->localizeNeigh(curNeigh)) = &__triangles.back();
last = curNeigh->base()[(k+1)%3];
}
__deleteTriangle(curNeigh);
return false;
}
return true;
}
void Triangulation::__resetTags() {
for(Triangles::iterator i = __triangles.begin(); i != __triangles.end(); ++i)
i->tag = 0;
}
bool Triangulation::__checkLine(const CurveVertex &L, const bool closed) {
ASSERT(L.size() > 1);
ASSERT(closed); // !closed not implemented
ASSERT(!closed || L.front() == L.back());
typedef TinyVector<2,PointIndex> Edge;
typedef std::map<Edge,TinyVector<2,TriangleIndex> > EdgeMapping;
typedef std::set<Edge> BorderEdges;
BorderEdges borderEdges;
PointIndex lastPoint = L.front();
for(CurveVertex::const_iterator i=L.begin(); ++i != L.end();) {
Edge edge;
edge[0]=lastPoint;
edge[1]=*i;
if(edge[0] > edge[1]) std::swap(edge[0],edge[1]); // arithmtique des pointeurs ou entiers
borderEdges.insert(edge);
lastPoint = *i;
}
EdgeMapping edgeMapping;
do {
edgeMapping.clear();
for(Triangles::iterator tr = __triangles.begin();tr !=__triangles.end(); ++tr) {
for(unsigned j=0; j<3; ++j) {
PointIndex p1 = tr->base()[(j+1)%3];
PointIndex p2 = tr->base()[(j+2)%3];
if(p1 > p2) std::swap(p1,p2); // arithmtique des pointeurs ou des entiers
Edge edge;
edge[0]=p1;
edge[1]=p2;
const TriangleIndex tr2 = tr->neigh(j);
if ((tr2 != NULL) && !(tr->tag & HSide[j]) && edgeMapping.find(edge) == edgeMapping.end()) {
Edge edgetemp;
edgetemp[0]=p1;
edgetemp[1]=p2;
TinyVector<2,TriangleIndex> trtemp;
trtemp[0]=&*tr;
trtemp[1]=tr2;
edgeMapping.insert(EdgeMapping::value_type(edgetemp,trtemp));
}
}
}
for(BorderEdges::const_iterator i = borderEdges.begin(); i != borderEdges.end();) {
EdgeMapping::iterator foundi = edgeMapping.find(*i);
if (foundi != edgeMapping.end()) {
// say() << "Edge " << *(*i)[0] << " " << *(*i)[1] << " found\n";
TriangleIndex tr0 = foundi->second[0];
TriangleIndex tr1 = foundi->second[1];
ASSERT(tr0 != NULL && tr1 != NULL);
tr0->tag |= HSide[tr0->localizeNeigh(tr1)];
tr1->tag |= HSide[tr1->localizeNeigh(tr0)];
edgeMapping.erase(foundi);
BorderEdges::const_iterator i_to_rm = i++;
borderEdges.erase(i_to_rm);
} else ++i;
}
if (!borderEdges.empty()) {
bool done = false;
while (!done) {
EdgeMapping::iterator e = edgeMapping.begin();
for(unsigned n = myrand(edgeMapping.size());n>0;--n) ++e;
const PointIndex P1 = e->first[0];
const PointIndex P2 = e->first[1];
TriangleIndex tr0 = e->second[0];
TriangleIndex tr1 = e->second[1];
ASSERT(tr0 != NULL && tr1 != NULL);
const unsigned l0 = tr0->localizeNeigh(tr1);
const unsigned l1 = tr1->localizeNeigh(tr0);
const PointIndex Q1 = tr0->base()[l0];
const PointIndex Q2 = tr1->base()[l1];
if (checkInter(*P1,*P2,*Q1,*Q2)) {
// say() << "Flip edge " << *P1 << " " << *P2 << "\n";
done = true;
__permutation(tr0,tr0->localizeNeigh(tr1));
}
}
}
} while (!borderEdges.empty());
return false;
}
bool Triangulation::__colorize(const Points & points) {
// les variables pour les diffrents parcours
std::stack<TriangleIndex> Stack;
// On parcours les triangles pour propager la couleur de l'extrieur : 1
// Recherche de l'lment initial : ie un triangle de coin de la box englobante
Stack.push(__find_P_in_elt(__corners[0]));
while(!Stack.empty()) {
TriangleIndex itr = Stack.top(); Stack.pop();
if (!(itr->tag & HMask)) { // pas de couleur
itr->tag |= 1; // 1 : Couleur de remplissage
// parcours des voisins (il faudrait peut tre pas de boucle pour tre plus rapide)
for(unsigned j=0;j<3;j++) {
TriangleIndex ktr = itr->neigh(j);
if (ktr != NULL && // hors domaine
!(itr->tag & HSide[j]) &&
!(ktr->tag & HMask) // Dj colori
) {
Stack.push(ktr);
}
}
}
}
Triangles::iterator llf = __triangles.begin();
const Triangles::iterator lle = __triangles.end();
while(llf != lle && !((llf->tag & HUMask) &&
!(llf->tag & HMask) &&
(llf->neigh(0)->tag & HMask ||
llf->neigh(1)->tag & HMask ||
llf->neigh(2)->tag & HMask))) {
++llf;
}
if (llf == lle) {
#ifndef STEP_DETAILS
ffout(4) << "Bad Domain\n"; // c'est pour cela qu'il faut commencer par colorier l'extrieur
#endif /* STEP_DETAILS */
return false;
}
// On parcours les triangles pour propager une couleur
// Initialisation
Stack.push(&*llf);
while(!Stack.empty()) {
TriangleIndex itr = Stack.top();
Stack.pop();
if (!(itr->tag & HMask)) { // pas de couleur
itr->tag |= 2; // colorisation du domaine avec la couleur 2
// parcours des voisins (il faudrait peut tre pas de boucle pour tre plus rapide)
for(unsigned j=0;j<3;j++) {
TriangleIndex ktr = itr->neigh(j);
if (ktr != NULL &&
!(itr->tag & HSide[j]) &&
!(ktr->tag & HMask)) {
Stack.push(ktr);
}
}
}
}
return true;
}
TriangleIndex Triangulation::__findVertex(const Vertex * const p) {
for(Triangles::iterator i = __triangles.begin(); i != __triangles.end(); ++i) {
for(unsigned j=0;j<3;++j)
if (i->base()[j] == p)
return &*i;
}
ASSERT(false);
return NULL;
}
bool Triangulation::triangulize(Points & points,
const CurveVertex & envVertex,
const std::vector<CurveVertex> & holes,
const std::vector<CurveVertex> & curves) {
ASSERT(points.size() > 0);
__points = &points;
const Vertex & firstPoint = points[0];
PointNumType
xmin=firstPoint[0],xmax=firstPoint[0],
ymin=firstPoint[1],ymax=firstPoint[1];
for(unsigned i=0;i<points.size();++i) {
const PointNumType x = points[i][0];
const PointNumType y = points[i][1];
xmin = std::min(xmin,x);
xmax = std::max(xmax,x);
ymin = std::min(ymin,y);
ymax = std::max(ymax,y);
}
PointNumType
dx=(xmax-xmin)/10,
dy=(ymax-ymin)/10;
setBox(xmin-2.*dx,ymin-dy,xmax+2.*dx,ymax+dy);
#ifdef STEP_DETAILS
ffout(4) << "Box : (" << xmin-dx << " , " << ymin-dy << " ) , ( " << xmax+dx << " , " << ymax+dy << " )" << "\n";
#endif /* STEP_DETAILS */
#ifdef STEP_DETAILS
//Timer timer(true);
#endif /* STEP_DETAILS */
for(unsigned i=0;i<points.size();++i) {
ASSERT(__checkAll());
__findInsert(points[i]);
}
#ifdef STEP_DETAILS
// ffout(4) << "Insertion " << points.size() << " points [" << __triangles.size() << " triangles] in " << timer.//stop() << "_s" << "\n";
#endif /* STEP_DETAILS */
#ifdef STEP_WRITE
{ ASSERT(__checkAll()); std::ofstream o("insert.mesh"); export_mesh(o); }
#endif /* STEP_WRITE */
// Important car c'est les CheckLine qui marque les cots frontires dans les triangles
__resetTags();
#ifdef STEP_DETAILS
ffout(4) << "Check Enveloppe" << "\n";
for(unsigned i=1;checkLine(envVertex,true);)
{
ffout(4) << "\tPasse " << ++i << "\n";
ASSERT(checkAll()); std::ofstream o("passe.mesh"); export_mesh(o);
}
#else /* STEP_DETAILS */
while(__checkLine(envVertex,true));
#endif /* STEP_DETAILS */
//#ifdef STEP_DETAILS
// ffout(4) << "Check Trou(s)" << "\n";
// for(std::vector<CurveVertex>::iterator lf = holes.begin(); lf != holes.end(); ++lf)
// for(unsigned i=1;checkLine(*lf,true);)
// ffout(4) << "\tPasse " << ++i << "\n";
//#else /* STEP_DETAILS */
for(std::vector<CurveVertex>::const_iterator lf = holes.begin(); lf != holes.end(); ++lf)
while(__checkLine(*lf,true));
//#endif /* STEP_DETAILS */
#ifdef STEP_DETAILS
ffout(4) << "Check Courbe(s)" << "\n";
for(std::vector<CurveVertex>::const_iterator lf = curves.begin(); lf != curves.end(); ++lf)
for(unsigned i=1;checkLine(*lf,false);)
ffout(4) << "\tPasse " << ++i << "\n";
#else /* STEP_DETAILS */
for(std::vector<CurveVertex>::const_iterator lf = curves.begin(); lf != curves.end(); ++lf)
while(__checkLine(*lf,false));
#endif /* STEP_DETAILS */
#ifdef STEP_WRITE
{ ASSERT(__checkAll()); std::ofstream o("check.mesh"); export_mesh(o); }
#endif /* STEP_WRITE */
#ifdef STEP_DETAILS
ffout(4) << "Colorize" << "\n";
#endif /* STEP_DETAILS */
__colorize(points);
#ifdef STEP_WRITE
{ ASSERT(__checkAll()); std::ofstream o("color.mesh"); export_mesh(o); }
#endif /* STEP_WRITE */
unsigned count_tri = 0;
const Triangles::iterator e = __triangles.end();
for(Triangles::iterator ff,f = __triangles.begin(); f != e;) {
ff = f++;
if ((ff->tag & HMask) != 2) {
__elimination_tri(&*ff);
++count_tri;
}
}
#ifdef STEP_DETAILS
ffout(4) << count_tri << " Triangles Deleted " << __triangles.size() << " triangles left" << "\n";
#endif /* STEP_DETAILS */
#ifdef STEP_WRITE
{ ASSERT(__checkAll()); std::ofstream o("final.mesh"); export_mesh(o); }
#endif /* STEP_WRITE */
return true;
}
void Triangulation::export_mesh(std::ostream &o) const {
const Points & points = *__points;
o << "MeshVersionFormatted 1\n"
"Dimension\n"
"\t2\n"
"Vertices\n"
"\t" << points.size()+4
<< "\n";
for(unsigned i=0;i<4;++i)
o << __corners[i][0] << " " << __corners[i][1] << " 1\n";
for(Points::const_iterator i = points.begin(); i != points.end(); ++i) {
o << i->x() << " " << i->y() << " 1\n";
}
class BestRef {
private:
public:
const Vertex *points, *__corners;
public:
// BestRef(const Points * _points, const Point * ___corners)
// : points(_points),
// __corners(___corners) { }
unsigned operator()(const Vertex & p) {
const unsigned n = &p-__corners;
if (n < 4)
return n+1;
else
return &p-points+5;
}
} ref; // (firstPointer(points),firstPointer(__corners));
ref.points = &points[0];
ref.__corners = __corners;
o << "Triangles\n\t" << __triangles.size() << "\n";
for(Triangles::const_iterator i = __triangles.begin(); i != __triangles.end(); ++i) {
o << " " << ref((*i)(0))
<< " " << ref((*i)(1))
<< " " << ref((*i)(2)) << " " << (i->tag & HMask) << "\n";
}
}
bool Triangulation::__checkAll() {
bool ok = true;
unsigned jj = 0;
for(Triangles::iterator i = __triangles.begin(); i != __triangles.end(); ++i,++jj) {
// Check Orientation
const double S = SDet((*i)(0),(*i)(1),(*i)(2));
const bool check_point =
((i->base())[0] != (i->base())[1]) &&
((i->base())[0] != (i->base())[2]) &&
((i->base())[1] != (i->base())[2]);
if (!check_point) {
ffout(4) << "Warning : Triangle " << jj << " Plat : " /*<< *i */<< "\n";
ok = false;
}
ASSERT(check_point);
if (S <= 0) {
ffout(4) << "< triangle " << jj << " " /*<< *i */<< " : " << S << "\n";
ok = false;
}
// Check Neigh
for(unsigned j=0;j<3;j++) {
const TriangleIndex k = i->neigh(j);
if (k != NULL) {
const unsigned l = k->localizeNeigh(&*i);
if (!(k->base()[(l+2)%3] == i->base()[(j+1)%3] &&
k->base()[(l+1)%3] == i->base()[(j+2)%3])) {
ffout(4) << "Bad edge " << j << " for triangle : " /*<< *i */<< "\n";
ok = false;
}
}
}
}
return ok;
}
|