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
|
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#ifndef itkFastMarchingQuadEdgeMeshFilterBase_hxx
#define itkFastMarchingQuadEdgeMeshFilterBase_hxx
#include "itkMath.h"
#include "itkVector.h"
#include "itkMatrix.h"
namespace itk
{
template <typename TInput, typename TOutput>
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::FastMarchingQuadEdgeMeshFilterBase()
: Superclass()
{
this->m_InputMesh = nullptr;
}
template <typename TInput, typename TOutput>
IdentifierType
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::GetTotalNumberOfNodes() const
{
return this->GetInput()->GetNumberOfPoints();
}
template <typename TInput, typename TOutput>
void
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::SetOutputValue(OutputMeshType * oMesh,
const NodeType & iNode,
const OutputPixelType & iValue)
{
oMesh->SetPointData(iNode, iValue);
}
template <typename TInput, typename TOutput>
const typename FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::OutputPixelType
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::GetOutputValue(OutputMeshType * oMesh,
const NodeType & iNode) const
{
OutputPixelType outputValue{};
oMesh->GetPointData(iNode, &outputValue);
return outputValue;
}
template <typename TInput, typename TOutput>
unsigned char
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::GetLabelValueForGivenNode(const NodeType & iNode) const
{
auto it = m_Label.find(iNode);
if (it != m_Label.end())
{
return it->second;
}
else
{
return Traits::Far;
}
}
template <typename TInput, typename TOutput>
void
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::SetLabelValueForGivenNode(const NodeType & iNode,
const LabelType & iLabel)
{
m_Label[iNode] = iLabel;
}
template <typename TInput, typename TOutput>
void
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::UpdateNeighbors(OutputMeshType * oMesh, const NodeType & iNode)
{
OutputPointType p;
oMesh->GetPoint(iNode, &p);
OutputQEType * qe = p.GetEdge();
if (qe)
{
OutputQEType * qe_it = qe;
this->m_InputMesh = this->GetInput();
do
{
if (qe_it)
{
OutputPointIdentifierType neigh_id = qe_it->GetDestination();
const char label = this->GetLabelValueForGivenNode(neigh_id);
if ((label != Traits::Alive) && (label != Traits::InitialTrial) && (label != Traits::Forbidden))
{
this->UpdateValue(oMesh, neigh_id);
}
}
else
{
itkGenericExceptionMacro("qe_it is nullptr");
}
qe_it = qe_it->GetOnext();
} while (qe_it != qe);
}
else
{
itkGenericExceptionMacro("qe is nullptr");
}
}
template <typename TInput, typename TOutput>
void
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::UpdateValue(OutputMeshType * oMesh, const NodeType & iNode)
{
OutputPointType p;
oMesh->GetPoint(iNode, &p);
InputPixelType F{};
this->m_InputMesh->GetPointData(iNode, &F);
if (F < 0.)
{
itkGenericExceptionMacro("F < 0.");
}
OutputQEType * qe = p.GetEdge();
OutputPixelType outputPixel = this->m_LargeValue;
if (qe)
{
OutputQEType * qe_it = qe;
qe_it = qe_it->GetOnext();
do
{
OutputQEType * qe_it2 = qe_it->GetOnext();
if (qe_it2)
{
if (qe_it->GetLeft() != OutputMeshType::m_NoFace)
{
OutputPointIdentifierType id1 = qe_it->GetDestination();
OutputPointIdentifierType id2 = qe_it2->GetDestination();
const auto label1 = static_cast<LabelType>(this->GetLabelValueForGivenNode(id1));
const auto label2 = static_cast<LabelType>(this->GetLabelValueForGivenNode(id2));
bool IsFar1 = (label1 != Traits::Far);
bool IsFar2 = (label2 != Traits::Far);
if (IsFar1 || IsFar2)
{
OutputPointType q1 = oMesh->GetPoint(id1);
OutputPointType q2 = oMesh->GetPoint(id2);
auto val1 = static_cast<OutputVectorRealType>(this->GetOutputValue(oMesh, id1));
auto val2 = static_cast<OutputVectorRealType>(this->GetOutputValue(oMesh, id2));
if (val1 > val2)
{
OutputPointType temp_pt(q1);
q1 = q2;
q2 = temp_pt;
std::swap(val1, val2);
std::swap(IsFar1, IsFar2);
std::swap(id1, id2);
}
const OutputVectorRealType temp =
this->Solve(oMesh, iNode, p, F, id1, q1, IsFar1, val1, id2, q2, IsFar2, val2);
outputPixel = std::min(outputPixel, static_cast<OutputPixelType>(temp));
}
}
qe_it = qe_it2;
}
else
{
// throw one exception here
itkGenericExceptionMacro("qe_it2 is nullptr");
}
} while (qe_it != qe);
if (outputPixel < this->m_LargeValue)
{
this->SetOutputValue(oMesh, iNode, outputPixel);
this->SetLabelValueForGivenNode(iNode, Traits::Trial);
this->m_Heap.push(NodePairType(iNode, outputPixel));
}
}
else
{
// throw one exception
itkGenericExceptionMacro("qe_it is nullptr");
}
}
template <typename TInput, typename TOutput>
const typename FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::OutputVectorRealType
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::Solve(OutputMeshType * oMesh,
const NodeType & iId,
const OutputPointType & iCurrentPoint,
const OutputVectorRealType & iF,
const NodeType & iId1,
const OutputPointType & iP1,
const bool iIsFar1,
const OutputVectorRealType iVal1,
const NodeType & iId2,
const OutputPointType & iP2,
const bool iIsFar2,
const OutputVectorRealType & iVal2) const
{
OutputVectorType Edge1 = iP1 - iCurrentPoint;
OutputVectorType Edge2 = iP2 - iCurrentPoint;
OutputVectorRealType sq_norm1 = Edge1.GetSquaredNorm();
OutputVectorRealType norm1 = 0.;
OutputVectorRealType epsilon = NumericTraits<OutputVectorRealType>::epsilon();
if (sq_norm1 > epsilon)
{
norm1 = std::sqrt(sq_norm1);
OutputVectorRealType inv_norm1 = 1. / norm1;
Edge1 *= inv_norm1;
}
OutputVectorRealType sq_norm2 = Edge2.GetSquaredNorm();
OutputVectorRealType norm2 = 0.;
if (sq_norm2 > epsilon)
{
norm2 = std::sqrt(sq_norm2);
OutputVectorRealType inv_norm2 = 1. / norm2;
Edge2 *= inv_norm2;
}
if (!iIsFar1 && iIsFar2)
{
// only one point is a contributor
return iVal2 + norm2 * iF;
}
if (iIsFar1 && !iIsFar2)
{
// only one point is a contributor
return iVal1 + norm1 * iF;
}
if (iIsFar1 && iIsFar2)
{
auto dot = static_cast<OutputVectorRealType>(Edge1 * Edge2);
if (dot >= 0.)
{
return ComputeUpdate(iVal1, iVal2, norm2, sq_norm2, norm1, sq_norm1, dot, iF);
}
else
{
OutputVectorRealType sq_norm3, norm3, dot1, dot2;
OutputPointIdentifierType new_id;
bool unfolded =
UnfoldTriangle(oMesh, iId, iCurrentPoint, iId1, iP1, iId2, iP2, norm3, sq_norm3, dot1, dot2, new_id);
if (unfolded)
{
OutputVectorRealType t_sq_norm3 = norm3 * norm3;
auto val3 = static_cast<OutputVectorRealType>(this->GetOutputValue(oMesh, new_id));
OutputVectorRealType t1 = ComputeUpdate(iVal1, val3, norm3, t_sq_norm3, norm1, sq_norm1, dot1, iF);
OutputVectorRealType t2 = ComputeUpdate(iVal2, val3, norm3, t_sq_norm3, norm2, sq_norm2, dot2, iF);
return std::min(t1, t2);
}
else
{
return ComputeUpdate(iVal1, iVal2, norm2, sq_norm2, norm1, sq_norm1, dot, iF);
}
}
}
return this->m_LargeValue;
}
template <typename TInput, typename TOutput>
const typename FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::OutputVectorRealType
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::ComputeUpdate(const OutputVectorRealType & iVal1,
const OutputVectorRealType & iVal2,
const OutputVectorRealType & iNorm1,
const OutputVectorRealType & iSqNorm1,
const OutputVectorRealType & iNorm2,
const OutputVectorRealType & iSqNorm2,
const OutputVectorRealType & iDot,
const OutputVectorRealType & iF) const
{
auto large_value = static_cast<OutputVectorRealType>(this->m_LargeValue);
OutputVectorRealType t = large_value;
OutputVectorRealType CosAngle = iDot;
OutputVectorRealType SinAngle = std::sqrt(1. - iDot * iDot);
OutputVectorRealType u = iVal2 - iVal1;
OutputVectorRealType sq_u = u * u;
OutputVectorRealType f2 = iSqNorm1 + iSqNorm2 - 2. * iNorm1 * iNorm2 * CosAngle;
OutputVectorRealType f1 = iNorm2 * u * (iNorm1 * CosAngle - iNorm2);
OutputVectorRealType f0 = iSqNorm2 * (sq_u - iF * iF * iSqNorm1 * SinAngle * SinAngle);
OutputVectorRealType delta = f1 * f1 - f0 * f2;
OutputVectorRealType epsilon = NumericTraits<OutputVectorRealType>::epsilon();
if (delta >= 0.)
{
if (itk::Math::abs(f2) > epsilon)
{
t = (-f1 - std::sqrt(delta)) / f2;
// test if we must must choose the other solution
if ((t < u) || (iNorm2 * (t - u) / t < iNorm1 * CosAngle) || (iNorm1 / CosAngle < iNorm2 * (t - u) / t))
{
t = (-f1 + std::sqrt(delta)) / f2;
}
}
else
{
if (itk::Math::abs(f1) > epsilon)
{
t = -f0 / f1;
}
else
{
t = -large_value;
}
}
}
else
{
t = -large_value;
}
// choose the update from the 2 vertex only if upwind criterion is met
if ((u < t) && (iNorm1 * CosAngle < iNorm2 * (t - u) / t) && (iNorm2 * (t - u) / t < iNorm1 / CosAngle))
{
return t + iVal1;
}
else
{
return std::min(iNorm2 * iF + iVal1, iNorm1 * iF + iVal2);
}
}
template <typename TInput, typename TOutput>
bool
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::UnfoldTriangle(OutputMeshType * oMesh,
const OutputPointIdentifierType & iId,
const OutputPointType & iP,
const OutputPointIdentifierType & iId1,
const OutputPointType & iP1,
const OutputPointIdentifierType & iId2,
const OutputPointType & iP2,
OutputVectorRealType & oNorm,
OutputVectorRealType & oSqNorm,
OutputVectorRealType & oDot1,
OutputVectorRealType & oDot2,
OutputPointIdentifierType & oId) const
{
(void)iId;
OutputVectorType Edge1 = iP1 - iP;
OutputVectorRealType Norm1 = Edge1.GetNorm();
OutputVectorRealType epsilon = NumericTraits<OutputVectorRealType>::epsilon();
if (Norm1 > epsilon)
{
OutputVectorRealType inv_Norm = 1. / Norm1;
Edge1 *= inv_Norm;
}
OutputVectorType Edge2 = iP2 - iP;
OutputVectorRealType Norm2 = Edge2.GetNorm();
if (Norm2 > epsilon)
{
OutputVectorRealType inv_Norm = 1. / Norm2;
Edge2 *= inv_Norm;
}
auto dot = static_cast<OutputVectorRealType>(Edge1 * Edge2);
// the equation of the lines defining the unfolding region
// [e.g. line 1 : {x ; <x,eq1>=0} ]
using Vector2DType = Vector<OutputVectorRealType, 2>;
using Matrix2DType = Matrix<OutputVectorRealType, 2, 2>;
Vector2DType v1;
v1[0] = dot;
v1[1] = std::sqrt(1. - dot * dot);
Vector2DType v2;
v2[0] = 1.;
v2[1] = 0.;
Vector2DType x1;
x1[0] = Norm1;
x1[1] = 0.;
Vector2DType x2 = Norm2 * v1;
// keep track of the starting point
Vector2DType x_start1(x1);
Vector2DType x_start2(x2);
OutputPointIdentifierType id1 = iId1;
OutputPointIdentifierType id2 = iId2;
OutputQEType * qe = oMesh->FindEdge(id1, id2);
qe = qe->GetSym();
OutputPointType t_pt1 = iP1;
OutputPointType t_pt2 = iP2;
OutputPointType t_pt = t_pt1;
unsigned int nNum = 0;
while (nNum < 50 && qe->GetLeft() != OutputMeshType::m_NoFace)
{
OutputQEType * qe_Lnext = qe->GetLnext();
OutputPointIdentifierType t_id = qe_Lnext->GetDestination();
oMesh->GetPoint(t_id, &t_pt);
Edge1 = t_pt2 - t_pt1;
Norm1 = Edge1.GetNorm();
if (Norm1 > epsilon)
{
OutputVectorRealType inv_Norm = 1. / Norm1;
Edge1 *= inv_Norm;
}
Edge2 = t_pt - t_pt1;
Norm2 = Edge2.GetNorm();
if (Norm2 > epsilon)
{
OutputVectorRealType inv_Norm = 1. / Norm2;
Edge2 *= inv_Norm;
}
/* compute the position of the new point x on the unfolding plane (via a rotation of -alpha on (x2-x1)/rNorm1 )
| cos(alpha) sin(alpha)|
x = |-sin(alpha) cos(alpha)| * [x2-x1]*rNorm2/rNorm1 + x1 where cos(alpha)=dot
*/
Vector2DType vv;
if (Norm1 > epsilon)
{
vv = (x2 - x1) * Norm2 / Norm1;
}
else
{
vv.Fill(0.);
}
dot = Edge1 * Edge2;
Matrix2DType rotation;
rotation[0][0] = dot;
rotation[0][1] = std::sqrt(1. - dot * dot);
rotation[1][0] = -rotation[0][1];
rotation[1][1] = dot;
Vector2DType x = rotation * vv + x1;
/* compute the intersection points.
We look for x=x1+lambda*(x-x1) or x=x2+lambda*(x-x2) with <x,eqi>=0, so */
OutputVectorRealType lambda11 = -(x1 * v1) / ((x - x1) * v1); // left most
OutputVectorRealType lambda12 = -(x1 * v2) / ((x - x1) * v2); // right most
OutputVectorRealType lambda21 = -(x2 * v1) / ((x - x2) * v1); // left most
OutputVectorRealType lambda22 = -(x2 * v2) / ((x - x2) * v2); // right most
bool bIntersect11 = (lambda11 >= 0.) && (lambda11 <= 1.);
bool bIntersect12 = (lambda12 >= 0.) && (lambda12 <= 1.);
bool bIntersect21 = (lambda21 >= 0.) && (lambda21 <= 1.);
bool bIntersect22 = (lambda22 >= 0.) && (lambda22 <= 1.);
if (bIntersect11 && bIntersect12)
{
qe = oMesh->FindEdge(id1, t_id);
qe = qe->GetSym();
id2 = t_id;
t_pt2 = t_pt;
x2 = x;
}
else
{
if (bIntersect21 && bIntersect22)
{
qe = oMesh->FindEdge(id2, t_id);
qe = qe->GetSym();
id1 = t_id;
t_pt1 = t_pt;
x1 = x;
}
else
{ // bIntersect11 && !bIntersect12 && !bIntersect21 && bIntersect22
// that's it, we have found the point
oSqNorm = x.GetSquaredNorm();
if (oSqNorm > epsilon)
{
oNorm = std::sqrt(oSqNorm);
OutputVectorRealType temp_norm = x_start1.GetNorm();
if (temp_norm > epsilon)
{
oDot1 = x * x_start1 / (oNorm * temp_norm);
}
else
{
oDot1 = 0.;
}
temp_norm = x_start2.GetNorm();
if (temp_norm > epsilon)
{
oDot2 = x * x_start2 / (oNorm * temp_norm);
}
else
{
oDot2 = 0.;
}
}
else
{
oNorm = 0.;
oDot1 = 0.;
oDot2 = 0.;
}
oId = t_id;
return true;
}
}
++nNum;
}
return false;
}
template <typename TInput, typename TOutput>
bool
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::CheckTopology(OutputMeshType * oMesh, const NodeType & iNode)
{
(void)oMesh;
(void)iNode;
return true;
}
template <typename TInput, typename TOutput>
void
FastMarchingQuadEdgeMeshFilterBase<TInput, TOutput>::InitializeOutput(OutputMeshType * oMesh)
{
this->CopyInputMeshToOutputMeshGeometry();
// Check that the input mesh is made of triangles
{
OutputCellsContainerPointer cells = oMesh->GetCells();
OutputCellsContainerConstIterator c_it = cells->Begin();
OutputCellsContainerConstIterator c_end = cells->End();
while (c_it != c_end)
{
OutputCellType * cell = c_it.Value();
if (cell->GetNumberOfPoints() != 3)
{
itkGenericExceptionMacro("Input mesh has non triangular faces");
}
++c_it;
}
}
OutputPointsContainerPointer points = oMesh->GetPoints();
OutputPointDataContainerPointer pointdata = OutputPointDataContainer::New();
pointdata->Reserve(points->Size());
OutputPointsContainerIterator p_it = points->Begin();
OutputPointsContainerIterator p_end = points->End();
while (p_it != p_end)
{
pointdata->SetElement(p_it->Index(), this->m_LargeValue);
++p_it;
}
oMesh->SetPointData(pointdata);
m_Label.clear();
if (this->m_AlivePoints)
{
NodePairContainerConstIterator pointsIter = this->m_AlivePoints->Begin();
NodePairContainerConstIterator pointsEnd = this->m_AlivePoints->End();
while (pointsIter != pointsEnd)
{
// get node from alive points container
NodeType idx = pointsIter->Value().GetNode();
if (points->IndexExists(idx))
{
OutputPixelType outputPixel = pointsIter->Value().GetValue();
this->SetLabelValueForGivenNode(idx, Traits::Alive);
this->SetOutputValue(oMesh, idx, outputPixel);
}
++pointsIter;
}
}
if (this->m_ForbiddenPoints)
{
NodePairContainerConstIterator pointsIter = this->m_ForbiddenPoints->Begin();
NodePairContainerConstIterator pointsEnd = this->m_ForbiddenPoints->End();
OutputPixelType zero{};
while (pointsIter != pointsEnd)
{
NodeType idx = pointsIter->Value().GetNode();
if (points->IndexExists(idx))
{
this->SetLabelValueForGivenNode(idx, Traits::Forbidden);
this->SetOutputValue(oMesh, idx, zero);
}
++pointsIter;
}
}
if (this->m_TrialPoints)
{
NodePairContainerConstIterator pointsIter = this->m_TrialPoints->Begin();
NodePairContainerConstIterator pointsEnd = this->m_TrialPoints->End();
while (pointsIter != pointsEnd)
{
NodeType idx = pointsIter->Value().GetNode();
if (points->IndexExists(idx))
{
OutputPixelType outputPixel = pointsIter->Value().GetValue();
this->SetLabelValueForGivenNode(idx, Traits::InitialTrial);
this->SetOutputValue(oMesh, idx, outputPixel);
this->m_Heap.push(pointsIter->Value());
}
++pointsIter;
}
}
}
} // namespace itk
#endif // itkFastMarchingQuadEdgeMeshFilterBase_hxx
|