/*
* Project: MoleCuilder
* Description: creates and alters molecular systems
* Copyright (C) 2010-2012 University of Bonn. All rights reserved.
*
*
* This file is part of MoleCuilder.
*
* MoleCuilder 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 of the License, or
* (at your option) any later version.
*
* MoleCuilder 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 MoleCuilder. If not, see .
*/
/*
* tesselation.cpp
*
* Created on: Aug 3, 2009
* Author: heber
*/
// include config.h
#ifdef HAVE_CONFIG_H
#include
#endif
#include "CodePatterns/MemDebug.hpp"
#include
#include
#include
#include
#include
#include
#include "tesselation.hpp"
#include "BoundaryPointSet.hpp"
#include "BoundaryLineSet.hpp"
#include "BoundaryTriangleSet.hpp"
#include "BoundaryPolygonSet.hpp"
#include "CandidateForTesselation.hpp"
#include "CodePatterns/Assert.hpp"
#include "CodePatterns/Info.hpp"
#include "CodePatterns/IteratorAdaptors.hpp"
#include "CodePatterns/Log.hpp"
#include "CodePatterns/Verbose.hpp"
#include "Helpers/helpers.hpp"
#include "LinearAlgebra/Exceptions.hpp"
#include "LinearAlgebra/Line.hpp"
#include "LinearAlgebra/Plane.hpp"
#include "LinearAlgebra/Vector.hpp"
#include "LinearAlgebra/vector_ops.hpp"
#include "LinkedCell/IPointCloud.hpp"
#include "LinkedCell/linkedcell.hpp"
#include "LinkedCell/PointCloudAdaptor.hpp"
#include "tesselationhelpers.hpp"
#include "Atom/TesselPoint.hpp"
#include "triangleintersectionlist.hpp"
class molecule;
const char *TecplotSuffix = ".dat";
const char *Raster3DSuffix = ".r3d";
const char *VRMLSUffix = ".wrl";
const double ParallelEpsilon = 1e-3;
const double Tesselation::HULLEPSILON = 1e-9;
/** Constructor of class Tesselation.
*/
Tesselation::Tesselation() :
PointsOnBoundaryCount(0), LinesOnBoundaryCount(0), TrianglesOnBoundaryCount(0), LastTriangle(NULL), TriangleFilesWritten(0), InternalPointer(PointsOnBoundary.begin())
{
//Info FunctionInfo(__func__);
}
;
/** Destructor of class Tesselation.
* We have to free all points, lines and triangles.
*/
Tesselation::~Tesselation()
{
//Info FunctionInfo(__func__);
LOG(2, "INFO: Free'ing TesselStruct ... ");
for (TriangleMap::iterator runner = TrianglesOnBoundary.begin(); runner != TrianglesOnBoundary.end(); runner++) {
if (runner->second != NULL) {
delete (runner->second);
runner->second = NULL;
} else
ELOG(1, "The triangle " << runner->first << " has already been free'd.");
}
LOG(1, "INFO: This envelope was written to file " << TriangleFilesWritten << " times(s).");
}
/** Performs tesselation of a given point \a cloud with rolling sphere of
* \a SPHERERADIUS.
*
* @param cloud point cloud to tesselate
* @param SPHERERADIUS radius of the rolling sphere
*/
void Tesselation::operator()(IPointCloud & cloud, const double SPHERERADIUS)
{
// create linkedcell
LinkedCell_deprecated *LinkedList = new LinkedCell_deprecated(cloud, 2. * SPHERERADIUS);
// check for at least three points
{
bool ThreePointsFound = true;
cloud.GoToFirst();
for (size_t i = 0; i < 3; ++i, cloud.GoToNext())
ThreePointsFound &= (!cloud.IsEnd());
cloud.GoToFirst();
if (ThreePointsFound == false) {
ELOG(2, "Less than 3 points in cloud, not enough for tesselation.");
return;
}
}
// find a starting triangle
FindStartingTriangle(SPHERERADIUS, LinkedList);
CandidateForTesselation *baseline = NULL;
BoundaryTriangleSet *T = NULL;
bool OneLoopWithoutSuccessFlag = true;
while ((!OpenLines.empty()) && (OneLoopWithoutSuccessFlag)) {
// 2a. fill all new OpenLines
for (CandidateMap::iterator Runner = OpenLines.begin(); Runner != OpenLines.end(); Runner++) {
baseline = Runner->second;
if (baseline->pointlist.empty()) {
T = (((baseline->BaseLine->triangles.begin()))->second);
//the line is there, so there is a triangle, but only one.
const bool TesselationFailFlag = FindNextSuitableTriangle(*baseline, *T, SPHERERADIUS, LinkedList);
ASSERT(TesselationFailFlag, "Tesselation::operator() - no suitable candidate triangle found.");
}
}
// 2b. search for smallest ShortestAngle among all candidates
double ShortestAngle = 4. * M_PI;
for (CandidateMap::iterator Runner = OpenLines.begin(); Runner != OpenLines.end(); Runner++) {
if (Runner->second->ShortestAngle < ShortestAngle) {
baseline = Runner->second;
ShortestAngle = baseline->ShortestAngle;
}
}
if ((ShortestAngle == 4. * M_PI) || (baseline->pointlist.empty()))
OneLoopWithoutSuccessFlag = false;
else {
AddCandidatePolygon(*baseline, SPHERERADIUS, LinkedList);
}
}
delete LinkedList;
}
/** Determines the volume of a tesselated convex envelope.
*
* @param IsAngstroem unit of length is angstroem or bohr radii
* \return determined volume of envelope assumed being convex
*/
double Tesselation::getVolumeOfConvexEnvelope(const bool IsAngstroem) const
{
// calculate center of gravity
Vector center;
if (!PointsOnBoundary.empty()) {
for (PointMap::const_iterator iter = PointsOnBoundary.begin();
iter != PointsOnBoundary.end(); ++iter)
center += iter->second->node->getPosition();
center *= 1./(double)PointsOnBoundary.size();
}
// 6a. Every triangle forms a pyramid with the center of gravity as its peak, sum up the volumes
double volume = 0.;
for (TriangleMap::const_iterator runner = TrianglesOnBoundary.begin(); runner != TrianglesOnBoundary.end(); runner++) { // go through every triangle, calculate volume of its pyramid with CoG as peak
const double TetrahedronVolume = CalculateVolumeofGeneralTetraeder(
runner->second->endpoints[0]->getPosition(),
runner->second->endpoints[1]->getPosition(),
runner->second->endpoints[2]->getPosition(),
center);
LOG(1, "INFO: volume of tetrahedron is " << setprecision(10) << TetrahedronVolume
<< " " << (IsAngstroem ? "angstrom" : "atomiclength") << "^3.");
volume += TetrahedronVolume;
}
LOG(0, "RESULT: The summed volume is " << setprecision(6) << volume
<< " " << (IsAngstroem ? "angstrom" : "atomiclength") << "^3.");
return volume;
}
/** Determines the area of a tesselated envelope.
*
* @param IsAngstroem unit of length is angstroem or bohr radii
* \return determined surface area of the envelope
*/
double Tesselation::getAreaOfEnvelope(const bool IsAngstroem) const
{
double surfacearea = 0.;
Vector x;
Vector y;
// 6a. Every triangle forms a pyramid with the center of gravity as its peak, sum up the volumes
for (TriangleMap::const_iterator runner = TrianglesOnBoundary.begin(); runner != TrianglesOnBoundary.end(); runner++) { // go through every triangle, calculate volume of its pyramid with CoG as peak
const double area = runner->second->getArea();
LOG(1, "INFO: Area of triangle is " << setprecision(10) << area << " " << (IsAngstroem ? "angstrom" : "atomiclength") << "^2.");
surfacearea += area;
}
LOG(0, "RESULT: The summed surface area is " << setprecision(6) << surfacearea << " " << (IsAngstroem ? "angstrom" : "atomiclength") << "^3.");
return surfacearea;
}
/** Gueses first starting triangle of the convex envelope.
* We guess the starting triangle by taking the smallest distance between two points and looking for a fitting third.
* \param *out output stream for debugging
* \param PointsOnBoundary set of boundary points defining the convex envelope of the cluster
*/
void Tesselation::GuessStartingTriangle()
{
//Info FunctionInfo(__func__);
// 4b. create a starting triangle
// 4b1. create all distances
DistanceMultiMap DistanceMMap;
double distance, tmp;
Vector PlaneVector, TrialVector;
PointMap::iterator A, B, C; // three nodes of the first triangle
A = PointsOnBoundary.begin(); // the first may be chosen arbitrarily
// with A chosen, take each pair B,C and sort
if (A != PointsOnBoundary.end()) {
B = A;
B++;
for (; B != PointsOnBoundary.end(); B++) {
C = B;
C++;
for (; C != PointsOnBoundary.end(); C++) {
tmp = A->second->node->DistanceSquared(B->second->node->getPosition());
distance = tmp * tmp;
tmp = A->second->node->DistanceSquared(C->second->node->getPosition());
distance += tmp * tmp;
tmp = B->second->node->DistanceSquared(C->second->node->getPosition());
distance += tmp * tmp;
DistanceMMap.insert(DistanceMultiMapPair(distance, pair(B, C)));
}
}
}
// // listing distances
// if (DoLog(1)) {
// std::stringstream output;
// output << "Listing DistanceMMap:";
// for(DistanceMultiMap::iterator runner = DistanceMMap.begin(); runner != DistanceMMap.end(); runner++) {
// output << " " << runner->first << "(" << *runner->second.first->second << ", " << *runner->second.second->second << ")";
// }
// LOG(1, output.str());
// }
// 4b2. pick three baselines forming a triangle
// 1. we take from the smallest sum of squared distance as the base line BC (with peak A) onward as the triangle candidate
DistanceMultiMap::iterator baseline = DistanceMMap.begin();
for (; baseline != DistanceMMap.end(); baseline++) {
// we take from the smallest sum of squared distance as the base line BC (with peak A) onward as the triangle candidate
// 2. next, we have to check whether all points reside on only one side of the triangle
// 3. construct plane vector
PlaneVector = Plane(A->second->node->getPosition(), baseline->second.first->second->node->getPosition(), baseline->second.second->second->node->getPosition()).getNormal();
LOG(2, "Plane vector of candidate triangle is " << PlaneVector);
// 4. loop over all points
double sign = 0.;
PointMap::iterator checker = PointsOnBoundary.begin();
for (; checker != PointsOnBoundary.end(); checker++) {
// (neglecting A,B,C)
if ((checker == A) || (checker == baseline->second.first) || (checker == baseline->second.second))
continue;
// 4a. project onto plane vector
TrialVector = (checker->second->node->getPosition() - A->second->node->getPosition());
distance = TrialVector.ScalarProduct(PlaneVector);
if (fabs(distance) < 1e-4) // we need to have a small epsilon around 0 which is still ok
continue;
LOG(2, "Projection of " << checker->second->node->getName() << " yields distance of " << distance << ".");
tmp = distance / fabs(distance);
// 4b. Any have different sign to than before? (i.e. would lie outside convex hull with this starting triangle)
if ((sign != 0) && (tmp != sign)) {
// 4c. If so, break 4. loop and continue with next candidate in 1. loop
LOG(2, "Current candidates: " << A->second->node->getName() << "," << baseline->second.first->second->node->getName() << "," << baseline->second.second->second->node->getName() << " leaves " << checker->second->node->getName() << " outside the convex hull.");
break;
} else { // note the sign for later
LOG(2, "Current candidates: " << A->second->node->getName() << "," << baseline->second.first->second->node->getName() << "," << baseline->second.second->second->node->getName() << " leave " << checker->second->node->getName() << " inside the convex hull.");
sign = tmp;
}
// 4d. Check whether the point is inside the triangle (check distance to each node
tmp = checker->second->node->DistanceSquared(A->second->node->getPosition());
int innerpoint = 0;
if ((tmp < A->second->node->DistanceSquared(baseline->second.first->second->node->getPosition())) && (tmp < A->second->node->DistanceSquared(baseline->second.second->second->node->getPosition())))
innerpoint++;
tmp = checker->second->node->DistanceSquared(baseline->second.first->second->node->getPosition());
if ((tmp < baseline->second.first->second->node->DistanceSquared(A->second->node->getPosition())) && (tmp < baseline->second.first->second->node->DistanceSquared(baseline->second.second->second->node->getPosition())))
innerpoint++;
tmp = checker->second->node->DistanceSquared(baseline->second.second->second->node->getPosition());
if ((tmp < baseline->second.second->second->node->DistanceSquared(baseline->second.first->second->node->getPosition())) && (tmp < baseline->second.second->second->node->DistanceSquared(A->second->node->getPosition())))
innerpoint++;
// 4e. If so, break 4. loop and continue with next candidate in 1. loop
if (innerpoint == 3)
break;
}
// 5. come this far, all on same side? Then break 1. loop and construct triangle
if (checker == PointsOnBoundary.end()) {
LOG(2, "Looks like we have a candidate!");
break;
}
}
if (baseline != DistanceMMap.end()) {
BPS[0] = baseline->second.first->second;
BPS[1] = baseline->second.second->second;
BLS[0] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
BPS[0] = A->second;
BPS[1] = baseline->second.second->second;
BLS[1] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
BPS[0] = baseline->second.first->second;
BPS[1] = A->second;
BLS[2] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
// 4b3. insert created triangle
BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
TrianglesOnBoundary.insert(TrianglePair(TrianglesOnBoundaryCount, BTS));
TrianglesOnBoundaryCount++;
for (int i = 0; i < NDIM; i++) {
LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount, BTS->lines[i]));
LinesOnBoundaryCount++;
}
LOG(1, "Starting triangle is " << *BTS << ".");
} else {
ELOG(0, "No starting triangle found.");
}
}
;
/** Tesselates the convex envelope of a cluster from a single starting triangle.
* The starting triangle is made out of three baselines. Each line in the final tesselated cluster may belong to at most
* 2 triangles. Hence, we go through all current lines:
* -# if the lines contains to only one triangle
* -# We search all points in the boundary
* -# if the triangle is in forward direction of the baseline (at most 90 degrees angle between vector orthogonal to
* baseline in triangle plane pointing out of the triangle and normal vector of new triangle)
* -# if the triangle with the baseline and the current point has the smallest of angles (comparison between normal vectors)
* -# then we have a new triangle, whose baselines we again add (or increase their TriangleCount)
* \param *out output stream for debugging
* \param *configuration for IsAngstroem
* \param *cloud cluster of points
*/
void Tesselation::TesselateOnBoundary(IPointCloud & cloud)
{
//Info FunctionInfo(__func__);
bool flag;
PointMap::iterator winner;
class BoundaryPointSet *peak = NULL;
double SmallestAngle, TempAngle;
Vector NormalVector, VirtualNormalVector, CenterVector, TempVector, helper, PropagationVector, *Center = NULL;
LineMap::iterator LineChecker[2];
Center = cloud.GetCenter();
// create a first tesselation with the given BoundaryPoints
do {
flag = false;
for (LineMap::iterator baseline = LinesOnBoundary.begin(); baseline != LinesOnBoundary.end(); baseline++)
if (baseline->second->triangles.size() == 1) {
// 5a. go through each boundary point if not _both_ edges between either endpoint of the current line and this point exist (and belong to 2 triangles)
SmallestAngle = M_PI;
// get peak point with respect to this base line's only triangle
BTS = baseline->second->triangles.begin()->second; // there is only one triangle so far
LOG(3, "DEBUG: Current baseline is between " << *(baseline->second) << ".");
for (int i = 0; i < 3; i++)
if ((BTS->endpoints[i] != baseline->second->endpoints[0]) && (BTS->endpoints[i] != baseline->second->endpoints[1]))
peak = BTS->endpoints[i];
LOG(3, "DEBUG: and has peak " << *peak << ".");
// prepare some auxiliary vectors
Vector BaseLineCenter, BaseLine;
BaseLineCenter = 0.5 * ((baseline->second->endpoints[0]->node->getPosition()) + (baseline->second->endpoints[1]->node->getPosition()));
BaseLine = (baseline->second->endpoints[0]->node->getPosition()) - (baseline->second->endpoints[1]->node->getPosition());
// offset to center of triangle
CenterVector.Zero();
for (int i = 0; i < 3; i++)
CenterVector += BTS->getEndpoint(i);
CenterVector.Scale(1. / 3.);
LOG(2, "CenterVector of base triangle is " << CenterVector);
// normal vector of triangle
NormalVector = (*Center) - CenterVector;
BTS->GetNormalVector(NormalVector);
NormalVector = BTS->NormalVector;
LOG(4, "DEBUG: NormalVector of base triangle is " << NormalVector);
// vector in propagation direction (out of triangle)
// project center vector onto triangle plane (points from intersection plane-NormalVector to plane-CenterVector intersection)
PropagationVector = Plane(BaseLine, NormalVector, 0).getNormal();
TempVector = CenterVector - (baseline->second->endpoints[0]->node->getPosition()); // TempVector is vector on triangle plane pointing from one baseline egde towards center!
//LOG(0, "Projection of propagation onto temp: " << PropagationVector.Projection(&TempVector) << ".");
if (PropagationVector.ScalarProduct(TempVector) > 0) // make sure normal propagation vector points outward from baseline
PropagationVector.Scale(-1.);
LOG(4, "DEBUG: PropagationVector of base triangle is " << PropagationVector);
winner = PointsOnBoundary.end();
// loop over all points and calculate angle between normal vector of new and present triangle
for (PointMap::iterator target = PointsOnBoundary.begin(); target != PointsOnBoundary.end(); target++) {
if ((target->second != baseline->second->endpoints[0]) && (target->second != baseline->second->endpoints[1])) { // don't take the same endpoints
LOG(4, "DEBUG: Target point is " << *(target->second) << ":");
// first check direction, so that triangles don't intersect
VirtualNormalVector = (target->second->node->getPosition()) - BaseLineCenter;
VirtualNormalVector.ProjectOntoPlane(NormalVector);
TempAngle = VirtualNormalVector.Angle(PropagationVector);
LOG(5, "DEBUG: VirtualNormalVector is " << VirtualNormalVector << " and PropagationVector is " << PropagationVector << ".");
if (TempAngle > (M_PI / 2.)) { // no bends bigger than Pi/2 (90 degrees)
LOG(5, "DEBUG: Angle on triangle plane between propagation direction and base line to " << *(target->second) << " is " << TempAngle << ", bad direction!");
continue;
} else
LOG(5, "DEBUG: Angle on triangle plane between propagation direction and base line to " << *(target->second) << " is " << TempAngle << ", good direction!");
// check first and second endpoint (if any connecting line goes to target has at least not more than 1 triangle)
LineChecker[0] = baseline->second->endpoints[0]->lines.find(target->first);
LineChecker[1] = baseline->second->endpoints[1]->lines.find(target->first);
if (((LineChecker[0] != baseline->second->endpoints[0]->lines.end()) && (LineChecker[0]->second->triangles.size() == 2))) {
LOG(5, "DEBUG: " << *(baseline->second->endpoints[0]) << " has line " << *(LineChecker[0]->second) << " to " << *(target->second) << " as endpoint with " << LineChecker[0]->second->triangles.size() << " triangles.");
continue;
}
if (((LineChecker[1] != baseline->second->endpoints[1]->lines.end()) && (LineChecker[1]->second->triangles.size() == 2))) {
LOG(5, "DEBUG: " << *(baseline->second->endpoints[1]) << " has line " << *(LineChecker[1]->second) << " to " << *(target->second) << " as endpoint with " << LineChecker[1]->second->triangles.size() << " triangles.");
continue;
}
// check whether the envisaged triangle does not already exist (if both lines exist and have same endpoint)
if ((((LineChecker[0] != baseline->second->endpoints[0]->lines.end()) && (LineChecker[1] != baseline->second->endpoints[1]->lines.end()) && (GetCommonEndpoint(LineChecker[0]->second, LineChecker[1]->second) == peak)))) {
LOG(6, "DEBUG: Current target is peak!");
continue;
}
// check for linear dependence
TempVector = (baseline->second->endpoints[0]->node->getPosition()) - (target->second->node->getPosition());
helper = (baseline->second->endpoints[1]->node->getPosition()) - (target->second->node->getPosition());
helper.ProjectOntoPlane(TempVector);
if (fabs(helper.NormSquared()) < MYEPSILON) {
LOG(2, "Chosen set of vectors is linear dependent.");
continue;
}
// in case NOT both were found, create virtually this triangle, get its normal vector, calculate angle
flag = true;
VirtualNormalVector = Plane((baseline->second->endpoints[0]->node->getPosition()), (baseline->second->endpoints[1]->node->getPosition()), (target->second->node->getPosition())).getNormal();
TempVector = (1. / 3.) * ((baseline->second->endpoints[0]->node->getPosition()) + (baseline->second->endpoints[1]->node->getPosition()) + (target->second->node->getPosition()));
TempVector -= (*Center);
// make it always point outward
if (VirtualNormalVector.ScalarProduct(TempVector) < 0)
VirtualNormalVector.Scale(-1.);
// calculate angle
TempAngle = NormalVector.Angle(VirtualNormalVector);
LOG(5, "DEBUG: NormalVector is " << VirtualNormalVector << " and the angle is " << TempAngle << ".");
if ((SmallestAngle - TempAngle) > MYEPSILON) { // set to new possible winner
SmallestAngle = TempAngle;
winner = target;
LOG(5, "DEBUG: New winner " << *winner->second->node << " due to smaller angle between normal vectors.");
} else
if (fabs(SmallestAngle - TempAngle) < MYEPSILON) { // check the angle to propagation, both possible targets are in one plane! (their normals have same angle)
// hence, check the angles to some normal direction from our base line but in this common plane of both targets...
helper = (target->second->node->getPosition()) - BaseLineCenter;
helper.ProjectOntoPlane(BaseLine);
// ...the one with the smaller angle is the better candidate
TempVector = (target->second->node->getPosition()) - BaseLineCenter;
TempVector.ProjectOntoPlane(VirtualNormalVector);
TempAngle = TempVector.Angle(helper);
TempVector = (winner->second->node->getPosition()) - BaseLineCenter;
TempVector.ProjectOntoPlane(VirtualNormalVector);
if (TempAngle < TempVector.Angle(helper)) {
TempAngle = NormalVector.Angle(VirtualNormalVector);
SmallestAngle = TempAngle;
winner = target;
LOG(5, "DEBUG: New winner " << *winner->second->node << " due to smaller angle " << TempAngle << " to propagation direction.");
} else
LOG(5, "DEBUG: Keeping old winner " << *winner->second->node << " due to smaller angle to propagation direction.");
} else
LOG(5, "DEBUG: Keeping old winner " << *winner->second->node << " due to smaller angle between normal vectors.");
}
} // end of loop over all boundary points
// 5b. The point of the above whose triangle has the greatest angle with the triangle the current line belongs to (it only belongs to one, remember!): New triangle
if (winner != PointsOnBoundary.end()) {
LOG(3, "DEBUG: Winning target point is " << *(winner->second) << " with angle " << SmallestAngle << ".");
// create the lins of not yet present
BLS[0] = baseline->second;
// 5c. add lines to the line set if those were new (not yet part of a triangle), delete lines that belong to two triangles)
LineChecker[0] = baseline->second->endpoints[0]->lines.find(winner->first);
LineChecker[1] = baseline->second->endpoints[1]->lines.find(winner->first);
if (LineChecker[0] == baseline->second->endpoints[0]->lines.end()) { // create
BPS[0] = baseline->second->endpoints[0];
BPS[1] = winner->second;
BLS[1] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount, BLS[1]));
LinesOnBoundaryCount++;
} else
BLS[1] = LineChecker[0]->second;
if (LineChecker[1] == baseline->second->endpoints[1]->lines.end()) { // create
BPS[0] = baseline->second->endpoints[1];
BPS[1] = winner->second;
BLS[2] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount, BLS[2]));
LinesOnBoundaryCount++;
} else
BLS[2] = LineChecker[1]->second;
BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
BTS->GetCenter(helper);
helper -= (*Center);
helper *= -1;
BTS->GetNormalVector(helper);
TrianglesOnBoundary.insert(TrianglePair(TrianglesOnBoundaryCount, BTS));
TrianglesOnBoundaryCount++;
} else {
ELOG(2, "I could not determine a winner for this baseline " << *(baseline->second) << ".");
}
// 5d. If the set of lines is not yet empty, go to 5. and continue
} else
LOG(3, "DEBUG: Baseline candidate " << *(baseline->second) << " has a triangle count of " << baseline->second->triangles.size() << ".");
} while (flag);
// exit
delete (Center);
}
;
/** Inserts all points outside of the tesselated surface into it by adding new triangles.
* \param *out output stream for debugging
* \param *cloud cluster of points
* \param *LC LinkedCell_deprecated structure to find nearest point quickly
* \return true - all straddling points insert, false - something went wrong
*/
bool Tesselation::InsertStraddlingPoints(IPointCloud & cloud, const LinkedCell_deprecated *LC)
{
//Info FunctionInfo(__func__);
Vector Intersection, Normal;
TesselPoint *Walker = NULL;
Vector *Center = cloud.GetCenter();
TriangleList *triangles = NULL;
bool AddFlag = false;
LinkedCell_deprecated *BoundaryPoints = NULL;
bool SuccessFlag = true;
cloud.GoToFirst();
PointCloudAdaptor > newcloud(this, cloud.GetName());
BoundaryPoints = new LinkedCell_deprecated(newcloud, 5.);
while (!cloud.IsEnd()) { // we only have to go once through all points, as boundary can become only bigger
if (AddFlag) {
delete (BoundaryPoints);
BoundaryPoints = new LinkedCell_deprecated(newcloud, 5.);
AddFlag = false;
}
Walker = cloud.GetPoint();
LOG(3, "DEBUG: Current point is " << *Walker << ".");
// get the next triangle
triangles = FindClosestTrianglesToVector(Walker->getPosition(), BoundaryPoints);
if (triangles != NULL)
BTS = triangles->front();
else
BTS = NULL;
delete triangles;
if ((BTS == NULL) || (BTS->ContainsBoundaryPoint(Walker))) {
LOG(3, "DEBUG: No triangles found, probably a tesselation point itself.");
cloud.GoToNext();
continue;
} else {
}
LOG(3, "DEBUG: Closest triangle is " << *BTS << ".");
// get the intersection point
if (BTS->GetIntersectionInsideTriangle(*Center, Walker->getPosition(), Intersection)) {
LOG(3, "DEBUG: We have an intersection at " << Intersection << ".");
// we have the intersection, check whether in- or outside of boundary
if ((Center->DistanceSquared(Walker->getPosition()) - Center->DistanceSquared(Intersection)) < -MYEPSILON) {
// inside, next!
LOG(3, "DEBUG: " << *Walker << " is inside wrt triangle " << *BTS << ".");
} else {
// outside!
LOG(3, "DEBUG: " << *Walker << " is outside wrt triangle " << *BTS << ".");
class BoundaryLineSet *OldLines[3], *NewLines[3];
class BoundaryPointSet *OldPoints[3], *NewPoint;
// store the three old lines and old points
for (int i = 0; i < 3; i++) {
OldLines[i] = BTS->lines[i];
OldPoints[i] = BTS->endpoints[i];
}
Normal = BTS->NormalVector;
// add Walker to boundary points
LOG(3, "DEBUG: Adding " << *Walker << " to BoundaryPoints.");
AddFlag = true;
if (AddBoundaryPoint(Walker, 0))
NewPoint = BPS[0];
else
continue;
// remove triangle
LOG(3, "DEBUG: Erasing triangle " << *BTS << ".");
TrianglesOnBoundary.erase(BTS->Nr);
delete (BTS);
// create three new boundary lines
for (int i = 0; i < 3; i++) {
BPS[0] = NewPoint;
BPS[1] = OldPoints[i];
NewLines[i] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
LOG(4, "DEBUG: Creating new line " << *NewLines[i] << ".");
LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount, NewLines[i])); // no need for check for unique insertion as BPS[0] is definitely a new one
LinesOnBoundaryCount++;
}
// create three new triangle with new point
for (int i = 0; i < 3; i++) { // find all baselines
BLS[0] = OldLines[i];
int n = 1;
for (int j = 0; j < 3; j++) {
if (NewLines[j]->IsConnectedTo(BLS[0])) {
if (n > 2) {
ELOG(2, BLS[0] << " connects to all of the new lines?!");
return false;
} else
BLS[n++] = NewLines[j];
}
}
// create the triangle
BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
Normal.Scale(-1.);
BTS->GetNormalVector(Normal);
Normal.Scale(-1.);
LOG(3, "DEBUG: Created new triangle " << *BTS << ".");
TrianglesOnBoundary.insert(TrianglePair(TrianglesOnBoundaryCount, BTS));
TrianglesOnBoundaryCount++;
}
}
} else { // something is wrong with FindClosestTriangleToPoint!
ELOG(1, "The closest triangle did not produce an intersection!");
SuccessFlag = false;
break;
}
cloud.GoToNext();
}
// exit
delete (Center);
delete (BoundaryPoints);
return SuccessFlag;
}
;
/** Adds a point to the tesselation::PointsOnBoundary list.
* \param *Walker point to add
* \param n TesselStruct::BPS index to put pointer into
* \return true - new point was added, false - point already present
*/
bool Tesselation::AddBoundaryPoint(TesselPoint * Walker, const int n)
{
//Info FunctionInfo(__func__);
PointTestPair InsertUnique;
BPS[n] = new class BoundaryPointSet(Walker);
InsertUnique = PointsOnBoundary.insert(PointPair(Walker->getNr(), BPS[n]));
if (InsertUnique.second) { // if new point was not present before, increase counter
PointsOnBoundaryCount++;
return true;
} else {
delete (BPS[n]);
BPS[n] = InsertUnique.first->second;
return false;
}
}
;
/** Adds point to Tesselation::PointsOnBoundary if not yet present.
* Tesselation::TPS is set to either this new BoundaryPointSet or to the existing one of not unique.
* @param Candidate point to add
* @param n index for this point in Tesselation::TPS array
*/
void Tesselation::AddTesselationPoint(TesselPoint* Candidate, const int n)
{
//Info FunctionInfo(__func__);
PointTestPair InsertUnique;
TPS[n] = new class BoundaryPointSet(Candidate);
InsertUnique = PointsOnBoundary.insert(PointPair(Candidate->getNr(), TPS[n]));
if (InsertUnique.second) { // if new point was not present before, increase counter
PointsOnBoundaryCount++;
} else {
delete TPS[n];
LOG(4, "DEBUG: Node " << *((InsertUnique.first)->second->node) << " is already present in PointsOnBoundary.");
TPS[n] = (InsertUnique.first)->second;
}
}
;
/** Sets point to a present Tesselation::PointsOnBoundary.
* Tesselation::TPS is set to the existing one or NULL if not found.
* @param Candidate point to set to
* @param n index for this point in Tesselation::TPS array
*/
void Tesselation::SetTesselationPoint(TesselPoint* Candidate, const int n) const
{
//Info FunctionInfo(__func__);
PointMap::const_iterator FindPoint = PointsOnBoundary.find(Candidate->getNr());
if (FindPoint != PointsOnBoundary.end())
TPS[n] = FindPoint->second;
else
TPS[n] = NULL;
}
;
/** Function tries to add line from current Points in BPS to BoundaryLineSet.
* If successful it raises the line count and inserts the new line into the BLS,
* if unsuccessful, it writes the line which had been present into the BLS, deleting the new constructed one.
* @param *OptCenter desired OptCenter if there are more than one candidate line
* @param *candidate third point of the triangle to be, for checking between multiple open line candidates
* @param *a first endpoint
* @param *b second endpoint
* @param n index of Tesselation::BLS giving the line with both endpoints
* @return true - inserted a new line, false - present line used
*/
bool Tesselation::AddTesselationLine(const Vector * const OptCenter, const BoundaryPointSet * const candidate, class BoundaryPointSet *a, class BoundaryPointSet *b, const int n)
{
bool insertNewLine = true;
LineMap::iterator FindLine = a->lines.find(b->node->getNr());
BoundaryLineSet *WinningLine = NULL;
if (FindLine != a->lines.end()) {
LOG(3, "DEBUG: There is at least one line between " << *a << " and " << *b << ": " << *(FindLine->second) << ".");
pair FindPair;
FindPair = a->lines.equal_range(b->node->getNr());
for (FindLine = FindPair.first; (FindLine != FindPair.second) && (insertNewLine); FindLine++) {
LOG(3, "DEBUG: Checking line " << *(FindLine->second) << " ...");
// If there is a line with less than two attached triangles, we don't need a new line.
if (FindLine->second->triangles.size() == 1) {
CandidateMap::iterator Finder = OpenLines.find(FindLine->second);
ASSERT(Finder != OpenLines.end(), "Tesselation::AddTesselationLine() - " + toString(*FindLine->second) + " is not a new line and not in OpenLines.");
if (Finder->second != NULL) {
if (!Finder->second->pointlist.empty())
LOG(4, "DEBUG: line " << *(FindLine->second) << " is open with candidate " << **(Finder->second->pointlist.begin()) << ".");
else {
LOG(4, "ACCEPT: line " << *(FindLine->second) << " is open with no candidate.");
insertNewLine = false;
WinningLine = FindLine->second;
}
// get open line
for (TesselPointList::const_iterator CandidateChecker = Finder->second->pointlist.begin(); CandidateChecker != Finder->second->pointlist.end(); ++CandidateChecker) {
if ((*(CandidateChecker) == candidate->node) && (OptCenter == NULL || OptCenter->DistanceSquared(Finder->second->OptCenter) < MYEPSILON)) { // stop searching if candidate matches
LOG(4, "ACCEPT: Candidate " << *(*CandidateChecker) << " has the right center " << Finder->second->OptCenter << ".");
insertNewLine = false;
WinningLine = FindLine->second;
break;
} else {
LOG(5, "REJECT: Candidate " << *(*CandidateChecker) << "'s center " << Finder->second->OptCenter << " does not match desired on " << *OptCenter << ".");
}
}
} else {
LOG(4, "ACCEPT: There are no candidates for present open line, use what we have.");
insertNewLine = false;
WinningLine = FindLine->second;
}
}
}
}
if (insertNewLine) {
AddNewTesselationTriangleLine(a, b, n);
} else {
AddExistingTesselationTriangleLine(WinningLine, n);
}
return insertNewLine;
}
;
/**
* Adds lines from each of the current points in the BPS to BoundaryLineSet.
* Raises the line count and inserts the new line into the BLS.
*
* @param *a first endpoint
* @param *b second endpoint
* @param n index of Tesselation::BLS giving the line with both endpoints
*/
void Tesselation::AddNewTesselationTriangleLine(class BoundaryPointSet *a, class BoundaryPointSet *b, const int n)
{
//Info FunctionInfo(__func__);
LOG(2, "DEBUG: Adding open line [" << LinesOnBoundaryCount << "|" << *(a->node) << " and " << *(b->node) << ".");
BPS[0] = a;
BPS[1] = b;
BLS[n] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount); // this also adds the line to the local maps
// add line to global map
LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount, BLS[n]));
// increase counter
LinesOnBoundaryCount++;
// also add to open lines
CandidateForTesselation *CFT = new CandidateForTesselation(BLS[n]);
OpenLines.insert(pair(BLS[n], CFT));
}
;
/** Uses an existing line for a new triangle.
* Sets Tesselation::BLS[\a n] and removes the lines from Tesselation::OpenLines.
* \param *FindLine the line to add
* \param n index of the line to set in Tesselation::BLS
*/
void Tesselation::AddExistingTesselationTriangleLine(class BoundaryLineSet *Line, int n)
{
//Info FunctionInfo(__func__);
LOG(5, "DEBUG: Using existing line " << *Line);
// set endpoints and line
BPS[0] = Line->endpoints[0];
BPS[1] = Line->endpoints[1];
BLS[n] = Line;
// remove existing line from OpenLines
CandidateMap::iterator CandidateLine = OpenLines.find(BLS[n]);
if (CandidateLine != OpenLines.end()) {
LOG(6, "DEBUG: Removing line from OpenLines.");
delete (CandidateLine->second);
OpenLines.erase(CandidateLine);
} else {
ELOG(1, "Line exists and is attached to less than two triangles, but not in OpenLines!");
}
}
;
/** Function adds triangle to global list.
* Furthermore, the triangle receives the next free id and id counter \a TrianglesOnBoundaryCount is increased.
*/
void Tesselation::AddTesselationTriangle()
{
//Info FunctionInfo(__func__);
LOG(4, "DEBUG: Adding triangle to global TrianglesOnBoundary map.");
// add triangle to global map
TrianglesOnBoundary.insert(TrianglePair(TrianglesOnBoundaryCount, BTS));
TrianglesOnBoundaryCount++;
// set as last new triangle
LastTriangle = BTS;
// NOTE: add triangle to local maps is done in constructor of BoundaryTriangleSet
}
;
/** Function adds triangle to global list.
* Furthermore, the triangle number is set to \a Nr.
* \param getNr() triangle number
*/
void Tesselation::AddTesselationTriangle(const int nr)
{
//Info FunctionInfo(__func__);
LOG(4, "DEBUG: Adding triangle to global TrianglesOnBoundary map.");
// add triangle to global map
TrianglesOnBoundary.insert(TrianglePair(nr, BTS));
// set as last new triangle
LastTriangle = BTS;
// NOTE: add triangle to local maps is done in constructor of BoundaryTriangleSet
}
;
/** Removes a triangle from the tesselation.
* Removes itself from the TriangleMap's of its lines, calls for them RemoveTriangleLine() if they are no more connected.
* Removes itself from memory.
* \param *triangle to remove
*/
void Tesselation::RemoveTesselationTriangle(class BoundaryTriangleSet *triangle)
{
//Info FunctionInfo(__func__);
if (triangle == NULL)
return;
for (int i = 0; i < 3; i++) {
if (triangle->lines[i] != NULL) {
LOG(4, "DEBUG: Removing triangle Nr." << triangle->Nr << " in line " << *triangle->lines[i] << ".");
triangle->lines[i]->triangles.erase(triangle->Nr);
std::stringstream output;
output << *triangle->lines[i] << " is ";
if (triangle->lines[i]->triangles.empty()) {
output << "no more attached to any triangle, erasing.";
RemoveTesselationLine(triangle->lines[i]);
} else {
output << "still attached to another triangle: ";
for (TriangleMap::iterator TriangleRunner = triangle->lines[i]->triangles.begin(); TriangleRunner != triangle->lines[i]->triangles.end(); TriangleRunner++)
output << "\t[" << (TriangleRunner->second)->Nr << "|" << *((TriangleRunner->second)->endpoints[0]) << ", " << *((TriangleRunner->second)->endpoints[1]) << ", " << *((TriangleRunner->second)->endpoints[2]) << "] \t";
OpenLines.insert(pair(triangle->lines[i], NULL));
}
LOG(3, "DEBUG: " << output.str());
triangle->lines[i] = NULL; // free'd or not: disconnect
} else
ELOG(1, "This line " << i << " has already been free'd.");
}
if (TrianglesOnBoundary.erase(triangle->Nr))
LOG(3, "DEBUG: Removing triangle Nr. " << triangle->Nr << ".");
delete (triangle);
}
;
/** Removes a line from the tesselation.
* Removes itself from each endpoints' LineMap, then removes itself from global LinesOnBoundary list and free's the line.
* \param *line line to remove
*/
void Tesselation::RemoveTesselationLine(class BoundaryLineSet *line)
{
//Info FunctionInfo(__func__);
int Numbers[2];
if (line == NULL)
return;
// get other endpoint number for finding copies of same line
if (line->endpoints[1] != NULL)
Numbers[0] = line->endpoints[1]->Nr;
else
Numbers[0] = -1;
if (line->endpoints[0] != NULL)
Numbers[1] = line->endpoints[0]->Nr;
else
Numbers[1] = -1;
// erase from OpenLines if present
{
CandidateMap::iterator openlineiter = OpenLines.find(line);
if (openlineiter != OpenLines.end())
OpenLines.erase(openlineiter);
}
for (int i = 0; i < 2; i++) {
if (line->endpoints[i] != NULL) {
if (Numbers[i] != -1) { // as there may be multiple lines with same endpoints, we have to go through each and find in the endpoint's line list this line set
pair erasor = line->endpoints[i]->lines.equal_range(Numbers[i]);
for (LineMap::iterator Runner = erasor.first; Runner != erasor.second; Runner++)
if ((*Runner).second == line) {
LOG(4, "DEBUG: Removing Line Nr. " << line->Nr << " in boundary point " << *line->endpoints[i] << ".");
line->endpoints[i]->lines.erase(Runner);
break;
}
} else { // there's just a single line left
if (line->endpoints[i]->lines.erase(line->Nr))
LOG(4, "DEBUG: Removing Line Nr. " << line->Nr << " in boundary point " << *line->endpoints[i] << ".");
}
if (line->endpoints[i]->lines.empty()) {
LOG(4, "DEBUG: " << *line->endpoints[i] << " has no more lines it's attached to, erasing.");
RemoveTesselationPoint(line->endpoints[i]);
} else
if (DoLog(0)) {
std::stringstream output;
output << "DEBUG: " << *line->endpoints[i] << " has still lines it's attached to: ";
for (LineMap::iterator LineRunner = line->endpoints[i]->lines.begin(); LineRunner != line->endpoints[i]->lines.end(); LineRunner++)
output << "[" << *(LineRunner->second) << "] \t";
LOG(4, output.str());
}
line->endpoints[i] = NULL; // free'd or not: disconnect
} else
ELOG(4, "DEBUG: Endpoint " << i << " has already been free'd.");
}
if (!line->triangles.empty())
ELOG(2, "Memory Leak! I " << *line << " am still connected to some triangles.");
if (LinesOnBoundary.erase(line->Nr))
LOG(4, "DEBUG: Removing line Nr. " << line->Nr << ".");
delete (line);
}
;
/** Removes a point from the tesselation.
* Checks whether there are still lines connected, removes from global PointsOnBoundary list, then free's the point.
* \note If a point should be removed, while keep the tesselated surface intact (i.e. closed), use RemovePointFromTesselatedSurface()
* \param *point point to remove
*/
void Tesselation::RemoveTesselationPoint(class BoundaryPointSet *point)
{
//Info FunctionInfo(__func__);
if (point == NULL)
return;
if (PointsOnBoundary.erase(point->Nr))
LOG(4, "DEBUG: Removing point Nr. " << point->Nr << ".");
delete (point);
}
;
/** Checks validity of a given sphere of a candidate line.
* \sa CandidateForTesselation::CheckValidity(), which is more evolved.
* We check CandidateForTesselation::OtherOptCenter
* \param &CandidateLine contains other degenerated candidates which we have to subtract as well
* \param RADIUS radius of sphere
* \param *LC LinkedCell_deprecated structure with other atoms
* \return true - candidate triangle is degenerated, false - candidate triangle is not degenerated
*/
bool Tesselation::CheckDegeneracy(CandidateForTesselation &CandidateLine, const double RADIUS, const LinkedCell_deprecated *LC) const
{
//Info FunctionInfo(__func__);
LOG(3, "DEBUG: Checking degeneracy by whether sphere contains no others points ...");
bool flag = true;
LOG(3, "DEBUG: Check by: draw sphere {" << CandidateLine.OtherOptCenter[0] << " " << CandidateLine.OtherOptCenter[1] << " " << CandidateLine.OtherOptCenter[2] << "} radius " << RADIUS << " resolution 30");
// get all points inside the sphere
TesselPointList *ListofPoints = LC->GetPointsInsideSphere(RADIUS, &CandidateLine.OtherOptCenter);
LOG(3, "DEBUG: CheckDegeneracy: There are " << ListofPoints->size() << " points inside the sphere.");
LOG(4, "DEBUG: The following atoms are inside sphere at " << CandidateLine.OtherOptCenter << ":");
for (TesselPointList::const_iterator Runner = ListofPoints->begin(); Runner != ListofPoints->end(); ++Runner)
LOG(4, "DEBUG: " << *(*Runner) << " with distance " << (*Runner)->distance(CandidateLine.OtherOptCenter) << ".");
// remove triangles's endpoints
for (int i = 0; i < 2; i++)
ListofPoints->remove(CandidateLine.BaseLine->endpoints[i]->node);
// remove other candidates
for (TesselPointList::const_iterator Runner = CandidateLine.pointlist.begin(); Runner != CandidateLine.pointlist.end(); ++Runner)
ListofPoints->remove(*Runner);
// check for other points
if (!ListofPoints->empty()) {
LOG(3, "DEBUG: CheckDegeneracy: There are still " << ListofPoints->size() << " points inside the sphere.");
flag = false;
LOG(4, "DEBUG: External atoms inside of sphere at " << CandidateLine.OtherOptCenter << ":");
for (TesselPointList::const_iterator Runner = ListofPoints->begin(); Runner != ListofPoints->end(); ++Runner)
LOG(4, "DEBUG: " << *(*Runner) << " with distance " << (*Runner)->distance(CandidateLine.OtherOptCenter) << ".");
}
delete ListofPoints;
return flag;
}
;
/** Checks whether the triangle consisting of the three points is already present.
* Searches for the points in Tesselation::PointsOnBoundary and checks their
* lines. If any of the three edges already has two triangles attached, false is
* returned.
* \param *out output stream for debugging
* \param *Candidates endpoints of the triangle candidate
* \return integer 0 if no triangle exists, 1 if one triangle exists, 2 if two
* triangles exist which is the maximum for three points
*/
int Tesselation::CheckPresenceOfTriangle(TesselPoint *Candidates[3]) const
{
//Info FunctionInfo(__func__);
int adjacentTriangleCount = 0;
class BoundaryPointSet *Points[3];
// builds a triangle point set (Points) of the end points
for (int i = 0; i < 3; i++) {
PointMap::const_iterator FindPoint = PointsOnBoundary.find(Candidates[i]->getNr());
if (FindPoint != PointsOnBoundary.end()) {
Points[i] = FindPoint->second;
} else {
Points[i] = NULL;
}
}
// checks lines between the points in the Points for their adjacent triangles
for (int i = 0; i < 3; i++) {
if (Points[i] != NULL) {
for (int j = i; j < 3; j++) {
if (Points[j] != NULL) {
LineMap::const_iterator FindLine = Points[i]->lines.find(Points[j]->node->getNr());
for (; (FindLine != Points[i]->lines.end()) && (FindLine->first == Points[j]->node->getNr()); FindLine++) {
TriangleMap *triangles = &FindLine->second->triangles;
LOG(5, "DEBUG: Current line is " << FindLine->first << ": " << *(FindLine->second) << " with triangles " << triangles << ".");
for (TriangleMap::const_iterator FindTriangle = triangles->begin(); FindTriangle != triangles->end(); FindTriangle++) {
if (FindTriangle->second->IsPresentTupel(Points)) {
adjacentTriangleCount++;
}
}
}
// Only one of the triangle lines must be considered for the triangle count.
//LOG(5, "DEBUG: Found " << adjacentTriangleCount << " adjacent triangles for the point set.");
//return adjacentTriangleCount;
}
}
}
}
LOG(3, "DEBUG: Found " << adjacentTriangleCount << " adjacent triangles for the point set.");
return adjacentTriangleCount;
}
;
/** Checks whether the triangle consisting of the three points is already present.
* Searches for the points in Tesselation::PointsOnBoundary and checks their
* lines. If any of the three edges already has two triangles attached, false is
* returned.
* \param *out output stream for debugging
* \param *Candidates endpoints of the triangle candidate
* \return NULL - none found or pointer to triangle
*/
class BoundaryTriangleSet * Tesselation::GetPresentTriangle(TesselPoint *Candidates[3])
{
//Info FunctionInfo(__func__);
class BoundaryTriangleSet *triangle = NULL;
class BoundaryPointSet *Points[3];
// builds a triangle point set (Points) of the end points
for (int i = 0; i < 3; i++) {
PointMap::iterator FindPoint = PointsOnBoundary.find(Candidates[i]->getNr());
if (FindPoint != PointsOnBoundary.end()) {
Points[i] = FindPoint->second;
} else {
Points[i] = NULL;
}
}
// checks lines between the points in the Points for their adjacent triangles
for (int i = 0; i < 3; i++) {
if (Points[i] != NULL) {
for (int j = i; j < 3; j++) {
if (Points[j] != NULL) {
LineMap::iterator FindLine = Points[i]->lines.find(Points[j]->node->getNr());
for (; (FindLine != Points[i]->lines.end()) && (FindLine->first == Points[j]->node->getNr()); FindLine++) {
TriangleMap *triangles = &FindLine->second->triangles;
for (TriangleMap::iterator FindTriangle = triangles->begin(); FindTriangle != triangles->end(); FindTriangle++) {
if (FindTriangle->second->IsPresentTupel(Points)) {
if ((triangle == NULL) || (triangle->Nr > FindTriangle->second->Nr))
triangle = FindTriangle->second;
}
}
}
// Only one of the triangle lines must be considered for the triangle count.
//LOG(5, "DEBUG: Found " << adjacentTriangleCount << " adjacent triangles for the point set.");
//return adjacentTriangleCount;
}
}
}
}
return triangle;
}
;
/** Finds the starting triangle for FindNonConvexBorder().
* Looks at the outermost point per axis, then FindSecondPointForTesselation()
* for the second and FindNextSuitablePointViaAngleOfSphere() for the third
* point are called.
* \param *out output stream for debugging
* \param RADIUS radius of virtual rolling sphere
* \param *LC LinkedCell_deprecated structure with neighbouring TesselPoint's
* \return true - a starting triangle has been created, false - no valid triple of points found
*/
bool Tesselation::FindStartingTriangle(const double RADIUS, const LinkedCell_deprecated *LC)
{
//Info FunctionInfo(__func__);
int i = 0;
TesselPoint* MaxPoint[NDIM];
TesselPoint* Temporary;
double maxCoordinate[NDIM];
BoundaryLineSet *BaseLine = NULL;
Vector helper;
Vector Chord;
Vector SearchDirection;
Vector CircleCenter; // center of the circle, i.e. of the band of sphere's centers
Vector CirclePlaneNormal; // normal vector defining the plane this circle lives in
Vector SphereCenter;
Vector NormalVector;
NormalVector.Zero();
for (i = 0; i < 3; i++) {
MaxPoint[i] = NULL;
maxCoordinate[i] = -10e30;
}
// 1. searching topmost point with respect to each axis
LOG(2, "INFO: Searching for topmost pointer from each axis.");
for (int i = 0; i < NDIM; i++) { // each axis
LC->n[i] = LC->N[i] - 1; // current axis is topmost cell
const int map[NDIM] = {i, (i + 1) % NDIM, (i + 2) % NDIM};
for (LC->n[map[1]] = 0; LC->n[map[1]] < LC->N[map[1]]; LC->n[map[1]]++)
for (LC->n[map[2]] = 0; LC->n[map[2]] < LC->N[map[2]]; LC->n[map[2]]++) {
const TesselPointSTLList *List = LC->GetCurrentCell();
//LOG(1, "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << ".");
if (List != NULL) {
for (TesselPointSTLList::const_iterator Runner = List->begin(); Runner != List->end(); Runner++) {
if ((*Runner)->at(map[0]) > maxCoordinate[map[0]]) {
LOG(4, "DEBUG: New maximal for axis " << map[0] << " node is " << *(*Runner) << " at " << (*Runner)->getPosition() << ".");
maxCoordinate[map[0]] = (*Runner)->at(map[0]);
MaxPoint[map[0]] = (*Runner);
}
}
} else {
ELOG(1, "The current cell " << LC->n[0] << "," << LC->n[1] << "," << LC->n[2] << " is invalid!");
}
}
}
if (DoLog(3)) {
std::stringstream output;
output << "Found maximum coordinates: ";
for (int i = 0; i < NDIM; i++)
output << i << ": " << *MaxPoint[i] << "\t";
LOG(3, "DEBUG: " << output.str());
}
BTS = NULL;
for (int k = 0; k < NDIM; k++) {
NormalVector.Zero();
NormalVector[k] = 1.;
BaseLine = new BoundaryLineSet();
BaseLine->endpoints[0] = new BoundaryPointSet(MaxPoint[k]);
LOG(2, "INFO: Coordinates of start node at " << *BaseLine->endpoints[0]->node << ".");
double ShortestAngle;
ShortestAngle = 999999.; // This will contain the angle, which will be always positive (when looking for second point), when looking for third point this will be the quadrant.
Temporary = NULL;
FindSecondPointForTesselation(BaseLine->endpoints[0]->node, NormalVector, Temporary, &ShortestAngle, RADIUS, LC); // we give same point as next candidate as its bonds are looked into in find_second_...
if (Temporary == NULL) {
// have we found a second point?
delete BaseLine;
continue;
}
BaseLine->endpoints[1] = new BoundaryPointSet(Temporary);
LOG(2, "INFO: Second node is at " << *Temporary << ".");
// construct center of circle
CircleCenter = 0.5 * ((BaseLine->endpoints[0]->node->getPosition()) + (BaseLine->endpoints[1]->node->getPosition()));
LOG(4, "DEBUG: CircleCenter is at " << CircleCenter << ".");
// construct normal vector of circle
CirclePlaneNormal = (BaseLine->endpoints[0]->node->getPosition()) - (BaseLine->endpoints[1]->node->getPosition());
LOG(4, "DEBUG: CirclePlaneNormal is at " << CirclePlaneNormal << ".");
double radius = CirclePlaneNormal.NormSquared();
double CircleRadius = sqrt(RADIUS * RADIUS - radius / 4.);
NormalVector.ProjectOntoPlane(CirclePlaneNormal);
NormalVector.Normalize();
LOG(4, "DEBUG: NormalVector is at " << NormalVector << ".");
ShortestAngle = 2. * M_PI; // This will indicate the quadrant.
SphereCenter = (CircleRadius * NormalVector) + CircleCenter;
// Now, NormalVector and SphereCenter are two orthonormalized vectors in the plane defined by CirclePlaneNormal (not normalized)
// look in one direction of baseline for initial candidate
try {
SearchDirection = Plane(CirclePlaneNormal, NormalVector, 0).getNormal(); // whether we look "left" first or "right" first is not important ...
} catch (LinearAlgebraException &e) {
ELOG(1, "Vectors are linear dependent: " << CirclePlaneNormal << ", " << NormalVector << ".");
delete BaseLine;
continue;
}
// adding point 1 and point 2 and add the line between them
LOG(3, "DEBUG: Found second point is at " << *BaseLine->endpoints[1]->node << ".");
//LOG(1, "INFO: OldSphereCenter is at " << helper << ".");
CandidateForTesselation OptCandidates(BaseLine);
FindThirdPointForTesselation(NormalVector, SearchDirection, SphereCenter, OptCandidates, NULL, RADIUS, LC);
{
std::stringstream output;
for (TesselPointList::iterator it = OptCandidates.pointlist.begin(); it != OptCandidates.pointlist.end(); it++)
output << *(*it);
LOG(3, "DEBUG: List of third Points is: " << output.str());
}
if (!OptCandidates.pointlist.empty()) {
BTS = NULL;
AddCandidatePolygon(OptCandidates, RADIUS, LC);
} else {
delete BaseLine;
continue;
}
if (BTS != NULL) { // we have created one starting triangle
delete BaseLine;
break;
} else {
// remove all candidates from the list and then the list itself
OptCandidates.pointlist.clear();
}
delete BaseLine;
}
return (BTS != NULL);
}
;
/** Checks for a given baseline and a third point candidate whether baselines of the found triangle don't have even better candidates.
* This is supposed to prevent early closing of the tesselation.
* \param CandidateLine CandidateForTesselation with baseline and shortestangle , i.e. not \a *OptCandidate
* \param *ThirdNode third point in triangle, not in BoundaryLineSet::endpoints
* \param RADIUS radius of sphere
* \param *LC LinkedCell_deprecated structure
* \return true - there is a better candidate (smaller angle than \a ShortestAngle), false - no better TesselPoint candidate found
*/
//bool Tesselation::HasOtherBaselineBetterCandidate(CandidateForTesselation &CandidateLine, const TesselPoint * const ThirdNode, double RADIUS, const LinkedCell_deprecated * const LC) const
//{
// //Info FunctionInfo(__func__);
// bool result = false;
// Vector CircleCenter;
// Vector CirclePlaneNormal;
// Vector OldSphereCenter;
// Vector SearchDirection;
// Vector helper;
// TesselPoint *OtherOptCandidate = NULL;
// double OtherShortestAngle = 2.*M_PI; // This will indicate the quadrant.
// double radius, CircleRadius;
// BoundaryLineSet *Line = NULL;
// BoundaryTriangleSet *T = NULL;
//
// // check both other lines
// PointMap::const_iterator FindPoint = PointsOnBoundary.find(ThirdNode->getNr());
// if (FindPoint != PointsOnBoundary.end()) {
// for (int i=0;i<2;i++) {
// LineMap::const_iterator FindLine = (FindPoint->second)->lines.find(BaseRay->endpoints[0]->node->getNr());
// if (FindLine != (FindPoint->second)->lines.end()) {
// Line = FindLine->second;
// LOG(0, "Found line " << *Line << ".");
// if (Line->triangles.size() == 1) {
// T = Line->triangles.begin()->second;
// // construct center of circle
// CircleCenter.CopyVector(Line->endpoints[0]->node->node);
// CircleCenter.AddVector(Line->endpoints[1]->node->node);
// CircleCenter.Scale(0.5);
//
// // construct normal vector of circle
// CirclePlaneNormal.CopyVector(Line->endpoints[0]->node->node);
// CirclePlaneNormal.SubtractVector(Line->endpoints[1]->node->node);
//
// // calculate squared radius of circle
// radius = CirclePlaneNormal.ScalarProduct(&CirclePlaneNormal);
// if (radius/4. < RADIUS*RADIUS) {
// CircleRadius = RADIUS*RADIUS - radius/4.;
// CirclePlaneNormal.Normalize();
// //LOG(1, "INFO: CircleCenter is at " << CircleCenter << ", CirclePlaneNormal is " << CirclePlaneNormal << " with circle radius " << sqrt(CircleRadius) << ".");
//
// // construct old center
// GetCenterofCircumcircle(&OldSphereCenter, *T->endpoints[0]->node->node, *T->endpoints[1]->node->node, *T->endpoints[2]->node->node);
// helper.CopyVector(&T->NormalVector); // normal vector ensures that this is correct center of the two possible ones
// radius = Line->endpoints[0]->node->node->DistanceSquared(&OldSphereCenter);
// helper.Scale(sqrt(RADIUS*RADIUS - radius));
// OldSphereCenter.AddVector(&helper);
// OldSphereCenter.SubtractVector(&CircleCenter);
// //LOG(1, "INFO: OldSphereCenter is at " << OldSphereCenter << ".");
//
// // construct SearchDirection
// SearchDirection.MakeNormalVector(&T->NormalVector, &CirclePlaneNormal);
// helper.CopyVector(Line->endpoints[0]->node->node);
// helper.SubtractVector(ThirdNode->node);
// if (helper.ScalarProduct(&SearchDirection) < -HULLEPSILON)// ohoh, SearchDirection points inwards!
// SearchDirection.Scale(-1.);
// SearchDirection.ProjectOntoPlane(&OldSphereCenter);
// SearchDirection.Normalize();
// LOG(1, "INFO: SearchDirection is " << SearchDirection << ".");
// if (fabs(OldSphereCenter.ScalarProduct(&SearchDirection)) > HULLEPSILON) {
// // rotated the wrong way!
// ELOG(1, "SearchDirection and RelativeOldSphereCenter are still not orthogonal!");
// }
//
// // add third point
// FindThirdPointForTesselation(T->NormalVector, SearchDirection, OldSphereCenter, OptCandidates, ThirdNode, RADIUS, LC);
// for (TesselPointList::iterator it = OptCandidates.pointlist.begin(); it != OptCandidates.pointlist.end(); ++it) {
// if (((*it) == BaseRay->endpoints[0]->node) || ((*it) == BaseRay->endpoints[1]->node)) // skip if it's the same triangle than suggested
// continue;
// LOG(1, "INFO: Third point candidate is " << (*it)
// << " with circumsphere's center at " << (*it)->OptCenter << ".");
// LOG(1, "INFO: Baseline is " << *BaseRay);
//
// // check whether all edges of the new triangle still have space for one more triangle (i.e. TriangleCount <2)
// TesselPoint *PointCandidates[3];
// PointCandidates[0] = (*it);
// PointCandidates[1] = BaseRay->endpoints[0]->node;
// PointCandidates[2] = BaseRay->endpoints[1]->node;
// bool check=false;
// int existentTrianglesCount = CheckPresenceOfTriangle(PointCandidates);
// // If there is no triangle, add it regularly.
// if (existentTrianglesCount == 0) {
// SetTesselationPoint((*it), 0);
// SetTesselationPoint(BaseRay->endpoints[0]->node, 1);
// SetTesselationPoint(BaseRay->endpoints[1]->node, 2);
//
// if (CheckLineCriteriaForDegeneratedTriangle((const BoundaryPointSet ** const )TPS)) {
// OtherOptCandidate = (*it);
// check = true;
// }
// } else if ((existentTrianglesCount >= 1) && (existentTrianglesCount <= 3)) { // If there is a planar region within the structure, we need this triangle a second time.
// SetTesselationPoint((*it), 0);
// SetTesselationPoint(BaseRay->endpoints[0]->node, 1);
// SetTesselationPoint(BaseRay->endpoints[1]->node, 2);
//
// // We demand that at most one new degenerate line is created and that this line also already exists (which has to be the case due to existentTrianglesCount == 1)
// // i.e. at least one of the three lines must be present with TriangleCount <= 1
// if (CheckLineCriteriaForDegeneratedTriangle((const BoundaryPointSet ** const)TPS)) {
// OtherOptCandidate = (*it);
// check = true;
// }
// }
//
// if (check) {
// if (ShortestAngle > OtherShortestAngle) {
// LOG(0, "There is a better candidate than " << *ThirdNode << " with " << ShortestAngle << " from baseline " << *Line << ": " << *OtherOptCandidate << " with " << OtherShortestAngle << ".");
// result = true;
// break;
// }
// }
// }
// delete(OptCandidates);
// if (result)
// break;
// } else {
// LOG(0, "Circumcircle for base line " << *Line << " and base triangle " << T << " is too big!");
// }
// } else {
// ELOG(2, "Baseline is connected to two triangles already?");
// }
// } else {
// LOG(1, "No present baseline between " << BaseRay->endpoints[0] << " and candidate " << *ThirdNode << ".");
// }
// }
// } else {
// ELOG(1, "Could not find the TesselPoint " << *ThirdNode << ".");
// }
//
// return result;
//};
/** This function finds a triangle to a line, adjacent to an existing one.
* @param out output stream for debugging
* @param CandidateLine current cadndiate baseline to search from
* @param T current triangle which \a Line is edge of
* @param RADIUS radius of the rolling ball
* @param N number of found triangles
* @param *LC LinkedCell_deprecated structure with neighbouring points
* @return false - no suitable candidate found
*/
bool Tesselation::FindNextSuitableTriangle(CandidateForTesselation &CandidateLine, const BoundaryTriangleSet &T, const double& RADIUS, const LinkedCell_deprecated *LC)
{
//Info FunctionInfo(__func__);
Vector CircleCenter;
Vector CirclePlaneNormal;
Vector RelativeSphereCenter;
Vector SearchDirection;
Vector helper;
BoundaryPointSet *ThirdPoint = NULL;
LineMap::iterator testline;
double radius, CircleRadius;
for (int i = 0; i < 3; i++)
if ((T.endpoints[i] != CandidateLine.BaseLine->endpoints[0]) && (T.endpoints[i] != CandidateLine.BaseLine->endpoints[1])) {
ThirdPoint = T.endpoints[i];
break;
}
LOG(3, "DEBUG: Current baseline is " << *CandidateLine.BaseLine << " with ThirdPoint " << *ThirdPoint << " of triangle " << T << ".");
CandidateLine.T = &T;
// construct center of circle
CircleCenter = 0.5 * ((CandidateLine.BaseLine->endpoints[0]->node->getPosition()) + (CandidateLine.BaseLine->endpoints[1]->node->getPosition()));
// construct normal vector of circle
CirclePlaneNormal = (CandidateLine.BaseLine->endpoints[0]->node->getPosition()) - (CandidateLine.BaseLine->endpoints[1]->node->getPosition());
// calculate squared radius of circle
radius = CirclePlaneNormal.ScalarProduct(CirclePlaneNormal);
if (radius / 4. < RADIUS * RADIUS) {
// construct relative sphere center with now known CircleCenter
RelativeSphereCenter = T.SphereCenter - CircleCenter;
CircleRadius = RADIUS * RADIUS - radius / 4.;
CirclePlaneNormal.Normalize();
LOG(4, "DEBUG: CircleCenter is at " << CircleCenter << ", CirclePlaneNormal is " << CirclePlaneNormal << " with circle radius " << sqrt(CircleRadius) << ".");
LOG(4, "DEBUG: OldSphereCenter is at " << T.SphereCenter << ".");
// construct SearchDirection and an "outward pointer"
SearchDirection = Plane(RelativeSphereCenter, CirclePlaneNormal, 0).getNormal();
helper = CircleCenter - (ThirdPoint->node->getPosition());
if (helper.ScalarProduct(SearchDirection) < -HULLEPSILON) // ohoh, SearchDirection points inwards!
SearchDirection.Scale(-1.);
LOG(4, "DEBUG: SearchDirection is " << SearchDirection << ".");
if (fabs(RelativeSphereCenter.ScalarProduct(SearchDirection)) > HULLEPSILON) {
// rotated the wrong way!
ELOG(3, "DEBUG: SearchDirection and RelativeOldSphereCenter are still not orthogonal!");
}
// add third point
FindThirdPointForTesselation(T.NormalVector, SearchDirection, T.SphereCenter, CandidateLine, ThirdPoint, RADIUS, LC);
} else {
LOG(4, "DEBUG: Circumcircle for base line " << *CandidateLine.BaseLine << " and base triangle " << T << " is too big!");
}
if (CandidateLine.pointlist.empty()) {
ELOG(4, "DEBUG: Could not find a suitable candidate.");
return false;
}
{
std::stringstream output;
for (TesselPointList::iterator it = CandidateLine.pointlist.begin(); it != CandidateLine.pointlist.end(); ++it)
output << " " << *(*it);
LOG(3, "DEBUG: Third Points are: " << output.str());
}
return true;
}
;
/** Walks through Tesselation::OpenLines() and finds candidates for newly created ones.
* \param *&LCList atoms in LinkedCell_deprecated list
* \param RADIUS radius of the virtual sphere
* \return true - for all open lines without candidates so far, a candidate has been found,
* false - at least one open line without candidate still
*/
bool Tesselation::FindCandidatesforOpenLines(const double RADIUS, const LinkedCell_deprecated *&LCList)
{
bool TesselationFailFlag = true;
CandidateForTesselation *baseline = NULL;
BoundaryTriangleSet *T = NULL;
for (CandidateMap::iterator Runner = OpenLines.begin(); Runner != OpenLines.end(); Runner++) {
baseline = Runner->second;
if (baseline->pointlist.empty()) {
ASSERT((baseline->BaseLine->triangles.size() == 1), "Open line without exactly one attached triangle");
T = (((baseline->BaseLine->triangles.begin()))->second);
LOG(4, "DEBUG: Finding best candidate for open line " << *baseline->BaseLine << " of triangle " << *T);
TesselationFailFlag = TesselationFailFlag && FindNextSuitableTriangle(*baseline, *T, RADIUS, LCList); //the line is there, so there is a triangle, but only one.
}
}
return TesselationFailFlag;
}
;
/** Adds the present line and candidate point from \a &CandidateLine to the Tesselation.
* \param CandidateLine triangle to add
* \param RADIUS Radius of sphere
* \param *LC LinkedCell_deprecated structure
* \NOTE we need the copy operator here as the original CandidateForTesselation is removed in
* AddTesselationLine() in AddCandidateTriangle()
*/
void Tesselation::AddCandidatePolygon(CandidateForTesselation CandidateLine, const double RADIUS, const LinkedCell_deprecated *LC)
{
//Info FunctionInfo(__func__);
Vector Center;
TesselPoint * const TurningPoint = CandidateLine.BaseLine->endpoints[0]->node;
TesselPointList::iterator Runner;
TesselPointList::iterator Sprinter;
// fill the set of neighbours
TesselPointSet SetOfNeighbours;
SetOfNeighbours.insert(CandidateLine.BaseLine->endpoints[1]->node);
for (TesselPointList::iterator Runner = CandidateLine.pointlist.begin(); Runner != CandidateLine.pointlist.end(); Runner++)
SetOfNeighbours.insert(*Runner);
TesselPointList *connectedClosestPoints = GetCircleOfSetOfPoints(&SetOfNeighbours, TurningPoint, CandidateLine.BaseLine->endpoints[1]->node->getPosition());
if (DoLog(3)) {
std::stringstream output;
for (TesselPointList::iterator TesselRunner = connectedClosestPoints->begin(); TesselRunner != connectedClosestPoints->end(); ++TesselRunner)
output << **TesselRunner;
LOG(3, "DEBUG: List of Candidates for Turning Point " << *TurningPoint << ":" << output.str());
}
// go through all angle-sorted candidates (in degenerate n-nodes case we may have to add multiple triangles)
Runner = connectedClosestPoints->begin();
Sprinter = Runner;
Sprinter++;
while (Sprinter != connectedClosestPoints->end()) {
LOG(3, "DEBUG: Current Runner is " << *(*Runner) << " and sprinter is " << *(*Sprinter) << ".");
AddTesselationPoint(TurningPoint, 0);
AddTesselationPoint(*Runner, 1);
AddTesselationPoint(*Sprinter, 2);
AddCandidateTriangle(CandidateLine, Opt);
Runner = Sprinter;
Sprinter++;
if (Sprinter != connectedClosestPoints->end()) {
// fill the internal open lines with its respective candidate (otherwise lines in degenerate case are not picked)
FindDegeneratedCandidatesforOpenLines(*Sprinter, &CandidateLine.OptCenter); // Assume BTS contains last triangle
LOG(2, "DEBUG: There are still more triangles to add.");
}
// pick candidates for other open lines as well
FindCandidatesforOpenLines(RADIUS, LC);
// check whether we add a degenerate or a normal triangle
if (CheckDegeneracy(CandidateLine, RADIUS, LC)) {
// add normal and degenerate triangles
LOG(3, "DEBUG: Triangle of endpoints " << *TPS[0] << "," << *TPS[1] << " and " << *TPS[2] << " is degenerated, adding both sides.");
AddCandidateTriangle(CandidateLine, OtherOpt);
if (Sprinter != connectedClosestPoints->end()) {
// fill the internal open lines with its respective candidate (otherwise lines in degenerate case are not picked)
FindDegeneratedCandidatesforOpenLines(*Sprinter, &CandidateLine.OtherOptCenter);
}
// pick candidates for other open lines as well
FindCandidatesforOpenLines(RADIUS, LC);
}
}
delete (connectedClosestPoints);
}
;
/** for polygons (multiple candidates for a baseline) sets internal edges to the correct next candidate.
* \param *Sprinter next candidate to which internal open lines are set
* \param *OptCenter OptCenter for this candidate
*/
void Tesselation::FindDegeneratedCandidatesforOpenLines(TesselPoint * const Sprinter, const Vector * const OptCenter)
{
//Info FunctionInfo(__func__);
pair FindPair = TPS[0]->lines.equal_range(TPS[2]->node->getNr());
for (LineMap::const_iterator FindLine = FindPair.first; FindLine != FindPair.second; FindLine++) {
LOG(4, "DEBUG: Checking line " << *(FindLine->second) << " ...");
// If there is a line with less than two attached triangles, we don't need a new line.
if (FindLine->second->triangles.size() == 1) {
CandidateMap::iterator Finder = OpenLines.find(FindLine->second);
if (!Finder->second->pointlist.empty())
LOG(4, "DEBUG: line " << *(FindLine->second) << " is open with candidate " << **(Finder->second->pointlist.begin()) << ".");
else {
LOG(4, "DEBUG: line " << *(FindLine->second) << " is open with no candidate, setting to next Sprinter" << (*Sprinter));
Finder->second->T = BTS; // is last triangle
Finder->second->pointlist.push_back(Sprinter);
Finder->second->ShortestAngle = 0.;
Finder->second->OptCenter = *OptCenter;
}
}
}
}
;
/** If a given \a *triangle is degenerated, this adds both sides.
* i.e. the triangle with same BoundaryPointSet's but NormalVector in opposite direction.
* Note that endpoints are stored in Tesselation::TPS
* \param CandidateLine CanddiateForTesselation structure for the desired BoundaryLine
* \param RADIUS radius of sphere
* \param *LC pointer to LinkedCell_deprecated structure
*/
void Tesselation::AddDegeneratedTriangle(CandidateForTesselation &CandidateLine, const double RADIUS, const LinkedCell_deprecated *LC)
{
//Info FunctionInfo(__func__);
Vector Center;
CandidateMap::const_iterator CandidateCheck = OpenLines.end();
BoundaryTriangleSet *triangle = NULL;
/// 1. Create or pick the lines for the first triangle
LOG(3, "DEBUG: Creating/Picking lines for first triangle ...");
for (int i = 0; i < 3; i++) {
BLS[i] = NULL;
LOG(3, "DEBUG: Current line is between " << *TPS[(i + 0) % 3] << " and " << *TPS[(i + 1) % 3] << ":");
AddTesselationLine(&CandidateLine.OptCenter, TPS[(i + 2) % 3], TPS[(i + 0) % 3], TPS[(i + 1) % 3], i);
}
/// 2. create the first triangle and NormalVector and so on
LOG(3, "DEBUG: Adding first triangle with center at " << CandidateLine.OptCenter << " ...");
BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
AddTesselationTriangle();
// create normal vector
BTS->GetCenter(Center);
Center -= CandidateLine.OptCenter;
BTS->SphereCenter = CandidateLine.OptCenter;
BTS->GetNormalVector(Center);
// give some verbose output about the whole procedure
if (CandidateLine.T != NULL)
LOG(2, "INFO: --> New triangle with " << *BTS << " and normal vector " << BTS->NormalVector << ", from " << *CandidateLine.T << " and angle " << CandidateLine.ShortestAngle << ".");
else
LOG(2, "INFO: --> New starting triangle with " << *BTS << " and normal vector " << BTS->NormalVector << " and no top triangle.");
triangle = BTS;
/// 3. Gather candidates for each new line
LOG(3, "DEBUG: Adding candidates to new lines ...");
for (int i = 0; i < 3; i++) {
LOG(4, "DEBUG: Current line is between " << *TPS[(i + 0) % 3] << " and " << *TPS[(i + 1) % 3] << ":");
CandidateCheck = OpenLines.find(BLS[i]);
if ((CandidateCheck != OpenLines.end()) && (CandidateCheck->second->pointlist.empty())) {
if (CandidateCheck->second->T == NULL)
CandidateCheck->second->T = triangle;
FindNextSuitableTriangle(*(CandidateCheck->second), *CandidateCheck->second->T, RADIUS, LC);
}
}
/// 4. Create or pick the lines for the second triangle
LOG(3, "DEBUG: Creating/Picking lines for second triangle ...");
for (int i = 0; i < 3; i++) {
LOG(4, "DEBUG: Current line is between " << *TPS[(i + 0) % 3] << " and " << *TPS[(i + 1) % 3] << ":");
AddTesselationLine(&CandidateLine.OtherOptCenter, TPS[(i + 2) % 3], TPS[(i + 0) % 3], TPS[(i + 1) % 3], i);
}
/// 5. create the second triangle and NormalVector and so on
LOG(3, "DEBUG: Adding second triangle with center at " << CandidateLine.OtherOptCenter << " ...");
BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
AddTesselationTriangle();
BTS->SphereCenter = CandidateLine.OtherOptCenter;
// create normal vector in other direction
BTS->GetNormalVector(triangle->NormalVector);
BTS->NormalVector.Scale(-1.);
// give some verbose output about the whole procedure
if (CandidateLine.T != NULL)
LOG(2, "DEBUG: --> New degenerate triangle with " << *BTS << " and normal vector " << BTS->NormalVector << ", from " << *CandidateLine.T << " and angle " << CandidateLine.ShortestAngle << ".");
else
LOG(2, "DEBUG: --> New degenerate starting triangle with " << *BTS << " and normal vector " << BTS->NormalVector << " and no top triangle.");
/// 6. Adding triangle to new lines
LOG(3, "DEBUG: Adding second triangles to new lines ...");
for (int i = 0; i < 3; i++) {
LOG(4, "DEBUG: Current line is between " << *TPS[(i + 0) % 3] << " and " << *TPS[(i + 1) % 3] << ":");
CandidateCheck = OpenLines.find(BLS[i]);
if ((CandidateCheck != OpenLines.end()) && (CandidateCheck->second->pointlist.empty())) {
if (CandidateCheck->second->T == NULL)
CandidateCheck->second->T = BTS;
}
}
}
;
/** Adds a triangle to the Tesselation structure from three given TesselPoint's.
* Note that endpoints are in Tesselation::TPS.
* \param CandidateLine CandidateForTesselation structure contains other information
* \param type which opt center to add (i.e. which side) and thus which NormalVector to take
*/
void Tesselation::AddCandidateTriangle(CandidateForTesselation &CandidateLine, enum centers type)
{
//Info FunctionInfo(__func__);
Vector Center;
Vector *OptCenter = (type == Opt) ? &CandidateLine.OptCenter : &CandidateLine.OtherOptCenter;
// add the lines
AddTesselationLine(OptCenter, TPS[2], TPS[0], TPS[1], 0);
AddTesselationLine(OptCenter, TPS[1], TPS[0], TPS[2], 1);
AddTesselationLine(OptCenter, TPS[0], TPS[1], TPS[2], 2);
// add the triangles
BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
AddTesselationTriangle();
// create normal vector
BTS->GetCenter(Center);
Center.SubtractVector(*OptCenter);
BTS->SphereCenter = *OptCenter;
BTS->GetNormalVector(Center);
// give some verbose output about the whole procedure
if (CandidateLine.T != NULL)
LOG(2, "INFO: --> New" << ((type == OtherOpt) ? " degenerate " : " ") << "triangle with " << *BTS << " and normal vector " << BTS->NormalVector << ", from " << *CandidateLine.T << " and angle " << CandidateLine.ShortestAngle << ".");
else
LOG(2, "INFO: --> New" << ((type == OtherOpt) ? " degenerate " : " ") << "starting triangle with " << *BTS << " and normal vector " << BTS->NormalVector << " and no top triangle.");
}
;
bool Tesselation::isConvex() const
{
bool status = true;
// go through all lines on boundary
for (LineMap::const_iterator lineiter = LinesOnBoundary.begin();
lineiter != LinesOnBoundary.end(); ++lineiter) {
if (!lineiter->second->CheckConvexityCriterion()) {
LOG(2, "DEBUG: Line " << *lineiter->second << " is not convex.");
status = false;
}
}
return status;
}
/** Checks whether the quadragon of the two triangles connect to \a *Base is convex.
* We look whether the closest point on \a *Base with respect to the other baseline is outside
* of the segment formed by both endpoints (concave) or not (convex).
* \param *out output stream for debugging
* \param *Base line to be flipped
* \return NULL - convex, otherwise endpoint that makes it concave
*/
class BoundaryPointSet *Tesselation::IsConvexRectangle(class BoundaryLineSet *Base)
{
//Info FunctionInfo(__func__);
class BoundaryPointSet *Spot = NULL;
class BoundaryLineSet *OtherBase;
Vector *ClosestPoint;
int m = 0;
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); runner++)
for (int j = 0; j < 3; j++) // all of their endpoints and baselines
if (!Base->ContainsBoundaryPoint(runner->second->endpoints[j])) // and neither of its endpoints
BPS[m++] = runner->second->endpoints[j];
OtherBase = new class BoundaryLineSet(BPS, -1);
LOG(3, "DEBUG: Current base line is " << *Base << ".");
LOG(3, "DEBUG: Other base line is " << *OtherBase << ".");
// get the closest point on each line to the other line
ClosestPoint = GetClosestPointBetweenLine(Base, OtherBase);
// delete the temporary other base line
delete (OtherBase);
// get the distance vector from Base line to OtherBase line
Vector DistanceToIntersection[2], BaseLine;
double distance[2];
BaseLine = (Base->endpoints[1]->node->getPosition()) - (Base->endpoints[0]->node->getPosition());
for (int i = 0; i < 2; i++) {
DistanceToIntersection[i] = (*ClosestPoint) - (Base->endpoints[i]->node->getPosition());
distance[i] = BaseLine.ScalarProduct(DistanceToIntersection[i]);
}
delete (ClosestPoint);
if ((distance[0] * distance[1]) > 0) { // have same sign?
LOG(4, "REJECT: Both SKPs have same sign: " << distance[0] << " and " << distance[1] << ". " << *Base << "' rectangle is concave.");
if (distance[0] < distance[1]) {
Spot = Base->endpoints[0];
} else {
Spot = Base->endpoints[1];
}
return Spot;
} else { // different sign, i.e. we are in between
LOG(3, "ACCEPT: Rectangle of triangles of base line " << *Base << " is convex.");
return NULL;
}
}
;
void Tesselation::PrintAllBoundaryPoints(ofstream *out) const
{
//Info FunctionInfo(__func__);
// print all lines
std::stringstream output;
for (PointMap::const_iterator PointRunner = PointsOnBoundary.begin(); PointRunner != PointsOnBoundary.end(); PointRunner++)
output << " " << *(PointRunner->second);
LOG(3, "DEBUG: Printing all boundary points for debugging:" << output.str());
}
;
void Tesselation::PrintAllBoundaryLines(ofstream *out) const
{
//Info FunctionInfo(__func__);
// print all lines
std::stringstream output;
for (LineMap::const_iterator LineRunner = LinesOnBoundary.begin(); LineRunner != LinesOnBoundary.end(); LineRunner++)
output << " " << *(LineRunner->second);
LOG(3, "DEBUG: Printing all boundary lines for debugging:" << output.str());
}
;
void Tesselation::PrintAllBoundaryTriangles(ofstream *out) const
{
//Info FunctionInfo(__func__);
// print all triangles
std::stringstream output;
for (TriangleMap::const_iterator TriangleRunner = TrianglesOnBoundary.begin(); TriangleRunner != TrianglesOnBoundary.end(); TriangleRunner++)
output << " " << *(TriangleRunner->second);
LOG(3, "DEBUG: Printing all boundary triangles for debugging:" << output.str());
}
;
/** For a given boundary line \a *Base and its two triangles, picks the central baseline that is "higher".
* \param *out output stream for debugging
* \param *Base line to be flipped
* \return volume change due to flipping (0 - then no flipped occured)
*/
double Tesselation::PickFarthestofTwoBaselines(class BoundaryLineSet *Base)
{
//Info FunctionInfo(__func__);
class BoundaryLineSet *OtherBase;
Vector *ClosestPoint[2];
double volume;
int m = 0;
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); runner++)
for (int j = 0; j < 3; j++) // all of their endpoints and baselines
if (!Base->ContainsBoundaryPoint(runner->second->endpoints[j])) // and neither of its endpoints
BPS[m++] = runner->second->endpoints[j];
OtherBase = new class BoundaryLineSet(BPS, -1);
LOG(3, "DEBUG: Current base line is " << *Base << ".");
LOG(3, "DEBUG: Other base line is " << *OtherBase << ".");
// get the closest point on each line to the other line
ClosestPoint[0] = GetClosestPointBetweenLine(Base, OtherBase);
ClosestPoint[1] = GetClosestPointBetweenLine(OtherBase, Base);
// get the distance vector from Base line to OtherBase line
Vector Distance = (*ClosestPoint[1]) - (*ClosestPoint[0]);
// calculate volume
volume = CalculateVolumeofGeneralTetraeder(Base->endpoints[1]->node->getPosition(), OtherBase->endpoints[0]->node->getPosition(), OtherBase->endpoints[1]->node->getPosition(), Base->endpoints[0]->node->getPosition());
// delete the temporary other base line and the closest points
delete (ClosestPoint[0]);
delete (ClosestPoint[1]);
delete (OtherBase);
if (Distance.NormSquared() < MYEPSILON) { // check for intersection
LOG(3, "REJECT: Both lines have an intersection: Nothing to do.");
return false;
} else { // check for sign against BaseLineNormal
Vector BaseLineNormal;
BaseLineNormal.Zero();
if (Base->triangles.size() < 2) {
ELOG(1, "Less than two triangles are attached to this baseline!");
return 0.;
}
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); runner++) {
LOG(4, "DEBUG: Adding NormalVector " << runner->second->NormalVector << " of triangle " << *(runner->second) << ".");
BaseLineNormal += (runner->second->NormalVector);
}
BaseLineNormal.Scale(1. / 2.);
if (Distance.ScalarProduct(BaseLineNormal) > MYEPSILON) { // Distance points outwards, hence OtherBase higher than Base -> flip
LOG(3, "ACCEPT: Other base line would be higher: Flipping baseline.");
// calculate volume summand as a general tetraeder
return volume;
} else { // Base higher than OtherBase -> do nothing
LOG(3, "REJECT: Base line is higher: Nothing to do.");
return 0.;
}
}
}
;
/** For a given baseline and its two connected triangles, flips the baseline.
* I.e. we create the new baseline between the other two endpoints of these four
* endpoints and reconstruct the two triangles accordingly.
* \param *out output stream for debugging
* \param *Base line to be flipped
* \return pointer to allocated new baseline - flipping successful, NULL - something went awry
*/
class BoundaryLineSet * Tesselation::FlipBaseline(class BoundaryLineSet *Base)
{
//Info FunctionInfo(__func__);
class BoundaryLineSet *OldLines[4], *NewLine;
class BoundaryPointSet *OldPoints[2];
Vector BaseLineNormal[2];
Vector OtherEndpoint[2]; // fourth point to either triangle
int OldTriangleNrs[2], OldBaseLineNr;
int i, m;
// calculate NormalVector for later use
if (Base->triangles.size() < 2) {
ELOG(1, "Less than two triangles are attached to this baseline!");
return NULL;
}
{
int i=0;
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); runner++) {
LOG(5, "INFO: Adding NormalVector " << runner->second->NormalVector << " of triangle " << *(runner->second) << ".");
OtherEndpoint[(i+1)%2] = runner->second->GetThirdEndpoint(Base)->node->getPosition();
BaseLineNormal[i++] = (runner->second->NormalVector);
ASSERT( i <= 2,
"Tesselation::FlipBaseline() - not exactly two triangles found.");
}
}
// get the two triangles
// gather four endpoints and four lines
for (int j = 0; j < 4; j++)
OldLines[j] = NULL;
for (int j = 0; j < 2; j++)
OldPoints[j] = NULL;
i = 0;
m = 0;
// print OldLines and OldPoints for debugging
if (DoLog(3)) {
std::stringstream output;
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); runner++)
for (int j = 0; j < 3; j++) // all of their endpoints and baselines
if (runner->second->lines[j] != Base) // pick not the central baseline
output << *runner->second->lines[j] << "\t";
LOG(3, "DEBUG: The four old lines are: " << output.str());
}
if (DoLog(3)) {
std::stringstream output;
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); runner++)
for (int j = 0; j < 3; j++) // all of their endpoints and baselines
if (!Base->ContainsBoundaryPoint(runner->second->endpoints[j])) // and neither of its endpoints
output << *runner->second->endpoints[j] << "\t";
LOG(3, "DEBUG: The two old points are: " << output.str());
}
// index OldLines and OldPoints
// note that oldlines[0,1] belong to first triangle and hence first normal
// vector and oldlines[2,3] belong to second triangle and its normal vector
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); runner++)
for (int j = 0; j < 3; j++) // all of their endpoints and baselines
if (runner->second->lines[j] != Base) // pick not the central baseline
OldLines[i++] = runner->second->lines[j];
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); runner++)
for (int j = 0; j < 3; j++) // all of their endpoints and baselines
if (!Base->ContainsBoundaryPoint(runner->second->endpoints[j])) // and neither of its endpoints
OldPoints[m++] = runner->second->endpoints[j];
Vector BasePoints[2]; // endpoints of Base
if (OldLines[0]->ContainsBoundaryPoint(Base->endpoints[0])) {
BasePoints[0] = Base->endpoints[0]->node->getPosition();
BasePoints[1] = Base->endpoints[1]->node->getPosition();
} else {
BasePoints[0] = Base->endpoints[1]->node->getPosition();
BasePoints[1] = Base->endpoints[0]->node->getPosition();
}
// check whether everything is in place to create new lines and triangles
if (i < 4) {
ELOG(1, "We have not gathered enough baselines!");
return NULL;
}
for (int j = 0; j < 4; j++)
if (OldLines[j] == NULL) {
ELOG(1, "We have not gathered enough baselines!");
return NULL;
}
for (int j = 0; j < 2; j++)
if (OldPoints[j] == NULL) {
ELOG(1, "We have not gathered enough endpoints!");
return NULL;
}
// construct plane of first triangle for calculating normal vectors later
const Plane triangleplane = Base->triangles.begin()->second->getPlane();
// get fourth point projected down onto this plane
const Vector Intersection = triangleplane.getClosestPoint(OtherEndpoint[0]);
// remove triangles and baseline removes itself
LOG(3, "DEBUG: Deleting baseline " << *Base << " from global list.");
OldBaseLineNr = Base->Nr;
m = 0;
// first obtain all triangle to delete ... (otherwise we pull the carpet (Base) from under the for-loop's feet)
list TrianglesOfBase;
for (TriangleMap::iterator runner = Base->triangles.begin(); runner != Base->triangles.end(); ++runner)
TrianglesOfBase.push_back(runner->second);
// .. then delete each triangle (which deletes the line as well)
for (list::iterator runner = TrianglesOfBase.begin(); !TrianglesOfBase.empty(); runner = TrianglesOfBase.begin()) {
LOG(3, "DEBUG: Deleting triangle " << *(*runner) << ".");
OldTriangleNrs[m++] = (*runner)->Nr;
RemoveTesselationTriangle((*runner));
TrianglesOfBase.erase(runner);
}
// construct new baseline (with same number as old one)
BPS[0] = OldPoints[0];
BPS[1] = OldPoints[1];
NewLine = new class BoundaryLineSet(BPS, OldBaseLineNr);
LinesOnBoundary.insert(LinePair(OldBaseLineNr, NewLine)); // no need for check for unique insertion as NewLine is definitely a new one
LOG(3, "DEBUG: Created new baseline " << *NewLine << ".");
// Explanation for normal vector choice:
// Decisive for the normal vector of the new triangle is whether the fourth
// endpoint is with respect to the joining boundary line on one side or on
// the other with respect to the endpoint of the plane triangle that is not
// contained in the joining boundary line.
// construct new triangles with flipped baseline
i = -1;
if (OldLines[0]->IsConnectedTo(OldLines[2]))
i = 2;
if (OldLines[0]->IsConnectedTo(OldLines[3]))
i = 3;
if (i != -1) {
BLS[0] = OldLines[0];
BLS[1] = OldLines[i];
BLS[2] = NewLine;
BTS = new class BoundaryTriangleSet(BLS, OldTriangleNrs[0]);
{
Line joiningline = makeLineThrough(
OldLines[0]->endpoints[0]->node->getPosition(),
OldLines[0]->endpoints[1]->node->getPosition());
// BasePoints[1] is not contained in OldLines[0], hence is third endpoint
// of plane triangle. BasePoints[0] is in the joining OldLines[0] and
// we check whether Intersection is on same side as BasePoints[1] or not.
const Vector pointinginside =
joiningline.getClosestPoint(BasePoints[1]) - BasePoints[1];
const Vector pointingtointersection =
joiningline.getClosestPoint(Intersection) - Intersection;
const double sign_of_intersection =
pointingtointersection.ScalarProduct(pointinginside);
LOG(4, "DEBUG: Sign of SKP between both lines w.r.t " << *OldLines[0]
<< " is " << sign_of_intersection << ".");
const double sign_of_normal = (sign_of_intersection >= 0) ? -1. : +1.;
LOG(4, "DEBUG: Opposite normal direction is "
<< sign_of_normal * BaseLineNormal[0] << ".");
BTS->GetNormalVector(sign_of_normal * BaseLineNormal[0]);
}
AddTesselationTriangle(OldTriangleNrs[0]);
LOG(3, "DEBUG: Created new triangle " << *BTS << ".");
BLS[0] = (i == 2 ? OldLines[3] : OldLines[2]);
BLS[1] = OldLines[1];
BLS[2] = NewLine;
BTS = new class BoundaryTriangleSet(BLS, OldTriangleNrs[1]);
{
Line joiningline = makeLineThrough(
OldLines[1]->endpoints[0]->node->getPosition(),
OldLines[1]->endpoints[1]->node->getPosition());
// BasePoints[0] is not contained in OldLines[1], hence is third endpoint
// of plane triangle. BasePoints[1] is in th1e joining OldLines[1] and
// we check whether Intersection is on same side as BasePoints[0] or not.
const Vector pointinginside =
joiningline.getClosestPoint(BasePoints[0]) - BasePoints[0];
const Vector pointingtointersection =
joiningline.getClosestPoint(Intersection) - Intersection;
const double sign_of_intersection =
pointingtointersection.ScalarProduct(pointinginside);
LOG(4, "DEBUG: Sign of SKP between both lines w.r.t " << *OldLines[1]
<< " is " << sign_of_intersection << ".");
const double sign_of_normal = (sign_of_intersection >= 0) ? -1. : +1.;
LOG(4, "DEBUG: Opposite normal direction is "
<< sign_of_normal * BaseLineNormal[0] << ".");
BTS->GetNormalVector(sign_of_normal * BaseLineNormal[0]);
}
AddTesselationTriangle(OldTriangleNrs[1]);
LOG(3, "DEBUG: Created new triangle " << *BTS << ".");
} else {
ELOG(0, "The four old lines do not connect, something's utterly wrong here!");
return NULL;
}
return NewLine;
}
;
/** Finds the second point of starting triangle.
* \param *a first node
* \param Oben vector indicating the outside
* \param OptCandidate reference to recommended candidate on return
* \param Storage[3] array storing angles and other candidate information
* \param RADIUS radius of virtual sphere
* \param *LC LinkedCell_deprecated structure with neighbouring points
*/
void Tesselation::FindSecondPointForTesselation(TesselPoint* a, Vector Oben, TesselPoint*& OptCandidate, double Storage[3], double RADIUS, const LinkedCell_deprecated *LC)
{
//Info FunctionInfo(__func__);
Vector AngleCheck;
class TesselPoint* Candidate = NULL;
double norm = -1.;
double angle = 0.;
int N[NDIM];
int Nlower[NDIM];
int Nupper[NDIM];
if (LC->SetIndexToNode(a)) { // get cell for the starting point
for (int i = 0; i < NDIM; i++) // store indices of this cell
N[i] = LC->n[i];
} else {
ELOG(1, "Point " << *a << " is not found in cell " << LC->index << ".");
return;
}
// then go through the current and all neighbouring cells and check the contained points for possible candidates
for (int i = 0; i < NDIM; i++) {
Nlower[i] = ((N[i] - 1) >= 0) ? N[i] - 1 : 0;
Nupper[i] = ((N[i] + 1) < LC->N[i]) ? N[i] + 1 : LC->N[i] - 1;
}
LOG(3, "DEBUG: LC Intervals from [" << N[0] << "<->" << LC->N[0] << ", " << N[1] << "<->" << LC->N[1] << ", " << N[2] << "<->" << LC->N[2] << "] :" << " [" << Nlower[0] << "," << Nupper[0] << "], " << " [" << Nlower[1] << "," << Nupper[1] << "], " << " [" << Nlower[2] << "," << Nupper[2] << "], ");
for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
const TesselPointSTLList *List = LC->GetCurrentCell();
//LOG(1, "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << ".");
if (List != NULL) {
for (TesselPointSTLList::const_iterator Runner = List->begin(); Runner != List->end(); Runner++) {
Candidate = (*Runner);
// check if we only have one unique point yet ...
if (a != Candidate) {
// Calculate center of the circle with radius RADIUS through points a and Candidate
Vector OrthogonalizedOben, aCandidate, Center;
double distance, scaleFactor;
OrthogonalizedOben = Oben;
aCandidate = (a->getPosition()) - (Candidate->getPosition());
OrthogonalizedOben.ProjectOntoPlane(aCandidate);
OrthogonalizedOben.Normalize();
distance = 0.5 * aCandidate.Norm();
scaleFactor = sqrt(((RADIUS * RADIUS) - (distance * distance)));
OrthogonalizedOben.Scale(scaleFactor);
Center = 0.5 * ((Candidate->getPosition()) + (a->getPosition()));
Center += OrthogonalizedOben;
AngleCheck = Center - (a->getPosition());
norm = aCandidate.Norm();
// second point shall have smallest angle with respect to Oben vector
if (norm < RADIUS * 2.) {
angle = AngleCheck.Angle(Oben);
if (angle < Storage[0]) {
//LOG(1, "INFO: Old values of Storage is " << Storage[0] << ", " << Storage[1]);
LOG(4, "DEBUG: Current candidate is " << *Candidate << ": Is a better candidate with distance " << norm << " and angle " << angle << " to oben " << Oben << ".");
OptCandidate = Candidate;
Storage[0] = angle;
//LOG(4, "DEBUG: Changing something in Storage is " << Storage[0] << ", " << Storage[1]);
} else {
//LOG(4, "DEBUG: Current candidate is " << *Candidate << ": Looses with angle " << angle << " to a better candidate " << *OptCandidate);
}
} else {
//LOG(4, "DEBUG: Current candidate is " << *Candidate << ": Refused due to Radius " << norm);
}
} else {
//LOG(4, "DEBUG: Current candidate is " << *Candidate << ": Candidate is equal to first endpoint." << *a << ".");
}
}
} else {
LOG(4, "DEBUG: Linked cell list is empty.");
}
}
}
;
/** This recursive function finds a third point, to form a triangle with two given ones.
* Note that this function is for the starting triangle.
* The idea is as follows: A sphere with fixed radius is (almost) uniquely defined in space by three points
* that sit on its boundary. Hence, when two points are given and we look for the (next) third point, then
* the center of the sphere is still fixed up to a single parameter. The band of possible values
* describes a circle in 3D-space. The old center of the sphere for the current base triangle gives
* us the "null" on this circle, the new center of the candidate point will be some way along this
* circle. The shorter the way the better is the candidate. Note that the direction is clearly given
* by the normal vector of the base triangle that always points outwards by construction.
* Hence, we construct a Center of this circle which sits right in the middle of the current base line.
* We construct the normal vector that defines the plane this circle lies in, it is just in the
* direction of the baseline. And finally, we need the radius of the circle, which is given by the rest
* with respect to the length of the baseline and the sphere's fixed \a RADIUS.
* Note that there is one difficulty: The circumcircle is uniquely defined, but for the circumsphere's center
* there are two possibilities which becomes clear from the construction as seen below. Hence, we must check
* both.
* Note also that the acos() function is not unique on [0, 2.*M_PI). Hence, we need an additional check
* to decide for one of the two possible angles. Therefore we need a SearchDirection and to make this check
* sensible we need OldSphereCenter to be orthogonal to it. Either we construct SearchDirection orthogonal
* right away, or -- what we do here -- we rotate the relative sphere centers such that this orthogonality
* holds. Then, the normalized projection onto the SearchDirection is either +1 or -1 and thus states whether
* the angle is uniquely in either (0,M_PI] or [M_PI, 2.*M_PI).
* @param NormalVector normal direction of the base triangle (here the unit axis vector, \sa FindStartingTriangle())
* @param SearchDirection general direction where to search for the next point, relative to center of BaseLine
* @param OldSphereCenter center of sphere for base triangle, relative to center of BaseLine, giving null angle for the parameter circle
* @param CandidateLine CandidateForTesselation with the current base line and list of candidates and ShortestAngle
* @param ThirdPoint third point to avoid in search
* @param RADIUS radius of sphere
* @param *LC LinkedCell_deprecated structure with neighbouring points
*/
void Tesselation::FindThirdPointForTesselation(const Vector &NormalVector, const Vector &SearchDirection, const Vector &OldSphereCenter, CandidateForTesselation &CandidateLine, const class BoundaryPointSet * const ThirdPoint, const double RADIUS, const LinkedCell_deprecated *LC) const
{
//Info FunctionInfo(__func__);
Vector CircleCenter; // center of the circle, i.e. of the band of sphere's centers
Vector CirclePlaneNormal; // normal vector defining the plane this circle lives in
Vector SphereCenter;
Vector NewSphereCenter; // center of the sphere defined by the two points of BaseLine and the one of Candidate, first possibility
Vector OtherNewSphereCenter; // center of the sphere defined by the two points of BaseLine and the one of Candidate, second possibility
Vector NewNormalVector; // normal vector of the Candidate's triangle
Vector helper, OptCandidateCenter, OtherOptCandidateCenter;
Vector RelativeOldSphereCenter;
Vector NewPlaneCenter;
double CircleRadius; // radius of this circle
double radius;
double otherradius;
double alpha, Otheralpha; // angles (i.e. parameter for the circle).
int N[NDIM], Nlower[NDIM], Nupper[NDIM];
TesselPoint *Candidate = NULL;
LOG(4, "DEBUG: NormalVector of BaseTriangle is " << NormalVector << ".");
// copy old center
CandidateLine.OldCenter = OldSphereCenter;
CandidateLine.ThirdPoint = ThirdPoint;
CandidateLine.pointlist.clear();
// construct center of circle
CircleCenter = 0.5 * ((CandidateLine.BaseLine->endpoints[0]->node->getPosition()) + (CandidateLine.BaseLine->endpoints[1]->node->getPosition()));
// construct normal vector of circle
CirclePlaneNormal = (CandidateLine.BaseLine->endpoints[0]->node->getPosition()) - (CandidateLine.BaseLine->endpoints[1]->node->getPosition());
RelativeOldSphereCenter = OldSphereCenter - CircleCenter;
// calculate squared radius TesselPoint *ThirdPoint,f circle
radius = CirclePlaneNormal.NormSquared() / 4.;
if (radius < RADIUS * RADIUS) {
CircleRadius = RADIUS * RADIUS - radius;
CirclePlaneNormal.Normalize();
LOG(4, "DEBUG: CircleCenter is at " << CircleCenter << ", CirclePlaneNormal is " << CirclePlaneNormal << " with circle radius " << sqrt(CircleRadius) << ".");
// test whether old center is on the band's plane
if (fabs(RelativeOldSphereCenter.ScalarProduct(CirclePlaneNormal)) > HULLEPSILON) {
ELOG(1, "Something's very wrong here: RelativeOldSphereCenter is not on the band's plane as desired by " << fabs(RelativeOldSphereCenter.ScalarProduct(CirclePlaneNormal)) << "!");
RelativeOldSphereCenter.ProjectOntoPlane(CirclePlaneNormal);
}
radius = RelativeOldSphereCenter.NormSquared();
if (fabs(radius - CircleRadius) < HULLEPSILON) {
LOG(4, "DEBUG: RelativeOldSphereCenter is at " << RelativeOldSphereCenter << ".");
// check SearchDirection
LOG(4, "DEBUG: SearchDirection is " << SearchDirection << ".");
if (fabs(RelativeOldSphereCenter.ScalarProduct(SearchDirection)) > HULLEPSILON) { // rotated the wrong way!
ELOG(1, "SearchDirection and RelativeOldSphereCenter are not orthogonal!");
}
// get cell for the starting point
if (LC->SetIndexToVector(CircleCenter)) {
for (int i = 0; i < NDIM; i++) // store indices of this cell
N[i] = LC->n[i];
//LOG(1, "INFO: Center cell is " << N[0] << ", " << N[1] << ", " << N[2] << " with No. " << LC->index << ".");
} else {
ELOG(1, "Vector " << CircleCenter << " is outside of LinkedCell's bounding box.");
return;
}
// then go through the current and all neighbouring cells and check the contained points for possible candidates
// if (DoLog(3)) {
// std::stringstream output;
// output << "LC Intervals:";
// for (int i = 0; i < NDIM; i++)
// output << " [" << Nlower[i] << "," << Nupper[i] << "] ";
// LOG(0, output.str());
// }
for (int i = 0; i < NDIM; i++) {
Nlower[i] = ((N[i] - 1) >= 0) ? N[i] - 1 : 0;
Nupper[i] = ((N[i] + 1) < LC->N[i]) ? N[i] + 1 : LC->N[i] - 1;
}
for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
const TesselPointSTLList *List = LC->GetCurrentCell();
//LOG(1, "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << ".");
if (List != NULL) {
for (TesselPointSTLList::const_iterator Runner = List->begin(); Runner != List->end(); Runner++) {
Candidate = (*Runner);
// check for three unique points
LOG(4, "DEBUG: Current Candidate is " << *Candidate << " for BaseLine " << *CandidateLine.BaseLine << " with OldSphereCenter " << OldSphereCenter << ".");
if ((Candidate != CandidateLine.BaseLine->endpoints[0]->node) && (Candidate != CandidateLine.BaseLine->endpoints[1]->node)) {
// find center on the plane
GetCenterofCircumcircle(NewPlaneCenter, CandidateLine.BaseLine->endpoints[0]->node->getPosition(), CandidateLine.BaseLine->endpoints[1]->node->getPosition(), Candidate->getPosition());
LOG(5, "DEBUG: NewPlaneCenter is " << NewPlaneCenter << ".");
try {
NewNormalVector = Plane((CandidateLine.BaseLine->endpoints[0]->node->getPosition()), (CandidateLine.BaseLine->endpoints[1]->node->getPosition()), (Candidate->getPosition())).getNormal();
LOG(5, "DEBUG: NewNormalVector is " << NewNormalVector << ".");
radius = CandidateLine.BaseLine->endpoints[0]->node->DistanceSquared(NewPlaneCenter);
LOG(5, "DEBUG: CircleCenter is at " << CircleCenter << ", CirclePlaneNormal is " << CirclePlaneNormal << " with circle radius " << sqrt(CircleRadius) << ".");
LOG(5, "DEBUG: SearchDirection is " << SearchDirection << ".");
LOG(5, "DEBUG: Radius of CircumCenterCircle is " << radius << ".");
if (radius < RADIUS * RADIUS) {
otherradius = CandidateLine.BaseLine->endpoints[1]->node->DistanceSquared(NewPlaneCenter);
if (fabs(radius - otherradius) < HULLEPSILON) {
// construct both new centers
NewSphereCenter = NewPlaneCenter;
OtherNewSphereCenter = NewPlaneCenter;
helper = NewNormalVector;
helper.Scale(sqrt(RADIUS * RADIUS - radius));
LOG(5, "DEBUG: Distance of NewPlaneCenter " << NewPlaneCenter << " to either NewSphereCenter is " << helper.Norm() << " of vector " << helper << " with sphere radius " << RADIUS << ".");
NewSphereCenter += helper;
LOG(5, "DEBUG: NewSphereCenter is at " << NewSphereCenter << ".");
// OtherNewSphereCenter is created by the same vector just in the other direction
helper.Scale(-1.);
OtherNewSphereCenter += helper;
LOG(5, "DEBUG: OtherNewSphereCenter is at " << OtherNewSphereCenter << ".");
alpha = GetPathLengthonCircumCircle(CircleCenter, CirclePlaneNormal, CircleRadius, NewSphereCenter, OldSphereCenter, NormalVector, SearchDirection, HULLEPSILON);
Otheralpha = GetPathLengthonCircumCircle(CircleCenter, CirclePlaneNormal, CircleRadius, OtherNewSphereCenter, OldSphereCenter, NormalVector, SearchDirection, HULLEPSILON);
if ((ThirdPoint != NULL) && (Candidate == ThirdPoint->node)) { // in that case only the other circlecenter is valid
if (OldSphereCenter.DistanceSquared(NewSphereCenter) < OldSphereCenter.DistanceSquared(OtherNewSphereCenter))
alpha = Otheralpha;
} else
alpha = min(alpha, Otheralpha);
// if there is a better candidate, drop the current list and add the new candidate
// otherwise ignore the new candidate and keep the list
if (CandidateLine.ShortestAngle > (alpha - HULLEPSILON)) {
if (fabs(alpha - Otheralpha) > MYEPSILON) {
CandidateLine.OptCenter = NewSphereCenter;
CandidateLine.OtherOptCenter = OtherNewSphereCenter;
} else {
CandidateLine.OptCenter = OtherNewSphereCenter;
CandidateLine.OtherOptCenter = NewSphereCenter;
}
// if there is an equal candidate, add it to the list without clearing the list
if ((CandidateLine.ShortestAngle - HULLEPSILON) < alpha) {
CandidateLine.pointlist.push_back(Candidate);
LOG(3, "ACCEPT: We have found an equally good candidate: " << *(Candidate) << " with " << alpha << " and circumsphere's center at " << CandidateLine.OptCenter << ".");
} else {
// remove all candidates from the list and then the list itself
CandidateLine.pointlist.clear();
CandidateLine.pointlist.push_back(Candidate);
LOG(3, "ACCEPT: We have found a better candidate: " << *(Candidate) << " with " << alpha << " and circumsphere's center at " << CandidateLine.OptCenter << ".");
}
CandidateLine.ShortestAngle = alpha;
LOG(4, "DEBUG: There are " << CandidateLine.pointlist.size() << " candidates in the list now.");
} else {
if ((Candidate != NULL) && (CandidateLine.pointlist.begin() != CandidateLine.pointlist.end())) {
LOG(4, "REJECT: Old candidate " << *(*CandidateLine.pointlist.begin()) << " with " << CandidateLine.ShortestAngle << " is better than new one " << *Candidate << " with " << alpha << " .");
} else {
LOG(4, "REJECT: Candidate " << *Candidate << " with " << alpha << " was rejected.");
}
}
} else {
ASSERT(0, std::string("FindThirdPointForTesselation() - ") + std::string("REJECT: Distance to center of circumcircle is not the same from each corner of the triangle: ") + toString(fabs(radius - otherradius)));
}
} else {
LOG(4, "REJECT: NewSphereCenter " << NewSphereCenter << " for " << *Candidate << " is too far away: " << radius << ".");
}
} catch (LinearDependenceException &excp) {
LOG(4, boost::diagnostic_information(excp));
LOG(4, "REJECT: Three points from " << *CandidateLine.BaseLine << " and Candidate " << *Candidate << " are linear-dependent.");
}
} else {
if (ThirdPoint != NULL) {
LOG(4, "REJECT: Base triangle " << *CandidateLine.BaseLine << " and " << *ThirdPoint << " contains Candidate " << *Candidate << ".");
} else {
LOG(4, "REJECT: Base triangle " << *CandidateLine.BaseLine << " contains Candidate " << *Candidate << ".");
}
}
}
}
}
} else {
ELOG(1, "The projected center of the old sphere has radius " << radius << " instead of " << CircleRadius << ".");
}
} else {
if (ThirdPoint != NULL)
LOG(4, "DEBUG: Circumcircle for base line " << *CandidateLine.BaseLine << " and third node " << *ThirdPoint << " is too big!");
else
LOG(4, "DEBUG: Circumcircle for base line " << *CandidateLine.BaseLine << " is too big!");
}
LOG(5, "DEBUG: Sorting candidate list ...");
if (CandidateLine.pointlist.size() > 1) {
CandidateLine.pointlist.unique();
CandidateLine.pointlist.sort(); //SortCandidates);
}
ASSERT(CandidateLine.pointlist.empty() || (CandidateLine.CheckValidity(RADIUS, LC)), std::string("Tesselation::FindThirdPointForTesselation()") + std::string("There were other points contained in the rolling sphere as well!"));
}
;
/** Finds the endpoint two lines are sharing.
* \param *line1 first line
* \param *line2 second line
* \return point which is shared or NULL if none
*/
class BoundaryPointSet *Tesselation::GetCommonEndpoint(const BoundaryLineSet * line1, const BoundaryLineSet * line2) const
{
//Info FunctionInfo(__func__);
const BoundaryLineSet * lines[2] = {line1, line2};
class BoundaryPointSet *node = NULL;
PointMap OrderMap;
PointTestPair OrderTest;
for (int i = 0; i < 2; i++)
// for both lines
for (int j = 0; j < 2; j++) { // for both endpoints
OrderTest = OrderMap.insert(pair(lines[i]->endpoints[j]->Nr, lines[i]->endpoints[j]));
if (!OrderTest.second) { // if insertion fails, we have common endpoint
node = OrderTest.first->second;
LOG(1, "Common endpoint of lines " << *line1 << " and " << *line2 << " is: " << *node << ".");
j = 2;
i = 2;
break;
}
}
return node;
}
;
/** Finds the boundary points that are closest to a given Vector \a *x.
* \param *out output stream for debugging
* \param *x Vector to look from
* \return map of BoundaryPointSet of closest points sorted by squared distance or NULL.
*/
DistanceToPointMap * Tesselation::FindClosestBoundaryPointsToVector(const Vector &x, const LinkedCell_deprecated* LC) const
{
//Info FunctionInfo(__func__);
PointMap::const_iterator FindPoint;
int N[NDIM], Nlower[NDIM], Nupper[NDIM];
if (LinesOnBoundary.empty()) {
ELOG(1, "There is no tesselation structure to compare the point with, please create one first.");
return NULL;
}
// gather all points close to the desired one
LC->SetIndexToVector(x); // ignore status as we calculate bounds below sensibly
for (int i = 0; i < NDIM; i++) // store indices of this cell
N[i] = LC->n[i];
LOG(2, "DEBUG: Center cell is " << N[0] << ", " << N[1] << ", " << N[2] << " with No. " << LC->index << ".");
DistanceToPointMap * points = new DistanceToPointMap;
LC->GetNeighbourBounds(Nlower, Nupper);
for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
const TesselPointSTLList *List = LC->GetCurrentCell();
//LOG(1, "The current cell " << LC->n[0] << "," << LC->n[1] << "," << LC->n[2]);
if (List != NULL) {
for (TesselPointSTLList::const_iterator Runner = List->begin(); Runner != List->end(); Runner++) {
FindPoint = PointsOnBoundary.find((*Runner)->getNr());
if (FindPoint != PointsOnBoundary.end()) {
// when the closest point is on the edge of a triangle (and hence
// we find two closes triangles due to it having an adjacent one)
// we should make sure that both triangles end up in the same entry
// in the distance multimap. Hence, we round to 6 digit precision.
const double distance = 1e-6 * floor(FindPoint->second->node->DistanceSquared(x) * 1e+6);
points->insert(DistanceToPointPair(distance, FindPoint->second));
LOG(3, "DEBUG: Putting " << *FindPoint->second << " into the list.");
}
}
} else {
ELOG(1, "The current cell " << LC->n[0] << "," << LC->n[1] << "," << LC->n[2] << " is invalid!");
}
}
// check whether we found some points
if (points->empty()) {
ELOG(1, "There is no nearest point: too far away from the surface.");
delete (points);
return NULL;
}
return points;
}
;
/** Finds the boundary line that is closest to a given Vector \a *x.
* \param *out output stream for debugging
* \param *x Vector to look from
* \return closest BoundaryLineSet or NULL in degenerate case.
*/
BoundaryLineSet * Tesselation::FindClosestBoundaryLineToVector(const Vector &x, const LinkedCell_deprecated* LC) const
{
//Info FunctionInfo(__func__);
// get closest points
DistanceToPointMap * points = FindClosestBoundaryPointsToVector(x, LC);
if (points == NULL) {
ELOG(1, "There is no nearest point: too far away from the surface.");
return NULL;
}
// for each point, check its lines, remember closest
LOG(1, "Finding closest BoundaryLine to " << x << " ... ");
BoundaryLineSet *ClosestLine = NULL;
double MinDistance = -1.;
Vector helper;
Vector Center;
Vector BaseLine;
for (DistanceToPointMap::iterator Runner = points->begin(); Runner != points->end(); Runner++) {
for (LineMap::iterator LineRunner = Runner->second->lines.begin(); LineRunner != Runner->second->lines.end(); LineRunner++) {
// calculate closest point on line to desired point
helper = 0.5 * (((LineRunner->second)->endpoints[0]->node->getPosition()) + ((LineRunner->second)->endpoints[1]->node->getPosition()));
Center = (x) - helper;
BaseLine = ((LineRunner->second)->endpoints[0]->node->getPosition()) - ((LineRunner->second)->endpoints[1]->node->getPosition());
Center.ProjectOntoPlane(BaseLine);
const double distance = Center.NormSquared();
if ((ClosestLine == NULL) || (distance < MinDistance)) {
// additionally calculate intersection on line (whether it's on the line section or not)
helper = (x) - ((LineRunner->second)->endpoints[0]->node->getPosition()) - Center;
const double lengthA = helper.ScalarProduct(BaseLine);
helper = (x) - ((LineRunner->second)->endpoints[1]->node->getPosition()) - Center;
const double lengthB = helper.ScalarProduct(BaseLine);
if (lengthB * lengthA < 0) { // if have different sign
ClosestLine = LineRunner->second;
MinDistance = distance;
LOG(1, "ACCEPT: New closest line is " << *ClosestLine << " with projected distance " << MinDistance << ".");
} else {
LOG(1, "REJECT: Intersection is outside of the line section: " << lengthA << " and " << lengthB << ".");
}
} else {
LOG(1, "REJECT: Point is too further away than present line: " << distance << " >> " << MinDistance << ".");
}
}
}
delete (points);
// check whether closest line is "too close" :), then it's inside
if (ClosestLine == NULL) {
LOG(2, "DEBUG: Is the only point, no one else is closeby.");
return NULL;
}
return ClosestLine;
}
;
/** Finds the triangle that is closest to a given Vector \a *x.
* \param *out output stream for debugging
* \param *x Vector to look from
* \return BoundaryTriangleSet of nearest triangle or NULL.
*/
TriangleList * Tesselation::FindClosestTrianglesToVector(const Vector &x, const LinkedCell_deprecated* LC) const
{
//Info FunctionInfo(__func__);
// get closest points
DistanceToPointMap * points = FindClosestBoundaryPointsToVector(x, LC);
if (points == NULL) {
ELOG(1, "There is no nearest point: too far away from the surface.");
return NULL;
}
// for each point, check its lines, remember closest
LOG(1, "Finding closest BoundaryTriangle to " << x << " ... ");
LineSet ClosestLines;
double MinDistance = 1e+16;
Vector BaseLineIntersection;
Vector Center;
Vector BaseLine;
Vector BaseLineCenter;
for (DistanceToPointMap::iterator Runner = points->begin(); Runner != points->end(); Runner++) {
for (LineMap::iterator LineRunner = Runner->second->lines.begin(); LineRunner != Runner->second->lines.end(); LineRunner++) {
BaseLine = ((LineRunner->second)->endpoints[0]->node->getPosition()) - ((LineRunner->second)->endpoints[1]->node->getPosition());
const double lengthBase = BaseLine.NormSquared();
BaseLineIntersection = (x) - ((LineRunner->second)->endpoints[0]->node->getPosition());
const double lengthEndA = BaseLineIntersection.NormSquared();
BaseLineIntersection = (x) - ((LineRunner->second)->endpoints[1]->node->getPosition());
const double lengthEndB = BaseLineIntersection.NormSquared();
if ((lengthEndA > lengthBase) || (lengthEndB > lengthBase) || ((lengthEndA < MYEPSILON) || (lengthEndB < MYEPSILON))) { // intersection would be outside, take closer endpoint
const double lengthEnd = std::min(lengthEndA, lengthEndB);
if (lengthEnd - MinDistance < -MYEPSILON) { // new best line
ClosestLines.clear();
ClosestLines.insert(LineRunner->second);
MinDistance = lengthEnd;
LOG(1, "ACCEPT: Line " << *LineRunner->second << " to endpoint " << *LineRunner->second->endpoints[0]->node << " is closer with " << lengthEnd << ".");
} else
if (fabs(lengthEnd - MinDistance) < MYEPSILON) { // additional best candidate
ClosestLines.insert(LineRunner->second);
LOG(1, "ACCEPT: Line " << *LineRunner->second << " to endpoint " << *LineRunner->second->endpoints[1]->node << " is equally good with " << lengthEnd << ".");
} else { // line is worse
LOG(1, "REJECT: Line " << *LineRunner->second << " to either endpoints is further away than present closest line candidate: " << lengthEndA << ", " << lengthEndB << ", and distance is longer than baseline:" << lengthBase << ".");
}
} else { // intersection is closer, calculate
// calculate closest point on line to desired point
BaseLineIntersection = (x) - ((LineRunner->second)->endpoints[1]->node->getPosition());
Center = BaseLineIntersection;
Center.ProjectOntoPlane(BaseLine);
BaseLineIntersection -= Center;
const double distance = BaseLineIntersection.NormSquared();
if (Center.NormSquared() > BaseLine.NormSquared()) {
ELOG(0, "Algorithmic error: In second case we have intersection outside of baseline!");
}
if ((ClosestLines.empty()) || (distance < MinDistance)) {
ClosestLines.insert(LineRunner->second);
MinDistance = distance;
LOG(1, "ACCEPT: Intersection in between endpoints, new closest line " << *LineRunner->second << " is " << *ClosestLines.begin() << " with projected distance " << MinDistance << ".");
} else {
LOG(2, "REJECT: Point is further away from line " << *LineRunner->second << " than present closest line: " << distance << " >> " << MinDistance << ".");
}
}
}
}
delete (points);
// check whether closest line is "too close" :), then it's inside
if (ClosestLines.empty()) {
LOG(2, "DEBUG: Is the only point, no one else is closeby.");
return NULL;
}
TriangleList * candidates = new TriangleList;
for (LineSet::iterator LineRunner = ClosestLines.begin(); LineRunner != ClosestLines.end(); LineRunner++)
for (TriangleMap::iterator Runner = (*LineRunner)->triangles.begin(); Runner != (*LineRunner)->triangles.end(); Runner++) {
candidates->push_back(Runner->second);
}
return candidates;
}
;
/** Finds closest triangle to a point.
* This basically just takes care of the degenerate case, which is not handled in FindClosestTrianglesToPoint().
* \param *out output stream for debugging
* \param *x Vector to look from
* \param &distance contains found distance on return
* \return list of BoundaryTriangleSet of nearest triangles or NULL.
*/
class BoundaryTriangleSet * Tesselation::FindClosestTriangleToVector(const Vector &x, const LinkedCell_deprecated* LC) const
{
//Info FunctionInfo(__func__);
class BoundaryTriangleSet *result = NULL;
TriangleList *triangles = FindClosestTrianglesToVector(x, LC);
TriangleList candidates;
Vector Center;
Vector helper;
if ((triangles == NULL) || (triangles->empty()))
return NULL;
// go through all and pick the one with the best alignment to x
double MinAlignment = 2. * M_PI;
for (TriangleList::iterator Runner = triangles->begin(); Runner != triangles->end(); Runner++) {
(*Runner)->GetCenter(Center);
helper = (x) - Center;
const double Alignment = helper.Angle((*Runner)->NormalVector);
if (Alignment < MinAlignment) {
result = *Runner;
MinAlignment = Alignment;
LOG(1, "ACCEPT: Triangle " << *result << " is better aligned with " << MinAlignment << ".");
} else {
LOG(1, "REJECT: Triangle " << *result << " is worse aligned with " << MinAlignment << ".");
}
}
delete (triangles);
return result;
}
;
/** Checks whether the provided Vector is within the Tesselation structure.
* Basically calls Tesselation::GetDistanceToSurface() and checks the sign of the return value.
* @param point of which to check the position
* @param *LC LinkedCell_deprecated structure
*
* @return true if the point is inside the Tesselation structure, false otherwise
*/
bool Tesselation::IsInnerPoint(const Vector &Point, const LinkedCell_deprecated* const LC) const
{
TriangleIntersectionList Intersections(Point, this, LC);
return Intersections.IsInside();
}
Vector Tesselation::getNormal(const Vector &Point, const LinkedCell_deprecated* const LC) const
{
TriangleIntersectionList Intersections(Point, this, LC);
BoundaryTriangleSet *triangle = Intersections.GetClosestTriangle();
if (triangle != NULL) {
return triangle->NormalVector;
} else
return zeroVec;
}
/** Returns the distance to the surface given by the tesselation.
* Calls FindClosestTriangleToVector() and checks whether the resulting triangle's BoundaryTriangleSet#NormalVector points
* towards or away from the given \a &Point. Additionally, we check whether it's normal to the normal vector, i.e. on the
* closest triangle's plane. Then, we have to check whether \a Point is inside the triangle or not to determine whether it's
* an inside or outside point. This is done by calling BoundaryTriangleSet::GetIntersectionInsideTriangle().
* In the end we additionally find the point on the triangle who was smallest distance to \a Point:
* -# Separate distance from point to center in vector in NormalDirection and on the triangle plane.
* -# Check whether vector on triangle plane points inside the triangle or crosses triangle bounds.
* -# If inside, take it to calculate closest distance
* -# If not, take intersection with BoundaryLine as distance
*
* @note distance is squared despite it still contains a sign to determine in-/outside!
*
* @param point of which to check the position
* @param *LC LinkedCell_deprecated structure
*
* @return >0 if outside, ==0 if on surface, <0 if inside
*/
double Tesselation::GetDistanceSquaredToTriangle(const Vector &Point, const BoundaryTriangleSet* const triangle) const
{
//Info FunctionInfo(__func__);
Vector Center;
Vector helper;
Vector DistanceToCenter;
Vector Intersection;
double distance = 0.;
if (triangle == NULL) { // is boundary point or only point in point cloud?
LOG(1, "No triangle given!");
return -1.;
} else {
LOG(1, "INFO: Closest triangle found is " << *triangle << " with normal vector " << triangle->NormalVector << ".");
}
triangle->GetCenter(Center);
LOG(2, "INFO: Central point of the triangle is " << Center << ".");
DistanceToCenter = Center - Point;
LOG(2, "INFO: Vector from point to test to center is " << DistanceToCenter << ".");
// check whether we are on boundary
if (fabs(DistanceToCenter.ScalarProduct(triangle->NormalVector)) < MYEPSILON) {
// calculate whether inside of triangle
DistanceToCenter = Point + triangle->NormalVector; // points outside
Center = Point - triangle->NormalVector; // points towards MolCenter
LOG(1, "INFO: Calling Intersection with " << Center << " and " << DistanceToCenter << ".");
if (triangle->GetIntersectionInsideTriangle(Center, DistanceToCenter, Intersection)) {
LOG(1, Point << " is inner point: sufficiently close to boundary, " << Intersection << ".");
return 0.;
} else {
LOG(1, Point << " is NOT an inner point: on triangle plane but outside of triangle bounds.");
return false;
}
} else {
// calculate smallest distance
distance = triangle->GetClosestPointInsideTriangle(Point, Intersection);
LOG(1, "Closest point on triangle is " << Intersection << ".");
// then check direction to boundary
if (DistanceToCenter.ScalarProduct(triangle->NormalVector) > MYEPSILON) {
LOG(1, Point << " is an inner point, " << distance << " below surface.");
return -distance;
} else {
LOG(1, Point << " is NOT an inner point, " << distance << " above surface.");
return +distance;
}
}
}
;
/** Calculates minimum distance from \a&Point to a tesselated surface.
* Combines \sa FindClosestTrianglesToVector() and \sa GetDistanceSquaredToTriangle().
* \param &Point point to calculate distance from
* \param *LC needed for finding closest points fast
* \return distance squared to closest point on surface
*/
double Tesselation::GetDistanceToSurface(const Vector &Point, const LinkedCell_deprecated* const LC) const
{
//Info FunctionInfo(__func__);
TriangleIntersectionList Intersections(Point, this, LC);
return Intersections.GetSmallestDistance();
}
;
/** Calculates minimum distance from \a&Point to a tesselated surface.
* Combines \sa FindClosestTrianglesToVector() and \sa GetDistanceSquaredToTriangle().
* \param &Point point to calculate distance from
* \param *LC needed for finding closest points fast
* \return distance squared to closest point on surface
*/
BoundaryTriangleSet * Tesselation::GetClosestTriangleOnSurface(const Vector &Point, const LinkedCell_deprecated* const LC) const
{
//Info FunctionInfo(__func__);
TriangleIntersectionList Intersections(Point, this, LC);
return Intersections.GetClosestTriangle();
}
;
/** Gets all points connected to the provided point by triangulation lines.
*
* @param *Point of which get all connected points
*
* @return set of the all points linked to the provided one
*/
TesselPointSet * Tesselation::GetAllConnectedPoints(const TesselPoint* const Point) const
{
//Info FunctionInfo(__func__);
TesselPointSet *connectedPoints = new TesselPointSet;
class BoundaryPointSet *ReferencePoint = NULL;
TesselPoint* current;
bool takePoint = false;
// find the respective boundary point
PointMap::const_iterator PointRunner = PointsOnBoundary.find(Point->getNr());
if (PointRunner != PointsOnBoundary.end()) {
ReferencePoint = PointRunner->second;
} else {
ELOG(2, "GetAllConnectedPoints() could not find the BoundaryPoint belonging to " << *Point << ".");
ReferencePoint = NULL;
}
// little trick so that we look just through lines connect to the BoundaryPoint
// OR fall-back to look through all lines if there is no such BoundaryPoint
const LineMap *Lines;
;
if (ReferencePoint != NULL)
Lines = &(ReferencePoint->lines);
else
Lines = &LinesOnBoundary;
LineMap::const_iterator findLines = Lines->begin();
while (findLines != Lines->end()) {
takePoint = false;
if (findLines->second->endpoints[0]->Nr == Point->getNr()) {
takePoint = true;
current = findLines->second->endpoints[1]->node;
} else
if (findLines->second->endpoints[1]->Nr == Point->getNr()) {
takePoint = true;
current = findLines->second->endpoints[0]->node;
}
if (takePoint) {
LOG(1, "INFO: Endpoint " << *current << " of line " << *(findLines->second) << " is enlisted.");
connectedPoints->insert(current);
}
findLines++;
}
if (connectedPoints->empty()) { // if have not found any points
ELOG(1, "We have not found any connected points to " << *Point << ".");
return NULL;
}
return connectedPoints;
}
;
/** Gets all points connected to the provided point by triangulation lines, ordered such that we have the circle round the point.
* Maps them down onto the plane designated by the axis \a *Point and \a *Reference. The center of all points
* connected in the tesselation to \a *Point is mapped to spherical coordinates with the zero angle being given
* by the mapped down \a *Reference. Hence, the biggest and the smallest angles are those of the two shanks of the
* triangle we are looking for.
*
* @param *out output stream for debugging
* @param *SetOfNeighbours all points for which the angle should be calculated
* @param *Point of which get all connected points
* @param *Reference Reference vector for zero angle or NULL for no preference
* @return list of the all points linked to the provided one
*/
TesselPointList * Tesselation::GetCircleOfConnectedTriangles(TesselPointSet *SetOfNeighbours, const TesselPoint* const Point, const Vector &Reference) const
{
//Info FunctionInfo(__func__);
map anglesOfPoints;
TesselPointList *connectedCircle = new TesselPointList;
Vector PlaneNormal;
Vector AngleZero;
Vector OrthogonalVector;
Vector helper;
const TesselPoint * const TrianglePoints[3] = {Point, NULL, NULL};
TriangleList *triangles = NULL;
if (SetOfNeighbours == NULL) {
ELOG(2, "Could not find any connected points!");
delete (connectedCircle);
return NULL;
}
// calculate central point
triangles = FindTriangles(TrianglePoints);
ASSERT((triangles == NULL) || (triangles->empty()), std::string("Tesselation::GetCircleOfConnectedTriangles()") + std::string("Could not find any triangles for point " + toString(*Point) + "."));
for (TriangleList::iterator Runner = triangles->begin(); Runner != triangles->end(); Runner++)
PlaneNormal += (*Runner)->NormalVector;
PlaneNormal.Scale(1.0 / triangles->size());
LOG(4, "DEBUG: Calculated PlaneNormal of all circle points is " << PlaneNormal << ".");
PlaneNormal.Normalize();
// construct one orthogonal vector
AngleZero = (Reference) - (Point->getPosition());
AngleZero.ProjectOntoPlane(PlaneNormal);
if ((AngleZero.NormSquared() < MYEPSILON)) {
LOG(4, "DEBUG: Using alternatively " << (*SetOfNeighbours->begin())->getPosition() << " as angle 0 referencer.");
AngleZero = ((*SetOfNeighbours->begin())->getPosition()) - (Point->getPosition());
AngleZero.ProjectOntoPlane(PlaneNormal);
ASSERT(AngleZero.NormSquared() > MYEPSILON, std::string("Tesselation::GetCircleOfConnectedTriangles() - ") + std::string("AngleZero is 0 even with alternative reference.") + std::string("The algorithm has to be changed here!"));
}
LOG(4, "DEBUG: Reference vector on this plane representing angle 0 is " << AngleZero << ".");
if (AngleZero.NormSquared() > MYEPSILON)
OrthogonalVector = Plane(PlaneNormal, AngleZero, 0).getNormal();
else
OrthogonalVector.MakeNormalTo(PlaneNormal);
LOG(4, "DEBUG: OrthogonalVector on plane is " << OrthogonalVector << ".");
// go through all connected points and calculate angle
for (TesselPointSet::iterator listRunner = SetOfNeighbours->begin(); listRunner != SetOfNeighbours->end(); listRunner++) {
helper = ((*listRunner)->getPosition()) - (Point->getPosition());
helper.ProjectOntoPlane(PlaneNormal);
double angle = GetAngle(helper, AngleZero, OrthogonalVector);
LOG(4, "DEBUG" << angle << " for point " << **listRunner << ".");
anglesOfPoints.insert(pair(angle, (*listRunner)));
}
for (map::iterator AngleRunner = anglesOfPoints.begin(); AngleRunner != anglesOfPoints.end(); AngleRunner++) {
connectedCircle->push_back(AngleRunner->second);
}
return connectedCircle;
}
/** Gets all points connected to the provided point by triangulation lines, ordered such that we have the circle round the point.
* Maps them down onto the plane designated by the axis \a *Point and \a *Reference. The center of all points
* connected in the tesselation to \a *Point is mapped to spherical coordinates with the zero angle being given
* by the mapped down \a *Reference. Hence, the biggest and the smallest angles are those of the two shanks of the
* triangle we are looking for.
*
* @param *SetOfNeighbours all points for which the angle should be calculated
* @param *Point of which get all connected points
* @param *Reference Reference vector for zero angle or (0,0,0) for no preference
* @return list of the all points linked to the provided one
*/
TesselPointList * Tesselation::GetCircleOfSetOfPoints(TesselPointSet *SetOfNeighbours, const TesselPoint* const Point, const Vector &Reference) const
{
//Info FunctionInfo(__func__);
map anglesOfPoints;
TesselPointList *connectedCircle = new TesselPointList;
Vector center;
Vector PlaneNormal;
Vector AngleZero;
Vector OrthogonalVector;
Vector helper;
if (SetOfNeighbours == NULL) {
ELOG(2, "Could not find any connected points!");
delete (connectedCircle);
return NULL;
}
// check whether there's something to do
if (SetOfNeighbours->size() < 3) {
for (TesselPointSet::iterator TesselRunner = SetOfNeighbours->begin(); TesselRunner != SetOfNeighbours->end(); TesselRunner++)
connectedCircle->push_back(*TesselRunner);
return connectedCircle;
}
LOG(1, "INFO: Point is " << *Point << " and Reference is " << Reference << ".");
// calculate central point
TesselPointSet::const_iterator TesselA = SetOfNeighbours->begin();
TesselPointSet::const_iterator TesselB = SetOfNeighbours->begin();
TesselPointSet::const_iterator TesselC = SetOfNeighbours->begin();
TesselB++;
TesselC++;
TesselC++;
int counter = 0;
while (TesselC != SetOfNeighbours->end()) {
helper = Plane(((*TesselA)->getPosition()), ((*TesselB)->getPosition()), ((*TesselC)->getPosition())).getNormal();
LOG(5, "DEBUG: Making normal vector out of " << *(*TesselA) << ", " << *(*TesselB) << " and " << *(*TesselC) << ":" << helper);
counter++;
TesselA++;
TesselB++;
TesselC++;
PlaneNormal += helper;
}
//LOG(0, "Summed vectors " << center << "; number of points " << connectedPoints.size() << "; scale factor " << counter);
PlaneNormal.Scale(1.0 / (double)counter);
// LOG(1, "INFO: Calculated center of all circle points is " << center << ".");
//
// // projection plane of the circle is at the closes Point and normal is pointing away from center of all circle points
// PlaneNormal.CopyVector(Point->node);
// PlaneNormal.SubtractVector(¢er);
// PlaneNormal.Normalize();
LOG(4, "DEBUG: Calculated plane normal of circle is " << PlaneNormal << ".");
// construct one orthogonal vector
if (!Reference.IsZero()) {
AngleZero = (Reference) - (Point->getPosition());
AngleZero.ProjectOntoPlane(PlaneNormal);
}
if ((Reference.IsZero()) || (AngleZero.NormSquared() < MYEPSILON)) {
LOG(4, "DEBUG: Using alternatively " << (*SetOfNeighbours->begin())->getPosition() << " as angle 0 referencer.");
AngleZero = ((*SetOfNeighbours->begin())->getPosition()) - (Point->getPosition());
AngleZero.ProjectOntoPlane(PlaneNormal);
ASSERT(AngleZero.NormSquared() > MYEPSILON, std::string("Tesselation::GetCircleOfSetOfPoints() - ") + std::string("AngleZero is 0 even with alternative reference.") + std::string("The algorithm has to be changed here!"));
}
LOG(4, "DEBUG: Reference vector on this plane representing angle 0 is " << AngleZero << ".");
if (AngleZero.NormSquared() > MYEPSILON)
OrthogonalVector = Plane(PlaneNormal, AngleZero, 0).getNormal();
else
OrthogonalVector.MakeNormalTo(PlaneNormal);
LOG(4, "DEBUG: OrthogonalVector on plane is " << OrthogonalVector << ".");
// go through all connected points and calculate angle
pair