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