[357fba] | 1 | /*
|
---|
| 2 | * TesselationHelpers.cpp
|
---|
| 3 | *
|
---|
| 4 | * Created on: Aug 3, 2009
|
---|
| 5 | * Author: heber
|
---|
| 6 | */
|
---|
| 7 |
|
---|
[112b09] | 8 | #include "Helpers/MemDebug.hpp"
|
---|
| 9 |
|
---|
[f66195] | 10 | #include <fstream>
|
---|
| 11 |
|
---|
[f67b6e] | 12 | #include "info.hpp"
|
---|
[f66195] | 13 | #include "linkedcell.hpp"
|
---|
[e138de] | 14 | #include "log.hpp"
|
---|
[f66195] | 15 | #include "tesselation.hpp"
|
---|
[357fba] | 16 | #include "tesselationhelpers.hpp"
|
---|
[f66195] | 17 | #include "vector.hpp"
|
---|
[643e76] | 18 | #include "Line.hpp"
|
---|
[0a4f7f] | 19 | #include "vector_ops.hpp"
|
---|
[f66195] | 20 | #include "verbose.hpp"
|
---|
[d4c9ae] | 21 | #include "Plane.hpp"
|
---|
[357fba] | 22 |
|
---|
[f67b6e] | 23 | double DetGet(gsl_matrix * const A, const int inPlace)
|
---|
| 24 | {
|
---|
| 25 | Info FunctionInfo(__func__);
|
---|
[357fba] | 26 | /*
|
---|
| 27 | inPlace = 1 => A is replaced with the LU decomposed copy.
|
---|
| 28 | inPlace = 0 => A is retained, and a copy is used for LU.
|
---|
| 29 | */
|
---|
| 30 |
|
---|
| 31 | double det;
|
---|
| 32 | int signum;
|
---|
| 33 | gsl_permutation *p = gsl_permutation_alloc(A->size1);
|
---|
[24a5e0] | 34 | gsl_matrix *tmpA=0;
|
---|
[357fba] | 35 |
|
---|
| 36 | if (inPlace)
|
---|
| 37 | tmpA = A;
|
---|
| 38 | else {
|
---|
| 39 | gsl_matrix *tmpA = gsl_matrix_alloc(A->size1, A->size2);
|
---|
| 40 | gsl_matrix_memcpy(tmpA , A);
|
---|
| 41 | }
|
---|
| 42 |
|
---|
| 43 |
|
---|
| 44 | gsl_linalg_LU_decomp(tmpA , p , &signum);
|
---|
| 45 | det = gsl_linalg_LU_det(tmpA , signum);
|
---|
| 46 | gsl_permutation_free(p);
|
---|
| 47 | if (! inPlace)
|
---|
| 48 | gsl_matrix_free(tmpA);
|
---|
| 49 |
|
---|
| 50 | return det;
|
---|
| 51 | };
|
---|
| 52 |
|
---|
[c0f6c6] | 53 | void GetSphere(Vector * const center, const Vector &a, const Vector &b, const Vector &c, const double RADIUS)
|
---|
[357fba] | 54 | {
|
---|
[f67b6e] | 55 | Info FunctionInfo(__func__);
|
---|
[357fba] | 56 | gsl_matrix *A = gsl_matrix_calloc(3,3);
|
---|
| 57 | double m11, m12, m13, m14;
|
---|
| 58 |
|
---|
| 59 | for(int i=0;i<3;i++) {
|
---|
[0a4f7f] | 60 | gsl_matrix_set(A, i, 0, a[i]);
|
---|
| 61 | gsl_matrix_set(A, i, 1, b[i]);
|
---|
| 62 | gsl_matrix_set(A, i, 2, c[i]);
|
---|
[357fba] | 63 | }
|
---|
[f1cccd] | 64 | m11 = DetGet(A, 1);
|
---|
[357fba] | 65 |
|
---|
| 66 | for(int i=0;i<3;i++) {
|
---|
[0a4f7f] | 67 | gsl_matrix_set(A, i, 0, a[i]*a[i] + b[i]*b[i] + c[i]*c[i]);
|
---|
| 68 | gsl_matrix_set(A, i, 1, b[i]);
|
---|
| 69 | gsl_matrix_set(A, i, 2, c[i]);
|
---|
[357fba] | 70 | }
|
---|
[f1cccd] | 71 | m12 = DetGet(A, 1);
|
---|
[357fba] | 72 |
|
---|
| 73 | for(int i=0;i<3;i++) {
|
---|
[0a4f7f] | 74 | gsl_matrix_set(A, i, 0, a[i]*a[i] + b[i]*b[i] + c[i]*c[i]);
|
---|
| 75 | gsl_matrix_set(A, i, 1, a[i]);
|
---|
| 76 | gsl_matrix_set(A, i, 2, c[i]);
|
---|
[357fba] | 77 | }
|
---|
[f1cccd] | 78 | m13 = DetGet(A, 1);
|
---|
[357fba] | 79 |
|
---|
| 80 | for(int i=0;i<3;i++) {
|
---|
[0a4f7f] | 81 | gsl_matrix_set(A, i, 0, a[i]*a[i] + b[i]*b[i] + c[i]*c[i]);
|
---|
| 82 | gsl_matrix_set(A, i, 1, a[i]);
|
---|
| 83 | gsl_matrix_set(A, i, 2, b[i]);
|
---|
[357fba] | 84 | }
|
---|
[f1cccd] | 85 | m14 = DetGet(A, 1);
|
---|
[357fba] | 86 |
|
---|
| 87 | if (fabs(m11) < MYEPSILON)
|
---|
[58ed4a] | 88 | DoeLog(1) && (eLog()<< Verbose(1) << "three points are colinear." << endl);
|
---|
[357fba] | 89 |
|
---|
[0a4f7f] | 90 | center->at(0) = 0.5 * m12/ m11;
|
---|
| 91 | center->at(1) = -0.5 * m13/ m11;
|
---|
| 92 | center->at(2) = 0.5 * m14/ m11;
|
---|
[357fba] | 93 |
|
---|
[1513a74] | 94 | if (fabs(a.distance(*center) - RADIUS) > MYEPSILON)
|
---|
| 95 | DoeLog(1) && (eLog()<< Verbose(1) << "The given center is further way by " << fabs(a.distance(*center) - RADIUS) << " from a than RADIUS." << endl);
|
---|
[357fba] | 96 |
|
---|
| 97 | gsl_matrix_free(A);
|
---|
| 98 | };
|
---|
| 99 |
|
---|
| 100 |
|
---|
| 101 |
|
---|
| 102 | /**
|
---|
| 103 | * Function returns center of sphere with RADIUS, which rests on points a, b, c
|
---|
| 104 | * @param Center this vector will be used for return
|
---|
| 105 | * @param a vector first point of triangle
|
---|
| 106 | * @param b vector second point of triangle
|
---|
| 107 | * @param c vector third point of triangle
|
---|
[c0f6c6] | 108 | * @param *Umkreismittelpunkt new center point of circumference
|
---|
[357fba] | 109 | * @param Direction vector indicates up/down
|
---|
[c0f6c6] | 110 | * @param AlternativeDirection Vector, needed in case the triangles have 90 deg angle
|
---|
[357fba] | 111 | * @param Halfplaneindicator double indicates whether Direction is up or down
|
---|
[c0f6c6] | 112 | * @param AlternativeIndicator double indicates in case of orthogonal triangles which direction of AlternativeDirection is suitable
|
---|
[357fba] | 113 | * @param alpha double angle at a
|
---|
| 114 | * @param beta double, angle at b
|
---|
| 115 | * @param gamma, double, angle at c
|
---|
| 116 | * @param Radius, double
|
---|
| 117 | * @param Umkreisradius double radius of circumscribing circle
|
---|
| 118 | */
|
---|
[c0f6c6] | 119 | void GetCenterOfSphere(Vector* const & Center, const Vector &a, const Vector &b, const Vector &c, Vector * const NewUmkreismittelpunkt, const Vector* const Direction, const Vector* const AlternativeDirection,
|
---|
| 120 | const double HalfplaneIndicator, const double AlternativeIndicator, const double alpha, const double beta, const double gamma, const double RADIUS, const double Umkreisradius)
|
---|
[357fba] | 121 | {
|
---|
[f67b6e] | 122 | Info FunctionInfo(__func__);
|
---|
[357fba] | 123 | Vector TempNormal, helper;
|
---|
| 124 | double Restradius;
|
---|
| 125 | Vector OtherCenter;
|
---|
| 126 | Center->Zero();
|
---|
[273382] | 127 | helper = sin(2.*alpha) * a;
|
---|
| 128 | (*Center) += helper;
|
---|
| 129 | helper = sin(2.*beta) * b;
|
---|
| 130 | (*Center) += helper;
|
---|
| 131 | helper = sin(2.*gamma) * c;
|
---|
| 132 | (*Center) += helper;
|
---|
[357fba] | 133 | //*Center = a * sin(2.*alpha) + b * sin(2.*beta) + c * sin(2.*gamma) ;
|
---|
| 134 | Center->Scale(1./(sin(2.*alpha) + sin(2.*beta) + sin(2.*gamma)));
|
---|
[273382] | 135 | (*NewUmkreismittelpunkt) = (*Center);
|
---|
[a67d19] | 136 | DoLog(1) && (Log() << Verbose(1) << "Center of new circumference is " << *NewUmkreismittelpunkt << ".\n");
|
---|
[357fba] | 137 | // Here we calculated center of circumscribing circle, using barycentric coordinates
|
---|
[a67d19] | 138 | DoLog(1) && (Log() << Verbose(1) << "Center of circumference is " << *Center << " in direction " << *Direction << ".\n");
|
---|
[357fba] | 139 |
|
---|
[273382] | 140 | TempNormal = a - b;
|
---|
| 141 | helper = a - c;
|
---|
| 142 | TempNormal.VectorProduct(helper);
|
---|
[357fba] | 143 | if (fabs(HalfplaneIndicator) < MYEPSILON)
|
---|
| 144 | {
|
---|
[273382] | 145 | if ((TempNormal.ScalarProduct(*AlternativeDirection) <0 && AlternativeIndicator >0) || (TempNormal.ScalarProduct(*AlternativeDirection) >0 && AlternativeIndicator <0))
|
---|
[357fba] | 146 | {
|
---|
[273382] | 147 | TempNormal *= -1;
|
---|
[357fba] | 148 | }
|
---|
| 149 | }
|
---|
| 150 | else
|
---|
| 151 | {
|
---|
[273382] | 152 | if (((TempNormal.ScalarProduct(*Direction)<0) && (HalfplaneIndicator >0)) || ((TempNormal.ScalarProduct(*Direction)>0) && (HalfplaneIndicator<0)))
|
---|
[357fba] | 153 | {
|
---|
[273382] | 154 | TempNormal *= -1;
|
---|
[357fba] | 155 | }
|
---|
| 156 | }
|
---|
| 157 |
|
---|
| 158 | TempNormal.Normalize();
|
---|
| 159 | Restradius = sqrt(RADIUS*RADIUS - Umkreisradius*Umkreisradius);
|
---|
[a67d19] | 160 | DoLog(1) && (Log() << Verbose(1) << "Height of center of circumference to center of sphere is " << Restradius << ".\n");
|
---|
[357fba] | 161 | TempNormal.Scale(Restradius);
|
---|
[a67d19] | 162 | DoLog(1) && (Log() << Verbose(1) << "Shift vector to sphere of circumference is " << TempNormal << ".\n");
|
---|
[273382] | 163 | (*Center) += TempNormal;
|
---|
[a67d19] | 164 | DoLog(1) && (Log() << Verbose(1) << "Center of sphere of circumference is " << *Center << ".\n");
|
---|
[f1cccd] | 165 | GetSphere(&OtherCenter, a, b, c, RADIUS);
|
---|
[a67d19] | 166 | DoLog(1) && (Log() << Verbose(1) << "OtherCenter of sphere of circumference is " << OtherCenter << ".\n");
|
---|
[357fba] | 167 | };
|
---|
| 168 |
|
---|
| 169 |
|
---|
| 170 | /** Constructs the center of the circumcircle defined by three points \a *a, \a *b and \a *c.
|
---|
| 171 | * \param *Center new center on return
|
---|
| 172 | * \param *a first point
|
---|
| 173 | * \param *b second point
|
---|
| 174 | * \param *c third point
|
---|
| 175 | */
|
---|
[c0f6c6] | 176 | void GetCenterofCircumcircle(Vector * const Center, const Vector &a, const Vector &b, const Vector &c)
|
---|
[357fba] | 177 | {
|
---|
[f67b6e] | 178 | Info FunctionInfo(__func__);
|
---|
[357fba] | 179 | Vector helper;
|
---|
| 180 | double alpha, beta, gamma;
|
---|
[273382] | 181 | Vector SideA = b - c;
|
---|
| 182 | Vector SideB = c - a;
|
---|
| 183 | Vector SideC = a - b;
|
---|
| 184 | alpha = M_PI - SideB.Angle(SideC);
|
---|
| 185 | beta = M_PI - SideC.Angle(SideA);
|
---|
| 186 | gamma = M_PI - SideA.Angle(SideB);
|
---|
[f67b6e] | 187 | //Log() << Verbose(1) << "INFO: alpha = " << alpha/M_PI*180. << ", beta = " << beta/M_PI*180. << ", gamma = " << gamma/M_PI*180. << "." << endl;
|
---|
[e359a8] | 188 | if (fabs(M_PI - alpha - beta - gamma) > HULLEPSILON) {
|
---|
[299554] | 189 | DoeLog(2) && (eLog()<< Verbose(2) << "GetCenterofCircumcircle: Sum of angles " << (alpha+beta+gamma)/M_PI*180. << " > 180 degrees by " << fabs(M_PI - alpha - beta - gamma)/M_PI*180. << "!" << endl);
|
---|
[e359a8] | 190 | }
|
---|
[357fba] | 191 |
|
---|
| 192 | Center->Zero();
|
---|
[273382] | 193 | helper = sin(2.*alpha) * a;
|
---|
| 194 | (*Center) += helper;
|
---|
| 195 | helper = sin(2.*beta) * b;
|
---|
| 196 | (*Center) += helper;
|
---|
| 197 | helper = sin(2.*gamma) * c;
|
---|
| 198 | (*Center) += helper;
|
---|
[357fba] | 199 | Center->Scale(1./(sin(2.*alpha) + sin(2.*beta) + sin(2.*gamma)));
|
---|
| 200 | };
|
---|
| 201 |
|
---|
| 202 | /** Returns the parameter "path length" for a given \a NewSphereCenter relative to \a OldSphereCenter on a circle on the plane \a CirclePlaneNormal with center \a CircleCenter and radius \a CircleRadius.
|
---|
| 203 | * Test whether the \a NewSphereCenter is really on the given plane and in distance \a CircleRadius from \a CircleCenter.
|
---|
| 204 | * It calculates the angle, making it unique on [0,2.*M_PI) by comparing to SearchDirection.
|
---|
| 205 | * Also the new center is invalid if it the same as the old one and does not lie right above (\a NormalVector) the base line (\a CircleCenter).
|
---|
| 206 | * \param CircleCenter Center of the parameter circle
|
---|
| 207 | * \param CirclePlaneNormal normal vector to plane of the parameter circle
|
---|
| 208 | * \param CircleRadius radius of the parameter circle
|
---|
| 209 | * \param NewSphereCenter new center of a circumcircle
|
---|
| 210 | * \param OldSphereCenter old center of a circumcircle, defining the zero "path length" on the parameter circle
|
---|
| 211 | * \param NormalVector normal vector
|
---|
| 212 | * \param SearchDirection search direction to make angle unique on return.
|
---|
| 213 | * \return Angle between \a NewSphereCenter and \a OldSphereCenter relative to \a CircleCenter, 2.*M_PI if one test fails
|
---|
| 214 | */
|
---|
[c0f6c6] | 215 | double GetPathLengthonCircumCircle(const Vector &CircleCenter, const Vector &CirclePlaneNormal, const double CircleRadius, const Vector &NewSphereCenter, const Vector &OldSphereCenter, const Vector &NormalVector, const Vector &SearchDirection)
|
---|
[357fba] | 216 | {
|
---|
[f67b6e] | 217 | Info FunctionInfo(__func__);
|
---|
[357fba] | 218 | Vector helper;
|
---|
| 219 | double radius, alpha;
|
---|
[273382] | 220 |
|
---|
| 221 | Vector RelativeOldSphereCenter = OldSphereCenter - CircleCenter;
|
---|
| 222 | Vector RelativeNewSphereCenter = NewSphereCenter - CircleCenter;
|
---|
| 223 | helper = RelativeNewSphereCenter;
|
---|
[357fba] | 224 | // test whether new center is on the parameter circle's plane
|
---|
[273382] | 225 | if (fabs(helper.ScalarProduct(CirclePlaneNormal)) > HULLEPSILON) {
|
---|
[8cbb97] | 226 | DoeLog(1) && (eLog()<< Verbose(1) << "Something's very wrong here: NewSphereCenter is not on the band's plane as desired by " <<fabs(helper.ScalarProduct(CirclePlaneNormal)) << "!" << endl);
|
---|
[273382] | 227 | helper.ProjectOntoPlane(CirclePlaneNormal);
|
---|
[357fba] | 228 | }
|
---|
[b998c3] | 229 | radius = helper.NormSquared();
|
---|
[357fba] | 230 | // test whether the new center vector has length of CircleRadius
|
---|
| 231 | if (fabs(radius - CircleRadius) > HULLEPSILON)
|
---|
[58ed4a] | 232 | DoeLog(1) && (eLog()<< Verbose(1) << "The projected center of the new sphere has radius " << radius << " instead of " << CircleRadius << "." << endl);
|
---|
[273382] | 233 | alpha = helper.Angle(RelativeOldSphereCenter);
|
---|
[357fba] | 234 | // make the angle unique by checking the halfplanes/search direction
|
---|
[273382] | 235 | if (helper.ScalarProduct(SearchDirection) < -HULLEPSILON) // acos is not unique on [0, 2.*M_PI), hence extra check to decide between two half intervals
|
---|
[357fba] | 236 | alpha = 2.*M_PI - alpha;
|
---|
[a67d19] | 237 | DoLog(1) && (Log() << Verbose(1) << "INFO: RelativeNewSphereCenter is " << helper << ", RelativeOldSphereCenter is " << RelativeOldSphereCenter << " and resulting angle is " << alpha << "." << endl);
|
---|
[1513a74] | 238 | radius = helper.distance(RelativeOldSphereCenter);
|
---|
[273382] | 239 | helper.ProjectOntoPlane(NormalVector);
|
---|
[357fba] | 240 | // check whether new center is somewhat away or at least right over the current baseline to prevent intersecting triangles
|
---|
| 241 | if ((radius > HULLEPSILON) || (helper.Norm() < HULLEPSILON)) {
|
---|
[a67d19] | 242 | DoLog(1) && (Log() << Verbose(1) << "INFO: Distance between old and new center is " << radius << " and between new center and baseline center is " << helper.Norm() << "." << endl);
|
---|
[357fba] | 243 | return alpha;
|
---|
| 244 | } else {
|
---|
[a67d19] | 245 | DoLog(1) && (Log() << Verbose(1) << "INFO: NewSphereCenter " << RelativeNewSphereCenter << " is too close to RelativeOldSphereCenter" << RelativeOldSphereCenter << "." << endl);
|
---|
[357fba] | 246 | return 2.*M_PI;
|
---|
| 247 | }
|
---|
| 248 | };
|
---|
| 249 |
|
---|
| 250 | struct Intersection {
|
---|
| 251 | Vector x1;
|
---|
| 252 | Vector x2;
|
---|
| 253 | Vector x3;
|
---|
| 254 | Vector x4;
|
---|
| 255 | };
|
---|
| 256 |
|
---|
| 257 | /**
|
---|
| 258 | * Intersection calculation function.
|
---|
| 259 | *
|
---|
| 260 | * @param x to find the result for
|
---|
| 261 | * @param function parameter
|
---|
| 262 | */
|
---|
| 263 | double MinIntersectDistance(const gsl_vector * x, void *params)
|
---|
| 264 | {
|
---|
[f67b6e] | 265 | Info FunctionInfo(__func__);
|
---|
[357fba] | 266 | double retval = 0;
|
---|
| 267 | struct Intersection *I = (struct Intersection *)params;
|
---|
| 268 | Vector intersection;
|
---|
| 269 | for (int i=0;i<NDIM;i++)
|
---|
[0a4f7f] | 270 | intersection[i] = gsl_vector_get(x, i);
|
---|
[357fba] | 271 |
|
---|
[273382] | 272 | Vector SideA = I->x1 -I->x2 ;
|
---|
| 273 | Vector HeightA = intersection - I->x1;
|
---|
| 274 | HeightA.ProjectOntoPlane(SideA);
|
---|
[357fba] | 275 |
|
---|
[273382] | 276 | Vector SideB = I->x3 - I->x4;
|
---|
| 277 | Vector HeightB = intersection - I->x3;
|
---|
| 278 | HeightB.ProjectOntoPlane(SideB);
|
---|
[357fba] | 279 |
|
---|
[273382] | 280 | retval = HeightA.ScalarProduct(HeightA) + HeightB.ScalarProduct(HeightB);
|
---|
[f67b6e] | 281 | //Log() << Verbose(1) << "MinIntersectDistance called, result: " << retval << endl;
|
---|
[357fba] | 282 |
|
---|
| 283 | return retval;
|
---|
| 284 | };
|
---|
| 285 |
|
---|
| 286 |
|
---|
| 287 | /**
|
---|
| 288 | * Calculates whether there is an intersection between two lines. The first line
|
---|
| 289 | * always goes through point 1 and point 2 and the second line is given by the
|
---|
| 290 | * connection between point 4 and point 5.
|
---|
| 291 | *
|
---|
| 292 | * @param point 1 of line 1
|
---|
| 293 | * @param point 2 of line 1
|
---|
| 294 | * @param point 1 of line 2
|
---|
| 295 | * @param point 2 of line 2
|
---|
| 296 | *
|
---|
| 297 | * @return true if there is an intersection between the given lines, false otherwise
|
---|
| 298 | */
|
---|
[c0f6c6] | 299 | bool existsIntersection(const Vector &point1, const Vector &point2, const Vector &point3, const Vector &point4)
|
---|
[357fba] | 300 | {
|
---|
[f67b6e] | 301 | Info FunctionInfo(__func__);
|
---|
[357fba] | 302 | bool result;
|
---|
| 303 |
|
---|
| 304 | struct Intersection par;
|
---|
[273382] | 305 | par.x1 = point1;
|
---|
| 306 | par.x2 = point2;
|
---|
| 307 | par.x3 = point3;
|
---|
| 308 | par.x4 = point4;
|
---|
[357fba] | 309 |
|
---|
| 310 | const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
|
---|
| 311 | gsl_multimin_fminimizer *s = NULL;
|
---|
| 312 | gsl_vector *ss, *x;
|
---|
[f1cccd] | 313 | gsl_multimin_function minexFunction;
|
---|
[357fba] | 314 |
|
---|
| 315 | size_t iter = 0;
|
---|
| 316 | int status;
|
---|
| 317 | double size;
|
---|
| 318 |
|
---|
| 319 | /* Starting point */
|
---|
| 320 | x = gsl_vector_alloc(NDIM);
|
---|
[0a4f7f] | 321 | gsl_vector_set(x, 0, point1[0]);
|
---|
| 322 | gsl_vector_set(x, 1, point1[1]);
|
---|
| 323 | gsl_vector_set(x, 2, point1[2]);
|
---|
[357fba] | 324 |
|
---|
| 325 | /* Set initial step sizes to 1 */
|
---|
| 326 | ss = gsl_vector_alloc(NDIM);
|
---|
| 327 | gsl_vector_set_all(ss, 1.0);
|
---|
| 328 |
|
---|
| 329 | /* Initialize method and iterate */
|
---|
[f1cccd] | 330 | minexFunction.n = NDIM;
|
---|
| 331 | minexFunction.f = &MinIntersectDistance;
|
---|
| 332 | minexFunction.params = (void *)∥
|
---|
[357fba] | 333 |
|
---|
| 334 | s = gsl_multimin_fminimizer_alloc(T, NDIM);
|
---|
[f1cccd] | 335 | gsl_multimin_fminimizer_set(s, &minexFunction, x, ss);
|
---|
[357fba] | 336 |
|
---|
| 337 | do {
|
---|
| 338 | iter++;
|
---|
| 339 | status = gsl_multimin_fminimizer_iterate(s);
|
---|
| 340 |
|
---|
| 341 | if (status) {
|
---|
| 342 | break;
|
---|
| 343 | }
|
---|
| 344 |
|
---|
| 345 | size = gsl_multimin_fminimizer_size(s);
|
---|
| 346 | status = gsl_multimin_test_size(size, 1e-2);
|
---|
| 347 |
|
---|
| 348 | if (status == GSL_SUCCESS) {
|
---|
[a67d19] | 349 | DoLog(1) && (Log() << Verbose(1) << "converged to minimum" << endl);
|
---|
[357fba] | 350 | }
|
---|
| 351 | } while (status == GSL_CONTINUE && iter < 100);
|
---|
| 352 |
|
---|
| 353 | // check whether intersection is in between or not
|
---|
[273382] | 354 | Vector intersection;
|
---|
[357fba] | 355 | double t1, t2;
|
---|
| 356 | for (int i = 0; i < NDIM; i++) {
|
---|
[0a4f7f] | 357 | intersection[i] = gsl_vector_get(s->x, i);
|
---|
[357fba] | 358 | }
|
---|
| 359 |
|
---|
[273382] | 360 | Vector SideA = par.x2 - par.x1;
|
---|
| 361 | Vector HeightA = intersection - par.x1;
|
---|
[357fba] | 362 |
|
---|
[273382] | 363 | t1 = HeightA.ScalarProduct(SideA)/SideA.ScalarProduct(SideA);
|
---|
[357fba] | 364 |
|
---|
[273382] | 365 | Vector SideB = par.x4 - par.x3;
|
---|
| 366 | Vector HeightB = intersection - par.x3;
|
---|
[357fba] | 367 |
|
---|
[273382] | 368 | t2 = HeightB.ScalarProduct(SideB)/SideB.ScalarProduct(SideB);
|
---|
[357fba] | 369 |
|
---|
[f67b6e] | 370 | Log() << Verbose(1) << "Intersection " << intersection << " is at "
|
---|
[357fba] | 371 | << t1 << " for (" << point1 << "," << point2 << ") and at "
|
---|
| 372 | << t2 << " for (" << point3 << "," << point4 << "): ";
|
---|
| 373 |
|
---|
| 374 | if (((t1 >= 0) && (t1 <= 1)) && ((t2 >= 0) && (t2 <= 1))) {
|
---|
[a67d19] | 375 | DoLog(1) && (Log() << Verbose(1) << "true intersection." << endl);
|
---|
[357fba] | 376 | result = true;
|
---|
| 377 | } else {
|
---|
[a67d19] | 378 | DoLog(1) && (Log() << Verbose(1) << "intersection out of region of interest." << endl);
|
---|
[357fba] | 379 | result = false;
|
---|
| 380 | }
|
---|
| 381 |
|
---|
| 382 | // free minimizer stuff
|
---|
| 383 | gsl_vector_free(x);
|
---|
| 384 | gsl_vector_free(ss);
|
---|
| 385 | gsl_multimin_fminimizer_free(s);
|
---|
| 386 |
|
---|
| 387 | return result;
|
---|
[91e7e4a] | 388 | };
|
---|
| 389 |
|
---|
[57066a] | 390 | /** Gets the angle between a point and a reference relative to the provided center.
|
---|
| 391 | * We have two shanks point and reference between which the angle is calculated
|
---|
| 392 | * and by scalar product with OrthogonalVector we decide the interval.
|
---|
| 393 | * @param point to calculate the angle for
|
---|
| 394 | * @param reference to which to calculate the angle
|
---|
| 395 | * @param OrthogonalVector points in direction of [pi,2pi] interval
|
---|
| 396 | *
|
---|
| 397 | * @return angle between point and reference
|
---|
| 398 | */
|
---|
[c0f6c6] | 399 | double GetAngle(const Vector &point, const Vector &reference, const Vector &OrthogonalVector)
|
---|
[57066a] | 400 | {
|
---|
[f67b6e] | 401 | Info FunctionInfo(__func__);
|
---|
[57066a] | 402 | if (reference.IsZero())
|
---|
| 403 | return M_PI;
|
---|
| 404 |
|
---|
| 405 | // calculate both angles and correct with in-plane vector
|
---|
| 406 | if (point.IsZero())
|
---|
| 407 | return M_PI;
|
---|
[273382] | 408 | double phi = point.Angle(reference);
|
---|
| 409 | if (OrthogonalVector.ScalarProduct(point) > 0) {
|
---|
[57066a] | 410 | phi = 2.*M_PI - phi;
|
---|
| 411 | }
|
---|
| 412 |
|
---|
[a67d19] | 413 | DoLog(1) && (Log() << Verbose(1) << "INFO: " << point << " has angle " << phi << " with respect to reference " << reference << "." << endl);
|
---|
[57066a] | 414 |
|
---|
| 415 | return phi;
|
---|
| 416 | }
|
---|
| 417 |
|
---|
[91e7e4a] | 418 |
|
---|
| 419 | /** Calculates the volume of a general tetraeder.
|
---|
| 420 | * \param *a first vector
|
---|
| 421 | * \param *a first vector
|
---|
| 422 | * \param *a first vector
|
---|
| 423 | * \param *a first vector
|
---|
| 424 | * \return \f$ \frac{1}{6} \cdot ((a-d) \times (a-c) \cdot (a-b)) \f$
|
---|
| 425 | */
|
---|
[c0f6c6] | 426 | double CalculateVolumeofGeneralTetraeder(const Vector &a, const Vector &b, const Vector &c, const Vector &d)
|
---|
[91e7e4a] | 427 | {
|
---|
[f67b6e] | 428 | Info FunctionInfo(__func__);
|
---|
[91e7e4a] | 429 | Vector Point, TetraederVector[3];
|
---|
| 430 | double volume;
|
---|
| 431 |
|
---|
[1bd79e] | 432 | TetraederVector[0] = a;
|
---|
| 433 | TetraederVector[1] = b;
|
---|
| 434 | TetraederVector[2] = c;
|
---|
[91e7e4a] | 435 | for (int j=0;j<3;j++)
|
---|
[273382] | 436 | TetraederVector[j].SubtractVector(d);
|
---|
[1bd79e] | 437 | Point = TetraederVector[0];
|
---|
[273382] | 438 | Point.VectorProduct(TetraederVector[1]);
|
---|
| 439 | volume = 1./6. * fabs(Point.ScalarProduct(TetraederVector[2]));
|
---|
[91e7e4a] | 440 | return volume;
|
---|
| 441 | };
|
---|
[357fba] | 442 |
|
---|
[57066a] | 443 |
|
---|
| 444 | /** Checks for a new special triangle whether one of its edges is already present with one one triangle connected.
|
---|
| 445 | * This enforces that special triangles (i.e. degenerated ones) should at last close the open-edge frontier and not
|
---|
| 446 | * make it bigger (i.e. closing one (the baseline) and opening two new ones).
|
---|
| 447 | * \param TPS[3] nodes of the triangle
|
---|
| 448 | * \return true - there is such a line (i.e. creation of degenerated triangle is valid), false - no such line (don't create)
|
---|
| 449 | */
|
---|
[c0f6c6] | 450 | bool CheckLineCriteriaForDegeneratedTriangle(const BoundaryPointSet * const nodes[3])
|
---|
[57066a] | 451 | {
|
---|
[f67b6e] | 452 | Info FunctionInfo(__func__);
|
---|
[57066a] | 453 | bool result = false;
|
---|
| 454 | int counter = 0;
|
---|
| 455 |
|
---|
| 456 | // check all three points
|
---|
| 457 | for (int i=0;i<3;i++)
|
---|
| 458 | for (int j=i+1; j<3; j++) {
|
---|
[f1ef60a] | 459 | if (nodes[i] == NULL) {
|
---|
[a67d19] | 460 | DoLog(1) && (Log() << Verbose(1) << "Node nr. " << i << " is not yet present." << endl);
|
---|
[f1ef60a] | 461 | result = true;
|
---|
| 462 | } else if (nodes[i]->lines.find(nodes[j]->node->nr) != nodes[i]->lines.end()) { // there already is a line
|
---|
[776b64] | 463 | LineMap::const_iterator FindLine;
|
---|
| 464 | pair<LineMap::const_iterator,LineMap::const_iterator> FindPair;
|
---|
[57066a] | 465 | FindPair = nodes[i]->lines.equal_range(nodes[j]->node->nr);
|
---|
| 466 | for (FindLine = FindPair.first; FindLine != FindPair.second; ++FindLine) {
|
---|
| 467 | // If there is a line with less than two attached triangles, we don't need a new line.
|
---|
| 468 | if (FindLine->second->triangles.size() < 2) {
|
---|
| 469 | counter++;
|
---|
| 470 | break; // increase counter only once per edge
|
---|
| 471 | }
|
---|
| 472 | }
|
---|
| 473 | } else { // no line
|
---|
[a67d19] | 474 | DoLog(1) && (Log() << Verbose(1) << "The line between " << *nodes[i] << " and " << *nodes[j] << " is not yet present, hence no need for a degenerate triangle." << endl);
|
---|
[57066a] | 475 | result = true;
|
---|
| 476 | }
|
---|
| 477 | }
|
---|
| 478 | if ((!result) && (counter > 1)) {
|
---|
[a67d19] | 479 | DoLog(1) && (Log() << Verbose(1) << "INFO: Degenerate triangle is ok, at least two, here " << counter << ", existing lines are used." << endl);
|
---|
[57066a] | 480 | result = true;
|
---|
| 481 | }
|
---|
| 482 | return result;
|
---|
| 483 | };
|
---|
| 484 |
|
---|
| 485 |
|
---|
[f67b6e] | 486 | ///** Sort function for the candidate list.
|
---|
| 487 | // */
|
---|
| 488 | //bool SortCandidates(const CandidateForTesselation* candidate1, const CandidateForTesselation* candidate2)
|
---|
| 489 | //{
|
---|
| 490 | // Info FunctionInfo(__func__);
|
---|
| 491 | // Vector BaseLineVector, OrthogonalVector, helper;
|
---|
| 492 | // if (candidate1->BaseLine != candidate2->BaseLine) { // sanity check
|
---|
[58ed4a] | 493 | // DoeLog(1) && (eLog()<< Verbose(1) << "sortCandidates was called for two different baselines: " << candidate1->BaseLine << " and " << candidate2->BaseLine << "." << endl);
|
---|
[f67b6e] | 494 | // //return false;
|
---|
| 495 | // exit(1);
|
---|
| 496 | // }
|
---|
| 497 | // // create baseline vector
|
---|
| 498 | // BaseLineVector.CopyVector(candidate1->BaseLine->endpoints[1]->node->node);
|
---|
| 499 | // BaseLineVector.SubtractVector(candidate1->BaseLine->endpoints[0]->node->node);
|
---|
| 500 | // BaseLineVector.Normalize();
|
---|
| 501 | //
|
---|
| 502 | // // create normal in-plane vector to cope with acos() non-uniqueness on [0,2pi] (note that is pointing in the "right" direction already, hence ">0" test!)
|
---|
| 503 | // helper.CopyVector(candidate1->BaseLine->endpoints[0]->node->node);
|
---|
| 504 | // helper.SubtractVector(candidate1->point->node);
|
---|
| 505 | // OrthogonalVector.CopyVector(&helper);
|
---|
| 506 | // helper.VectorProduct(&BaseLineVector);
|
---|
| 507 | // OrthogonalVector.SubtractVector(&helper);
|
---|
| 508 | // OrthogonalVector.Normalize();
|
---|
| 509 | //
|
---|
| 510 | // // calculate both angles and correct with in-plane vector
|
---|
| 511 | // helper.CopyVector(candidate1->point->node);
|
---|
| 512 | // helper.SubtractVector(candidate1->BaseLine->endpoints[0]->node->node);
|
---|
| 513 | // double phi = BaseLineVector.Angle(&helper);
|
---|
| 514 | // if (OrthogonalVector.ScalarProduct(&helper) > 0) {
|
---|
| 515 | // phi = 2.*M_PI - phi;
|
---|
| 516 | // }
|
---|
| 517 | // helper.CopyVector(candidate2->point->node);
|
---|
| 518 | // helper.SubtractVector(candidate1->BaseLine->endpoints[0]->node->node);
|
---|
| 519 | // double psi = BaseLineVector.Angle(&helper);
|
---|
| 520 | // if (OrthogonalVector.ScalarProduct(&helper) > 0) {
|
---|
| 521 | // psi = 2.*M_PI - psi;
|
---|
| 522 | // }
|
---|
| 523 | //
|
---|
| 524 | // Log() << Verbose(1) << *candidate1->point << " has angle " << phi << endl;
|
---|
| 525 | // Log() << Verbose(1) << *candidate2->point << " has angle " << psi << endl;
|
---|
| 526 | //
|
---|
| 527 | // // return comparison
|
---|
| 528 | // return phi < psi;
|
---|
| 529 | //};
|
---|
[57066a] | 530 |
|
---|
| 531 | /**
|
---|
| 532 | * Finds the point which is second closest to the provided one.
|
---|
| 533 | *
|
---|
| 534 | * @param Point to which to find the second closest other point
|
---|
| 535 | * @param linked cell structure
|
---|
| 536 | *
|
---|
| 537 | * @return point which is second closest to the provided one
|
---|
| 538 | */
|
---|
[71b20e] | 539 | TesselPoint* FindSecondClosestTesselPoint(const Vector* Point, const LinkedCell* const LC)
|
---|
[57066a] | 540 | {
|
---|
[f67b6e] | 541 | Info FunctionInfo(__func__);
|
---|
[57066a] | 542 | TesselPoint* closestPoint = NULL;
|
---|
| 543 | TesselPoint* secondClosestPoint = NULL;
|
---|
| 544 | double distance = 1e16;
|
---|
| 545 | double secondDistance = 1e16;
|
---|
| 546 | Vector helper;
|
---|
| 547 | int N[NDIM], Nlower[NDIM], Nupper[NDIM];
|
---|
| 548 |
|
---|
| 549 | LC->SetIndexToVector(Point); // ignore status as we calculate bounds below sensibly
|
---|
| 550 | for(int i=0;i<NDIM;i++) // store indices of this cell
|
---|
| 551 | N[i] = LC->n[i];
|
---|
[a67d19] | 552 | DoLog(1) && (Log() << Verbose(1) << "INFO: Center cell is " << N[0] << ", " << N[1] << ", " << N[2] << " with No. " << LC->index << "." << endl);
|
---|
[57066a] | 553 |
|
---|
| 554 | LC->GetNeighbourBounds(Nlower, Nupper);
|
---|
[f67b6e] | 555 | //Log() << Verbose(1) << endl;
|
---|
[57066a] | 556 | for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
|
---|
| 557 | for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
|
---|
| 558 | for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
|
---|
[734816] | 559 | const LinkedCell::LinkedNodes *List = LC->GetCurrentCell();
|
---|
[f67b6e] | 560 | //Log() << Verbose(1) << "The current cell " << LC->n[0] << "," << LC->n[1] << "," << LC->n[2] << endl;
|
---|
[57066a] | 561 | if (List != NULL) {
|
---|
[734816] | 562 | for (LinkedCell::LinkedNodes::const_iterator Runner = List->begin(); Runner != List->end(); Runner++) {
|
---|
[273382] | 563 | helper = (*Point) - (*(*Runner)->node);
|
---|
[57066a] | 564 | double currentNorm = helper. Norm();
|
---|
| 565 | if (currentNorm < distance) {
|
---|
| 566 | // remember second point
|
---|
| 567 | secondDistance = distance;
|
---|
| 568 | secondClosestPoint = closestPoint;
|
---|
| 569 | // mark down new closest point
|
---|
| 570 | distance = currentNorm;
|
---|
| 571 | closestPoint = (*Runner);
|
---|
[e138de] | 572 | //Log() << Verbose(2) << "INFO: New Second Nearest Neighbour is " << *secondClosestPoint << "." << endl;
|
---|
[57066a] | 573 | }
|
---|
| 574 | }
|
---|
| 575 | } else {
|
---|
[717e0c] | 576 | eLog() << Verbose(1) << "The current cell " << LC->n[0] << "," << LC->n[1] << ","
|
---|
[57066a] | 577 | << LC->n[2] << " is invalid!" << endl;
|
---|
| 578 | }
|
---|
| 579 | }
|
---|
| 580 |
|
---|
| 581 | return secondClosestPoint;
|
---|
| 582 | };
|
---|
| 583 |
|
---|
| 584 | /**
|
---|
| 585 | * Finds the point which is closest to the provided one.
|
---|
| 586 | *
|
---|
| 587 | * @param Point to which to find the closest other point
|
---|
| 588 | * @param SecondPoint the second closest other point on return, NULL if none found
|
---|
| 589 | * @param linked cell structure
|
---|
| 590 | *
|
---|
| 591 | * @return point which is closest to the provided one, NULL if none found
|
---|
| 592 | */
|
---|
[71b20e] | 593 | TesselPoint* FindClosestTesselPoint(const Vector* Point, TesselPoint *&SecondPoint, const LinkedCell* const LC)
|
---|
[57066a] | 594 | {
|
---|
[f67b6e] | 595 | Info FunctionInfo(__func__);
|
---|
[57066a] | 596 | TesselPoint* closestPoint = NULL;
|
---|
| 597 | SecondPoint = NULL;
|
---|
| 598 | double distance = 1e16;
|
---|
| 599 | double secondDistance = 1e16;
|
---|
| 600 | Vector helper;
|
---|
| 601 | int N[NDIM], Nlower[NDIM], Nupper[NDIM];
|
---|
| 602 |
|
---|
| 603 | LC->SetIndexToVector(Point); // ignore status as we calculate bounds below sensibly
|
---|
| 604 | for(int i=0;i<NDIM;i++) // store indices of this cell
|
---|
| 605 | N[i] = LC->n[i];
|
---|
[a67d19] | 606 | DoLog(1) && (Log() << Verbose(1) << "INFO: Center cell is " << N[0] << ", " << N[1] << ", " << N[2] << " with No. " << LC->index << "." << endl);
|
---|
[57066a] | 607 |
|
---|
| 608 | LC->GetNeighbourBounds(Nlower, Nupper);
|
---|
[f67b6e] | 609 | //Log() << Verbose(1) << endl;
|
---|
[57066a] | 610 | for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
|
---|
| 611 | for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
|
---|
| 612 | for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
|
---|
[734816] | 613 | const LinkedCell::LinkedNodes *List = LC->GetCurrentCell();
|
---|
[f67b6e] | 614 | //Log() << Verbose(1) << "The current cell " << LC->n[0] << "," << LC->n[1] << "," << LC->n[2] << endl;
|
---|
[57066a] | 615 | if (List != NULL) {
|
---|
[734816] | 616 | for (LinkedCell::LinkedNodes::const_iterator Runner = List->begin(); Runner != List->end(); Runner++) {
|
---|
[273382] | 617 | helper = (*Point) - (*(*Runner)->node);
|
---|
[71b20e] | 618 | double currentNorm = helper.NormSquared();
|
---|
[57066a] | 619 | if (currentNorm < distance) {
|
---|
| 620 | secondDistance = distance;
|
---|
| 621 | SecondPoint = closestPoint;
|
---|
| 622 | distance = currentNorm;
|
---|
| 623 | closestPoint = (*Runner);
|
---|
[f67b6e] | 624 | //Log() << Verbose(1) << "INFO: New Nearest Neighbour is " << *closestPoint << "." << endl;
|
---|
[57066a] | 625 | } else if (currentNorm < secondDistance) {
|
---|
| 626 | secondDistance = currentNorm;
|
---|
| 627 | SecondPoint = (*Runner);
|
---|
[f67b6e] | 628 | //Log() << Verbose(1) << "INFO: New Second Nearest Neighbour is " << *SecondPoint << "." << endl;
|
---|
[57066a] | 629 | }
|
---|
| 630 | }
|
---|
| 631 | } else {
|
---|
[717e0c] | 632 | eLog() << Verbose(1) << "The current cell " << LC->n[0] << "," << LC->n[1] << ","
|
---|
[57066a] | 633 | << LC->n[2] << " is invalid!" << endl;
|
---|
| 634 | }
|
---|
| 635 | }
|
---|
[a2028e] | 636 | // output
|
---|
| 637 | if (closestPoint != NULL) {
|
---|
[a67d19] | 638 | DoLog(1) && (Log() << Verbose(1) << "Closest point is " << *closestPoint);
|
---|
[a2028e] | 639 | if (SecondPoint != NULL)
|
---|
[a67d19] | 640 | DoLog(0) && (Log() << Verbose(0) << " and second closest is " << *SecondPoint);
|
---|
| 641 | DoLog(0) && (Log() << Verbose(0) << "." << endl);
|
---|
[a2028e] | 642 | }
|
---|
[57066a] | 643 | return closestPoint;
|
---|
| 644 | };
|
---|
| 645 |
|
---|
| 646 | /** Returns the closest point on \a *Base with respect to \a *OtherBase.
|
---|
| 647 | * \param *out output stream for debugging
|
---|
| 648 | * \param *Base reference line
|
---|
| 649 | * \param *OtherBase other base line
|
---|
| 650 | * \return Vector on reference line that has closest distance
|
---|
| 651 | */
|
---|
[e138de] | 652 | Vector * GetClosestPointBetweenLine(const BoundaryLineSet * const Base, const BoundaryLineSet * const OtherBase)
|
---|
[57066a] | 653 | {
|
---|
[f67b6e] | 654 | Info FunctionInfo(__func__);
|
---|
[57066a] | 655 | // construct the plane of the two baselines (i.e. take both their directional vectors)
|
---|
[273382] | 656 | Vector Baseline = (*Base->endpoints[1]->node->node) - (*Base->endpoints[0]->node->node);
|
---|
| 657 | Vector OtherBaseline = (*OtherBase->endpoints[1]->node->node) - (*OtherBase->endpoints[0]->node->node);
|
---|
| 658 | Vector Normal = Baseline;
|
---|
| 659 | Normal.VectorProduct(OtherBaseline);
|
---|
[57066a] | 660 | Normal.Normalize();
|
---|
[a67d19] | 661 | DoLog(1) && (Log() << Verbose(1) << "First direction is " << Baseline << ", second direction is " << OtherBaseline << ", normal of intersection plane is " << Normal << "." << endl);
|
---|
[57066a] | 662 |
|
---|
| 663 | // project one offset point of OtherBase onto this plane (and add plane offset vector)
|
---|
[273382] | 664 | Vector NewOffset = (*OtherBase->endpoints[0]->node->node) - (*Base->endpoints[0]->node->node);
|
---|
| 665 | NewOffset.ProjectOntoPlane(Normal);
|
---|
| 666 | NewOffset += (*Base->endpoints[0]->node->node);
|
---|
| 667 | Vector NewDirection = NewOffset + OtherBaseline;
|
---|
[57066a] | 668 |
|
---|
| 669 | // calculate the intersection between this projected baseline and Base
|
---|
| 670 | Vector *Intersection = new Vector;
|
---|
[643e76] | 671 | Line line1 = makeLineThrough(*(Base->endpoints[0]->node->node),*(Base->endpoints[1]->node->node));
|
---|
| 672 | Line line2 = makeLineThrough(NewOffset, NewDirection);
|
---|
| 673 | *Intersection = line1.getIntersection(line2);
|
---|
[273382] | 674 | Normal = (*Intersection) - (*Base->endpoints[0]->node->node);
|
---|
[8cbb97] | 675 | DoLog(1) && (Log() << Verbose(1) << "Found closest point on " << *Base << " at " << *Intersection << ", factor in line is " << fabs(Normal.ScalarProduct(Baseline)/Baseline.NormSquared()) << "." << endl);
|
---|
[57066a] | 676 |
|
---|
| 677 | return Intersection;
|
---|
| 678 | };
|
---|
| 679 |
|
---|
[c4d4df] | 680 | /** Returns the distance to the plane defined by \a *triangle
|
---|
| 681 | * \param *out output stream for debugging
|
---|
| 682 | * \param *x Vector to calculate distance to
|
---|
| 683 | * \param *triangle triangle defining plane
|
---|
| 684 | * \return distance between \a *x and plane defined by \a *triangle, -1 - if something went wrong
|
---|
| 685 | */
|
---|
[e138de] | 686 | double DistanceToTrianglePlane(const Vector *x, const BoundaryTriangleSet * const triangle)
|
---|
[c4d4df] | 687 | {
|
---|
[f67b6e] | 688 | Info FunctionInfo(__func__);
|
---|
[c4d4df] | 689 | double distance = 0.;
|
---|
| 690 | if (x == NULL) {
|
---|
| 691 | return -1;
|
---|
| 692 | }
|
---|
[d4c9ae] | 693 | distance = x->DistanceToSpace(triangle->getPlane());
|
---|
[c4d4df] | 694 | return distance;
|
---|
| 695 | };
|
---|
[57066a] | 696 |
|
---|
| 697 | /** Creates the objects in a VRML file.
|
---|
| 698 | * \param *out output stream for debugging
|
---|
| 699 | * \param *vrmlfile output stream for tecplot data
|
---|
| 700 | * \param *Tess Tesselation structure with constructed triangles
|
---|
| 701 | * \param *mol molecule structure with atom positions
|
---|
| 702 | */
|
---|
[e138de] | 703 | void WriteVrmlFile(ofstream * const vrmlfile, const Tesselation * const Tess, const PointCloud * const cloud)
|
---|
[57066a] | 704 | {
|
---|
[f67b6e] | 705 | Info FunctionInfo(__func__);
|
---|
[57066a] | 706 | TesselPoint *Walker = NULL;
|
---|
| 707 | int i;
|
---|
[e138de] | 708 | Vector *center = cloud->GetCenter();
|
---|
[57066a] | 709 | if (vrmlfile != NULL) {
|
---|
[e138de] | 710 | //Log() << Verbose(1) << "Writing Raster3D file ... ";
|
---|
[57066a] | 711 | *vrmlfile << "#VRML V2.0 utf8" << endl;
|
---|
| 712 | *vrmlfile << "#Created by molecuilder" << endl;
|
---|
| 713 | *vrmlfile << "#All atoms as spheres" << endl;
|
---|
| 714 | cloud->GoToFirst();
|
---|
| 715 | while (!cloud->IsEnd()) {
|
---|
| 716 | Walker = cloud->GetPoint();
|
---|
| 717 | *vrmlfile << "Sphere {" << endl << " "; // 2 is sphere type
|
---|
| 718 | for (i=0;i<NDIM;i++)
|
---|
[0a4f7f] | 719 | *vrmlfile << Walker->node->at(i)-center->at(i) << " ";
|
---|
[57066a] | 720 | *vrmlfile << "\t0.1\t1. 1. 1." << endl; // radius 0.05 and white as colour
|
---|
| 721 | cloud->GoToNext();
|
---|
| 722 | }
|
---|
| 723 |
|
---|
| 724 | *vrmlfile << "# All tesselation triangles" << endl;
|
---|
[776b64] | 725 | for (TriangleMap::const_iterator TriangleRunner = Tess->TrianglesOnBoundary.begin(); TriangleRunner != Tess->TrianglesOnBoundary.end(); TriangleRunner++) {
|
---|
[57066a] | 726 | *vrmlfile << "1" << endl << " "; // 1 is triangle type
|
---|
| 727 | for (i=0;i<3;i++) { // print each node
|
---|
| 728 | for (int j=0;j<NDIM;j++) // and for each node all NDIM coordinates
|
---|
[0a4f7f] | 729 | *vrmlfile << TriangleRunner->second->endpoints[i]->node->node->at(j)-center->at(j) << " ";
|
---|
[57066a] | 730 | *vrmlfile << "\t";
|
---|
| 731 | }
|
---|
| 732 | *vrmlfile << "1. 0. 0." << endl; // red as colour
|
---|
| 733 | *vrmlfile << "18" << endl << " 0.5 0.5 0.5" << endl; // 18 is transparency type for previous object
|
---|
| 734 | }
|
---|
| 735 | } else {
|
---|
[58ed4a] | 736 | DoeLog(1) && (eLog()<< Verbose(1) << "Given vrmlfile is " << vrmlfile << "." << endl);
|
---|
[57066a] | 737 | }
|
---|
| 738 | delete(center);
|
---|
| 739 | };
|
---|
| 740 |
|
---|
| 741 | /** Writes additionally the current sphere (i.e. the last triangle to file).
|
---|
| 742 | * \param *out output stream for debugging
|
---|
| 743 | * \param *rasterfile output stream for tecplot data
|
---|
| 744 | * \param *Tess Tesselation structure with constructed triangles
|
---|
| 745 | * \param *mol molecule structure with atom positions
|
---|
| 746 | */
|
---|
[e138de] | 747 | void IncludeSphereinRaster3D(ofstream * const rasterfile, const Tesselation * const Tess, const PointCloud * const cloud)
|
---|
[57066a] | 748 | {
|
---|
[f67b6e] | 749 | Info FunctionInfo(__func__);
|
---|
[57066a] | 750 | Vector helper;
|
---|
[6a7f78c] | 751 |
|
---|
| 752 | if (Tess->LastTriangle != NULL) {
|
---|
| 753 | // include the current position of the virtual sphere in the temporary raster3d file
|
---|
| 754 | Vector *center = cloud->GetCenter();
|
---|
| 755 | // make the circumsphere's center absolute again
|
---|
[273382] | 756 | Vector helper = (1./3.) * ((*Tess->LastTriangle->endpoints[0]->node->node) +
|
---|
| 757 | (*Tess->LastTriangle->endpoints[1]->node->node) +
|
---|
| 758 | (*Tess->LastTriangle->endpoints[2]->node->node));
|
---|
| 759 | helper -= (*center);
|
---|
[6a7f78c] | 760 | // and add to file plus translucency object
|
---|
| 761 | *rasterfile << "# current virtual sphere\n";
|
---|
| 762 | *rasterfile << "8\n 25.0 0.6 -1.0 -1.0 -1.0 0.2 0 0 0 0\n";
|
---|
[0a4f7f] | 763 | *rasterfile << "2\n " << helper[0] << " " << helper[1] << " " << helper[2] << "\t" << 5. << "\t1 0 0\n";
|
---|
[6a7f78c] | 764 | *rasterfile << "9\n terminating special property\n";
|
---|
| 765 | delete(center);
|
---|
| 766 | }
|
---|
[57066a] | 767 | };
|
---|
| 768 |
|
---|
| 769 | /** Creates the objects in a raster3d file (renderable with a header.r3d).
|
---|
| 770 | * \param *out output stream for debugging
|
---|
| 771 | * \param *rasterfile output stream for tecplot data
|
---|
| 772 | * \param *Tess Tesselation structure with constructed triangles
|
---|
| 773 | * \param *mol molecule structure with atom positions
|
---|
| 774 | */
|
---|
[e138de] | 775 | void WriteRaster3dFile(ofstream * const rasterfile, const Tesselation * const Tess, const PointCloud * const cloud)
|
---|
[57066a] | 776 | {
|
---|
[f67b6e] | 777 | Info FunctionInfo(__func__);
|
---|
[57066a] | 778 | TesselPoint *Walker = NULL;
|
---|
| 779 | int i;
|
---|
[fc9992] | 780 | Vector *center = cloud->GetCenter();
|
---|
[57066a] | 781 | if (rasterfile != NULL) {
|
---|
[e138de] | 782 | //Log() << Verbose(1) << "Writing Raster3D file ... ";
|
---|
[57066a] | 783 | *rasterfile << "# Raster3D object description, created by MoleCuilder" << endl;
|
---|
| 784 | *rasterfile << "@header.r3d" << endl;
|
---|
| 785 | *rasterfile << "# All atoms as spheres" << endl;
|
---|
| 786 | cloud->GoToFirst();
|
---|
| 787 | while (!cloud->IsEnd()) {
|
---|
| 788 | Walker = cloud->GetPoint();
|
---|
| 789 | *rasterfile << "2" << endl << " "; // 2 is sphere type
|
---|
[15b670] | 790 | for (int j=0;j<NDIM;j++) { // and for each node all NDIM coordinates
|
---|
| 791 | const double tmp = Walker->node->at(j)-center->at(j);
|
---|
| 792 | *rasterfile << ((fabs(tmp) < MYEPSILON) ? 0 : tmp) << " ";
|
---|
| 793 | }
|
---|
[57066a] | 794 | *rasterfile << "\t0.1\t1. 1. 1." << endl; // radius 0.05 and white as colour
|
---|
| 795 | cloud->GoToNext();
|
---|
| 796 | }
|
---|
| 797 |
|
---|
| 798 | *rasterfile << "# All tesselation triangles" << endl;
|
---|
| 799 | *rasterfile << "8\n 25. -1. 1. 1. 1. 0.0 0 0 0 2\n SOLID 1.0 0.0 0.0\n BACKFACE 0.3 0.3 1.0 0 0\n";
|
---|
[776b64] | 800 | for (TriangleMap::const_iterator TriangleRunner = Tess->TrianglesOnBoundary.begin(); TriangleRunner != Tess->TrianglesOnBoundary.end(); TriangleRunner++) {
|
---|
[57066a] | 801 | *rasterfile << "1" << endl << " "; // 1 is triangle type
|
---|
| 802 | for (i=0;i<3;i++) { // print each node
|
---|
[15b670] | 803 | for (int j=0;j<NDIM;j++) { // and for each node all NDIM coordinates
|
---|
| 804 | const double tmp = TriangleRunner->second->endpoints[i]->node->node->at(j)-center->at(j);
|
---|
| 805 | *rasterfile << ((fabs(tmp) < MYEPSILON) ? 0 : tmp) << " ";
|
---|
| 806 | }
|
---|
[57066a] | 807 | *rasterfile << "\t";
|
---|
| 808 | }
|
---|
| 809 | *rasterfile << "1. 0. 0." << endl; // red as colour
|
---|
| 810 | //*rasterfile << "18" << endl << " 0.5 0.5 0.5" << endl; // 18 is transparency type for previous object
|
---|
| 811 | }
|
---|
| 812 | *rasterfile << "9\n# terminating special property\n";
|
---|
| 813 | } else {
|
---|
[58ed4a] | 814 | DoeLog(1) && (eLog()<< Verbose(1) << "Given rasterfile is " << rasterfile << "." << endl);
|
---|
[57066a] | 815 | }
|
---|
[e138de] | 816 | IncludeSphereinRaster3D(rasterfile, Tess, cloud);
|
---|
[57066a] | 817 | delete(center);
|
---|
| 818 | };
|
---|
| 819 |
|
---|
| 820 | /** This function creates the tecplot file, displaying the tesselation of the hull.
|
---|
| 821 | * \param *out output stream for debugging
|
---|
| 822 | * \param *tecplot output stream for tecplot data
|
---|
| 823 | * \param N arbitrary number to differentiate various zones in the tecplot format
|
---|
| 824 | */
|
---|
[e138de] | 825 | void WriteTecplotFile(ofstream * const tecplot, const Tesselation * const TesselStruct, const PointCloud * const cloud, const int N)
|
---|
[57066a] | 826 | {
|
---|
[f67b6e] | 827 | Info FunctionInfo(__func__);
|
---|
[57066a] | 828 | if ((tecplot != NULL) && (TesselStruct != NULL)) {
|
---|
| 829 | // write header
|
---|
| 830 | *tecplot << "TITLE = \"3D CONVEX SHELL\"" << endl;
|
---|
| 831 | *tecplot << "VARIABLES = \"X\" \"Y\" \"Z\" \"U\"" << endl;
|
---|
[6a7f78c] | 832 | *tecplot << "ZONE T=\"";
|
---|
| 833 | if (N < 0) {
|
---|
| 834 | *tecplot << cloud->GetName();
|
---|
| 835 | } else {
|
---|
| 836 | *tecplot << N << "-";
|
---|
[b60a29] | 837 | if (TesselStruct->LastTriangle != NULL) {
|
---|
| 838 | for (int i=0;i<3;i++)
|
---|
[68f03d] | 839 | *tecplot << (i==0 ? "" : "_") << TesselStruct->LastTriangle->endpoints[i]->node->getName();
|
---|
[b60a29] | 840 | } else {
|
---|
| 841 | *tecplot << "none";
|
---|
| 842 | }
|
---|
[6a7f78c] | 843 | }
|
---|
[57066a] | 844 | *tecplot << "\", N=" << TesselStruct->PointsOnBoundary.size() << ", E=" << TesselStruct->TrianglesOnBoundary.size() << ", DATAPACKING=POINT, ZONETYPE=FETRIANGLE" << endl;
|
---|
[15b670] | 845 | const int MaxId=cloud->GetMaxId();
|
---|
| 846 | int *LookupList = new int[MaxId];
|
---|
| 847 | for (int i=0; i< MaxId ; i++){
|
---|
[57066a] | 848 | LookupList[i] = -1;
|
---|
[c72112] | 849 | }
|
---|
[57066a] | 850 |
|
---|
| 851 | // print atom coordinates
|
---|
| 852 | int Counter = 1;
|
---|
| 853 | TesselPoint *Walker = NULL;
|
---|
[c72112] | 854 | for (PointMap::const_iterator target = TesselStruct->PointsOnBoundary.begin(); target != TesselStruct->PointsOnBoundary.end(); ++target) {
|
---|
[57066a] | 855 | Walker = target->second->node;
|
---|
| 856 | LookupList[Walker->nr] = Counter++;
|
---|
[15b670] | 857 | for (int i=0;i<NDIM;i++) {
|
---|
| 858 | const double tmp = Walker->node->at(i);
|
---|
| 859 | *tecplot << ((fabs(tmp) < MYEPSILON) ? 0 : tmp) << " ";
|
---|
| 860 | }
|
---|
| 861 | *tecplot << target->second->value << endl;
|
---|
[57066a] | 862 | }
|
---|
| 863 | *tecplot << endl;
|
---|
| 864 | // print connectivity
|
---|
[a67d19] | 865 | DoLog(1) && (Log() << Verbose(1) << "The following triangles were created:" << endl);
|
---|
[776b64] | 866 | for (TriangleMap::const_iterator runner = TesselStruct->TrianglesOnBoundary.begin(); runner != TesselStruct->TrianglesOnBoundary.end(); runner++) {
|
---|
[68f03d] | 867 | DoLog(1) && (Log() << Verbose(1) << " " << runner->second->endpoints[0]->node->getName() << "<->" << runner->second->endpoints[1]->node->getName() << "<->" << runner->second->endpoints[2]->node->getName() << endl);
|
---|
[57066a] | 868 | *tecplot << LookupList[runner->second->endpoints[0]->node->nr] << " " << LookupList[runner->second->endpoints[1]->node->nr] << " " << LookupList[runner->second->endpoints[2]->node->nr] << endl;
|
---|
| 869 | }
|
---|
| 870 | delete[] (LookupList);
|
---|
| 871 | }
|
---|
| 872 | };
|
---|
[7dea7c] | 873 |
|
---|
| 874 | /** Calculates the concavity for each of the BoundaryPointSet's in a Tesselation.
|
---|
| 875 | * Sets BoundaryPointSet::value equal to the number of connected lines that are not convex.
|
---|
| 876 | * \param *out output stream for debugging
|
---|
| 877 | * \param *TesselStruct pointer to Tesselation structure
|
---|
| 878 | */
|
---|
[e138de] | 879 | void CalculateConcavityPerBoundaryPoint(const Tesselation * const TesselStruct)
|
---|
[7dea7c] | 880 | {
|
---|
[f67b6e] | 881 | Info FunctionInfo(__func__);
|
---|
[7dea7c] | 882 | class BoundaryPointSet *point = NULL;
|
---|
| 883 | class BoundaryLineSet *line = NULL;
|
---|
| 884 |
|
---|
| 885 | // calculate remaining concavity
|
---|
[776b64] | 886 | for (PointMap::const_iterator PointRunner = TesselStruct->PointsOnBoundary.begin(); PointRunner != TesselStruct->PointsOnBoundary.end(); PointRunner++) {
|
---|
[7dea7c] | 887 | point = PointRunner->second;
|
---|
[a67d19] | 888 | DoLog(1) && (Log() << Verbose(1) << "INFO: Current point is " << *point << "." << endl);
|
---|
[7dea7c] | 889 | point->value = 0;
|
---|
| 890 | for (LineMap::iterator LineRunner = point->lines.begin(); LineRunner != point->lines.end(); LineRunner++) {
|
---|
| 891 | line = LineRunner->second;
|
---|
[f67b6e] | 892 | //Log() << Verbose(1) << "INFO: Current line of point " << *point << " is " << *line << "." << endl;
|
---|
[e138de] | 893 | if (!line->CheckConvexityCriterion())
|
---|
[7dea7c] | 894 | point->value += 1;
|
---|
| 895 | }
|
---|
| 896 | }
|
---|
| 897 | };
|
---|
| 898 |
|
---|
| 899 |
|
---|
| 900 | /** Checks whether each BoundaryLineSet in the Tesselation has two triangles.
|
---|
| 901 | * \param *out output stream for debugging
|
---|
| 902 | * \param *TesselStruct
|
---|
| 903 | * \return true - all have exactly two triangles, false - some not, list is printed to screen
|
---|
| 904 | */
|
---|
[e138de] | 905 | bool CheckListOfBaselines(const Tesselation * const TesselStruct)
|
---|
[7dea7c] | 906 | {
|
---|
[f67b6e] | 907 | Info FunctionInfo(__func__);
|
---|
[776b64] | 908 | LineMap::const_iterator testline;
|
---|
[7dea7c] | 909 | bool result = false;
|
---|
| 910 | int counter = 0;
|
---|
| 911 |
|
---|
[a67d19] | 912 | DoLog(1) && (Log() << Verbose(1) << "Check: List of Baselines with not two connected triangles:" << endl);
|
---|
[7dea7c] | 913 | for (testline = TesselStruct->LinesOnBoundary.begin(); testline != TesselStruct->LinesOnBoundary.end(); testline++) {
|
---|
| 914 | if (testline->second->triangles.size() != 2) {
|
---|
[a67d19] | 915 | DoLog(2) && (Log() << Verbose(2) << *testline->second << "\t" << testline->second->triangles.size() << endl);
|
---|
[7dea7c] | 916 | counter++;
|
---|
| 917 | }
|
---|
| 918 | }
|
---|
| 919 | if (counter == 0) {
|
---|
[a67d19] | 920 | DoLog(1) && (Log() << Verbose(1) << "None." << endl);
|
---|
[7dea7c] | 921 | result = true;
|
---|
| 922 | }
|
---|
| 923 | return result;
|
---|
| 924 | }
|
---|
| 925 |
|
---|
[262bae] | 926 | /** Counts the number of triangle pairs that contain the given polygon.
|
---|
| 927 | * \param *P polygon with endpoints to look for
|
---|
| 928 | * \param *T set of triangles to create pairs from containing \a *P
|
---|
| 929 | */
|
---|
| 930 | int CountTrianglePairContainingPolygon(const BoundaryPolygonSet * const P, const TriangleSet * const T)
|
---|
| 931 | {
|
---|
| 932 | Info FunctionInfo(__func__);
|
---|
| 933 | // check number of endpoints in *P
|
---|
| 934 | if (P->endpoints.size() != 4) {
|
---|
[58ed4a] | 935 | DoeLog(1) && (eLog()<< Verbose(1) << "CountTrianglePairContainingPolygon works only on polygons with 4 nodes!" << endl);
|
---|
[262bae] | 936 | return 0;
|
---|
| 937 | }
|
---|
| 938 |
|
---|
| 939 | // check number of triangles in *T
|
---|
| 940 | if (T->size() < 2) {
|
---|
[58ed4a] | 941 | DoeLog(1) && (eLog()<< Verbose(1) << "Not enough triangles to have pairs!" << endl);
|
---|
[262bae] | 942 | return 0;
|
---|
| 943 | }
|
---|
| 944 |
|
---|
[a67d19] | 945 | DoLog(0) && (Log() << Verbose(0) << "Polygon is " << *P << endl);
|
---|
[262bae] | 946 | // create each pair, get the endpoints and check whether *P is contained.
|
---|
| 947 | int counter = 0;
|
---|
| 948 | PointSet Trianglenodes;
|
---|
| 949 | class BoundaryPolygonSet PairTrianglenodes;
|
---|
| 950 | for(TriangleSet::iterator Walker = T->begin(); Walker != T->end(); Walker++) {
|
---|
| 951 | for (int i=0;i<3;i++)
|
---|
| 952 | Trianglenodes.insert((*Walker)->endpoints[i]);
|
---|
| 953 |
|
---|
| 954 | for(TriangleSet::iterator PairWalker = Walker; PairWalker != T->end(); PairWalker++) {
|
---|
| 955 | if (Walker != PairWalker) { // skip first
|
---|
| 956 | PairTrianglenodes.endpoints = Trianglenodes;
|
---|
| 957 | for (int i=0;i<3;i++)
|
---|
| 958 | PairTrianglenodes.endpoints.insert((*PairWalker)->endpoints[i]);
|
---|
[856098] | 959 | const int size = PairTrianglenodes.endpoints.size();
|
---|
| 960 | if (size == 4) {
|
---|
[a67d19] | 961 | DoLog(0) && (Log() << Verbose(0) << " Current pair of triangles: " << **Walker << "," << **PairWalker << " with " << size << " distinct endpoints:" << PairTrianglenodes << endl);
|
---|
[856098] | 962 | // now check
|
---|
| 963 | if (PairTrianglenodes.ContainsPresentTupel(P)) {
|
---|
| 964 | counter++;
|
---|
[a67d19] | 965 | DoLog(0) && (Log() << Verbose(0) << " ACCEPT: Matches with " << *P << endl);
|
---|
[856098] | 966 | } else {
|
---|
[a67d19] | 967 | DoLog(0) && (Log() << Verbose(0) << " REJECT: No match with " << *P << endl);
|
---|
[856098] | 968 | }
|
---|
[262bae] | 969 | } else {
|
---|
[a67d19] | 970 | DoLog(0) && (Log() << Verbose(0) << " REJECT: Less than four endpoints." << endl);
|
---|
[262bae] | 971 | }
|
---|
| 972 | }
|
---|
| 973 | }
|
---|
[856098] | 974 | Trianglenodes.clear();
|
---|
[262bae] | 975 | }
|
---|
| 976 | return counter;
|
---|
| 977 | };
|
---|
| 978 |
|
---|
| 979 | /** Checks whether two give polygons have two or more points in common.
|
---|
| 980 | * \param *P1 first polygon
|
---|
| 981 | * \param *P2 second polygon
|
---|
| 982 | * \return true - are connected, false = are note
|
---|
| 983 | */
|
---|
| 984 | bool ArePolygonsEdgeConnected(const BoundaryPolygonSet * const P1, const BoundaryPolygonSet * const P2)
|
---|
| 985 | {
|
---|
| 986 | Info FunctionInfo(__func__);
|
---|
| 987 | int counter = 0;
|
---|
| 988 | for(PointSet::const_iterator Runner = P1->endpoints.begin(); Runner != P1->endpoints.end(); Runner++) {
|
---|
| 989 | if (P2->ContainsBoundaryPoint((*Runner))) {
|
---|
| 990 | counter++;
|
---|
[a67d19] | 991 | DoLog(1) && (Log() << Verbose(1) << *(*Runner) << " of second polygon is found in the first one." << endl);
|
---|
[262bae] | 992 | return true;
|
---|
| 993 | }
|
---|
| 994 | }
|
---|
| 995 | return false;
|
---|
| 996 | };
|
---|
| 997 |
|
---|
| 998 | /** Combines second into the first and deletes the second.
|
---|
| 999 | * \param *P1 first polygon, contains all nodes on return
|
---|
| 1000 | * \param *&P2 second polygon, is deleted.
|
---|
| 1001 | */
|
---|
| 1002 | void CombinePolygons(BoundaryPolygonSet * const P1, BoundaryPolygonSet * &P2)
|
---|
| 1003 | {
|
---|
| 1004 | Info FunctionInfo(__func__);
|
---|
[856098] | 1005 | pair <PointSet::iterator, bool> Tester;
|
---|
| 1006 | for(PointSet::iterator Runner = P2->endpoints.begin(); Runner != P2->endpoints.end(); Runner++) {
|
---|
| 1007 | Tester = P1->endpoints.insert((*Runner));
|
---|
| 1008 | if (Tester.second)
|
---|
[a67d19] | 1009 | DoLog(0) && (Log() << Verbose(0) << "Inserting endpoint " << *(*Runner) << " into first polygon." << endl);
|
---|
[262bae] | 1010 | }
|
---|
| 1011 | P2->endpoints.clear();
|
---|
| 1012 | delete(P2);
|
---|
| 1013 | };
|
---|
| 1014 |
|
---|