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