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