source: molecuilder/src/boundary.cpp@ dd69cc

Last change on this file since dd69cc was dd69cc, checked in by Frederik Heber <heber@…>, 16 years ago

FIX: DetermineCenterOfAll() did point in wrong direction, convex envelope working again, but algorithm is not stable

  • molecule::DetermineCenterOfAll() center was the sum of all scaled by -1/Number instead of 1/Number. This was corrected for in the code elsewhere. Now, it returns the correct center and where it's elsewhere called we subtract instead of add
  • Tesselation::TesselateOnBoundary() does now something meaningful. The NormalVector for the starting triangle could not be constructed, as we lacked the initial normal vector. It is however easy to construct from the center of all and the center of the starting triangle. However, the algorithm is faulty for 1_2_dimethylethane.
  • GetBoundaryPoints() - as we do not translate to center of gravity anymore, we have to subtract the center of all. Otherwise the radial method for determining boundary points does not work.

Is not yet working.

  • Property mode set to 100755
File size: 123.4 KB
Line 
1#include "boundary.hpp"
2#include "linkedcell.hpp"
3#include "molecules.hpp"
4#include <gsl/gsl_matrix.h>
5#include <gsl/gsl_linalg.h>
6#include <gsl/gsl_multimin.h>
7#include <gsl/gsl_permutation.h>
8
9#define DEBUG 1
10#define DoSingleStepOutput 0
11#define DoTecplotOutput 1
12#define DoRaster3DOutput 1
13#define DoVRMLOutput 1
14#define TecplotSuffix ".dat"
15#define Raster3DSuffix ".r3d"
16#define VRMLSUffix ".wrl"
17#define HULLEPSILON 1e-7
18
19// ======================================== Points on Boundary =================================
20
21BoundaryPointSet::BoundaryPointSet()
22{
23 LinesCount = 0;
24 Nr = -1;
25}
26;
27
28BoundaryPointSet::BoundaryPointSet(atom *Walker)
29{
30 node = Walker;
31 LinesCount = 0;
32 Nr = Walker->nr;
33}
34;
35
36BoundaryPointSet::~BoundaryPointSet()
37{
38 cout << Verbose(5) << "Erasing point nr. " << Nr << "." << endl;
39 if (!lines.empty())
40 cerr << "WARNING: Memory Leak! I " << *this << " am still connected to some lines." << endl;
41 node = NULL;
42}
43;
44
45void BoundaryPointSet::AddLine(class BoundaryLineSet *line)
46{
47 cout << Verbose(6) << "Adding " << *this << " to line " << *line << "."
48 << endl;
49 if (line->endpoints[0] == this)
50 {
51 lines.insert(LinePair(line->endpoints[1]->Nr, line));
52 }
53 else
54 {
55 lines.insert(LinePair(line->endpoints[0]->Nr, line));
56 }
57 LinesCount++;
58}
59;
60
61ostream &
62operator <<(ostream &ost, BoundaryPointSet &a)
63{
64 ost << "[" << a.Nr << "|" << a.node->Name << "]";
65 return ost;
66}
67;
68
69// ======================================== Lines on Boundary =================================
70
71BoundaryLineSet::BoundaryLineSet()
72{
73 for (int i = 0; i < 2; i++)
74 endpoints[i] = NULL;
75 TrianglesCount = 0;
76 Nr = -1;
77}
78;
79
80BoundaryLineSet::BoundaryLineSet(class BoundaryPointSet *Point[2], int number)
81{
82 // set number
83 Nr = number;
84 // set endpoints in ascending order
85 SetEndpointsOrdered(endpoints, Point[0], Point[1]);
86 // add this line to the hash maps of both endpoints
87 Point[0]->AddLine(this); //Taken out, to check whether we can avoid unwanted double adding.
88 Point[1]->AddLine(this); //
89 // clear triangles list
90 TrianglesCount = 0;
91 cout << Verbose(5) << "New Line with endpoints " << *this << "." << endl;
92}
93;
94
95BoundaryLineSet::~BoundaryLineSet()
96{
97 int Numbers[2];
98 Numbers[0] = endpoints[1]->Nr;
99 Numbers[1] = endpoints[0]->Nr;
100 for (int i = 0; i < 2; i++) {
101 cout << Verbose(5) << "Erasing Line Nr. " << Nr << " in boundary point " << *endpoints[i] << "." << endl;
102 // as there may be multiple lines with same endpoints, we have to go through each and find in the endpoint's line list this line set
103 pair<LineMap::iterator, LineMap::iterator> erasor = endpoints[i]->lines.equal_range(Numbers[i]);
104 for (LineMap::iterator Runner = erasor.first; Runner != erasor.second; Runner++)
105 if ((*Runner).second == this) {
106 endpoints[i]->lines.erase(Runner);
107 break;
108 }
109 if (endpoints[i]->lines.empty()) {
110 cout << Verbose(5) << *endpoints[i] << " has no more lines it's attached to, erasing." << endl;
111 if (endpoints[i] != NULL) {
112 delete(endpoints[i]);
113 endpoints[i] = NULL;
114 } else
115 cerr << "ERROR: Endpoint " << i << " has already been free'd." << endl;
116 } else
117 cout << Verbose(5) << *endpoints[i] << " has still lines it's attached to." << endl;
118 }
119 if (!triangles.empty())
120 cerr << "WARNING: Memory Leak! I " << *this << " am still connected to some triangles." << endl;
121}
122;
123
124void
125BoundaryLineSet::AddTriangle(class BoundaryTriangleSet *triangle)
126{
127 cout << Verbose(6) << "Add " << triangle->Nr << " to line " << *this << "."
128 << endl;
129 triangles.insert(TrianglePair(triangle->Nr, triangle));
130 TrianglesCount++;
131}
132;
133
134ostream &
135operator <<(ostream &ost, BoundaryLineSet &a)
136{
137 ost << "[" << a.Nr << "|" << a.endpoints[0]->node->Name << ","
138 << a.endpoints[1]->node->Name << "]";
139 return ost;
140}
141;
142
143// ======================================== Triangles on Boundary =================================
144
145
146BoundaryTriangleSet::BoundaryTriangleSet()
147{
148 for (int i = 0; i < 3; i++)
149 {
150 endpoints[i] = NULL;
151 lines[i] = NULL;
152 }
153 Nr = -1;
154}
155;
156
157BoundaryTriangleSet::BoundaryTriangleSet(class BoundaryLineSet *line[3], int number)
158{
159 // set number
160 Nr = number;
161 // set lines
162 cout << Verbose(5) << "New triangle " << Nr << ":" << endl;
163 for (int i = 0; i < 3; i++)
164 {
165 lines[i] = line[i];
166 lines[i]->AddTriangle(this);
167 }
168 // get ascending order of endpoints
169 map<int, class BoundaryPointSet *> OrderMap;
170 for (int i = 0; i < 3; i++)
171 // for all three lines
172 for (int j = 0; j < 2; j++)
173 { // for both endpoints
174 OrderMap.insert(pair<int, class BoundaryPointSet *> (
175 line[i]->endpoints[j]->Nr, line[i]->endpoints[j]));
176 // and we don't care whether insertion fails
177 }
178 // set endpoints
179 int Counter = 0;
180 cout << Verbose(6) << " with end points ";
181 for (map<int, class BoundaryPointSet *>::iterator runner = OrderMap.begin(); runner
182 != OrderMap.end(); runner++)
183 {
184 endpoints[Counter] = runner->second;
185 cout << " " << *endpoints[Counter];
186 Counter++;
187 }
188 if (Counter < 3)
189 {
190 cerr << "ERROR! We have a triangle with only two distinct endpoints!"
191 << endl;
192 //exit(1);
193 }
194 cout << "." << endl;
195}
196;
197
198BoundaryTriangleSet::~BoundaryTriangleSet()
199{
200 for (int i = 0; i < 3; i++) {
201 cout << Verbose(5) << "Erasing triangle Nr." << Nr << endl;
202 lines[i]->triangles.erase(Nr);
203 if (lines[i]->triangles.empty()) {
204 if (lines[i] != NULL) {
205 cout << Verbose(5) << *lines[i] << " is no more attached to any triangle, erasing." << endl;
206 delete (lines[i]);
207 lines[i] = NULL;
208 } else
209 cerr << "ERROR: This line " << i << " has already been free'd." << endl;
210 } else
211 cout << Verbose(5) << *lines[i] << " is still attached to another triangle." << endl;
212 }
213}
214;
215
216void
217BoundaryTriangleSet::GetNormalVector(Vector &OtherVector)
218{
219 // get normal vector
220 NormalVector.MakeNormalVector(&endpoints[0]->node->x, &endpoints[1]->node->x,
221 &endpoints[2]->node->x);
222
223 // make it always point inward (any offset vector onto plane projected onto normal vector suffices)
224 if (NormalVector.Projection(&OtherVector) > 0)
225 NormalVector.Scale(-1.);
226}
227;
228
229ostream &
230operator <<(ostream &ost, BoundaryTriangleSet &a)
231{
232 ost << "[" << a.Nr << "|" << a.endpoints[0]->node->Name << ","
233 << a.endpoints[1]->node->Name << "," << a.endpoints[2]->node->Name << "]";
234 return ost;
235}
236;
237
238
239// ============================ CandidateForTesselation =============================
240
241CandidateForTesselation::CandidateForTesselation(
242 atom *candidate, BoundaryLineSet* line, Vector OptCandidateCenter, Vector OtherOptCandidateCenter
243) {
244 point = candidate;
245 BaseLine = line;
246 OptCenter.CopyVector(&OptCandidateCenter);
247 OtherOptCenter.CopyVector(&OtherOptCandidateCenter);
248}
249
250CandidateForTesselation::~CandidateForTesselation() {
251 point = NULL;
252 BaseLine = NULL;
253}
254
255// ========================================== F U N C T I O N S =================================
256
257/** Finds the endpoint two lines are sharing.
258 * \param *line1 first line
259 * \param *line2 second line
260 * \return point which is shared or NULL if none
261 */
262class BoundaryPointSet *
263GetCommonEndpoint(class BoundaryLineSet * line1, class BoundaryLineSet * line2)
264{
265 class BoundaryLineSet * lines[2] =
266 { line1, line2 };
267 class BoundaryPointSet *node = NULL;
268 map<int, class BoundaryPointSet *> OrderMap;
269 pair<map<int, class BoundaryPointSet *>::iterator, bool> OrderTest;
270 for (int i = 0; i < 2; i++)
271 // for both lines
272 for (int j = 0; j < 2; j++)
273 { // for both endpoints
274 OrderTest = OrderMap.insert(pair<int, class BoundaryPointSet *> (
275 lines[i]->endpoints[j]->Nr, lines[i]->endpoints[j]));
276 if (!OrderTest.second)
277 { // if insertion fails, we have common endpoint
278 node = OrderTest.first->second;
279 cout << Verbose(5) << "Common endpoint of lines " << *line1
280 << " and " << *line2 << " is: " << *node << "." << endl;
281 j = 2;
282 i = 2;
283 break;
284 }
285 }
286 return node;
287}
288;
289
290/** Determines the boundary points of a cluster.
291 * Does a projection per axis onto the orthogonal plane, transforms into spherical coordinates, sorts them by the angle
292 * and looks at triples: if the middle has less a distance than the allowed maximum height of the triangle formed by the plane's
293 * center and first and last point in the triple, it is thrown out.
294 * \param *out output stream for debugging
295 * \param *mol molecule structure representing the cluster
296 */
297Boundaries *
298GetBoundaryPoints(ofstream *out, molecule *mol)
299{
300 atom *Walker = NULL;
301 PointMap PointsOnBoundary;
302 LineMap LinesOnBoundary;
303 TriangleMap TrianglesOnBoundary;
304 Vector *MolCenter = mol->DetermineCenterOfAll(out);
305 Vector helper;
306
307 *out << Verbose(1) << "Finding all boundary points." << endl;
308 Boundaries *BoundaryPoints = new Boundaries[NDIM]; // first is alpha, second is (r, nr)
309 BoundariesTestPair BoundaryTestPair;
310 Vector AxisVector, AngleReferenceVector, AngleReferenceNormalVector;
311 double radius, angle;
312 // 3a. Go through every axis
313 for (int axis = 0; axis < NDIM; axis++) {
314 AxisVector.Zero();
315 AngleReferenceVector.Zero();
316 AngleReferenceNormalVector.Zero();
317 AxisVector.x[axis] = 1.;
318 AngleReferenceVector.x[(axis + 1) % NDIM] = 1.;
319 AngleReferenceNormalVector.x[(axis + 2) % NDIM] = 1.;
320
321 //*out << Verbose(1) << "Axisvector is " << AxisVector << " and AngleReferenceVector is " << AngleReferenceVector << ", and AngleReferenceNormalVector is " << AngleReferenceNormalVector << "." << endl;
322
323 // 3b. construct set of all points, transformed into cylindrical system and with left and right neighbours
324 Walker = mol->start;
325 while (Walker->next != mol->end) {
326 Walker = Walker->next;
327 Vector ProjectedVector;
328 ProjectedVector.CopyVector(&Walker->x);
329 ProjectedVector.SubtractVector(MolCenter);
330 ProjectedVector.ProjectOntoPlane(&AxisVector);
331
332 // correct for negative side
333 radius = ProjectedVector.Norm();
334 if (fabs(radius) > MYEPSILON)
335 angle = ProjectedVector.Angle(&AngleReferenceVector);
336 else
337 angle = 0.; // otherwise it's a vector in Axis Direction and unimportant for boundary issues
338
339 //*out << "Checking sign in quadrant : " << ProjectedVector.Projection(&AngleReferenceNormalVector) << "." << endl;
340 if (ProjectedVector.Projection(&AngleReferenceNormalVector) > 0) {
341 angle = 2. * M_PI - angle;
342 }
343 //*out << Verbose(2) << "Inserting " << *Walker << ": (r, alpha) = (" << radius << "," << angle << "): " << ProjectedVector << endl;
344 BoundaryTestPair = BoundaryPoints[axis].insert(BoundariesPair(angle,
345 DistancePair (radius, Walker)));
346 if (BoundaryTestPair.second) { // successfully inserted
347 } else { // same point exists, check first r, then distance of original vectors to center of gravity
348 *out << Verbose(2) << "Encountered two vectors whose projection onto axis " << axis << " is equal: " << endl;
349 *out << Verbose(2) << "Present vector: " << BoundaryTestPair.first->second.second->x << endl;
350 *out << Verbose(2) << "New vector: " << Walker << endl;
351 double tmp = ProjectedVector.Norm();
352 if (tmp > BoundaryTestPair.first->second.first) {
353 BoundaryTestPair.first->second.first = tmp;
354 BoundaryTestPair.first->second.second = Walker;
355 *out << Verbose(2) << "Keeping new vector." << endl;
356 } else if (tmp == BoundaryTestPair.first->second.first) {
357 helper.CopyVector(&Walker->x);
358 helper.SubtractVector(MolCenter);
359 if (BoundaryTestPair.first->second.second->x.ScalarProduct(&BoundaryTestPair.first->second.second->x) < helper.ScalarProduct(&helper)) { // Norm() does a sqrt, which makes it a lot slower
360 BoundaryTestPair.first->second.second = Walker;
361 *out << Verbose(2) << "Keeping new vector." << endl;
362 } else {
363 *out << Verbose(2) << "Keeping present vector." << endl;
364 }
365 } else {
366 *out << Verbose(2) << "Keeping present vector." << endl;
367 }
368 }
369 }
370 // printing all inserted for debugging
371 // {
372 // *out << Verbose(2) << "Printing list of candidates for axis " << axis << " which we have inserted so far." << endl;
373 // int i=0;
374 // for(Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
375 // if (runner != BoundaryPoints[axis].begin())
376 // *out << ", " << i << ": " << *runner->second.second;
377 // else
378 // *out << i << ": " << *runner->second.second;
379 // i++;
380 // }
381 // *out << endl;
382 // }
383 // 3c. throw out points whose distance is less than the mean of left and right neighbours
384 bool flag = false;
385 *out << Verbose(1) << "Looking for candidates to kick out by convex condition ... " << endl;
386 do { // do as long as we still throw one out per round
387 flag = false;
388 Boundaries::iterator left = BoundaryPoints[axis].end();
389 Boundaries::iterator right = BoundaryPoints[axis].end();
390 for (Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
391 // set neighbours correctly
392 if (runner == BoundaryPoints[axis].begin()) {
393 left = BoundaryPoints[axis].end();
394 } else {
395 left = runner;
396 }
397 left--;
398 right = runner;
399 right++;
400 if (right == BoundaryPoints[axis].end()) {
401 right = BoundaryPoints[axis].begin();
402 }
403 // check distance
404
405 // construct the vector of each side of the triangle on the projected plane (defined by normal vector AxisVector)
406 {
407 Vector SideA, SideB, SideC, SideH;
408 SideA.CopyVector(&left->second.second->x);
409 SideA.SubtractVector(MolCenter);
410 SideA.ProjectOntoPlane(&AxisVector);
411 // *out << "SideA: ";
412 // SideA.Output(out);
413 // *out << endl;
414
415 SideB.CopyVector(&right->second.second->x);
416 SideB.SubtractVector(MolCenter);
417 SideB.ProjectOntoPlane(&AxisVector);
418 // *out << "SideB: ";
419 // SideB.Output(out);
420 // *out << endl;
421
422 SideC.CopyVector(&left->second.second->x);
423 SideC.SubtractVector(&right->second.second->x);
424 SideC.ProjectOntoPlane(&AxisVector);
425 // *out << "SideC: ";
426 // SideC.Output(out);
427 // *out << endl;
428
429 SideH.CopyVector(&runner->second.second->x);
430 SideH.SubtractVector(MolCenter);
431 SideH.ProjectOntoPlane(&AxisVector);
432 // *out << "SideH: ";
433 // SideH.Output(out);
434 // *out << endl;
435
436 // calculate each length
437 double a = SideA.Norm();
438 //double b = SideB.Norm();
439 //double c = SideC.Norm();
440 double h = SideH.Norm();
441 // calculate the angles
442 double alpha = SideA.Angle(&SideH);
443 double beta = SideA.Angle(&SideC);
444 double gamma = SideB.Angle(&SideH);
445 double delta = SideC.Angle(&SideH);
446 double MinDistance = a * sin(beta) / (sin(delta)) * (((alpha < M_PI / 2.) || (gamma < M_PI / 2.)) ? 1. : -1.);
447 //*out << Verbose(2) << " I calculated: a = " << a << ", h = " << h << ", beta(" << left->second.second->Name << "," << left->second.second->Name << "-" << right->second.second->Name << ") = " << beta << ", delta(" << left->second.second->Name << "," << runner->second.second->Name << ") = " << delta << ", Min = " << MinDistance << "." << endl;
448 //*out << Verbose(1) << "Checking CoG distance of runner " << *runner->second.second << " " << h << " against triangle's side length spanned by (" << *left->second.second << "," << *right->second.second << ") of " << MinDistance << "." << endl;
449 if ((fabs(h / fabs(h) - MinDistance / fabs(MinDistance)) < MYEPSILON) && (h < MinDistance)) {
450 // throw out point
451 //*out << Verbose(1) << "Throwing out " << *runner->second.second << "." << endl;
452 BoundaryPoints[axis].erase(runner);
453 flag = true;
454 }
455 }
456 }
457 } while (flag);
458 }
459 delete(MolCenter);
460 return BoundaryPoints;
461}
462;
463
464/** Determines greatest diameters of a cluster defined by its convex envelope.
465 * Looks at lines parallel to one axis and where they intersect on the projected planes
466 * \param *out output stream for debugging
467 * \param *BoundaryPoints NDIM set of boundary points defining the convex envelope on each projected plane
468 * \param *mol molecule structure representing the cluster
469 * \param IsAngstroem whether we have angstroem or atomic units
470 * \return NDIM array of the diameters
471 */
472double *
473GetDiametersOfCluster(ofstream *out, Boundaries *BoundaryPtr, molecule *mol,
474 bool IsAngstroem)
475{
476 // get points on boundary of NULL was given as parameter
477 bool BoundaryFreeFlag = false;
478 Boundaries *BoundaryPoints = BoundaryPtr;
479 if (BoundaryPoints == NULL)
480 {
481 BoundaryFreeFlag = true;
482 BoundaryPoints = GetBoundaryPoints(out, mol);
483 }
484 else
485 {
486 *out << Verbose(1) << "Using given boundary points set." << endl;
487 }
488 // determine biggest "diameter" of cluster for each axis
489 Boundaries::iterator Neighbour, OtherNeighbour;
490 double *GreatestDiameter = new double[NDIM];
491 for (int i = 0; i < NDIM; i++)
492 GreatestDiameter[i] = 0.;
493 double OldComponent, tmp, w1, w2;
494 Vector DistanceVector, OtherVector;
495 int component, Othercomponent;
496 for (int axis = 0; axis < NDIM; axis++)
497 { // regard each projected plane
498 //*out << Verbose(1) << "Current axis is " << axis << "." << endl;
499 for (int j = 0; j < 2; j++)
500 { // and for both axis on the current plane
501 component = (axis + j + 1) % NDIM;
502 Othercomponent = (axis + 1 + ((j + 1) & 1)) % NDIM;
503 //*out << Verbose(1) << "Current component is " << component << ", Othercomponent is " << Othercomponent << "." << endl;
504 for (Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner
505 != BoundaryPoints[axis].end(); runner++)
506 {
507 //*out << Verbose(2) << "Current runner is " << *(runner->second.second) << "." << endl;
508 // seek for the neighbours pair where the Othercomponent sign flips
509 Neighbour = runner;
510 Neighbour++;
511 if (Neighbour == BoundaryPoints[axis].end()) // make it wrap around
512 Neighbour = BoundaryPoints[axis].begin();
513 DistanceVector.CopyVector(&runner->second.second->x);
514 DistanceVector.SubtractVector(&Neighbour->second.second->x);
515 do
516 { // seek for neighbour pair where it flips
517 OldComponent = DistanceVector.x[Othercomponent];
518 Neighbour++;
519 if (Neighbour == BoundaryPoints[axis].end()) // make it wrap around
520 Neighbour = BoundaryPoints[axis].begin();
521 DistanceVector.CopyVector(&runner->second.second->x);
522 DistanceVector.SubtractVector(&Neighbour->second.second->x);
523 //*out << Verbose(3) << "OldComponent is " << OldComponent << ", new one is " << DistanceVector.x[Othercomponent] << "." << endl;
524 }
525 while ((runner != Neighbour) && (fabs(OldComponent / fabs(
526 OldComponent) - DistanceVector.x[Othercomponent] / fabs(
527 DistanceVector.x[Othercomponent])) < MYEPSILON)); // as long as sign does not flip
528 if (runner != Neighbour)
529 {
530 OtherNeighbour = Neighbour;
531 if (OtherNeighbour == BoundaryPoints[axis].begin()) // make it wrap around
532 OtherNeighbour = BoundaryPoints[axis].end();
533 OtherNeighbour--;
534 //*out << Verbose(2) << "The pair, where the sign of OtherComponent flips, is: " << *(Neighbour->second.second) << " and " << *(OtherNeighbour->second.second) << "." << endl;
535 // now we have found the pair: Neighbour and OtherNeighbour
536 OtherVector.CopyVector(&runner->second.second->x);
537 OtherVector.SubtractVector(&OtherNeighbour->second.second->x);
538 //*out << Verbose(2) << "Distances to Neighbour and OtherNeighbour are " << DistanceVector.x[component] << " and " << OtherVector.x[component] << "." << endl;
539 //*out << Verbose(2) << "OtherComponents to Neighbour and OtherNeighbour are " << DistanceVector.x[Othercomponent] << " and " << OtherVector.x[Othercomponent] << "." << endl;
540 // do linear interpolation between points (is exact) to extract exact intersection between Neighbour and OtherNeighbour
541 w1 = fabs(OtherVector.x[Othercomponent]);
542 w2 = fabs(DistanceVector.x[Othercomponent]);
543 tmp = fabs((w1 * DistanceVector.x[component] + w2
544 * OtherVector.x[component]) / (w1 + w2));
545 // mark if it has greater diameter
546 //*out << Verbose(2) << "Comparing current greatest " << GreatestDiameter[component] << " to new " << tmp << "." << endl;
547 GreatestDiameter[component] = (GreatestDiameter[component]
548 > tmp) ? GreatestDiameter[component] : tmp;
549 } //else
550 //*out << Verbose(2) << "Saw no sign flip, probably top or bottom node." << endl;
551 }
552 }
553 }
554 *out << Verbose(0) << "RESULT: The biggest diameters are "
555 << GreatestDiameter[0] << " and " << GreatestDiameter[1] << " and "
556 << GreatestDiameter[2] << " " << (IsAngstroem ? "angstrom"
557 : "atomiclength") << "." << endl;
558
559 // free reference lists
560 if (BoundaryFreeFlag)
561 delete[] (BoundaryPoints);
562
563 return GreatestDiameter;
564}
565;
566
567/** Creates the objects in a VRML file.
568 * \param *out output stream for debugging
569 * \param *vrmlfile output stream for tecplot data
570 * \param *Tess Tesselation structure with constructed triangles
571 * \param *mol molecule structure with atom positions
572 */
573void write_vrml_file(ofstream *out, ofstream *vrmlfile, class Tesselation *Tess, class molecule *mol)
574{
575 atom *Walker = mol->start;
576 bond *Binder = mol->first;
577 int i;
578 Vector *center = mol->DetermineCenterOfAll(out);
579 if (vrmlfile != NULL) {
580 //cout << Verbose(1) << "Writing Raster3D file ... ";
581 *vrmlfile << "#VRML V2.0 utf8" << endl;
582 *vrmlfile << "#Created by molecuilder" << endl;
583 *vrmlfile << "#All atoms as spheres" << endl;
584 while (Walker->next != mol->end) {
585 Walker = Walker->next;
586 *vrmlfile << "Sphere {" << endl << " "; // 2 is sphere type
587 for (i=0;i<NDIM;i++)
588 *vrmlfile << Walker->x.x[i]-center->x[i] << " ";
589 *vrmlfile << "\t0.1\t1. 1. 1." << endl; // radius 0.05 and white as colour
590 }
591
592 *vrmlfile << "# All bonds as vertices" << endl;
593 while (Binder->next != mol->last) {
594 Binder = Binder->next;
595 *vrmlfile << "3" << endl << " "; // 2 is round-ended cylinder type
596 for (i=0;i<NDIM;i++)
597 *vrmlfile << Binder->leftatom->x.x[i]-center->x[i] << " ";
598 *vrmlfile << "\t0.03\t";
599 for (i=0;i<NDIM;i++)
600 *vrmlfile << Binder->rightatom->x.x[i]-center->x[i] << " ";
601 *vrmlfile << "\t0.03\t0. 0. 1." << endl; // radius 0.05 and blue as colour
602 }
603
604 *vrmlfile << "# All tesselation triangles" << endl;
605 for (TriangleMap::iterator TriangleRunner = Tess->TrianglesOnBoundary.begin(); TriangleRunner != Tess->TrianglesOnBoundary.end(); TriangleRunner++) {
606 *vrmlfile << "1" << endl << " "; // 1 is triangle type
607 for (i=0;i<3;i++) { // print each node
608 for (int j=0;j<NDIM;j++) // and for each node all NDIM coordinates
609 *vrmlfile << TriangleRunner->second->endpoints[i]->node->x.x[j]-center->x[j] << " ";
610 *vrmlfile << "\t";
611 }
612 *vrmlfile << "1. 0. 0." << endl; // red as colour
613 *vrmlfile << "18" << endl << " 0.5 0.5 0.5" << endl; // 18 is transparency type for previous object
614 }
615 } else {
616 cerr << "ERROR: Given vrmlfile is " << vrmlfile << "." << endl;
617 }
618 delete(center);
619};
620
621/** Creates the objects in a raster3d file (renderable with a header.r3d).
622 * \param *out output stream for debugging
623 * \param *rasterfile output stream for tecplot data
624 * \param *Tess Tesselation structure with constructed triangles
625 * \param *mol molecule structure with atom positions
626 */
627void write_raster3d_file(ofstream *out, ofstream *rasterfile, class Tesselation *Tess, class molecule *mol)
628{
629 atom *Walker = mol->start;
630 bond *Binder = mol->first;
631 int i;
632 Vector *center = mol->DetermineCenterOfAll(out);
633 if (rasterfile != NULL) {
634 //cout << Verbose(1) << "Writing Raster3D file ... ";
635 *rasterfile << "# Raster3D object description, created by MoleCuilder" << endl;
636 *rasterfile << "@header.r3d" << endl;
637 *rasterfile << "# All atoms as spheres" << endl;
638 while (Walker->next != mol->end) {
639 Walker = Walker->next;
640 *rasterfile << "2" << endl << " "; // 2 is sphere type
641 for (i=0;i<NDIM;i++)
642 *rasterfile << Walker->x.x[i]-center->x[i] << " ";
643 *rasterfile << "\t0.1\t1. 1. 1." << endl; // radius 0.05 and white as colour
644 }
645
646 *rasterfile << "# All bonds as vertices" << endl;
647 while (Binder->next != mol->last) {
648 Binder = Binder->next;
649 *rasterfile << "3" << endl << " "; // 2 is round-ended cylinder type
650 for (i=0;i<NDIM;i++)
651 *rasterfile << Binder->leftatom->x.x[i]-center->x[i] << " ";
652 *rasterfile << "\t0.03\t";
653 for (i=0;i<NDIM;i++)
654 *rasterfile << Binder->rightatom->x.x[i]-center->x[i] << " ";
655 *rasterfile << "\t0.03\t0. 0. 1." << endl; // radius 0.05 and blue as colour
656 }
657
658 *rasterfile << "# All tesselation triangles" << endl;
659 *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";
660 for (TriangleMap::iterator TriangleRunner = Tess->TrianglesOnBoundary.begin(); TriangleRunner != Tess->TrianglesOnBoundary.end(); TriangleRunner++) {
661 *rasterfile << "1" << endl << " "; // 1 is triangle type
662 for (i=0;i<3;i++) { // print each node
663 for (int j=0;j<NDIM;j++) // and for each node all NDIM coordinates
664 *rasterfile << TriangleRunner->second->endpoints[i]->node->x.x[j]-center->x[j] << " ";
665 *rasterfile << "\t";
666 }
667 *rasterfile << "1. 0. 0." << endl; // red as colour
668 //*rasterfile << "18" << endl << " 0.5 0.5 0.5" << endl; // 18 is transparency type for previous object
669 }
670 *rasterfile << "9\n terminating special property\n";
671 } else {
672 cerr << "ERROR: Given rasterfile is " << rasterfile << "." << endl;
673 }
674 delete(center);
675};
676
677/** This function creates the tecplot file, displaying the tesselation of the hull.
678 * \param *out output stream for debugging
679 * \param *tecplot output stream for tecplot data
680 * \param N arbitrary number to differentiate various zones in the tecplot format
681 */
682void write_tecplot_file(ofstream *out, ofstream *tecplot, class Tesselation *TesselStruct, class molecule *mol, int N)
683{
684 if (tecplot != NULL)
685 {
686 *tecplot << "TITLE = \"3D CONVEX SHELL\"" << endl;
687 *tecplot << "VARIABLES = \"X\" \"Y\" \"Z\"" << endl;
688 *tecplot << "ZONE T=\"TRIANGLES" << N << "\", N="
689 << TesselStruct->PointsOnBoundaryCount << ", E="
690 << TesselStruct->TrianglesOnBoundaryCount
691 << ", DATAPACKING=POINT, ZONETYPE=FETRIANGLE" << endl;
692 int *LookupList = new int[mol->AtomCount];
693 for (int i = 0; i < mol->AtomCount; i++)
694 LookupList[i] = -1;
695
696 // print atom coordinates
697 *out << Verbose(2) << "The following triangles were created:";
698 int Counter = 1;
699 atom *Walker = NULL;
700 for (PointMap::iterator target = TesselStruct->PointsOnBoundary.begin(); target
701 != TesselStruct->PointsOnBoundary.end(); target++)
702 {
703 Walker = target->second->node;
704 LookupList[Walker->nr] = Counter++;
705 *tecplot << Walker->x.x[0] << " " << Walker->x.x[1] << " "
706 << Walker->x.x[2] << " " << endl;
707 }
708 *tecplot << endl;
709 // print connectivity
710 for (TriangleMap::iterator runner =
711 TesselStruct->TrianglesOnBoundary.begin(); runner
712 != TesselStruct->TrianglesOnBoundary.end(); runner++)
713 {
714 *out << " " << runner->second->endpoints[0]->node->Name << "<->"
715 << runner->second->endpoints[1]->node->Name << "<->"
716 << runner->second->endpoints[2]->node->Name;
717 *tecplot << LookupList[runner->second->endpoints[0]->node->nr] << " "
718 << LookupList[runner->second->endpoints[1]->node->nr] << " "
719 << LookupList[runner->second->endpoints[2]->node->nr] << endl;
720 }
721 delete[] (LookupList);
722 *out << endl;
723 }
724}
725
726/** Tesselates the convex boundary by finding all boundary points.
727 * \param *out output stream for debugging
728 * \param *mol molecule structure with Atom's and Bond's
729 * \param *TesselStruct Tesselation filled with points, lines and triangles on boundary on return
730 * \param *LCList atoms in LinkedCell list
731 * \param *filename filename prefix for output of vertex data
732 * \return *TesselStruct is filled with convex boundary and tesselation is stored under \a *filename.
733 */
734void Find_convex_border(ofstream *out, molecule* mol, class Tesselation *&TesselStruct, class LinkedCell *LCList, const char *filename)
735{
736 bool BoundaryFreeFlag = false;
737 Boundaries *BoundaryPoints = NULL;
738
739 cout << Verbose(1) << "Begin of find_convex_border" << endl;
740
741 if (TesselStruct != NULL) // free if allocated
742 delete(TesselStruct);
743 TesselStruct = new class Tesselation;
744
745 // 1. Find all points on the boundary
746 if (BoundaryPoints == NULL) {
747 BoundaryFreeFlag = true;
748 BoundaryPoints = GetBoundaryPoints(out, mol);
749 } else {
750 *out << Verbose(1) << "Using given boundary points set." << endl;
751 }
752
753// printing all inserted for debugging
754 for (int axis=0; axis < NDIM; axis++)
755 {
756 *out << Verbose(2) << "Printing list of candidates for axis " << axis << " which we have inserted so far." << endl;
757 int i=0;
758 for(Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
759 if (runner != BoundaryPoints[axis].begin())
760 *out << ", " << i << ": " << *runner->second.second;
761 else
762 *out << i << ": " << *runner->second.second;
763 i++;
764 }
765 *out << endl;
766 }
767
768 // 2. fill the boundary point list
769 for (int axis = 0; axis < NDIM; axis++)
770 for (Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++)
771 TesselStruct->AddPoint(runner->second.second);
772
773 *out << Verbose(2) << "I found " << TesselStruct->PointsOnBoundaryCount << " points on the convex boundary." << endl;
774 // now we have the whole set of edge points in the BoundaryList
775
776 // listing for debugging
777 // *out << Verbose(1) << "Listing PointsOnBoundary:";
778 // for(PointMap::iterator runner = PointsOnBoundary.begin(); runner != PointsOnBoundary.end(); runner++) {
779 // *out << " " << *runner->second;
780 // }
781 // *out << endl;
782
783 // 3a. guess starting triangle
784 TesselStruct->GuessStartingTriangle(out);
785
786 // 3b. go through all lines, that are not yet part of two triangles (only of one so far)
787 TesselStruct->TesselateOnBoundary(out, mol);
788
789 *out << Verbose(2) << "I created " << TesselStruct->TrianglesOnBoundaryCount << " triangles with " << TesselStruct->LinesOnBoundaryCount << " lines and " << TesselStruct->PointsOnBoundaryCount << " points." << endl;
790
791 // 4. Store triangles in tecplot file
792 if (filename != NULL) {
793 string OutputName(filename);
794 OutputName.append(TecplotSuffix);
795 ofstream *tecplot = new ofstream(OutputName.c_str());
796 write_tecplot_file(out, tecplot, TesselStruct, mol, 0);
797 tecplot->close();
798 delete(tecplot);
799 ofstream *rasterplot = new ofstream(OutputName.c_str());
800 write_raster3d_file(out, rasterplot, TesselStruct, mol);
801 rasterplot->close();
802 delete(rasterplot);
803 }
804
805 // free reference lists
806 if (BoundaryFreeFlag)
807 delete[] (BoundaryPoints);
808
809 cout << Verbose(1) << "End of find_convex_border" << endl;
810};
811
812
813/** Determines the volume of a cluster.
814 * Determines first the convex envelope, then tesselates it and calculates its volume.
815 * \param *out output stream for debugging
816 * \param *TesselStruct Tesselation filled with points, lines and triangles on boundary on return
817 * \param *configuration needed for path to store convex envelope file
818 * \return determined volume of the cluster in cubed config:GetIsAngstroem()
819 */
820double VolumeOfConvexEnvelope(ofstream *out, class Tesselation *TesselStruct, class config *configuration)
821{
822 bool IsAngstroem = configuration->GetIsAngstroem();
823 double volume = 0.;
824 double PyramidVolume = 0.;
825 double G, h;
826 Vector x, y;
827 double a, b, c;
828
829 // 6a. Every triangle forms a pyramid with the center of gravity as its peak, sum up the volumes
830 *out << Verbose(1)
831 << "Calculating the volume of the pyramids formed out of triangles and center of gravity."
832 << endl;
833 for (TriangleMap::iterator runner = TesselStruct->TrianglesOnBoundary.begin(); runner
834 != TesselStruct->TrianglesOnBoundary.end(); runner++)
835 { // go through every triangle, calculate volume of its pyramid with CoG as peak
836 x.CopyVector(&runner->second->endpoints[0]->node->x);
837 x.SubtractVector(&runner->second->endpoints[1]->node->x);
838 y.CopyVector(&runner->second->endpoints[0]->node->x);
839 y.SubtractVector(&runner->second->endpoints[2]->node->x);
840 a = sqrt(runner->second->endpoints[0]->node->x.DistanceSquared(
841 &runner->second->endpoints[1]->node->x));
842 b = sqrt(runner->second->endpoints[0]->node->x.DistanceSquared(
843 &runner->second->endpoints[2]->node->x));
844 c = sqrt(runner->second->endpoints[2]->node->x.DistanceSquared(
845 &runner->second->endpoints[1]->node->x));
846 G = sqrt(((a + b + c) * (a + b + c) - 2 * (a * a + b * b + c * c)) / 16.); // area of tesselated triangle
847 x.MakeNormalVector(&runner->second->endpoints[0]->node->x,
848 &runner->second->endpoints[1]->node->x,
849 &runner->second->endpoints[2]->node->x);
850 x.Scale(runner->second->endpoints[1]->node->x.Projection(&x));
851 h = x.Norm(); // distance of CoG to triangle
852 PyramidVolume = (1. / 3.) * G * h; // this formula holds for _all_ pyramids (independent of n-edge base or (not) centered peak)
853 *out << Verbose(2) << "Area of triangle is " << G << " "
854 << (IsAngstroem ? "angstrom" : "atomiclength") << "^2, height is "
855 << h << " and the volume is " << PyramidVolume << " "
856 << (IsAngstroem ? "angstrom" : "atomiclength") << "^3." << endl;
857 volume += PyramidVolume;
858 }
859 *out << Verbose(0) << "RESULT: The summed volume is " << setprecision(10)
860 << volume << " " << (IsAngstroem ? "angstrom" : "atomiclength") << "^3."
861 << endl;
862
863 return volume;
864}
865;
866
867/** Creates multiples of the by \a *mol given cluster and suspends them in water with a given final density.
868 * We get cluster volume by VolumeOfConvexEnvelope() and its diameters by GetDiametersOfCluster()
869 * \param *out output stream for debugging
870 * \param *configuration needed for path to store convex envelope file
871 * \param *mol molecule structure representing the cluster
872 * \param ClusterVolume guesstimated cluster volume, if equal 0 we used VolumeOfConvexEnvelope() instead.
873 * \param celldensity desired average density in final cell
874 */
875void
876PrepareClustersinWater(ofstream *out, config *configuration, molecule *mol,
877 double ClusterVolume, double celldensity)
878{
879 // transform to PAS
880 mol->PrincipalAxisSystem(out, true);
881
882 // some preparations beforehand
883 bool IsAngstroem = configuration->GetIsAngstroem();
884 Boundaries *BoundaryPoints = GetBoundaryPoints(out, mol);
885 class Tesselation *TesselStruct = NULL;
886 LinkedCell LCList(mol, 10.);
887 Find_convex_border(out, mol, TesselStruct, &LCList, NULL);
888 double clustervolume;
889 if (ClusterVolume == 0)
890 clustervolume = VolumeOfConvexEnvelope(out, TesselStruct, configuration);
891 else
892 clustervolume = ClusterVolume;
893 double *GreatestDiameter = GetDiametersOfCluster(out, BoundaryPoints, mol, IsAngstroem);
894 Vector BoxLengths;
895 int repetition[NDIM] =
896 { 1, 1, 1 };
897 int TotalNoClusters = 1;
898 for (int i = 0; i < NDIM; i++)
899 TotalNoClusters *= repetition[i];
900
901 // sum up the atomic masses
902 double totalmass = 0.;
903 atom *Walker = mol->start;
904 while (Walker->next != mol->end)
905 {
906 Walker = Walker->next;
907 totalmass += Walker->type->mass;
908 }
909 *out << Verbose(0) << "RESULT: The summed mass is " << setprecision(10)
910 << totalmass << " atomicmassunit." << endl;
911
912 *out << Verbose(0) << "RESULT: The average density is " << setprecision(10)
913 << totalmass / clustervolume << " atomicmassunit/"
914 << (IsAngstroem ? "angstrom" : "atomiclength") << "^3." << endl;
915
916 // solve cubic polynomial
917 *out << Verbose(1) << "Solving equidistant suspension in water problem ..."
918 << endl;
919 double cellvolume;
920 if (IsAngstroem)
921 cellvolume = (TotalNoClusters * totalmass / SOLVENTDENSITY_A - (totalmass
922 / clustervolume)) / (celldensity - 1);
923 else
924 cellvolume = (TotalNoClusters * totalmass / SOLVENTDENSITY_a0 - (totalmass
925 / clustervolume)) / (celldensity - 1);
926 *out << Verbose(1) << "Cellvolume needed for a density of " << celldensity
927 << " g/cm^3 is " << cellvolume << " " << (IsAngstroem ? "angstrom"
928 : "atomiclength") << "^3." << endl;
929
930 double minimumvolume = TotalNoClusters * (GreatestDiameter[0]
931 * GreatestDiameter[1] * GreatestDiameter[2]);
932 *out << Verbose(1)
933 << "Minimum volume of the convex envelope contained in a rectangular box is "
934 << minimumvolume << " atomicmassunit/" << (IsAngstroem ? "angstrom"
935 : "atomiclength") << "^3." << endl;
936 if (minimumvolume > cellvolume)
937 {
938 cerr << Verbose(0)
939 << "ERROR: the containing box already has a greater volume than the envisaged cell volume!"
940 << endl;
941 cout << Verbose(0)
942 << "Setting Box dimensions to minimum possible, the greatest diameters."
943 << endl;
944 for (int i = 0; i < NDIM; i++)
945 BoxLengths.x[i] = GreatestDiameter[i];
946 mol->CenterEdge(out, &BoxLengths);
947 }
948 else
949 {
950 BoxLengths.x[0] = (repetition[0] * GreatestDiameter[0] + repetition[1]
951 * GreatestDiameter[1] + repetition[2] * GreatestDiameter[2]);
952 BoxLengths.x[1] = (repetition[0] * repetition[1] * GreatestDiameter[0]
953 * GreatestDiameter[1] + repetition[0] * repetition[2]
954 * GreatestDiameter[0] * GreatestDiameter[2] + repetition[1]
955 * repetition[2] * GreatestDiameter[1] * GreatestDiameter[2]);
956 BoxLengths.x[2] = minimumvolume - cellvolume;
957 double x0 = 0., x1 = 0., x2 = 0.;
958 if (gsl_poly_solve_cubic(BoxLengths.x[0], BoxLengths.x[1],
959 BoxLengths.x[2], &x0, &x1, &x2) == 1) // either 1 or 3 on return
960 *out << Verbose(0) << "RESULT: The resulting spacing is: " << x0
961 << " ." << endl;
962 else
963 {
964 *out << Verbose(0) << "RESULT: The resulting spacings are: " << x0
965 << " and " << x1 << " and " << x2 << " ." << endl;
966 x0 = x2; // sorted in ascending order
967 }
968
969 cellvolume = 1;
970 for (int i = 0; i < NDIM; i++)
971 {
972 BoxLengths.x[i] = repetition[i] * (x0 + GreatestDiameter[i]);
973 cellvolume *= BoxLengths.x[i];
974 }
975
976 // set new box dimensions
977 *out << Verbose(0) << "Translating to box with these boundaries." << endl;
978 mol->SetBoxDimension(&BoxLengths);
979 mol->CenterInBox((ofstream *) &cout);
980 }
981 // update Box of atoms by boundary
982 mol->SetBoxDimension(&BoxLengths);
983 *out << Verbose(0) << "RESULT: The resulting cell dimensions are: "
984 << BoxLengths.x[0] << " and " << BoxLengths.x[1] << " and "
985 << BoxLengths.x[2] << " with total volume of " << cellvolume << " "
986 << (IsAngstroem ? "angstrom" : "atomiclength") << "^3." << endl;
987}
988;
989
990// =========================================================== class TESSELATION ===========================================
991
992/** Constructor of class Tesselation.
993 */
994Tesselation::Tesselation()
995{
996 PointsOnBoundaryCount = 0;
997 LinesOnBoundaryCount = 0;
998 TrianglesOnBoundaryCount = 0;
999 TriangleFilesWritten = 0;
1000}
1001;
1002
1003/** Constructor of class Tesselation.
1004 * We have to free all points, lines and triangles.
1005 */
1006Tesselation::~Tesselation()
1007{
1008 cout << Verbose(1) << "Free'ing TesselStruct ... " << endl;
1009 for (TriangleMap::iterator runner = TrianglesOnBoundary.begin(); runner != TrianglesOnBoundary.end(); runner++) {
1010 if (runner->second != NULL) {
1011 delete (runner->second);
1012 runner->second = NULL;
1013 } else
1014 cerr << "ERROR: The triangle " << runner->first << " has already been free'd." << endl;
1015 }
1016}
1017;
1018
1019/** Gueses first starting triangle of the convex envelope.
1020 * We guess the starting triangle by taking the smallest distance between two points and looking for a fitting third.
1021 * \param *out output stream for debugging
1022 * \param PointsOnBoundary set of boundary points defining the convex envelope of the cluster
1023 */
1024void
1025Tesselation::GuessStartingTriangle(ofstream *out)
1026{
1027 // 4b. create a starting triangle
1028 // 4b1. create all distances
1029 DistanceMultiMap DistanceMMap;
1030 double distance, tmp;
1031 Vector PlaneVector, TrialVector;
1032 PointMap::iterator A, B, C; // three nodes of the first triangle
1033 A = PointsOnBoundary.begin(); // the first may be chosen arbitrarily
1034
1035 // with A chosen, take each pair B,C and sort
1036 if (A != PointsOnBoundary.end())
1037 {
1038 B = A;
1039 B++;
1040 for (; B != PointsOnBoundary.end(); B++)
1041 {
1042 C = B;
1043 C++;
1044 for (; C != PointsOnBoundary.end(); C++)
1045 {
1046 tmp = A->second->node->x.DistanceSquared(&B->second->node->x);
1047 distance = tmp * tmp;
1048 tmp = A->second->node->x.DistanceSquared(&C->second->node->x);
1049 distance += tmp * tmp;
1050 tmp = B->second->node->x.DistanceSquared(&C->second->node->x);
1051 distance += tmp * tmp;
1052 DistanceMMap.insert(DistanceMultiMapPair(distance, pair<
1053 PointMap::iterator, PointMap::iterator> (B, C)));
1054 }
1055 }
1056 }
1057 // // listing distances
1058 // *out << Verbose(1) << "Listing DistanceMMap:";
1059 // for(DistanceMultiMap::iterator runner = DistanceMMap.begin(); runner != DistanceMMap.end(); runner++) {
1060 // *out << " " << runner->first << "(" << *runner->second.first->second << ", " << *runner->second.second->second << ")";
1061 // }
1062 // *out << endl;
1063 // 4b2. pick three baselines forming a triangle
1064 // 1. we take from the smallest sum of squared distance as the base line BC (with peak A) onward as the triangle candidate
1065 DistanceMultiMap::iterator baseline = DistanceMMap.begin();
1066 for (; baseline != DistanceMMap.end(); baseline++)
1067 {
1068 // we take from the smallest sum of squared distance as the base line BC (with peak A) onward as the triangle candidate
1069 // 2. next, we have to check whether all points reside on only one side of the triangle
1070 // 3. construct plane vector
1071 PlaneVector.MakeNormalVector(&A->second->node->x,
1072 &baseline->second.first->second->node->x,
1073 &baseline->second.second->second->node->x);
1074 *out << Verbose(2) << "Plane vector of candidate triangle is ";
1075 PlaneVector.Output(out);
1076 *out << endl;
1077 // 4. loop over all points
1078 double sign = 0.;
1079 PointMap::iterator checker = PointsOnBoundary.begin();
1080 for (; checker != PointsOnBoundary.end(); checker++)
1081 {
1082 // (neglecting A,B,C)
1083 if ((checker == A) || (checker == baseline->second.first) || (checker
1084 == baseline->second.second))
1085 continue;
1086 // 4a. project onto plane vector
1087 TrialVector.CopyVector(&checker->second->node->x);
1088 TrialVector.SubtractVector(&A->second->node->x);
1089 distance = TrialVector.Projection(&PlaneVector);
1090 if (fabs(distance) < 1e-4) // we need to have a small epsilon around 0 which is still ok
1091 continue;
1092 *out << Verbose(3) << "Projection of " << checker->second->node->Name
1093 << " yields distance of " << distance << "." << endl;
1094 tmp = distance / fabs(distance);
1095 // 4b. Any have different sign to than before? (i.e. would lie outside convex hull with this starting triangle)
1096 if ((sign != 0) && (tmp != sign))
1097 {
1098 // 4c. If so, break 4. loop and continue with next candidate in 1. loop
1099 *out << Verbose(2) << "Current candidates: "
1100 << A->second->node->Name << ","
1101 << baseline->second.first->second->node->Name << ","
1102 << baseline->second.second->second->node->Name << " leave "
1103 << checker->second->node->Name << " outside the convex hull."
1104 << endl;
1105 break;
1106 }
1107 else
1108 { // note the sign for later
1109 *out << Verbose(2) << "Current candidates: "
1110 << A->second->node->Name << ","
1111 << baseline->second.first->second->node->Name << ","
1112 << baseline->second.second->second->node->Name << " leave "
1113 << checker->second->node->Name << " inside the convex hull."
1114 << endl;
1115 sign = tmp;
1116 }
1117 // 4d. Check whether the point is inside the triangle (check distance to each node
1118 tmp = checker->second->node->x.DistanceSquared(&A->second->node->x);
1119 int innerpoint = 0;
1120 if ((tmp < A->second->node->x.DistanceSquared(
1121 &baseline->second.first->second->node->x)) && (tmp
1122 < A->second->node->x.DistanceSquared(
1123 &baseline->second.second->second->node->x)))
1124 innerpoint++;
1125 tmp = checker->second->node->x.DistanceSquared(
1126 &baseline->second.first->second->node->x);
1127 if ((tmp < baseline->second.first->second->node->x.DistanceSquared(
1128 &A->second->node->x)) && (tmp
1129 < baseline->second.first->second->node->x.DistanceSquared(
1130 &baseline->second.second->second->node->x)))
1131 innerpoint++;
1132 tmp = checker->second->node->x.DistanceSquared(
1133 &baseline->second.second->second->node->x);
1134 if ((tmp < baseline->second.second->second->node->x.DistanceSquared(
1135 &baseline->second.first->second->node->x)) && (tmp
1136 < baseline->second.second->second->node->x.DistanceSquared(
1137 &A->second->node->x)))
1138 innerpoint++;
1139 // 4e. If so, break 4. loop and continue with next candidate in 1. loop
1140 if (innerpoint == 3)
1141 break;
1142 }
1143 // 5. come this far, all on same side? Then break 1. loop and construct triangle
1144 if (checker == PointsOnBoundary.end())
1145 {
1146 *out << "Looks like we have a candidate!" << endl;
1147 break;
1148 }
1149 }
1150 if (baseline != DistanceMMap.end())
1151 {
1152 BPS[0] = baseline->second.first->second;
1153 BPS[1] = baseline->second.second->second;
1154 BLS[0] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
1155 BPS[0] = A->second;
1156 BPS[1] = baseline->second.second->second;
1157 BLS[1] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
1158 BPS[0] = baseline->second.first->second;
1159 BPS[1] = A->second;
1160 BLS[2] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount);
1161
1162 // 4b3. insert created triangle
1163 BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
1164 TrianglesOnBoundary.insert(TrianglePair(TrianglesOnBoundaryCount, BTS));
1165 TrianglesOnBoundaryCount++;
1166 for (int i = 0; i < NDIM; i++)
1167 {
1168 LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount, BTS->lines[i]));
1169 LinesOnBoundaryCount++;
1170 }
1171
1172 *out << Verbose(1) << "Starting triangle is " << *BTS << "." << endl;
1173 }
1174 else
1175 {
1176 *out << Verbose(1) << "No starting triangle found." << endl;
1177 exit(255);
1178 }
1179}
1180;
1181
1182/** Tesselates the convex envelope of a cluster from a single starting triangle.
1183 * The starting triangle is made out of three baselines. Each line in the final tesselated cluster may belong to at most
1184 * 2 triangles. Hence, we go through all current lines:
1185 * -# if the lines contains to only one triangle
1186 * -# We search all points in the boundary
1187 * -# if the triangle with the baseline and the current point has the smallest of angles (comparison between normal vectors
1188 * -# if the triangle is in forward direction of the baseline (at most 90 degrees angle between vector orthogonal to
1189 * baseline in triangle plane pointing out of the triangle and normal vector of new triangle)
1190 * -# then we have a new triangle, whose baselines we again add (or increase their TriangleCount)
1191 * \param *out output stream for debugging
1192 * \param *configuration for IsAngstroem
1193 * \param *mol the cluster as a molecule structure
1194 */
1195void Tesselation::TesselateOnBoundary(ofstream *out, molecule *mol)
1196{
1197 bool flag;
1198 PointMap::iterator winner;
1199 class BoundaryPointSet *peak = NULL;
1200 double SmallestAngle, TempAngle;
1201 Vector NormalVector, VirtualNormalVector, CenterVector, TempVector, PropagationVector, *MolCenter = NULL;
1202 LineMap::iterator LineChecker[2];
1203
1204 MolCenter = mol->DetermineCenterOfAll(out);
1205 do {
1206 flag = false;
1207 for (LineMap::iterator baseline = LinesOnBoundary.begin(); baseline != LinesOnBoundary.end(); baseline++)
1208 if (baseline->second->TrianglesCount == 1) {
1209 *out << Verbose(2) << "Current baseline is between " << *(baseline->second) << "." << endl;
1210 // 5a. go through each boundary point if not _both_ edges between either endpoint of the current line and this point exist (and belong to 2 triangles)
1211 SmallestAngle = M_PI;
1212 BTS = baseline->second->triangles.begin()->second; // there is only one triangle so far
1213 // get peak point with respect to this base line's only triangle
1214 for (int i = 0; i < 3; i++)
1215 if ((BTS->endpoints[i] != baseline->second->endpoints[0]) && (BTS->endpoints[i] != baseline->second->endpoints[1]))
1216 peak = BTS->endpoints[i];
1217 *out << Verbose(3) << " and has peak " << *peak << "." << endl;
1218 // offset to center of triangle
1219 CenterVector.Zero();
1220 for (int i = 0; i < 3; i++)
1221 CenterVector.AddVector(&BTS->endpoints[i]->node->x);
1222 CenterVector.Scale(1. / 3.);
1223 *out << Verbose(4) << "CenterVector of base triangle is " << CenterVector << endl;
1224 NormalVector.CopyVector(MolCenter);
1225 NormalVector.SubtractVector(&CenterVector);
1226 // normal vector of triangle
1227 BTS->GetNormalVector(NormalVector);
1228 *out << Verbose(4) << "NormalVector of base triangle is " << NormalVector << endl;
1229 // vector in propagation direction (out of triangle)
1230 // project center vector onto triangle plane (points from intersection plane-NormalVector to plane-CenterVector intersection)
1231 TempVector.CopyVector(&baseline->second->endpoints[0]->node->x);
1232 TempVector.SubtractVector(&baseline->second->endpoints[1]->node->x);
1233 PropagationVector.MakeNormalVector(&TempVector, &NormalVector);
1234 TempVector.CopyVector(&CenterVector);
1235 TempVector.SubtractVector(&baseline->second->endpoints[0]->node->x); // TempVector is vector on triangle plane pointing from one baseline egde towards center!
1236 //*out << Verbose(2) << "Projection of propagation onto temp: " << PropagationVector.Projection(&TempVector) << "." << endl;
1237 if (PropagationVector.Projection(&TempVector) > 0) // make sure normal propagation vector points outward from baseline
1238 PropagationVector.Scale(-1.);
1239 *out << Verbose(4) << "PropagationVector of base triangle is ";
1240 PropagationVector.Output(out);
1241 *out << endl;
1242 winner = PointsOnBoundary.end();
1243 for (PointMap::iterator target = PointsOnBoundary.begin(); target
1244 != PointsOnBoundary.end(); target++)
1245 if ((target->second != baseline->second->endpoints[0])
1246 && (target->second != baseline->second->endpoints[1]))
1247 { // don't take the same endpoints
1248 *out << Verbose(3) << "Target point is " << *(target->second)
1249 << ":";
1250 bool continueflag = true;
1251
1252 VirtualNormalVector.CopyVector(
1253 &baseline->second->endpoints[0]->node->x);
1254 VirtualNormalVector.AddVector(
1255 &baseline->second->endpoints[0]->node->x);
1256 VirtualNormalVector.Scale(-1. / 2.); // points now to center of base line
1257 VirtualNormalVector.AddVector(&target->second->node->x); // points from center of base line to target
1258 TempAngle = VirtualNormalVector.Angle(&PropagationVector);
1259 continueflag = continueflag && (TempAngle < (M_PI/2.)); // no bends bigger than Pi/2 (90 degrees)
1260 if (!continueflag)
1261 {
1262 *out << Verbose(4)
1263 << "Angle between propagation direction and base line to "
1264 << *(target->second) << " is " << TempAngle
1265 << ", bad direction!" << endl;
1266 continue;
1267 }
1268 else
1269 *out << Verbose(4)
1270 << "Angle between propagation direction and base line to "
1271 << *(target->second) << " is " << TempAngle
1272 << ", good direction!" << endl;
1273 LineChecker[0] = baseline->second->endpoints[0]->lines.find(
1274 target->first);
1275 LineChecker[1] = baseline->second->endpoints[1]->lines.find(
1276 target->first);
1277 // if (LineChecker[0] != baseline->second->endpoints[0]->lines.end())
1278 // *out << Verbose(4) << *(baseline->second->endpoints[0]) << " has line " << *(LineChecker[0]->second) << " to " << *(target->second) << " as endpoint with " << LineChecker[0]->second->TrianglesCount << " triangles." << endl;
1279 // else
1280 // *out << Verbose(4) << *(baseline->second->endpoints[0]) << " has no line to " << *(target->second) << " as endpoint." << endl;
1281 // if (LineChecker[1] != baseline->second->endpoints[1]->lines.end())
1282 // *out << Verbose(4) << *(baseline->second->endpoints[1]) << " has line " << *(LineChecker[1]->second) << " to " << *(target->second) << " as endpoint with " << LineChecker[1]->second->TrianglesCount << " triangles." << endl;
1283 // else
1284 // *out << Verbose(4) << *(baseline->second->endpoints[1]) << " has no line to " << *(target->second) << " as endpoint." << endl;
1285 // check first endpoint (if any connecting line goes to target or at least not more than 1)
1286 continueflag = continueflag && (((LineChecker[0]
1287 == baseline->second->endpoints[0]->lines.end())
1288 || (LineChecker[0]->second->TrianglesCount == 1)));
1289 if (!continueflag)
1290 {
1291 *out << Verbose(4) << *(baseline->second->endpoints[0])
1292 << " has line " << *(LineChecker[0]->second)
1293 << " to " << *(target->second)
1294 << " as endpoint with "
1295 << LineChecker[0]->second->TrianglesCount
1296 << " triangles." << endl;
1297 continue;
1298 }
1299 // check second endpoint (if any connecting line goes to target or at least not more than 1)
1300 continueflag = continueflag && (((LineChecker[1]
1301 == baseline->second->endpoints[1]->lines.end())
1302 || (LineChecker[1]->second->TrianglesCount == 1)));
1303 if (!continueflag)
1304 {
1305 *out << Verbose(4) << *(baseline->second->endpoints[1])
1306 << " has line " << *(LineChecker[1]->second)
1307 << " to " << *(target->second)
1308 << " as endpoint with "
1309 << LineChecker[1]->second->TrianglesCount
1310 << " triangles." << endl;
1311 continue;
1312 }
1313 // check whether the envisaged triangle does not already exist (if both lines exist and have same endpoint)
1314 continueflag = continueflag && (!(((LineChecker[0]
1315 != baseline->second->endpoints[0]->lines.end())
1316 && (LineChecker[1]
1317 != baseline->second->endpoints[1]->lines.end())
1318 && (GetCommonEndpoint(LineChecker[0]->second,
1319 LineChecker[1]->second) == peak))));
1320 if (!continueflag)
1321 {
1322 *out << Verbose(4) << "Current target is peak!" << endl;
1323 continue;
1324 }
1325 // in case NOT both were found
1326 if (continueflag)
1327 { // create virtually this triangle, get its normal vector, calculate angle
1328 flag = true;
1329 VirtualNormalVector.MakeNormalVector(
1330 &baseline->second->endpoints[0]->node->x,
1331 &baseline->second->endpoints[1]->node->x,
1332 &target->second->node->x);
1333 // make it always point inward
1334 if (baseline->second->endpoints[0]->node->x.Projection(
1335 &VirtualNormalVector) > 0)
1336 VirtualNormalVector.Scale(-1.);
1337 // calculate angle
1338 TempAngle = NormalVector.Angle(&VirtualNormalVector);
1339 *out << Verbose(4) << "NormalVector is ";
1340 VirtualNormalVector.Output(out);
1341 *out << " and the angle is " << TempAngle << "." << endl;
1342 if (SmallestAngle > TempAngle)
1343 { // set to new possible winner
1344 SmallestAngle = TempAngle;
1345 winner = target;
1346 }
1347 }
1348 }
1349 // 5b. The point of the above whose triangle has the greatest angle with the triangle the current line belongs to (it only belongs to one, remember!): New triangle
1350 if (winner != PointsOnBoundary.end())
1351 {
1352 *out << Verbose(2) << "Winning target point is "
1353 << *(winner->second) << " with angle " << SmallestAngle
1354 << "." << endl;
1355 // create the lins of not yet present
1356 BLS[0] = baseline->second;
1357 // 5c. add lines to the line set if those were new (not yet part of a triangle), delete lines that belong to two triangles)
1358 LineChecker[0] = baseline->second->endpoints[0]->lines.find(
1359 winner->first);
1360 LineChecker[1] = baseline->second->endpoints[1]->lines.find(
1361 winner->first);
1362 if (LineChecker[0]
1363 == baseline->second->endpoints[0]->lines.end())
1364 { // create
1365 BPS[0] = baseline->second->endpoints[0];
1366 BPS[1] = winner->second;
1367 BLS[1] = new class BoundaryLineSet(BPS,
1368 LinesOnBoundaryCount);
1369 LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount,
1370 BLS[1]));
1371 LinesOnBoundaryCount++;
1372 }
1373 else
1374 BLS[1] = LineChecker[0]->second;
1375 if (LineChecker[1]
1376 == baseline->second->endpoints[1]->lines.end())
1377 { // create
1378 BPS[0] = baseline->second->endpoints[1];
1379 BPS[1] = winner->second;
1380 BLS[2] = new class BoundaryLineSet(BPS,
1381 LinesOnBoundaryCount);
1382 LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount,
1383 BLS[2]));
1384 LinesOnBoundaryCount++;
1385 }
1386 else
1387 BLS[2] = LineChecker[1]->second;
1388 BTS = new class BoundaryTriangleSet(BLS,
1389 TrianglesOnBoundaryCount);
1390 TrianglesOnBoundary.insert(TrianglePair(
1391 TrianglesOnBoundaryCount, BTS));
1392 TrianglesOnBoundaryCount++;
1393 }
1394 else
1395 {
1396 *out << Verbose(1)
1397 << "I could not determine a winner for this baseline "
1398 << *(baseline->second) << "." << endl;
1399 }
1400
1401 // 5d. If the set of lines is not yet empty, go to 5. and continue
1402 }
1403 else
1404 *out << Verbose(2) << "Baseline candidate " << *(baseline->second)
1405 << " has a triangle count of "
1406 << baseline->second->TrianglesCount << "." << endl;
1407 } while (flag);
1408 delete(MolCenter);
1409};
1410
1411/** Adds an atom to the tesselation::PointsOnBoundary list.
1412 * \param *Walker atom to add
1413 */
1414void
1415Tesselation::AddPoint(atom *Walker)
1416{
1417 PointTestPair InsertUnique;
1418 BPS[0] = new class BoundaryPointSet(Walker);
1419 InsertUnique = PointsOnBoundary.insert(PointPair(Walker->nr, BPS[0]));
1420 if (InsertUnique.second) // if new point was not present before, increase counter
1421 PointsOnBoundaryCount++;
1422}
1423;
1424
1425/** Adds point to Tesselation::PointsOnBoundary if not yet present.
1426 * Tesselation::TPS is set to either this new BoundaryPointSet or to the existing one of not unique.
1427 * @param Candidate point to add
1428 * @param n index for this point in Tesselation::TPS array
1429 */
1430void
1431Tesselation::AddTrianglePoint(atom* Candidate, int n)
1432{
1433 PointTestPair InsertUnique;
1434 TPS[n] = new class BoundaryPointSet(Candidate);
1435 InsertUnique = PointsOnBoundary.insert(PointPair(Candidate->nr, TPS[n]));
1436 if (InsertUnique.second) { // if new point was not present before, increase counter
1437 PointsOnBoundaryCount++;
1438 } else {
1439 delete TPS[n];
1440 cout << Verbose(3) << "Atom " << *((InsertUnique.first)->second->node) << " is already present in PointsOnBoundary." << endl;
1441 TPS[n] = (InsertUnique.first)->second;
1442 }
1443}
1444;
1445
1446/** Function tries to add line from current Points in BPS to BoundaryLineSet.
1447 * If successful it raises the line count and inserts the new line into the BLS,
1448 * if unsuccessful, it writes the line which had been present into the BLS, deleting the new constructed one.
1449 * @param *a first endpoint
1450 * @param *b second endpoint
1451 * @param n index of Tesselation::BLS giving the line with both endpoints
1452 */
1453void Tesselation::AddTriangleLine(class BoundaryPointSet *a, class BoundaryPointSet *b, int n) {
1454 bool insertNewLine = true;
1455
1456 if (a->lines.find(b->node->nr) != a->lines.end()) {
1457 LineMap::iterator FindLine;
1458 pair<LineMap::iterator,LineMap::iterator> FindPair;
1459 FindPair = a->lines.equal_range(b->node->nr);
1460
1461 for (FindLine = FindPair.first; FindLine != FindPair.second; ++FindLine) {
1462 // If there is a line with less than two attached triangles, we don't need a new line.
1463 if (FindLine->second->TrianglesCount < 2) {
1464 insertNewLine = false;
1465 cout << Verbose(3) << "Using existing line " << *FindLine->second << endl;
1466
1467 BPS[0] = FindLine->second->endpoints[0];
1468 BPS[1] = FindLine->second->endpoints[1];
1469 BLS[n] = FindLine->second;
1470
1471 break;
1472 }
1473 }
1474 }
1475
1476 if (insertNewLine) {
1477 AlwaysAddTriangleLine(a, b, n);
1478 }
1479}
1480;
1481
1482/**
1483 * Adds lines from each of the current points in the BPS to BoundaryLineSet.
1484 * Raises the line count and inserts the new line into the BLS.
1485 *
1486 * @param *a first endpoint
1487 * @param *b second endpoint
1488 * @param n index of Tesselation::BLS giving the line with both endpoints
1489 */
1490void Tesselation::AlwaysAddTriangleLine(class BoundaryPointSet *a, class BoundaryPointSet *b, int n)
1491{
1492 cout << Verbose(3) << "Adding line between " << *(a->node) << " and " << *(b->node) << "." << endl;
1493 BPS[0] = a;
1494 BPS[1] = b;
1495 BLS[n] = new class BoundaryLineSet(BPS, LinesOnBoundaryCount); // this also adds the line to the local maps
1496 // add line to global map
1497 LinesOnBoundary.insert(LinePair(LinesOnBoundaryCount, BLS[n]));
1498 // increase counter
1499 LinesOnBoundaryCount++;
1500};
1501
1502/** Function tries to add Triangle just created to Triangle and remarks if already existent (Failure of algorithm).
1503 * Furthermore it adds the triangle to all of its lines, in order to recognize those which are saturated later.
1504 */
1505void
1506Tesselation::AddTriangle()
1507{
1508 cout << Verbose(1) << "Adding triangle to global TrianglesOnBoundary map." << endl;
1509
1510 // add triangle to global map
1511 TrianglesOnBoundary.insert(TrianglePair(TrianglesOnBoundaryCount, BTS));
1512 TrianglesOnBoundaryCount++;
1513
1514 // NOTE: add triangle to local maps is done in constructor of BoundaryTriangleSet
1515}
1516;
1517
1518
1519double det_get(gsl_matrix *A, int inPlace) {
1520 /*
1521 inPlace = 1 => A is replaced with the LU decomposed copy.
1522 inPlace = 0 => A is retained, and a copy is used for LU.
1523 */
1524
1525 double det;
1526 int signum;
1527 gsl_permutation *p = gsl_permutation_alloc(A->size1);
1528 gsl_matrix *tmpA;
1529
1530 if (inPlace)
1531 tmpA = A;
1532 else {
1533 gsl_matrix *tmpA = gsl_matrix_alloc(A->size1, A->size2);
1534 gsl_matrix_memcpy(tmpA , A);
1535 }
1536
1537
1538 gsl_linalg_LU_decomp(tmpA , p , &signum);
1539 det = gsl_linalg_LU_det(tmpA , signum);
1540 gsl_permutation_free(p);
1541 if (! inPlace)
1542 gsl_matrix_free(tmpA);
1543
1544 return det;
1545};
1546
1547void get_sphere(Vector *center, Vector &a, Vector &b, Vector &c, double RADIUS)
1548{
1549 gsl_matrix *A = gsl_matrix_calloc(3,3);
1550 double m11, m12, m13, m14;
1551
1552 for(int i=0;i<3;i++) {
1553 gsl_matrix_set(A, i, 0, a.x[i]);
1554 gsl_matrix_set(A, i, 1, b.x[i]);
1555 gsl_matrix_set(A, i, 2, c.x[i]);
1556 }
1557 m11 = det_get(A, 1);
1558
1559 for(int i=0;i<3;i++) {
1560 gsl_matrix_set(A, i, 0, a.x[i]*a.x[i] + b.x[i]*b.x[i] + c.x[i]*c.x[i]);
1561 gsl_matrix_set(A, i, 1, b.x[i]);
1562 gsl_matrix_set(A, i, 2, c.x[i]);
1563 }
1564 m12 = det_get(A, 1);
1565
1566 for(int i=0;i<3;i++) {
1567 gsl_matrix_set(A, i, 0, a.x[i]*a.x[i] + b.x[i]*b.x[i] + c.x[i]*c.x[i]);
1568 gsl_matrix_set(A, i, 1, a.x[i]);
1569 gsl_matrix_set(A, i, 2, c.x[i]);
1570 }
1571 m13 = det_get(A, 1);
1572
1573 for(int i=0;i<3;i++) {
1574 gsl_matrix_set(A, i, 0, a.x[i]*a.x[i] + b.x[i]*b.x[i] + c.x[i]*c.x[i]);
1575 gsl_matrix_set(A, i, 1, a.x[i]);
1576 gsl_matrix_set(A, i, 2, b.x[i]);
1577 }
1578 m14 = det_get(A, 1);
1579
1580 if (fabs(m11) < MYEPSILON)
1581 cerr << "ERROR: three points are colinear." << endl;
1582
1583 center->x[0] = 0.5 * m12/ m11;
1584 center->x[1] = -0.5 * m13/ m11;
1585 center->x[2] = 0.5 * m14/ m11;
1586
1587 if (fabs(a.Distance(center) - RADIUS) > MYEPSILON)
1588 cerr << "ERROR: The given center is further way by " << fabs(a.Distance(center) - RADIUS) << " from a than RADIUS." << endl;
1589
1590 gsl_matrix_free(A);
1591};
1592
1593
1594
1595/**
1596 * Function returns center of sphere with RADIUS, which rests on points a, b, c
1597 * @param Center this vector will be used for return
1598 * @param a vector first point of triangle
1599 * @param b vector second point of triangle
1600 * @param c vector third point of triangle
1601 * @param *Umkreismittelpunkt new cneter point of circumference
1602 * @param Direction vector indicates up/down
1603 * @param AlternativeDirection vecotr, needed in case the triangles have 90 deg angle
1604 * @param Halfplaneindicator double indicates whether Direction is up or down
1605 * @param AlternativeIndicator doube indicates in case of orthogonal triangles which direction of AlternativeDirection is suitable
1606 * @param alpha double angle at a
1607 * @param beta double, angle at b
1608 * @param gamma, double, angle at c
1609 * @param Radius, double
1610 * @param Umkreisradius double radius of circumscribing circle
1611 */
1612void Get_center_of_sphere(Vector* Center, Vector a, Vector b, Vector c, Vector *NewUmkreismittelpunkt, Vector* Direction, Vector* AlternativeDirection,
1613 double HalfplaneIndicator, double AlternativeIndicator, double alpha, double beta, double gamma, double RADIUS, double Umkreisradius)
1614{
1615 Vector TempNormal, helper;
1616 double Restradius;
1617 Vector OtherCenter;
1618 cout << Verbose(3) << "Begin of Get_center_of_sphere.\n";
1619 Center->Zero();
1620 helper.CopyVector(&a);
1621 helper.Scale(sin(2.*alpha));
1622 Center->AddVector(&helper);
1623 helper.CopyVector(&b);
1624 helper.Scale(sin(2.*beta));
1625 Center->AddVector(&helper);
1626 helper.CopyVector(&c);
1627 helper.Scale(sin(2.*gamma));
1628 Center->AddVector(&helper);
1629 //*Center = a * sin(2.*alpha) + b * sin(2.*beta) + c * sin(2.*gamma) ;
1630 Center->Scale(1./(sin(2.*alpha) + sin(2.*beta) + sin(2.*gamma)));
1631 NewUmkreismittelpunkt->CopyVector(Center);
1632 cout << Verbose(4) << "Center of new circumference is " << *NewUmkreismittelpunkt << ".\n";
1633 // Here we calculated center of circumscribing circle, using barycentric coordinates
1634 cout << Verbose(4) << "Center of circumference is " << *Center << " in direction " << *Direction << ".\n";
1635
1636 TempNormal.CopyVector(&a);
1637 TempNormal.SubtractVector(&b);
1638 helper.CopyVector(&a);
1639 helper.SubtractVector(&c);
1640 TempNormal.VectorProduct(&helper);
1641 if (fabs(HalfplaneIndicator) < MYEPSILON)
1642 {
1643 if ((TempNormal.ScalarProduct(AlternativeDirection) <0 and AlternativeIndicator >0) or (TempNormal.ScalarProduct(AlternativeDirection) >0 and AlternativeIndicator <0))
1644 {
1645 TempNormal.Scale(-1);
1646 }
1647 }
1648 else
1649 {
1650 if (TempNormal.ScalarProduct(Direction)<0 && HalfplaneIndicator >0 || TempNormal.ScalarProduct(Direction)>0 && HalfplaneIndicator<0)
1651 {
1652 TempNormal.Scale(-1);
1653 }
1654 }
1655
1656 TempNormal.Normalize();
1657 Restradius = sqrt(RADIUS*RADIUS - Umkreisradius*Umkreisradius);
1658 cout << Verbose(4) << "Height of center of circumference to center of sphere is " << Restradius << ".\n";
1659 TempNormal.Scale(Restradius);
1660 cout << Verbose(4) << "Shift vector to sphere of circumference is " << TempNormal << ".\n";
1661
1662 Center->AddVector(&TempNormal);
1663 cout << Verbose(0) << "Center of sphere of circumference is " << *Center << ".\n";
1664 get_sphere(&OtherCenter, a, b, c, RADIUS);
1665 cout << Verbose(0) << "OtherCenter of sphere of circumference is " << OtherCenter << ".\n";
1666 cout << Verbose(3) << "End of Get_center_of_sphere.\n";
1667};
1668
1669
1670/** Constructs the center of the circumcircle defined by three points \a *a, \a *b and \a *c.
1671 * \param *Center new center on return
1672 * \param *a first point
1673 * \param *b second point
1674 * \param *c third point
1675 */
1676void GetCenterofCircumcircle(Vector *Center, Vector *a, Vector *b, Vector *c)
1677{
1678 Vector helper;
1679 double alpha, beta, gamma;
1680 Vector SideA, SideB, SideC;
1681 SideA.CopyVector(b);
1682 SideA.SubtractVector(c);
1683 SideB.CopyVector(c);
1684 SideB.SubtractVector(a);
1685 SideC.CopyVector(a);
1686 SideC.SubtractVector(b);
1687 alpha = M_PI - SideB.Angle(&SideC);
1688 beta = M_PI - SideC.Angle(&SideA);
1689 gamma = M_PI - SideA.Angle(&SideB);
1690 //cout << Verbose(3) << "INFO: alpha = " << alpha/M_PI*180. << ", beta = " << beta/M_PI*180. << ", gamma = " << gamma/M_PI*180. << "." << endl;
1691 if (fabs(M_PI - alpha - beta - gamma) > HULLEPSILON)
1692 cerr << "GetCenterofCircumcircle: Sum of angles " << (alpha+beta+gamma)/M_PI*180. << " > 180 degrees by " << fabs(M_PI - alpha - beta - gamma)/M_PI*180. << "!" << endl;
1693
1694 Center->Zero();
1695 helper.CopyVector(a);
1696 helper.Scale(sin(2.*alpha));
1697 Center->AddVector(&helper);
1698 helper.CopyVector(b);
1699 helper.Scale(sin(2.*beta));
1700 Center->AddVector(&helper);
1701 helper.CopyVector(c);
1702 helper.Scale(sin(2.*gamma));
1703 Center->AddVector(&helper);
1704 Center->Scale(1./(sin(2.*alpha) + sin(2.*beta) + sin(2.*gamma)));
1705};
1706
1707/** 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.
1708 * Test whether the \a NewSphereCenter is really on the given plane and in distance \a CircleRadius from \a CircleCenter.
1709 * It calculates the angle, making it unique on [0,2.*M_PI) by comparing to SearchDirection.
1710 * 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).
1711 * \param CircleCenter Center of the parameter circle
1712 * \param CirclePlaneNormal normal vector to plane of the parameter circle
1713 * \param CircleRadius radius of the parameter circle
1714 * \param NewSphereCenter new center of a circumcircle
1715 * \param OldSphereCenter old center of a circumcircle, defining the zero "path length" on the parameter circle
1716 * \param NormalVector normal vector
1717 * \param SearchDirection search direction to make angle unique on return.
1718 * \return Angle between \a NewSphereCenter and \a OldSphereCenter relative to \a CircleCenter, 2.*M_PI if one test fails
1719 */
1720double GetPathLengthonCircumCircle(Vector &CircleCenter, Vector &CirclePlaneNormal, double CircleRadius, Vector &NewSphereCenter, Vector &OldSphereCenter, Vector &NormalVector, Vector &SearchDirection)
1721{
1722 Vector helper;
1723 double radius, alpha;
1724
1725 helper.CopyVector(&NewSphereCenter);
1726 // test whether new center is on the parameter circle's plane
1727 if (fabs(helper.ScalarProduct(&CirclePlaneNormal)) > HULLEPSILON) {
1728 cerr << "ERROR: Something's very wrong here: NewSphereCenter is not on the band's plane as desired by " <<fabs(helper.ScalarProduct(&CirclePlaneNormal)) << "!" << endl;
1729 helper.ProjectOntoPlane(&CirclePlaneNormal);
1730 }
1731 radius = helper.ScalarProduct(&helper);
1732 // test whether the new center vector has length of CircleRadius
1733 if (fabs(radius - CircleRadius) > HULLEPSILON)
1734 cerr << Verbose(1) << "ERROR: The projected center of the new sphere has radius " << radius << " instead of " << CircleRadius << "." << endl;
1735 alpha = helper.Angle(&OldSphereCenter);
1736 // make the angle unique by checking the halfplanes/search direction
1737 if (helper.ScalarProduct(&SearchDirection) < -HULLEPSILON) // acos is not unique on [0, 2.*M_PI), hence extra check to decide between two half intervals
1738 alpha = 2.*M_PI - alpha;
1739 //cout << Verbose(2) << "INFO: RelativeNewSphereCenter is " << helper << ", RelativeOldSphereCenter is " << OldSphereCenter << " and resulting angle is " << alpha << "." << endl;
1740 radius = helper.Distance(&OldSphereCenter);
1741 helper.ProjectOntoPlane(&NormalVector);
1742 // check whether new center is somewhat away or at least right over the current baseline to prevent intersecting triangles
1743 if ((radius > HULLEPSILON) || (helper.Norm() < HULLEPSILON)) {
1744 //cout << Verbose(2) << "INFO: Distance between old and new center is " << radius << " and between new center and baseline center is " << helper.Norm() << "." << endl;
1745 return alpha;
1746 } else {
1747 //cout << Verbose(1) << "INFO: NewSphereCenter " << helper << " is too close to OldSphereCenter" << OldSphereCenter << "." << endl;
1748 return 2.*M_PI;
1749 }
1750};
1751
1752
1753/** Checks whether the triangle consisting of the three atoms is already present.
1754 * Searches for the points in Tesselation::PointsOnBoundary and checks their
1755 * lines. If any of the three edges already has two triangles attached, false is
1756 * returned.
1757 * \param *out output stream for debugging
1758 * \param *Candidates endpoints of the triangle candidate
1759 * \return integer 0 if no triangle exists, 1 if one triangle exists, 2 if two
1760 * triangles exist which is the maximum for three points
1761 */
1762int Tesselation::CheckPresenceOfTriangle(ofstream *out, atom *Candidates[3]) {
1763 LineMap::iterator FindLine;
1764 PointMap::iterator FindPoint;
1765 TriangleMap::iterator FindTriangle;
1766 int adjacentTriangleCount = 0;
1767 class BoundaryPointSet *Points[3];
1768
1769 //*out << Verbose(2) << "Begin of CheckPresenceOfTriangle" << endl;
1770 // builds a triangle point set (Points) of the end points
1771 for (int i = 0; i < 3; i++) {
1772 FindPoint = PointsOnBoundary.find(Candidates[i]->nr);
1773 if (FindPoint != PointsOnBoundary.end()) {
1774 Points[i] = FindPoint->second;
1775 } else {
1776 Points[i] = NULL;
1777 }
1778 }
1779
1780 // checks lines between the points in the Points for their adjacent triangles
1781 for (int i = 0; i < 3; i++) {
1782 if (Points[i] != NULL) {
1783 for (int j = i; j < 3; j++) {
1784 if (Points[j] != NULL) {
1785 FindLine = Points[i]->lines.find(Points[j]->node->nr);
1786 if (FindLine != Points[i]->lines.end()) {
1787 for (; FindLine->first == Points[j]->node->nr; FindLine++) {
1788 FindTriangle = FindLine->second->triangles.begin();
1789 for (; FindTriangle != FindLine->second->triangles.end(); FindTriangle++) {
1790 if ((
1791 (FindTriangle->second->endpoints[0] == Points[0])
1792 || (FindTriangle->second->endpoints[0] == Points[1])
1793 || (FindTriangle->second->endpoints[0] == Points[2])
1794 ) && (
1795 (FindTriangle->second->endpoints[1] == Points[0])
1796 || (FindTriangle->second->endpoints[1] == Points[1])
1797 || (FindTriangle->second->endpoints[1] == Points[2])
1798 ) && (
1799 (FindTriangle->second->endpoints[2] == Points[0])
1800 || (FindTriangle->second->endpoints[2] == Points[1])
1801 || (FindTriangle->second->endpoints[2] == Points[2])
1802 )
1803 ) {
1804 adjacentTriangleCount++;
1805 }
1806 }
1807 }
1808 // Only one of the triangle lines must be considered for the triangle count.
1809 *out << Verbose(2) << "Found " << adjacentTriangleCount << " adjacent triangles for the point set." << endl;
1810 return adjacentTriangleCount;
1811
1812 }
1813 }
1814 }
1815 }
1816 }
1817
1818 *out << Verbose(2) << "Found " << adjacentTriangleCount << " adjacent triangles for the point set." << endl;
1819 return adjacentTriangleCount;
1820};
1821
1822/** This recursive function finds a third point, to form a triangle with two given ones.
1823 * Note that this function is for the starting triangle.
1824 * The idea is as follows: A sphere with fixed radius is (almost) uniquely defined in space by three points
1825 * that sit on its boundary. Hence, when two points are given and we look for the (next) third point, then
1826 * the center of the sphere is still fixed up to a single parameter. The band of possible values
1827 * describes a circle in 3D-space. The old center of the sphere for the current base triangle gives
1828 * us the "null" on this circle, the new center of the candidate point will be some way along this
1829 * circle. The shorter the way the better is the candidate. Note that the direction is clearly given
1830 * by the normal vector of the base triangle that always points outwards by construction.
1831 * Hence, we construct a Center of this circle which sits right in the middle of the current base line.
1832 * We construct the normal vector that defines the plane this circle lies in, it is just in the
1833 * direction of the baseline. And finally, we need the radius of the circle, which is given by the rest
1834 * with respect to the length of the baseline and the sphere's fixed \a RADIUS.
1835 * Note that there is one difficulty: The circumcircle is uniquely defined, but for the circumsphere's center
1836 * there are two possibilities which becomes clear from the construction as seen below. Hence, we must check
1837 * both.
1838 * Note also that the acos() function is not unique on [0, 2.*M_PI). Hence, we need an additional check
1839 * to decide for one of the two possible angles. Therefore we need a SearchDirection and to make this check
1840 * sensible we need OldSphereCenter to be orthogonal to it. Either we construct SearchDirection orthogonal
1841 * right away, or -- what we do here -- we rotate the relative sphere centers such that this orthogonality
1842 * holds. Then, the normalized projection onto the SearchDirection is either +1 or -1 and thus states whether
1843 * the angle is uniquely in either (0,M_PI] or [M_PI, 2.*M_PI).
1844 * @param NormalVector normal direction of the base triangle (here the unit axis vector, \sa Find_starting_triangle())
1845 * @param SearchDirection general direction where to search for the next point, relative to center of BaseLine
1846 * @param OldSphereCenter center of sphere for base triangle, relative to center of BaseLine, giving null angle for the parameter circle
1847 * @param BaseLine BoundaryLineSet with the current base line
1848 * @param ThirdNode third atom to avoid in search
1849 * @param candidates list of equally good candidates to return
1850 * @param ShortestAngle the current path length on this circle band for the current Opt_Candidate
1851 * @param RADIUS radius of sphere
1852 * @param *LC LinkedCell structure with neighbouring atoms
1853 */
1854void Find_third_point_for_Tesselation(
1855 Vector NormalVector, Vector SearchDirection, Vector OldSphereCenter,
1856 class BoundaryLineSet *BaseLine, atom *ThirdNode, CandidateList* &candidates,
1857 double *ShortestAngle, const double RADIUS, LinkedCell *LC
1858) {
1859 Vector CircleCenter; // center of the circle, i.e. of the band of sphere's centers
1860 Vector CirclePlaneNormal; // normal vector defining the plane this circle lives in
1861 Vector SphereCenter;
1862 Vector NewSphereCenter; // center of the sphere defined by the two points of BaseLine and the one of Candidate, first possibility
1863 Vector OtherNewSphereCenter; // center of the sphere defined by the two points of BaseLine and the one of Candidate, second possibility
1864 Vector NewNormalVector; // normal vector of the Candidate's triangle
1865 Vector helper, OptCandidateCenter, OtherOptCandidateCenter;
1866 LinkedAtoms *List = NULL;
1867 double CircleRadius; // radius of this circle
1868 double radius;
1869 double alpha, Otheralpha; // angles (i.e. parameter for the circle).
1870 int N[NDIM], Nlower[NDIM], Nupper[NDIM];
1871 atom *Candidate = NULL;
1872 CandidateForTesselation *optCandidate = NULL;
1873
1874 cout << Verbose(1) << "Begin of Find_third_point_for_Tesselation" << endl;
1875
1876 //cout << Verbose(2) << "INFO: NormalVector of BaseTriangle is " << NormalVector << "." << endl;
1877
1878 // construct center of circle
1879 CircleCenter.CopyVector(&(BaseLine->endpoints[0]->node->x));
1880 CircleCenter.AddVector(&BaseLine->endpoints[1]->node->x);
1881 CircleCenter.Scale(0.5);
1882
1883 // construct normal vector of circle
1884 CirclePlaneNormal.CopyVector(&BaseLine->endpoints[0]->node->x);
1885 CirclePlaneNormal.SubtractVector(&BaseLine->endpoints[1]->node->x);
1886
1887 // calculate squared radius atom *ThirdNode,f circle
1888 radius = CirclePlaneNormal.ScalarProduct(&CirclePlaneNormal);
1889 if (radius/4. < RADIUS*RADIUS) {
1890 CircleRadius = RADIUS*RADIUS - radius/4.;
1891 CirclePlaneNormal.Normalize();
1892 //cout << Verbose(2) << "INFO: CircleCenter is at " << CircleCenter << ", CirclePlaneNormal is " << CirclePlaneNormal << " with circle radius " << sqrt(CircleRadius) << "." << endl;
1893
1894 // test whether old center is on the band's plane
1895 if (fabs(OldSphereCenter.ScalarProduct(&CirclePlaneNormal)) > HULLEPSILON) {
1896 cerr << "ERROR: Something's very wrong here: OldSphereCenter is not on the band's plane as desired by " << fabs(OldSphereCenter.ScalarProduct(&CirclePlaneNormal)) << "!" << endl;
1897 OldSphereCenter.ProjectOntoPlane(&CirclePlaneNormal);
1898 }
1899 radius = OldSphereCenter.ScalarProduct(&OldSphereCenter);
1900 if (fabs(radius - CircleRadius) < HULLEPSILON) {
1901
1902 // check SearchDirection
1903 //cout << Verbose(2) << "INFO: SearchDirection is " << SearchDirection << "." << endl;
1904 if (fabs(OldSphereCenter.ScalarProduct(&SearchDirection)) > HULLEPSILON) { // rotated the wrong way!
1905 cerr << "ERROR: SearchDirection and RelativeOldSphereCenter are not orthogonal!" << endl;
1906 }
1907
1908 // get cell for the starting atom
1909 if (LC->SetIndexToVector(&CircleCenter)) {
1910 for(int i=0;i<NDIM;i++) // store indices of this cell
1911 N[i] = LC->n[i];
1912 //cout << Verbose(2) << "INFO: Center cell is " << N[0] << ", " << N[1] << ", " << N[2] << " with No. " << LC->index << "." << endl;
1913 } else {
1914 cerr << "ERROR: Vector " << CircleCenter << " is outside of LinkedCell's bounding box." << endl;
1915 return;
1916 }
1917 // then go through the current and all neighbouring cells and check the contained atoms for possible candidates
1918 //cout << Verbose(2) << "LC Intervals:";
1919 for (int i=0;i<NDIM;i++) {
1920 Nlower[i] = ((N[i]-1) >= 0) ? N[i]-1 : 0;
1921 Nupper[i] = ((N[i]+1) < LC->N[i]) ? N[i]+1 : LC->N[i]-1;
1922 //cout << " [" << Nlower[i] << "," << Nupper[i] << "] ";
1923 }
1924 //cout << endl;
1925 for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
1926 for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
1927 for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
1928 List = LC->GetCurrentCell();
1929 //cout << Verbose(2) << "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << "." << endl;
1930 if (List != NULL) {
1931 for (LinkedAtoms::iterator Runner = List->begin(); Runner != List->end(); Runner++) {
1932 Candidate = (*Runner);
1933
1934 // check for three unique points
1935 //cout << Verbose(2) << "INFO: Current Candidate is " << *Candidate << " at " << Candidate->x << "." << endl;
1936 if ((Candidate != BaseLine->endpoints[0]->node) && (Candidate != BaseLine->endpoints[1]->node) ){
1937
1938 // construct both new centers
1939 GetCenterofCircumcircle(&NewSphereCenter, &(BaseLine->endpoints[0]->node->x), &(BaseLine->endpoints[1]->node->x), &(Candidate->x));
1940 OtherNewSphereCenter.CopyVector(&NewSphereCenter);
1941
1942 if ((NewNormalVector.MakeNormalVector(&(BaseLine->endpoints[0]->node->x), &(BaseLine->endpoints[1]->node->x), &(Candidate->x)))
1943 && (fabs(NewNormalVector.ScalarProduct(&NewNormalVector)) > HULLEPSILON)
1944 ) {
1945 helper.CopyVector(&NewNormalVector);
1946 //cout << Verbose(2) << "INFO: NewNormalVector is " << NewNormalVector << "." << endl;
1947 radius = BaseLine->endpoints[0]->node->x.DistanceSquared(&NewSphereCenter);
1948 if (radius < RADIUS*RADIUS) {
1949 helper.Scale(sqrt(RADIUS*RADIUS - radius));
1950 //cout << Verbose(2) << "INFO: Distance of NewCircleCenter to NewSphereCenter is " << helper.Norm() << " with sphere radius " << RADIUS << "." << endl;
1951 NewSphereCenter.AddVector(&helper);
1952 NewSphereCenter.SubtractVector(&CircleCenter);
1953 //cout << Verbose(2) << "INFO: NewSphereCenter is at " << NewSphereCenter << "." << endl;
1954
1955 // OtherNewSphereCenter is created by the same vector just in the other direction
1956 helper.Scale(-1.);
1957 OtherNewSphereCenter.AddVector(&helper);
1958 OtherNewSphereCenter.SubtractVector(&CircleCenter);
1959 //cout << Verbose(2) << "INFO: OtherNewSphereCenter is at " << OtherNewSphereCenter << "." << endl;
1960
1961 alpha = GetPathLengthonCircumCircle(CircleCenter, CirclePlaneNormal, CircleRadius, NewSphereCenter, OldSphereCenter, NormalVector, SearchDirection);
1962 Otheralpha = GetPathLengthonCircumCircle(CircleCenter, CirclePlaneNormal, CircleRadius, OtherNewSphereCenter, OldSphereCenter, NormalVector, SearchDirection);
1963 alpha = min(alpha, Otheralpha);
1964 // if there is a better candidate, drop the current list and add the new candidate
1965 // otherwise ignore the new candidate and keep the list
1966 if (*ShortestAngle > (alpha - HULLEPSILON)) {
1967 optCandidate = new CandidateForTesselation(Candidate, BaseLine, OptCandidateCenter, OtherOptCandidateCenter);
1968 if (fabs(alpha - Otheralpha) > MYEPSILON) {
1969 optCandidate->OptCenter.CopyVector(&NewSphereCenter);
1970 optCandidate->OtherOptCenter.CopyVector(&OtherNewSphereCenter);
1971 } else {
1972 optCandidate->OptCenter.CopyVector(&OtherNewSphereCenter);
1973 optCandidate->OtherOptCenter.CopyVector(&NewSphereCenter);
1974 }
1975 // if there is an equal candidate, add it to the list without clearing the list
1976 if ((*ShortestAngle - HULLEPSILON) < alpha) {
1977 candidates->push_back(optCandidate);
1978 cout << Verbose(2) << "ACCEPT: We have found an equally good candidate: " << *(optCandidate->point) << " with "
1979 << alpha << " and circumsphere's center at " << optCandidate->OptCenter << "." << endl;
1980 } else {
1981 // remove all candidates from the list and then the list itself
1982 class CandidateForTesselation *remover = NULL;
1983 for (CandidateList::iterator it = candidates->begin(); it != candidates->end(); ++it) {
1984 remover = *it;
1985 delete(remover);
1986 }
1987 candidates->clear();
1988 candidates->push_back(optCandidate);
1989 cout << Verbose(2) << "ACCEPT: We have found a better candidate: " << *(optCandidate->point) << " with "
1990 << alpha << " and circumsphere's center at " << optCandidate->OptCenter << "." << endl;
1991 }
1992 *ShortestAngle = alpha;
1993 //cout << Verbose(2) << "INFO: There are " << candidates->size() << " candidates in the list now." << endl;
1994 } else {
1995 if ((optCandidate != NULL) && (optCandidate->point != NULL)) {
1996 //cout << Verbose(2) << "REJECT: Old candidate: " << *(optCandidate->point) << " is better than " << alpha << " with " << *ShortestAngle << "." << endl;
1997 } else {
1998 //cout << Verbose(2) << "REJECT: Candidate " << *Candidate << " with " << alpha << " was rejected." << endl;
1999 }
2000 }
2001
2002 } else {
2003 //cout << Verbose(2) << "REJECT: NewSphereCenter " << NewSphereCenter << " is too far away: " << radius << "." << endl;
2004 }
2005 } else {
2006 //cout << Verbose(2) << "REJECT: Three points from " << *BaseLine << " and Candidate " << *Candidate << " are linear-dependent." << endl;
2007 }
2008 } else {
2009 if (ThirdNode != NULL) {
2010 //cout << Verbose(2) << "REJECT: Base triangle " << *BaseLine << " and " << *ThirdNode << " contains Candidate " << *Candidate << "." << endl;
2011 } else {
2012 //cout << Verbose(2) << "REJECT: Base triangle " << *BaseLine << " contains Candidate " << *Candidate << "." << endl;
2013 }
2014 }
2015 }
2016 }
2017 }
2018 } else {
2019 cerr << Verbose(2) << "ERROR: The projected center of the old sphere has radius " << radius << " instead of " << CircleRadius << "." << endl;
2020 }
2021 } else {
2022 if (ThirdNode != NULL)
2023 cout << Verbose(2) << "Circumcircle for base line " << *BaseLine << " and third node " << *ThirdNode << " is too big!" << endl;
2024 else
2025 cout << Verbose(2) << "Circumcircle for base line " << *BaseLine << " is too big!" << endl;
2026 }
2027
2028 //cout << Verbose(2) << "INFO: Sorting candidate list ..." << endl;
2029 if (candidates->size() > 1) {
2030 candidates->unique();
2031 candidates->sort(sortCandidates);
2032 }
2033
2034 cout << Verbose(1) << "End of Find_third_point_for_Tesselation" << endl;
2035};
2036
2037struct Intersection {
2038 Vector x1;
2039 Vector x2;
2040 Vector x3;
2041 Vector x4;
2042};
2043
2044/**
2045 * Intersection calculation function.
2046 *
2047 * @param x to find the result for
2048 * @param function parameter
2049 */
2050double MinIntersectDistance(const gsl_vector * x, void *params) {
2051 double retval = 0;
2052 struct Intersection *I = (struct Intersection *)params;
2053 Vector intersection;
2054 Vector SideA,SideB,HeightA, HeightB;
2055 for (int i=0;i<NDIM;i++)
2056 intersection.x[i] = gsl_vector_get(x, i);
2057
2058 SideA.CopyVector(&(I->x1));
2059 SideA.SubtractVector(&I->x2);
2060 HeightA.CopyVector(&intersection);
2061 HeightA.SubtractVector(&I->x1);
2062 HeightA.ProjectOntoPlane(&SideA);
2063
2064 SideB.CopyVector(&I->x3);
2065 SideB.SubtractVector(&I->x4);
2066 HeightB.CopyVector(&intersection);
2067 HeightB.SubtractVector(&I->x3);
2068 HeightB.ProjectOntoPlane(&SideB);
2069
2070 retval = HeightA.ScalarProduct(&HeightA) + HeightB.ScalarProduct(&HeightB);
2071 //cout << Verbose(2) << "MinIntersectDistance called, result: " << retval << endl;
2072
2073 return retval;
2074};
2075
2076
2077/**
2078 * Calculates whether there is an intersection between two lines. The first line
2079 * always goes through point 1 and point 2 and the second line is given by the
2080 * connection between point 4 and point 5.
2081 *
2082 * @param point 1 of line 1
2083 * @param point 2 of line 1
2084 * @param point 1 of line 2
2085 * @param point 2 of line 2
2086 *
2087 * @return true if there is an intersection between the given lines, false otherwise
2088 */
2089bool existsIntersection(Vector point1, Vector point2, Vector point3, Vector point4) {
2090 bool result;
2091
2092 struct Intersection par;
2093 par.x1.CopyVector(&point1);
2094 par.x2.CopyVector(&point2);
2095 par.x3.CopyVector(&point3);
2096 par.x4.CopyVector(&point4);
2097
2098 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
2099 gsl_multimin_fminimizer *s = NULL;
2100 gsl_vector *ss, *x;
2101 gsl_multimin_function minex_func;
2102
2103 size_t iter = 0;
2104 int status;
2105 double size;
2106
2107 /* Starting point */
2108 x = gsl_vector_alloc(NDIM);
2109 gsl_vector_set(x, 0, point1.x[0]);
2110 gsl_vector_set(x, 1, point1.x[1]);
2111 gsl_vector_set(x, 2, point1.x[2]);
2112
2113 /* Set initial step sizes to 1 */
2114 ss = gsl_vector_alloc(NDIM);
2115 gsl_vector_set_all(ss, 1.0);
2116
2117 /* Initialize method and iterate */
2118 minex_func.n = NDIM;
2119 minex_func.f = &MinIntersectDistance;
2120 minex_func.params = (void *)&par;
2121
2122 s = gsl_multimin_fminimizer_alloc(T, NDIM);
2123 gsl_multimin_fminimizer_set(s, &minex_func, x, ss);
2124
2125 do {
2126 iter++;
2127 status = gsl_multimin_fminimizer_iterate(s);
2128
2129 if (status) {
2130 break;
2131 }
2132
2133 size = gsl_multimin_fminimizer_size(s);
2134 status = gsl_multimin_test_size(size, 1e-2);
2135
2136 if (status == GSL_SUCCESS) {
2137 cout << Verbose(2) << "converged to minimum" << endl;
2138 }
2139 } while (status == GSL_CONTINUE && iter < 100);
2140
2141 // check whether intersection is in between or not
2142 Vector intersection, SideA, SideB, HeightA, HeightB;
2143 double t1, t2;
2144 for (int i = 0; i < NDIM; i++) {
2145 intersection.x[i] = gsl_vector_get(s->x, i);
2146 }
2147
2148 SideA.CopyVector(&par.x2);
2149 SideA.SubtractVector(&par.x1);
2150 HeightA.CopyVector(&intersection);
2151 HeightA.SubtractVector(&par.x1);
2152
2153 t1 = HeightA.Projection(&SideA)/SideA.ScalarProduct(&SideA);
2154
2155 SideB.CopyVector(&par.x4);
2156 SideB.SubtractVector(&par.x3);
2157 HeightB.CopyVector(&intersection);
2158 HeightB.SubtractVector(&par.x3);
2159
2160 t2 = HeightB.Projection(&SideB)/SideB.ScalarProduct(&SideB);
2161
2162 cout << Verbose(2) << "Intersection " << intersection << " is at "
2163 << t1 << " for (" << point1 << "," << point2 << ") and at "
2164 << t2 << " for (" << point3 << "," << point4 << "): ";
2165
2166 if (((t1 >= 0) && (t1 <= 1)) && ((t2 >= 0) && (t2 <= 1))) {
2167 cout << "true intersection." << endl;
2168 result = true;
2169 } else {
2170 cout << "intersection out of region of interest." << endl;
2171 result = false;
2172 }
2173
2174 // free minimizer stuff
2175 gsl_vector_free(x);
2176 gsl_vector_free(ss);
2177 gsl_multimin_fminimizer_free(s);
2178
2179 return result;
2180}
2181
2182/** Finds the second point of starting triangle.
2183 * \param *a first atom
2184 * \param *Candidate pointer to candidate atom on return
2185 * \param Oben vector indicating the outside
2186 * \param Opt_Candidate reference to recommended candidate on return
2187 * \param Storage[3] array storing angles and other candidate information
2188 * \param RADIUS radius of virtual sphere
2189 * \param *LC LinkedCell structure with neighbouring atoms
2190 */
2191void Find_second_point_for_Tesselation(atom* a, atom* Candidate, Vector Oben, atom*& Opt_Candidate, double Storage[3], double RADIUS, LinkedCell *LC)
2192{
2193 cout << Verbose(2) << "Begin of Find_second_point_for_Tesselation" << endl;
2194 Vector AngleCheck;
2195 double norm = -1., angle;
2196 LinkedAtoms *List = NULL;
2197 int N[NDIM], Nlower[NDIM], Nupper[NDIM];
2198
2199 if (LC->SetIndexToAtom(a)) { // get cell for the starting atom
2200 for(int i=0;i<NDIM;i++) // store indices of this cell
2201 N[i] = LC->n[i];
2202 } else {
2203 cerr << "ERROR: Atom " << *a << " is not found in cell " << LC->index << "." << endl;
2204 return;
2205 }
2206 // then go through the current and all neighbouring cells and check the contained atoms for possible candidates
2207 cout << Verbose(3) << "LC Intervals from [";
2208 for (int i=0;i<NDIM;i++) {
2209 cout << " " << N[i] << "<->" << LC->N[i];
2210 }
2211 cout << "] :";
2212 for (int i=0;i<NDIM;i++) {
2213 Nlower[i] = ((N[i]-1) >= 0) ? N[i]-1 : 0;
2214 Nupper[i] = ((N[i]+1) < LC->N[i]) ? N[i]+1 : LC->N[i]-1;
2215 cout << " [" << Nlower[i] << "," << Nupper[i] << "] ";
2216 }
2217 cout << endl;
2218
2219
2220 for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
2221 for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
2222 for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
2223 List = LC->GetCurrentCell();
2224 //cout << Verbose(2) << "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << "." << endl;
2225 if (List != NULL) {
2226 for (LinkedAtoms::iterator Runner = List->begin(); Runner != List->end(); Runner++) {
2227 Candidate = (*Runner);
2228 // check if we only have one unique point yet ...
2229 if (a != Candidate) {
2230 // Calculate center of the circle with radius RADIUS through points a and Candidate
2231 Vector OrthogonalizedOben, a_Candidate, Center;
2232 double distance, scaleFactor;
2233
2234 OrthogonalizedOben.CopyVector(&Oben);
2235 a_Candidate.CopyVector(&(a->x));
2236 a_Candidate.SubtractVector(&(Candidate->x));
2237 OrthogonalizedOben.ProjectOntoPlane(&a_Candidate);
2238 OrthogonalizedOben.Normalize();
2239 distance = 0.5 * a_Candidate.Norm();
2240 scaleFactor = sqrt(((RADIUS * RADIUS) - (distance * distance)));
2241 OrthogonalizedOben.Scale(scaleFactor);
2242
2243 Center.CopyVector(&(Candidate->x));
2244 Center.AddVector(&(a->x));
2245 Center.Scale(0.5);
2246 Center.AddVector(&OrthogonalizedOben);
2247
2248 AngleCheck.CopyVector(&Center);
2249 AngleCheck.SubtractVector(&(a->x));
2250 norm = a_Candidate.Norm();
2251 // second point shall have smallest angle with respect to Oben vector
2252 if (norm < RADIUS*2.) {
2253 angle = AngleCheck.Angle(&Oben);
2254 if (angle < Storage[0]) {
2255 //cout << Verbose(3) << "Old values of Storage: %lf %lf \n", Storage[0], Storage[1]);
2256 cout << Verbose(3) << "Current candidate is " << *Candidate << ": Is a better candidate with distance " << norm << " and angle " << angle << " to oben " << Oben << ".\n";
2257 Opt_Candidate = Candidate;
2258 Storage[0] = angle;
2259 //cout << Verbose(3) << "Changing something in Storage: %lf %lf. \n", Storage[0], Storage[2]);
2260 } else {
2261 //cout << Verbose(3) << "Current candidate is " << *Candidate << ": Looses with angle " << angle << " to a better candidate " << *Opt_Candidate << endl;
2262 }
2263 } else {
2264 //cout << Verbose(3) << "Current candidate is " << *Candidate << ": Refused due to Radius " << norm << endl;
2265 }
2266 } else {
2267 //cout << Verbose(3) << "Current candidate is " << *Candidate << ": Candidate is equal to first endpoint." << *a << "." << endl;
2268 }
2269 }
2270 } else {
2271 cout << Verbose(3) << "Linked cell list is empty." << endl;
2272 }
2273 }
2274 cout << Verbose(2) << "End of Find_second_point_for_Tesselation" << endl;
2275};
2276
2277/** Finds the starting triangle for find_non_convex_border().
2278 * Looks at the outermost atom per axis, then Find_second_point_for_Tesselation()
2279 * for the second and Find_next_suitable_point_via_Angle_of_Sphere() for the third
2280 * point are called.
2281 * \param RADIUS radius of virtual rolling sphere
2282 * \param *LC LinkedCell structure with neighbouring atoms
2283 */
2284void Tesselation::Find_starting_triangle(ofstream *out, molecule *mol, const double RADIUS, LinkedCell *LC)
2285{
2286 cout << Verbose(1) << "Begin of Find_starting_triangle\n";
2287 int i = 0;
2288 LinkedAtoms *List = NULL;
2289 atom* FirstPoint = NULL;
2290 atom* SecondPoint = NULL;
2291 atom* MaxAtom[NDIM];
2292 double max_coordinate[NDIM];
2293 Vector Oben;
2294 Vector helper;
2295 Vector Chord;
2296 Vector SearchDirection;
2297
2298 Oben.Zero();
2299
2300 for (i = 0; i < 3; i++) {
2301 MaxAtom[i] = NULL;
2302 max_coordinate[i] = -1;
2303 }
2304
2305 // 1. searching topmost atom with respect to each axis
2306 for (int i=0;i<NDIM;i++) { // each axis
2307 LC->n[i] = LC->N[i]-1; // current axis is topmost cell
2308 for (LC->n[(i+1)%NDIM]=0;LC->n[(i+1)%NDIM]<LC->N[(i+1)%NDIM];LC->n[(i+1)%NDIM]++)
2309 for (LC->n[(i+2)%NDIM]=0;LC->n[(i+2)%NDIM]<LC->N[(i+2)%NDIM];LC->n[(i+2)%NDIM]++) {
2310 List = LC->GetCurrentCell();
2311 //cout << Verbose(2) << "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << "." << endl;
2312 if (List != NULL) {
2313 for (LinkedAtoms::iterator Runner = List->begin();Runner != List->end();Runner++) {
2314 if ((*Runner)->x.x[i] > max_coordinate[i]) {
2315 cout << Verbose(2) << "New maximal for axis " << i << " atom is " << *(*Runner) << " at " << (*Runner)->x << "." << endl;
2316 max_coordinate[i] = (*Runner)->x.x[i];
2317 MaxAtom[i] = (*Runner);
2318 }
2319 }
2320 } else {
2321 cerr << "ERROR: The current cell " << LC->n[0] << "," << LC->n[1] << "," << LC->n[2] << " is invalid!" << endl;
2322 }
2323 }
2324 }
2325
2326 cout << Verbose(2) << "Found maximum coordinates: ";
2327 for (int i=0;i<NDIM;i++)
2328 cout << i << ": " << *MaxAtom[i] << "\t";
2329 cout << endl;
2330
2331 BTS = NULL;
2332 CandidateList *Opt_Candidates = new CandidateList();
2333 for (int k=0;k<NDIM;k++) {
2334 Oben.x[k] = 1.;
2335 FirstPoint = MaxAtom[k];
2336 cout << Verbose(1) << "Coordinates of start atom " << *FirstPoint << " at " << FirstPoint->x << "." << endl;
2337
2338 double ShortestAngle;
2339 atom* Opt_Candidate = NULL;
2340 ShortestAngle = 999999.; // This will contain the angle, which will be always positive (when looking for second point), when looking for third point this will be the quadrant.
2341
2342 Find_second_point_for_Tesselation(FirstPoint, NULL, Oben, Opt_Candidate, &ShortestAngle, RADIUS, LC); // we give same point as next candidate as its bonds are looked into in find_second_...
2343 SecondPoint = Opt_Candidate;
2344 if (SecondPoint == NULL) // have we found a second point?
2345 continue;
2346 else
2347 cout << Verbose(1) << "Found second point is " << *SecondPoint << " at " << SecondPoint->x << ".\n";
2348
2349 helper.CopyVector(&(FirstPoint->x));
2350 helper.SubtractVector(&(SecondPoint->x));
2351 helper.Normalize();
2352 Oben.ProjectOntoPlane(&helper);
2353 Oben.Normalize();
2354 helper.VectorProduct(&Oben);
2355 ShortestAngle = 2.*M_PI; // This will indicate the quadrant.
2356
2357 Chord.CopyVector(&(FirstPoint->x)); // bring into calling function
2358 Chord.SubtractVector(&(SecondPoint->x));
2359 double radius = Chord.ScalarProduct(&Chord);
2360 double CircleRadius = sqrt(RADIUS*RADIUS - radius/4.);
2361 helper.CopyVector(&Oben);
2362 helper.Scale(CircleRadius);
2363 // Now, oben and helper are two orthonormalized vectors in the plane defined by Chord (not normalized)
2364
2365 // look in one direction of baseline for initial candidate
2366 SearchDirection.MakeNormalVector(&Chord, &Oben); // whether we look "left" first or "right" first is not important ...
2367
2368 // adding point 1 and point 2 and the line between them
2369 AddTrianglePoint(FirstPoint, 0);
2370 AddTrianglePoint(SecondPoint, 1);
2371 AddTriangleLine(TPS[0], TPS[1], 0);
2372
2373 //cout << Verbose(2) << "INFO: OldSphereCenter is at " << helper << ".\n";
2374 Find_third_point_for_Tesselation(
2375 Oben, SearchDirection, helper, BLS[0], NULL, *&Opt_Candidates, &ShortestAngle, RADIUS, LC
2376 );
2377 cout << Verbose(1) << "List of third Points is ";
2378 for (CandidateList::iterator it = Opt_Candidates->begin(); it != Opt_Candidates->end(); ++it) {
2379 cout << " " << *(*it)->point;
2380 }
2381 cout << endl;
2382
2383 for (CandidateList::iterator it = Opt_Candidates->begin(); it != Opt_Candidates->end(); ++it) {
2384 // add third triangle point
2385 AddTrianglePoint((*it)->point, 2);
2386 // add the second and third line
2387 AddTriangleLine(TPS[1], TPS[2], 1);
2388 AddTriangleLine(TPS[0], TPS[2], 2);
2389 // ... and triangles to the Maps of the Tesselation class
2390 BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
2391 AddTriangle();
2392 // ... and calculate its normal vector (with correct orientation)
2393 (*it)->OptCenter.Scale(-1.);
2394 cout << Verbose(2) << "Anti-Oben is currently " << (*it)->OptCenter << "." << endl;
2395 BTS->GetNormalVector((*it)->OptCenter); // vector to compare with should point inwards
2396 cout << Verbose(0) << "==> Found starting triangle consists of " << *FirstPoint << ", " << *SecondPoint << " and "
2397 << *(*it)->point << " with normal vector " << BTS->NormalVector << ".\n";
2398
2399 // if we do not reach the end with the next step of iteration, we need to setup a new first line
2400 if (it != Opt_Candidates->end()--) {
2401 FirstPoint = (*it)->BaseLine->endpoints[0]->node;
2402 SecondPoint = (*it)->point;
2403 // adding point 1 and point 2 and the line between them
2404 AddTrianglePoint(FirstPoint, 0);
2405 AddTrianglePoint(SecondPoint, 1);
2406 AddTriangleLine(TPS[0], TPS[1], 0);
2407 }
2408 cout << Verbose(2) << "Projection is " << BTS->NormalVector.Projection(&Oben) << "." << endl;
2409 }
2410 if (BTS != NULL) // we have created one starting triangle
2411 break;
2412 else {
2413 // remove all candidates from the list and then the list itself
2414 class CandidateForTesselation *remover = NULL;
2415 for (CandidateList::iterator it = Opt_Candidates->begin(); it != Opt_Candidates->end(); ++it) {
2416 remover = *it;
2417 delete(remover);
2418 }
2419 Opt_Candidates->clear();
2420 }
2421 }
2422
2423 // remove all candidates from the list and then the list itself
2424 class CandidateForTesselation *remover = NULL;
2425 for (CandidateList::iterator it = Opt_Candidates->begin(); it != Opt_Candidates->end(); ++it) {
2426 remover = *it;
2427 delete(remover);
2428 }
2429 delete(Opt_Candidates);
2430 cout << Verbose(1) << "End of Find_starting_triangle\n";
2431};
2432
2433/** Checks for a new special triangle whether one of its edges is already present with one one triangle connected.
2434 * This enforces that special triangles (i.e. degenerated ones) should at last close the open-edge frontier and not
2435 * make it bigger (i.e. closing one (the baseline) and opening two new ones).
2436 * \param TPS[3] nodes of the triangle
2437 * \return true - there is such a line (i.e. creation of degenerated triangle is valid), false - no such line (don't create)
2438 */
2439bool CheckLineCriteriaforDegeneratedTriangle(class BoundaryPointSet *nodes[3])
2440{
2441 bool result = false;
2442 int counter = 0;
2443
2444 // check all three points
2445 for (int i=0;i<3;i++)
2446 for (int j=i+1; j<3; j++) {
2447 if (nodes[i]->lines.find(nodes[j]->node->nr) != nodes[i]->lines.end()) { // there already is a line
2448 LineMap::iterator FindLine;
2449 pair<LineMap::iterator,LineMap::iterator> FindPair;
2450 FindPair = nodes[i]->lines.equal_range(nodes[j]->node->nr);
2451 for (FindLine = FindPair.first; FindLine != FindPair.second; ++FindLine) {
2452 // If there is a line with less than two attached triangles, we don't need a new line.
2453 if (FindLine->second->TrianglesCount < 2) {
2454 counter++;
2455 break; // increase counter only once per edge
2456 }
2457 }
2458 } else { // no line
2459 cout << Verbose(1) << "ERROR: The line between " << nodes[i] << " and " << nodes[j] << " is not yet present, hence no need for a degenerate triangle!" << endl;
2460 result = true;
2461 }
2462 }
2463 if (counter > 1) {
2464 cout << Verbose(2) << "INFO: Degenerate triangle is ok, at least two, here " << counter << ", existing lines are used." << endl;
2465 result = true;
2466 }
2467 return result;
2468};
2469
2470
2471/** This function finds a triangle to a line, adjacent to an existing one.
2472 * @param out output stream for debugging
2473 * @param *mol molecule with Atom's and Bond's
2474 * @param Line current baseline to search from
2475 * @param T current triangle which \a Line is edge of
2476 * @param RADIUS radius of the rolling ball
2477 * @param N number of found triangles
2478 * @param *filename filename base for intermediate envelopes
2479 * @param *LC LinkedCell structure with neighbouring atoms
2480 */
2481bool Tesselation::Find_next_suitable_triangle(ofstream *out,
2482 molecule *mol, BoundaryLineSet &Line, BoundaryTriangleSet &T,
2483 const double& RADIUS, int N, const char *tempbasename, LinkedCell *LC)
2484{
2485 cout << Verbose(0) << "Begin of Find_next_suitable_triangle\n";
2486 ofstream *tempstream = NULL;
2487 char NumberName[255];
2488 bool result = true;
2489 CandidateList *Opt_Candidates = new CandidateList();
2490
2491 Vector CircleCenter;
2492 Vector CirclePlaneNormal;
2493 Vector OldSphereCenter;
2494 Vector SearchDirection;
2495 Vector helper;
2496 atom *ThirdNode = NULL;
2497 LineMap::iterator testline;
2498 double ShortestAngle = 2.*M_PI; // This will indicate the quadrant.
2499 double radius, CircleRadius;
2500
2501 cout << Verbose(1) << "Current baseline is " << Line << " of triangle " << T << "." << endl;
2502 for (int i=0;i<3;i++)
2503 if ((T.endpoints[i]->node != Line.endpoints[0]->node) && (T.endpoints[i]->node != Line.endpoints[1]->node))
2504 ThirdNode = T.endpoints[i]->node;
2505
2506 // construct center of circle
2507 CircleCenter.CopyVector(&Line.endpoints[0]->node->x);
2508 CircleCenter.AddVector(&Line.endpoints[1]->node->x);
2509 CircleCenter.Scale(0.5);
2510
2511 // construct normal vector of circle
2512 CirclePlaneNormal.CopyVector(&Line.endpoints[0]->node->x);
2513 CirclePlaneNormal.SubtractVector(&Line.endpoints[1]->node->x);
2514
2515 // calculate squared radius of circle
2516 radius = CirclePlaneNormal.ScalarProduct(&CirclePlaneNormal);
2517 if (radius/4. < RADIUS*RADIUS) {
2518 CircleRadius = RADIUS*RADIUS - radius/4.;
2519 CirclePlaneNormal.Normalize();
2520 cout << Verbose(2) << "INFO: CircleCenter is at " << CircleCenter << ", CirclePlaneNormal is " << CirclePlaneNormal << " with circle radius " << sqrt(CircleRadius) << "." << endl;
2521
2522 // construct old center
2523 GetCenterofCircumcircle(&OldSphereCenter, &(T.endpoints[0]->node->x), &(T.endpoints[1]->node->x), &(T.endpoints[2]->node->x));
2524 helper.CopyVector(&T.NormalVector); // normal vector ensures that this is correct center of the two possible ones
2525 radius = Line.endpoints[0]->node->x.DistanceSquared(&OldSphereCenter);
2526 helper.Scale(sqrt(RADIUS*RADIUS - radius));
2527 OldSphereCenter.AddVector(&helper);
2528 OldSphereCenter.SubtractVector(&CircleCenter);
2529 //cout << Verbose(2) << "INFO: OldSphereCenter is at " << OldSphereCenter << "." << endl;
2530
2531 // construct SearchDirection
2532 SearchDirection.MakeNormalVector(&T.NormalVector, &CirclePlaneNormal);
2533 helper.CopyVector(&Line.endpoints[0]->node->x);
2534 helper.SubtractVector(&ThirdNode->x);
2535 if (helper.ScalarProduct(&SearchDirection) < -HULLEPSILON)// ohoh, SearchDirection points inwards!
2536 SearchDirection.Scale(-1.);
2537 SearchDirection.ProjectOntoPlane(&OldSphereCenter);
2538 SearchDirection.Normalize();
2539 cout << Verbose(2) << "INFO: SearchDirection is " << SearchDirection << "." << endl;
2540 if (fabs(OldSphereCenter.ScalarProduct(&SearchDirection)) > HULLEPSILON) {
2541 // rotated the wrong way!
2542 cerr << "ERROR: SearchDirection and RelativeOldSphereCenter are still not orthogonal!" << endl;
2543 }
2544
2545 // add third point
2546 Find_third_point_for_Tesselation(
2547 T.NormalVector, SearchDirection, OldSphereCenter, &Line, ThirdNode, Opt_Candidates,
2548 &ShortestAngle, RADIUS, LC
2549 );
2550
2551 } else {
2552 cout << Verbose(1) << "Circumcircle for base line " << Line << " and base triangle " << T << " is too big!" << endl;
2553 }
2554
2555 if (Opt_Candidates->begin() == Opt_Candidates->end()) {
2556 cerr << "WARNING: Could not find a suitable candidate." << endl;
2557 return false;
2558 }
2559 cout << Verbose(1) << "Third Points are ";
2560 for (CandidateList::iterator it = Opt_Candidates->begin(); it != Opt_Candidates->end(); ++it) {
2561 cout << " " << *(*it)->point;
2562 }
2563 cout << endl;
2564
2565 BoundaryLineSet *BaseRay = &Line;
2566 for (CandidateList::iterator it = Opt_Candidates->begin(); it != Opt_Candidates->end(); ++it) {
2567 cout << Verbose(1) << " Third point candidate is " << *(*it)->point
2568 << " with circumsphere's center at " << (*it)->OptCenter << "." << endl;
2569 cout << Verbose(1) << " Baseline is " << *BaseRay << endl;
2570
2571 // check whether all edges of the new triangle still have space for one more triangle (i.e. TriangleCount <2)
2572 atom *AtomCandidates[3];
2573 AtomCandidates[0] = (*it)->point;
2574 AtomCandidates[1] = BaseRay->endpoints[0]->node;
2575 AtomCandidates[2] = BaseRay->endpoints[1]->node;
2576 int existentTrianglesCount = CheckPresenceOfTriangle(out, AtomCandidates);
2577
2578 BTS = NULL;
2579 // If there is no triangle, add it regularly.
2580 if (existentTrianglesCount == 0) {
2581 AddTrianglePoint((*it)->point, 0);
2582 AddTrianglePoint(BaseRay->endpoints[0]->node, 1);
2583 AddTrianglePoint(BaseRay->endpoints[1]->node, 2);
2584
2585 AddTriangleLine(TPS[0], TPS[1], 0);
2586 AddTriangleLine(TPS[0], TPS[2], 1);
2587 AddTriangleLine(TPS[1], TPS[2], 2);
2588
2589 BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
2590 AddTriangle();
2591 (*it)->OptCenter.Scale(-1.);
2592 BTS->GetNormalVector((*it)->OptCenter);
2593 (*it)->OptCenter.Scale(-1.);
2594
2595 cout << "--> New triangle with " << *BTS << " and normal vector " << BTS->NormalVector
2596 << " for this triangle ... " << endl;
2597 //cout << Verbose(1) << "We have "<< TrianglesOnBoundaryCount << " for line " << *BaseRay << "." << endl;
2598 } else if (existentTrianglesCount == 1) { // If there is a planar region within the structure, we need this triangle a second time.
2599 AddTrianglePoint((*it)->point, 0);
2600 AddTrianglePoint(BaseRay->endpoints[0]->node, 1);
2601 AddTrianglePoint(BaseRay->endpoints[1]->node, 2);
2602
2603 // We demand that at most one new degenerate line is created and that this line also already exists (which has to be the case due to existentTrianglesCount == 1)
2604 // i.e. at least one of the three lines must be present with TriangleCount <= 1
2605 if (CheckLineCriteriaforDegeneratedTriangle(TPS)) {
2606 AddTriangleLine(TPS[0], TPS[1], 0);
2607 AddTriangleLine(TPS[0], TPS[2], 1);
2608 AddTriangleLine(TPS[1], TPS[2], 2);
2609
2610 BTS = new class BoundaryTriangleSet(BLS, TrianglesOnBoundaryCount);
2611 AddTriangle();
2612
2613 (*it)->OtherOptCenter.Scale(-1.);
2614 BTS->GetNormalVector((*it)->OtherOptCenter);
2615 (*it)->OtherOptCenter.Scale(-1.);
2616
2617 cout << "--> WARNING: Special new triangle with " << *BTS << " and normal vector " << BTS->NormalVector
2618 << " for this triangle ... " << endl;
2619 cout << Verbose(1) << "We have "<< BaseRay->TrianglesCount << " for line " << BaseRay << "." << endl;
2620 } else {
2621 cout << Verbose(1) << "WARNING: This triangle consisting of ";
2622 cout << *(*it)->point << ", ";
2623 cout << *BaseRay->endpoints[0]->node << " and ";
2624 cout << *BaseRay->endpoints[1]->node << " ";
2625 cout << "exists and is not added, as it does not seem helpful!" << endl;
2626 result = false;
2627 }
2628 } else {
2629 cout << Verbose(1) << "This triangle consisting of ";
2630 cout << *(*it)->point << ", ";
2631 cout << *BaseRay->endpoints[0]->node << " and ";
2632 cout << *BaseRay->endpoints[1]->node << " ";
2633 cout << "is invalid!" << endl;
2634 result = false;
2635 }
2636
2637 if ((result) && (existentTrianglesCount < 2) && (DoSingleStepOutput && (TrianglesOnBoundaryCount % 1 == 0))) { // if we have a new triangle and want to output each new triangle configuration
2638 sprintf(NumberName, "-%04d-%s_%s_%s", TriangleFilesWritten, BTS->endpoints[0]->node->Name, BTS->endpoints[1]->node->Name, BTS->endpoints[2]->node->Name);
2639 if (DoTecplotOutput) {
2640 string NameofTempFile(tempbasename);
2641 NameofTempFile.append(NumberName);
2642 for(size_t npos = NameofTempFile.find_first_of(' '); npos != string::npos; npos = NameofTempFile.find(' ', npos))
2643 NameofTempFile.erase(npos, 1);
2644 NameofTempFile.append(TecplotSuffix);
2645 cout << Verbose(1) << "Writing temporary non convex hull to file " << NameofTempFile << ".\n";
2646 tempstream = new ofstream(NameofTempFile.c_str(), ios::trunc);
2647 write_tecplot_file(out, tempstream, this, mol, TriangleFilesWritten);
2648 tempstream->close();
2649 tempstream->flush();
2650 delete(tempstream);
2651 }
2652
2653 if (DoRaster3DOutput) {
2654 string NameofTempFile(tempbasename);
2655 NameofTempFile.append(NumberName);
2656 for(size_t npos = NameofTempFile.find_first_of(' '); npos != string::npos; npos = NameofTempFile.find(' ', npos))
2657 NameofTempFile.erase(npos, 1);
2658 NameofTempFile.append(Raster3DSuffix);
2659 cout << Verbose(1) << "Writing temporary non convex hull to file " << NameofTempFile << ".\n";
2660 tempstream = new ofstream(NameofTempFile.c_str(), ios::trunc);
2661 write_raster3d_file(out, tempstream, this, mol);
2662 // include the current position of the virtual sphere in the temporary raster3d file
2663 // make the circumsphere's center absolute again
2664 helper.CopyVector(&BaseRay->endpoints[0]->node->x);
2665 helper.AddVector(&BaseRay->endpoints[1]->node->x);
2666 helper.Scale(0.5);
2667 (*it)->OptCenter.AddVector(&helper);
2668 Vector *center = mol->DetermineCenterOfAll(out);
2669 (*it)->OptCenter.SubtractVector(center);
2670 delete(center);
2671 // and add to file plus translucency object
2672 *tempstream << "# current virtual sphere\n";
2673 *tempstream << "8\n 25.0 0.6 -1.0 -1.0 -1.0 0.2 0 0 0 0\n";
2674 *tempstream << "2\n " << (*it)->OptCenter.x[0] << " "
2675 << (*it)->OptCenter.x[1] << " " << (*it)->OptCenter.x[2]
2676 << "\t" << RADIUS << "\t1 0 0\n";
2677 *tempstream << "9\n terminating special property\n";
2678 tempstream->close();
2679 tempstream->flush();
2680 delete(tempstream);
2681 }
2682 if (DoTecplotOutput || DoRaster3DOutput)
2683 TriangleFilesWritten++;
2684 }
2685
2686 // set baseline to new ray from ref point (here endpoints[0]->node) to current candidate (here (*it)->point))
2687 BaseRay = BLS[0];
2688 }
2689
2690 // remove all candidates from the list and then the list itself
2691 class CandidateForTesselation *remover = NULL;
2692 for (CandidateList::iterator it = Opt_Candidates->begin(); it != Opt_Candidates->end(); ++it) {
2693 remover = *it;
2694 delete(remover);
2695 }
2696 delete(Opt_Candidates);
2697 cout << Verbose(0) << "End of Find_next_suitable_triangle\n";
2698 return result;
2699};
2700
2701/**
2702 * Sort function for the candidate list.
2703 */
2704bool sortCandidates(CandidateForTesselation* candidate1, CandidateForTesselation* candidate2) {
2705 Vector BaseLineVector, OrthogonalVector, helper;
2706 if (candidate1->BaseLine != candidate2->BaseLine) { // sanity check
2707 cout << Verbose(0) << "ERROR: sortCandidates was called for two different baselines: " << candidate1->BaseLine << " and " << candidate2->BaseLine << "." << endl;
2708 //return false;
2709 exit(1);
2710 }
2711 // create baseline vector
2712 BaseLineVector.CopyVector(&(candidate1->BaseLine->endpoints[1]->node->x));
2713 BaseLineVector.SubtractVector(&(candidate1->BaseLine->endpoints[0]->node->x));
2714 BaseLineVector.Normalize();
2715
2716 // 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!)
2717 helper.CopyVector(&(candidate1->BaseLine->endpoints[0]->node->x));
2718 helper.SubtractVector(&(candidate1->point->x));
2719 OrthogonalVector.CopyVector(&helper);
2720 helper.VectorProduct(&BaseLineVector);
2721 OrthogonalVector.SubtractVector(&helper);
2722 OrthogonalVector.Normalize();
2723
2724 // calculate both angles and correct with in-plane vector
2725 helper.CopyVector(&(candidate1->point->x));
2726 helper.SubtractVector(&(candidate1->BaseLine->endpoints[0]->node->x));
2727 double phi = BaseLineVector.Angle(&helper);
2728 if (OrthogonalVector.ScalarProduct(&helper) > 0) {
2729 phi = 2.*M_PI - phi;
2730 }
2731 helper.CopyVector(&(candidate2->point->x));
2732 helper.SubtractVector(&(candidate1->BaseLine->endpoints[0]->node->x));
2733 double psi = BaseLineVector.Angle(&helper);
2734 if (OrthogonalVector.ScalarProduct(&helper) > 0) {
2735 psi = 2.*M_PI - psi;
2736 }
2737
2738 cout << Verbose(2) << *candidate1->point << " has angle " << phi << endl;
2739 cout << Verbose(2) << *candidate2->point << " has angle " << psi << endl;
2740
2741 // return comparison
2742 return phi < psi;
2743}
2744
2745/** Tesselates the non convex boundary by rolling a virtual sphere along the surface of the molecule.
2746 * \param *out output stream for debugging
2747 * \param *mol molecule structure with Atom's and Bond's
2748 * \param *Tess Tesselation filled with points, lines and triangles on boundary on return
2749 * \param *LCList atoms in LinkedCell list
2750 * \param *filename filename prefix for output of vertex data
2751 * \para RADIUS radius of the virtual sphere
2752 */
2753void Find_non_convex_border(ofstream *out, molecule* mol, class Tesselation *Tess, class LinkedCell *LCList, const char *filename, const double RADIUS)
2754{
2755 int N = 0;
2756 bool freeTess = false;
2757 bool freeLC = false;
2758 *out << Verbose(1) << "Entering search for non convex hull. " << endl;
2759 if (Tess == NULL) {
2760 *out << Verbose(1) << "Allocating Tesselation struct ..." << endl;
2761 Tess = new Tesselation;
2762 freeTess = true;
2763 }
2764 LineMap::iterator baseline;
2765 LineMap::iterator testline;
2766 *out << Verbose(0) << "Begin of Find_non_convex_border\n";
2767 bool flag = false; // marks whether we went once through all baselines without finding any without two triangles
2768 bool failflag = false;
2769
2770 if (LCList == NULL) {
2771 LCList = new LinkedCell(mol, 2.*RADIUS);
2772 freeLC = true;
2773 }
2774
2775 Tess->Find_starting_triangle(out, mol, RADIUS, LCList);
2776
2777 baseline = Tess->LinesOnBoundary.begin();
2778 while ((baseline != Tess->LinesOnBoundary.end()) || (flag)) {
2779 if (baseline->second->TrianglesCount == 1) {
2780 failflag = Tess->Find_next_suitable_triangle(out, mol, *(baseline->second), *(((baseline->second->triangles.begin()))->second), RADIUS, N, filename, LCList); //the line is there, so there is a triangle, but only one.
2781 flag = flag || failflag;
2782 if (!failflag)
2783 cerr << "WARNING: Find_next_suitable_triangle failed." << endl;
2784 } else {
2785 //cout << Verbose(1) << "Line " << *baseline->second << " has " << baseline->second->TrianglesCount << " triangles adjacent" << endl;
2786 if (baseline->second->TrianglesCount != 2)
2787 cout << Verbose(1) << "ERROR: TESSELATION FINISHED WITH INVALID TRIANGLE COUNT!" << endl;
2788 }
2789
2790 N++;
2791 baseline++;
2792 if ((baseline == Tess->LinesOnBoundary.end()) && (flag)) {
2793 baseline = Tess->LinesOnBoundary.begin(); // restart if we reach end due to newly inserted lines
2794 flag = false;
2795 }
2796 }
2797 if (1) { //failflag) {
2798 *out << Verbose(1) << "Writing final tecplot file\n";
2799 if (DoTecplotOutput) {
2800 string OutputName(filename);
2801 OutputName.append(TecplotSuffix);
2802 ofstream *tecplot = new ofstream(OutputName.c_str());
2803 write_tecplot_file(out, tecplot, Tess, mol, -1);
2804 tecplot->close();
2805 delete(tecplot);
2806 }
2807 if (DoRaster3DOutput) {
2808 string OutputName(filename);
2809 OutputName.append(Raster3DSuffix);
2810 ofstream *raster = new ofstream(OutputName.c_str());
2811 write_raster3d_file(out, raster, Tess, mol);
2812 raster->close();
2813 delete(raster);
2814 }
2815 } else {
2816 cerr << "ERROR: Could definitively not find all necessary triangles!" << endl;
2817 }
2818
2819 cout << Verbose(2) << "Check: List of Baselines with not two connected triangles:" << endl;
2820 int counter = 0;
2821 for (testline = Tess->LinesOnBoundary.begin(); testline != Tess->LinesOnBoundary.end(); testline++) {
2822 if (testline->second->TrianglesCount != 2) {
2823 cout << Verbose(2) << *testline->second << "\t" << testline->second->TrianglesCount << endl;
2824 counter++;
2825 }
2826 }
2827 if (counter == 0)
2828 cout << Verbose(2) << "None." << endl;
2829
2830 if (freeTess)
2831 delete(Tess);
2832 if (freeLC)
2833 delete(LCList);
2834 *out << Verbose(0) << "End of Find_non_convex_border\n";
2835};
2836
2837/** Finds a hole of sufficient size in \a this molecule to embed \a *srcmol into it.
2838 * \param *out output stream for debugging
2839 * \param *srcmol molecule to embed into
2840 * \return *Vector new center of \a *srcmol for embedding relative to \a this
2841 */
2842Vector* molecule::FindEmbeddingHole(ofstream *out, molecule *srcmol)
2843{
2844 Vector *Center = new Vector;
2845 Center->Zero();
2846 // calculate volume/shape of \a *srcmol
2847
2848 // find embedding holes
2849
2850 // if more than one, let user choose
2851
2852 // return embedding center
2853 return Center;
2854};
2855
Note: See TracBrowser for help on using the repository browser.