source: src/Tesselation/boundary.cpp@ 9b56f0a

Action_Thermostats Add_AtomRandomPerturbation Add_FitFragmentPartialChargesAction Add_RotateAroundBondAction Add_SelectAtomByNameAction Added_ParseSaveFragmentResults AddingActions_SaveParseParticleParameters Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_ParticleName_to_Atom Adding_StructOpt_integration_tests AtomFragments Automaking_mpqc_open AutomationFragmentation_failures Candidate_v1.5.4 Candidate_v1.6.0 Candidate_v1.6.1 Candidate_v1.7.0 ChangeBugEmailaddress ChangingTestPorts ChemicalSpaceEvaluator CombiningParticlePotentialParsing Combining_Subpackages Debian_Package_split Debian_package_split_molecuildergui_only Disabling_MemDebug Docu_Python_wait EmpiricalPotential_contain_HomologyGraph EmpiricalPotential_contain_HomologyGraph_documentation Enable_parallel_make_install Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph FitPartialCharges_GlobalError Fix_BoundInBox_CenterInBox_MoleculeActions Fix_ChargeSampling_PBC Fix_ChronosMutex Fix_FitPartialCharges Fix_FitPotential_needs_atomicnumbers Fix_ForceAnnealing Fix_IndependentFragmentGrids Fix_ParseParticles Fix_ParseParticles_split_forward_backward_Actions Fix_PopActions Fix_QtFragmentList_sorted_selection Fix_Restrictedkeyset_FragmentMolecule Fix_StatusMsg Fix_StepWorldTime_single_argument Fix_Verbose_Codepatterns Fix_fitting_potentials Fixes ForceAnnealing_goodresults ForceAnnealing_oldresults ForceAnnealing_tocheck ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion FragmentAction_writes_AtomFragments FragmentMolecule_checks_bonddegrees GeometryObjects Gui_Fixes Gui_displays_atomic_force_velocity ImplicitCharges IndependentFragmentGrids IndependentFragmentGrids_IndividualZeroInstances IndependentFragmentGrids_IntegrationTest IndependentFragmentGrids_Sole_NN_Calculation JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix MoreRobust_FragmentAutomation ODR_violation_mpqc_open PartialCharges_OrthogonalSummation PdbParser_setsAtomName PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks Rewrite_FitPartialCharges RotateToPrincipalAxisSystem_UndoRedo SaturateAtoms_findBestMatching SaturateAtoms_singleDegree StoppableMakroAction Subpackage_CodePatterns Subpackage_JobMarket Subpackage_LinearAlgebra Subpackage_levmar Subpackage_mpqc_open Subpackage_vmg Switchable_LogView ThirdParty_MPQC_rebuilt_buildsystem TrajectoryDependenant_MaxOrder TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps TremoloParser_setsAtomName Ubuntu_1604_changes stable
Last change on this file since 9b56f0a was f01769, checked in by Frederik Heber <heber@…>, 10 years ago

Replaced World::getAtom() wherever possible by const version.

  • some AtomSet member functions now have const atom ptr instead of atom ptr.
  • molecule can return const and non-const AtomSet.
  • added FromIdToConstAtom to allow iterate through atoms in molecule (which are stored by id, not by ptr) in const fashion.
  • in molecule::isInMolecule() is now const, ::CopyMolecule..() is non-const (because copying involves father atom who is stored non-const).
  • Property mode set to 100644
File size: 40.5 KB
RevLine 
[bcf653]1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
[0aa122]4 * Copyright (C) 2010-2012 University of Bonn. All rights reserved.
[5aaa43]5 * Copyright (C) 2013 Frederik Heber. All rights reserved.
[94d5ac6]6 *
7 *
8 * This file is part of MoleCuilder.
9 *
10 * MoleCuilder is free software: you can redistribute it and/or modify
11 * it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation, either version 2 of the License, or
13 * (at your option) any later version.
14 *
15 * MoleCuilder is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with MoleCuilder. If not, see <http://www.gnu.org/licenses/>.
[bcf653]22 */
23
[f66195]24/** \file boundary.cpp
[edb93c]25 *
26 * Implementations and super-function for envelopes
[2319ed]27 */
28
[bf3817]29// include config.h
30#ifdef HAVE_CONFIG_H
31#include <config.h>
32#endif
33
[ad011c]34#include "CodePatterns/MemDebug.hpp"
[112b09]35
[6f0841]36#include "Atom/atom.hpp"
[129204]37#include "Bond/bond.hpp"
[8eb17a]38#include "boundary.hpp"
[34c43a]39#include "BoundaryLineSet.hpp"
40#include "BoundaryPointSet.hpp"
41#include "BoundaryTriangleSet.hpp"
42#include "Box.hpp"
43#include "CandidateForTesselation.hpp"
44#include "CodePatterns/Info.hpp"
45#include "CodePatterns/Log.hpp"
46#include "CodePatterns/Verbose.hpp"
[f66195]47#include "config.hpp"
[3bdb6d]48#include "Element/element.hpp"
[34c43a]49#include "LinearAlgebra/Plane.hpp"
50#include "LinearAlgebra/RealSpaceMatrix.hpp"
[53c7fc]51#include "LinkedCell/linkedcell.hpp"
52#include "LinkedCell/PointCloudAdaptor.hpp"
[f66195]53#include "molecule.hpp"
[42127c]54#include "MoleculeListClass.hpp"
[34c43a]55#include "RandomNumbers/RandomNumberGeneratorFactory.hpp"
56#include "RandomNumbers/RandomNumberGenerator.hpp"
[f66195]57#include "tesselation.hpp"
58#include "tesselationhelpers.hpp"
[b34306]59#include "World.hpp"
[2319ed]60
[986ed3]61#include <iostream>
[36166d]62#include <iomanip>
63
[357fba]64#include<gsl/gsl_poly.h>
[d6eb80]65#include<time.h>
[2319ed]66
[357fba]67// ========================================== F U N C T I O N S =================================
[2319ed]68
69
[357fba]70/** Determines greatest diameters of a cluster defined by its convex envelope.
71 * Looks at lines parallel to one axis and where they intersect on the projected planes
[2319ed]72 * \param *out output stream for debugging
[357fba]73 * \param *BoundaryPoints NDIM set of boundary points defining the convex envelope on each projected plane
74 * \param *mol molecule structure representing the cluster
[776b64]75 * \param *&TesselStruct Tesselation structure with triangles
[357fba]76 * \param IsAngstroem whether we have angstroem or atomic units
77 * \return NDIM array of the diameters
[e4ea46]78 */
[e138de]79double *GetDiametersOfCluster(const Boundaries *BoundaryPtr, const molecule *mol, Tesselation *&TesselStruct, const bool IsAngstroem)
[caf5d6]80{
[ce7bfd]81 //Info FunctionInfo(__func__);
[357fba]82 // get points on boundary of NULL was given as parameter
83 bool BoundaryFreeFlag = false;
[ad37ab]84 double OldComponent = 0.;
85 double tmp = 0.;
86 double w1 = 0.;
87 double w2 = 0.;
88 Vector DistanceVector;
89 Vector OtherVector;
90 int component = 0;
91 int Othercomponent = 0;
[776b64]92 Boundaries::const_iterator Neighbour;
93 Boundaries::const_iterator OtherNeighbour;
[ad37ab]94 double *GreatestDiameter = new double[NDIM];
95
[776b64]96 const Boundaries *BoundaryPoints;
97 if (BoundaryPtr == NULL) {
[357fba]98 BoundaryFreeFlag = true;
[e138de]99 BoundaryPoints = GetBoundaryPoints(mol, TesselStruct);
[86234b]100 } else {
[776b64]101 BoundaryPoints = BoundaryPtr;
[47d041]102 LOG(0, "Using given boundary points set.");
[86234b]103 }
[357fba]104 // determine biggest "diameter" of cluster for each axis
105 for (int i = 0; i < NDIM; i++)
106 GreatestDiameter[i] = 0.;
107 for (int axis = 0; axis < NDIM; axis++)
108 { // regard each projected plane
[47d041]109 //LOG(1, "Current axis is " << axis << ".");
[357fba]110 for (int j = 0; j < 2; j++)
111 { // and for both axis on the current plane
112 component = (axis + j + 1) % NDIM;
113 Othercomponent = (axis + 1 + ((j + 1) & 1)) % NDIM;
[47d041]114 //LOG(1, "Current component is " << component << ", Othercomponent is " << Othercomponent << ".");
[776b64]115 for (Boundaries::const_iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
[47d041]116 //LOG(1, "Current runner is " << *(runner->second.second) << ".");
[357fba]117 // seek for the neighbours pair where the Othercomponent sign flips
118 Neighbour = runner;
119 Neighbour++;
120 if (Neighbour == BoundaryPoints[axis].end()) // make it wrap around
121 Neighbour = BoundaryPoints[axis].begin();
[d74077]122 DistanceVector = (runner->second.second->getPosition()) - (Neighbour->second.second->getPosition());
[776b64]123 do { // seek for neighbour pair where it flips
[0a4f7f]124 OldComponent = DistanceVector[Othercomponent];
[357fba]125 Neighbour++;
126 if (Neighbour == BoundaryPoints[axis].end()) // make it wrap around
127 Neighbour = BoundaryPoints[axis].begin();
[d74077]128 DistanceVector = (runner->second.second->getPosition()) - (Neighbour->second.second->getPosition());
[47d041]129 //LOG(2, "OldComponent is " << OldComponent << ", new one is " << DistanceVector.x[Othercomponent] << ".");
[776b64]130 } while ((runner != Neighbour) && (fabs(OldComponent / fabs(
[0a4f7f]131 OldComponent) - DistanceVector[Othercomponent] / fabs(
132 DistanceVector[Othercomponent])) < MYEPSILON)); // as long as sign does not flip
[776b64]133 if (runner != Neighbour) {
[357fba]134 OtherNeighbour = Neighbour;
135 if (OtherNeighbour == BoundaryPoints[axis].begin()) // make it wrap around
136 OtherNeighbour = BoundaryPoints[axis].end();
137 OtherNeighbour--;
[47d041]138 //LOG(1, "The pair, where the sign of OtherComponent flips, is: " << *(Neighbour->second.second) << " and " << *(OtherNeighbour->second.second) << ".");
[357fba]139 // now we have found the pair: Neighbour and OtherNeighbour
[d74077]140 OtherVector = (runner->second.second->getPosition()) - (OtherNeighbour->second.second->getPosition());
[47d041]141 //LOG(1, "Distances to Neighbour and OtherNeighbour are " << DistanceVector.x[component] << " and " << OtherVector.x[component] << ".");
142 //LOG(1, "OtherComponents to Neighbour and OtherNeighbour are " << DistanceVector.x[Othercomponent] << " and " << OtherVector.x[Othercomponent] << ".");
[357fba]143 // do linear interpolation between points (is exact) to extract exact intersection between Neighbour and OtherNeighbour
[0a4f7f]144 w1 = fabs(OtherVector[Othercomponent]);
145 w2 = fabs(DistanceVector[Othercomponent]);
146 tmp = fabs((w1 * DistanceVector[component] + w2
147 * OtherVector[component]) / (w1 + w2));
[357fba]148 // mark if it has greater diameter
[47d041]149 //LOG(1, "Comparing current greatest " << GreatestDiameter[component] << " to new " << tmp << ".");
[357fba]150 GreatestDiameter[component] = (GreatestDiameter[component]
151 > tmp) ? GreatestDiameter[component] : tmp;
152 } //else
[47d041]153 //LOG(1, "Saw no sign flip, probably top or bottom node.");
[3d919e]154 }
155 }
156 }
[47d041]157 LOG(0, "RESULT: The biggest diameters are "
[357fba]158 << GreatestDiameter[0] << " and " << GreatestDiameter[1] << " and "
159 << GreatestDiameter[2] << " " << (IsAngstroem ? "angstrom"
[47d041]160 : "atomiclength") << ".");
[03648b]161
[357fba]162 // free reference lists
163 if (BoundaryFreeFlag)
164 delete[] (BoundaryPoints);
[e4ea46]165
[357fba]166 return GreatestDiameter;
[e4ea46]167}
168;
[03648b]169
[042f82]170
[357fba]171/** Determines the boundary points of a cluster.
172 * Does a projection per axis onto the orthogonal plane, transforms into spherical coordinates, sorts them by the angle
173 * and looks at triples: if the middle has less a distance than the allowed maximum height of the triangle formed by the plane's
174 * center and first and last point in the triple, it is thrown out.
[f01769]175 *
176 * \todo When storing const ptrs in tesselation structures, remove const_cast
177 *
[357fba]178 * \param *out output stream for debugging
179 * \param *mol molecule structure representing the cluster
[776b64]180 * \param *&TesselStruct pointer to Tesselation structure
[e4ea46]181 */
[e138de]182Boundaries *GetBoundaryPoints(const molecule *mol, Tesselation *&TesselStruct)
[caf5d6]183{
[ce7bfd]184 //Info FunctionInfo(__func__);
[357fba]185 PointMap PointsOnBoundary;
186 LineMap LinesOnBoundary;
187 TriangleMap TrianglesOnBoundary;
[833b15]188 Vector MolCenter = mol->DetermineCenterOfAll();
[357fba]189 Vector helper;
[ad37ab]190 BoundariesTestPair BoundaryTestPair;
191 Vector AxisVector;
192 Vector AngleReferenceVector;
193 Vector AngleReferenceNormalVector;
194 Vector ProjectedVector;
[5309ba]195 Boundaries *BoundaryPoints = new Boundaries[NDIM]; // first is alpha, second is (r, Nr)
[ad37ab]196 double angle = 0.;
[042f82]197
[357fba]198 // 3a. Go through every axis
199 for (int axis = 0; axis < NDIM; axis++) {
200 AxisVector.Zero();
201 AngleReferenceVector.Zero();
202 AngleReferenceNormalVector.Zero();
[0a4f7f]203 AxisVector[axis] = 1.;
204 AngleReferenceVector[(axis + 1) % NDIM] = 1.;
205 AngleReferenceNormalVector[(axis + 2) % NDIM] = 1.;
[042f82]206
[47d041]207 LOG(1, "Axisvector is " << AxisVector << " and AngleReferenceVector is " << AngleReferenceVector << ", and AngleReferenceNormalVector is " << AngleReferenceNormalVector << ".");
[042f82]208
[357fba]209 // 3b. construct set of all points, transformed into cylindrical system and with left and right neighbours
[59fff1]210 // Boundaries stores non-const TesselPoint ref, hence we need iterator here
[f01769]211 for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
[833b15]212 ProjectedVector = (*iter)->getPosition() - (MolCenter);
[273382]213 ProjectedVector.ProjectOntoPlane(AxisVector);
[042f82]214
[357fba]215 // correct for negative side
[776b64]216 const double radius = ProjectedVector.NormSquared();
[357fba]217 if (fabs(radius) > MYEPSILON)
[273382]218 angle = ProjectedVector.Angle(AngleReferenceVector);
[357fba]219 else
220 angle = 0.; // otherwise it's a vector in Axis Direction and unimportant for boundary issues
[042f82]221
[47d041]222 //LOG(1, "Checking sign in quadrant : " << ProjectedVector.Projection(&AngleReferenceNormalVector) << ".");
[273382]223 if (ProjectedVector.ScalarProduct(AngleReferenceNormalVector) > 0) {
[357fba]224 angle = 2. * M_PI - angle;
225 }
[47d041]226 LOG(1, "Inserting " << **iter << ": (r, alpha) = (" << radius << "," << angle << "): " << ProjectedVector);
[f01769]227 BoundaryTestPair = BoundaryPoints[axis].insert(
228 BoundariesPair(angle, TesselPointDistancePair (radius, const_cast<atom *>(*iter))));
[357fba]229 if (!BoundaryTestPair.second) { // same point exists, check first r, then distance of original vectors to center of gravity
[47d041]230 LOG(2, "Encountered two vectors whose projection onto axis " << axis << " is equal: ");
231 LOG(2, "Present vector: " << *BoundaryTestPair.first->second.second);
232 LOG(2, "New vector: " << **iter);
[776b64]233 const double ProjectedVectorNorm = ProjectedVector.NormSquared();
234 if ((ProjectedVectorNorm - BoundaryTestPair.first->second.first) > MYEPSILON) {
235 BoundaryTestPair.first->second.first = ProjectedVectorNorm;
[f01769]236 BoundaryTestPair.first->second.second = const_cast<atom *>(*iter);
[47d041]237 LOG(2, "Keeping new vector due to larger projected distance " << ProjectedVectorNorm << ".");
[776b64]238 } else if (fabs(ProjectedVectorNorm - BoundaryTestPair.first->second.first) < MYEPSILON) {
[833b15]239 helper = (*iter)->getPosition() - (MolCenter);
[776b64]240 const double oldhelperNorm = helper.NormSquared();
[833b15]241 helper = BoundaryTestPair.first->second.second->getPosition() - (MolCenter);
[776b64]242 if (helper.NormSquared() < oldhelperNorm) {
[f01769]243 BoundaryTestPair.first->second.second = const_cast<atom *>(*iter);
[47d041]244 LOG(2, "Keeping new vector due to larger distance to molecule center " << helper.NormSquared() << ".");
[357fba]245 } else {
[47d041]246 LOG(2, "Keeping present vector due to larger distance to molecule center " << oldhelperNorm << ".");
[357fba]247 }
248 } else {
[47d041]249 LOG(2, "Keeping present vector due to larger projected distance " << ProjectedVectorNorm << ".");
[357fba]250 }
[018741]251 }
[3d919e]252 }
[357fba]253 // printing all inserted for debugging
254 // {
[47d041]255 // std::stringstream output;
256 // output << "Printing list of candidates for axis " << axis << " which we have inserted so far: ";
[357fba]257 // int i=0;
258 // for(Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
259 // if (runner != BoundaryPoints[axis].begin())
[47d041]260 // output << ", " << i << ": " << *runner->second.second;
[357fba]261 // else
[47d041]262 // output << i << ": " << *runner->second.second;
[357fba]263 // i++;
264 // }
[47d041]265 // LOG(1, output.str());
[357fba]266 // }
267 // 3c. throw out points whose distance is less than the mean of left and right neighbours
268 bool flag = false;
[47d041]269 LOG(1, "Looking for candidates to kick out by convex condition ... ");
[357fba]270 do { // do as long as we still throw one out per round
271 flag = false;
[accebe]272 Boundaries::iterator left = BoundaryPoints[axis].begin();
273 Boundaries::iterator right = BoundaryPoints[axis].begin();
274 Boundaries::iterator runner = BoundaryPoints[axis].begin();
275 bool LoopOnceDone = false;
276 while (!LoopOnceDone) {
277 runner = right;
278 right++;
[357fba]279 // set neighbours correctly
280 if (runner == BoundaryPoints[axis].begin()) {
281 left = BoundaryPoints[axis].end();
282 } else {
283 left = runner;
284 }
285 left--;
286 if (right == BoundaryPoints[axis].end()) {
287 right = BoundaryPoints[axis].begin();
[accebe]288 LoopOnceDone = true;
[357fba]289 }
290 // check distance
[3d919e]291
[357fba]292 // construct the vector of each side of the triangle on the projected plane (defined by normal vector AxisVector)
293 {
294 Vector SideA, SideB, SideC, SideH;
[833b15]295 SideA = left->second.second->getPosition() - (MolCenter);
[273382]296 SideA.ProjectOntoPlane(AxisVector);
[47d041]297 // LOG(1, "SideA: " << SideA);
[3d919e]298
[833b15]299 SideB = right->second.second->getPosition() -(MolCenter);
[273382]300 SideB.ProjectOntoPlane(AxisVector);
[47d041]301 // LOG(1, "SideB: " << SideB);
[3d919e]302
[d74077]303 SideC = left->second.second->getPosition() - right->second.second->getPosition();
[273382]304 SideC.ProjectOntoPlane(AxisVector);
[47d041]305 // LOG(1, "SideC: " << SideC);
[3d919e]306
[833b15]307 SideH = runner->second.second->getPosition() -(MolCenter);
[273382]308 SideH.ProjectOntoPlane(AxisVector);
[47d041]309 // LOG(1, "SideH: " << SideH);
[3d919e]310
[357fba]311 // calculate each length
[ad37ab]312 const double a = SideA.Norm();
313 //const double b = SideB.Norm();
314 //const double c = SideC.Norm();
315 const double h = SideH.Norm();
[357fba]316 // calculate the angles
[273382]317 const double alpha = SideA.Angle(SideH);
318 const double beta = SideA.Angle(SideC);
319 const double gamma = SideB.Angle(SideH);
320 const double delta = SideC.Angle(SideH);
[ad37ab]321 const double MinDistance = a * sin(beta) / (sin(delta)) * (((alpha < M_PI / 2.) || (gamma < M_PI / 2.)) ? 1. : -1.);
[47d041]322 //LOG(1, " 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 << ".");
323 LOG(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 << ".");
[357fba]324 if ((fabs(h / fabs(h) - MinDistance / fabs(MinDistance)) < MYEPSILON) && ((h - MinDistance)) < -MYEPSILON) {
325 // throw out point
[47d041]326 LOG(1, "Throwing out " << *runner->second.second << ".");
[357fba]327 BoundaryPoints[axis].erase(runner);
[accebe]328 runner = right;
[357fba]329 flag = true;
[3d919e]330 }
331 }
332 }
[357fba]333 } while (flag);
[3d919e]334 }
[357fba]335 return BoundaryPoints;
[6ac7ee]336};
337
[357fba]338/** Tesselates the convex boundary by finding all boundary points.
339 * \param *out output stream for debugging
[776b64]340 * \param *mol molecule structure with Atom's and Bond's.
[bdc91e]341 * \param *BoundaryPts set of boundary points to use or NULL
[357fba]342 * \param *TesselStruct Tesselation filled with points, lines and triangles on boundary on return
[6bd7e0]343 * \param *LCList atoms in LinkedCell_deprecated list
[357fba]344 * \param *filename filename prefix for output of vertex data
345 * \return *TesselStruct is filled with convex boundary and tesselation is stored under \a *filename.
[6ac7ee]346 */
[6bd7e0]347void FindConvexBorder(const molecule* mol, Boundaries *BoundaryPts, Tesselation *&TesselStruct, const LinkedCell_deprecated *LCList, const char *filename)
[6ac7ee]348{
[ce7bfd]349 //Info FunctionInfo(__func__);
[357fba]350 bool BoundaryFreeFlag = false;
351 Boundaries *BoundaryPoints = NULL;
[3d919e]352
[776b64]353 if (TesselStruct != NULL) // free if allocated
354 delete(TesselStruct);
355 TesselStruct = new class Tesselation;
[3d919e]356
[357fba]357 // 1. Find all points on the boundary
[bdc91e]358 if (BoundaryPts == NULL) {
359 BoundaryFreeFlag = true;
360 BoundaryPoints = GetBoundaryPoints(mol, TesselStruct);
[357fba]361 } else {
[bdc91e]362 BoundaryPoints = BoundaryPts;
[47d041]363 LOG(0, "Using given boundary points set.");
[3d919e]364 }
365
[357fba]366// printing all inserted for debugging
[47d041]367 if (DoLog(1)) {
368 for (int axis=0; axis < NDIM; axis++) {
369 std::stringstream output;
370 output << "Printing list of candidates for axis " << axis << " which we have inserted so far: ";
371 int i=0;
372 for(Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
373 if (runner != BoundaryPoints[axis].begin())
374 output << ", " << i << ": " << *runner->second.second;
375 else
376 output << i << ": " << *runner->second.second;
377 i++;
378 }
379 LOG(1, output.str());
[a37350]380 }
[bdc91e]381 }
[3d919e]382
[357fba]383 // 2. fill the boundary point list
384 for (int axis = 0; axis < NDIM; axis++)
385 for (Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++)
[776b64]386 if (!TesselStruct->AddBoundaryPoint(runner->second.second, 0))
[47d041]387 LOG(2, "Point " << *(runner->second.second) << " is already present.");
[e4ea46]388
[47d041]389 LOG(0, "I found " << TesselStruct->PointsOnBoundaryCount << " points on the convex boundary.");
[357fba]390 // now we have the whole set of edge points in the BoundaryList
[018741]391
[357fba]392 // listing for debugging
[47d041]393 //if (DoLog(1)) {
394 // std::stringstream output;
395 // output << "Listing PointsOnBoundary:";
[357fba]396 // for(PointMap::iterator runner = PointsOnBoundary.begin(); runner != PointsOnBoundary.end(); runner++) {
[47d041]397 // output << " " << *runner->second;
[357fba]398 // }
[47d041]399 // LOG(1, output.str());
400 //}
[018741]401
[357fba]402 // 3a. guess starting triangle
[e138de]403 TesselStruct->GuessStartingTriangle();
[018741]404
[357fba]405 // 3b. go through all lines, that are not yet part of two triangles (only of one so far)
[caa06ef]406 PointCloudAdaptor< molecule > cloud(const_cast<molecule *>(mol), mol->name);
[34c43a]407 TesselStruct->TesselateOnBoundary(cloud);
[3d919e]408
[357fba]409 // 3c. check whether all atoms lay inside the boundary, if not, add to boundary points, segment triangle into three with the new point
[34c43a]410 if (!TesselStruct->InsertStraddlingPoints(cloud, LCList))
[47d041]411 ELOG(1, "Insertion of straddling points failed!");
[3d919e]412
[47d041]413 LOG(0, "I created " << TesselStruct->TrianglesOnBoundary.size() << " intermediate triangles with " << TesselStruct->LinesOnBoundary.size() << " lines and " << TesselStruct->PointsOnBoundary.size() << " points.");
[ef0e6d]414
415 // 4. Store triangles in tecplot file
[d51975]416 StoreTrianglesinFile(mol, TesselStruct, filename, "_intermed");
[ef0e6d]417
[357fba]418 // 3d. check all baselines whether the peaks of the two adjacent triangles with respect to center of baseline are convex, if not, make the baseline between the two peaks and baseline endpoints become the new peaks
[ad37ab]419 bool AllConvex = true;
[093645]420 class BoundaryLineSet *line = NULL;
421 do {
422 AllConvex = true;
[776b64]423 for (LineMap::iterator LineRunner = TesselStruct->LinesOnBoundary.begin(); LineRunner != TesselStruct->LinesOnBoundary.end(); LineRunner++) {
[093645]424 line = LineRunner->second;
[47d041]425 LOG(1, "INFO: Current line is " << *line << ".");
[e138de]426 if (!line->CheckConvexityCriterion()) {
[47d041]427 LOG(1, "... line " << *line << " is concave, flipping it.");
[093645]428
429 // flip the line
[e138de]430 if (TesselStruct->PickFarthestofTwoBaselines(line) == 0.)
[47d041]431 ELOG(1, "Correction of concave baselines failed!");
[57066a]432 else {
[e138de]433 TesselStruct->FlipBaseline(line);
[47d041]434 LOG(1, "INFO: Correction of concave baselines worked.");
[accebe]435 LineRunner = TesselStruct->LinesOnBoundary.begin(); // LineRunner may have been erase if line was deleted from LinesOnBoundary
[57066a]436 }
[093645]437 }
438 }
439 } while (!AllConvex);
[3d919e]440
[ef0e6d]441 // 3e. we need another correction here, for TesselPoints that are below the surface (i.e. have an odd number of concave triangles surrounding it)
[776b64]442// if (!TesselStruct->CorrectConcaveTesselPoints(out))
[47d041]443// ELOG(1, "Correction of concave tesselpoints failed!");
[ef0e6d]444
[47d041]445 LOG(0, "I created " << TesselStruct->TrianglesOnBoundary.size() << " triangles with " << TesselStruct->LinesOnBoundary.size() << " lines and " << TesselStruct->PointsOnBoundary.size() << " points.");
[3d919e]446
[357fba]447 // 4. Store triangles in tecplot file
[d51975]448 StoreTrianglesinFile(mol, TesselStruct, filename, "");
[ef0e6d]449
[357fba]450 // free reference lists
451 if (BoundaryFreeFlag)
452 delete[] (BoundaryPoints);
[3d919e]453};
[6ac7ee]454
[d4fa23]455/** For testing removes one boundary point after another to check for leaks.
456 * \param *out output stream for debugging
457 * \param *TesselStruct Tesselation containing envelope with boundary points
458 * \param *mol molecule
459 * \param *filename name of file
460 * \return true - all removed, false - something went wrong
461 */
[e138de]462bool RemoveAllBoundaryPoints(class Tesselation *&TesselStruct, const molecule * const mol, const char * const filename)
[d4fa23]463{
[ce7bfd]464 //Info FunctionInfo(__func__);
[d4fa23]465 int i=0;
466 char number[MAXSTRINGSIZE];
467
468 if ((TesselStruct == NULL) || (TesselStruct->PointsOnBoundary.empty())) {
[47d041]469 ELOG(1, "TesselStruct is empty.");
[d4fa23]470 return false;
471 }
472
473 PointMap::iterator PointRunner;
474 while (!TesselStruct->PointsOnBoundary.empty()) {
[47d041]475 if (DoLog(1)) {
476 std::stringstream output;
477 output << "Remaining points are: ";
478 for (PointMap::iterator PointSprinter = TesselStruct->PointsOnBoundary.begin(); PointSprinter != TesselStruct->PointsOnBoundary.end(); PointSprinter++)
479 output << *(PointSprinter->second) << "\t";
480 LOG(1, output.str());
481 }
[d4fa23]482
483 PointRunner = TesselStruct->PointsOnBoundary.begin();
484 // remove point
[e138de]485 TesselStruct->RemovePointFromTesselatedSurface(PointRunner->second);
[d4fa23]486
487 // store envelope
488 sprintf(number, "-%04d", i++);
[e138de]489 StoreTrianglesinFile(mol, (const Tesselation *&)TesselStruct, filename, number);
[d4fa23]490 }
491
492 return true;
493};
494
[08ef35]495/** Creates a convex envelope from a given non-convex one.
[093645]496 * -# First step, remove concave spots, i.e. singular "dents"
497 * -# We go through all PointsOnBoundary.
498 * -# We CheckConvexityCriterion() for all its lines.
499 * -# If all its lines are concave, it cannot be on the convex envelope.
500 * -# Hence, we remove it and re-create all its triangles from its getCircleOfConnectedPoints()
501 * -# We calculate the additional volume.
502 * -# We go over all lines until none yields a concavity anymore.
503 * -# Second step, remove concave lines, i.e. line-shape "dents"
504 * -# We go through all LinesOnBoundary
505 * -# We CheckConvexityCriterion()
506 * -# If it returns concave, we flip the line in this quadrupel of points (abusing the degeneracy of the tesselation)
507 * -# We CheckConvexityCriterion(),
508 * -# if it's concave, we continue
509 * -# if not, we mark an error and stop
[08ef35]510 * Note: This routine - for free - calculates the difference in volume between convex and
[ee0032]511 * non-convex envelope, as the former is easy to calculate - Tesselation::getVolumeOfConvexEnvelope() - it
[08ef35]512 * can be used to compute volumes of arbitrary shapes.
513 * \param *out output stream for debugging
514 * \param *TesselStruct non-convex envelope, is changed in return!
[093645]515 * \param *mol molecule
516 * \param *filename name of file
[08ef35]517 * \return volume difference between the non- and the created convex envelope
518 */
[5f7b95]519double ConvexizeNonconvexEnvelope(
520 Tesselation *&TesselStruct,
521 const molecule * const mol,
522 const char * const filename,
523 bool DebugOutputEveryStep)
[08ef35]524{
[ce7bfd]525 //Info FunctionInfo(__func__);
[08ef35]526 double volume = 0;
527 class BoundaryPointSet *point = NULL;
528 class BoundaryLineSet *line = NULL;
[ad37ab]529 bool Concavity = false;
[57066a]530 char dummy[MAXSTRINGSIZE];
[ad37ab]531 PointMap::iterator PointRunner;
532 PointMap::iterator PointAdvance;
533 LineMap::iterator LineRunner;
534 LineMap::iterator LineAdvance;
535 TriangleMap::iterator TriangleRunner;
536 TriangleMap::iterator TriangleAdvance;
537 int run = 0;
[093645]538
539 // check whether there is something to work on
[08ef35]540 if (TesselStruct == NULL) {
[47d041]541 ELOG(1, "TesselStruct is empty!");
[08ef35]542 return volume;
543 }
544
[b8d215]545 LOG(1, "INFO: Making tesselated surface with " << TesselStruct->TrianglesOnBoundaryCount
546 << " convex ...");
547
[2ef2cc]548 // first purge all degenerate triangles
549 TesselStruct->RemoveDegeneratedTriangles();
550
[1d9b7aa]551 do {
552 Concavity = false;
[d289c6]553
[5f7b95]554 if (DebugOutputEveryStep) {
555 sprintf(dummy, "-%d", run++);
556 //CalculateConcavityPerBoundaryPoint(TesselStruct);
557 LOG(1, "INFO: Writing " << run << "th tesselation file.");
558 StoreTrianglesinFile(mol, (const Tesselation *&)TesselStruct, filename, dummy);
559 }
[57066a]560
[d289c6]561 // first step: remove all full-concave point
[1d9b7aa]562 PointRunner = TesselStruct->PointsOnBoundary.begin();
563 PointAdvance = PointRunner; // we need an advanced point, as the PointRunner might get removed
564 while (PointRunner != TesselStruct->PointsOnBoundary.end()) {
565 PointAdvance++;
566 point = PointRunner->second;
[b8d215]567 LOG(2, "INFO: Current point is " << *point << ".");
[aac19f]568 // check that at least a single line is concave
569 LineMap::iterator LineRunner = point->lines.begin();
570 for (; LineRunner != point->lines.end(); LineRunner++) {
571 const BoundaryLineSet * line = LineRunner->second;
[b8d215]572 LOG(3, "INFO: Current line of point " << *point << " is " << *line << ".");
[aac19f]573 if (!line->CheckConvexityCriterion())
574 break;
[d289c6]575 }
[aac19f]576 // remove the point if needed
577 if (LineRunner != point->lines.end()) {
578 const double tmp = TesselStruct->RemoveFullConcavePointFromTesselatedSurface(point);
579 if (tmp > 0.) {
580 volume += tmp;
581 Concavity = true;
582 }
[1d9b7aa]583 }
584 PointRunner = PointAdvance;
[093645]585 }
586
[5f7b95]587 if (DebugOutputEveryStep) {
588 sprintf(dummy, "-%d", run++);
589 //CalculateConcavityPerBoundaryPoint(TesselStruct);
590 LOG(1, "INFO: Writing " << run << "th tesselation file.");
591 StoreTrianglesinFile(mol, (const Tesselation *&)TesselStruct, filename, dummy);
592 }
[093645]593
[d289c6]594 // second step: flip baselines, i.e. add general tetraeder at concave lines
595 // when the tetraeder does not intersect with other already present triangles
[1d9b7aa]596 LineRunner = TesselStruct->LinesOnBoundary.begin();
597 LineAdvance = LineRunner; // we need an advanced line, as the LineRunner might get removed
[81ffd8]598 std::map<double, std::pair<BoundaryLineSet *, double> > GainMap;
[1d9b7aa]599 while (LineRunner != TesselStruct->LinesOnBoundary.end()) {
600 LineAdvance++;
601 line = LineRunner->second;
[d289c6]602 if (!line->CheckConvexityCriterion()) {
603 LOG(2, "INFO: concave line is " << *line << ".");
604 // gather the other points
605 BoundaryPointSet *BPS[4];
606 int m = 0;
607 {
608 for (TriangleMap::iterator runner = line->triangles.begin(); runner != line->triangles.end(); runner++)
609 for (int j = 0; j < 3; j++) // all of their endpoints and baselines
610 if (!line->ContainsBoundaryPoint(runner->second->endpoints[j])) // and neither of its endpoints
611 BPS[m++] = runner->second->endpoints[j];
612 }
613 BPS[2] = line->endpoints[0];
614 BPS[3] = line->endpoints[1];
615 LOG(3, "DEBUG: other line would consist of " << *BPS[0] << " and "
616 << *BPS[1] << ".");
617
618 // check for already present (third) side of the tetraeder as we then
619 // would create a degenerate triangle
620 bool TetraederSidePresent = false;
621 {
622 class TesselPoint *TriangleCandidates[3];
623 TriangleCandidates[0] = BPS[0]->node;
624 TriangleCandidates[1] = BPS[1]->node;
625 TriangleCandidates[2] = BPS[2]->node;
626 if ((TesselStruct->GetPresentTriangle(TriangleCandidates) != NULL)) {
627 LOG(2, "REJECT: Triangle side " << *TriangleCandidates[0] << ","
628 << *TriangleCandidates[1] << "," << *TriangleCandidates[2] << " present.");
629 TetraederSidePresent = true;
630 }
631 TriangleCandidates[2] = BPS[3]->node;
632 if ((TesselStruct->GetPresentTriangle(TriangleCandidates) != NULL)) {
633 LOG(2, "REJECT: Triangle side " << *TriangleCandidates[0] << ","
634 << *TriangleCandidates[1] << "," << *TriangleCandidates[2] << " present.");
635 TetraederSidePresent = true;
636 }
637 }
638
639 if ((BPS[0] != BPS[1]) && (m == 2) && (!TetraederSidePresent)) {
640 // check whether all adjacent triangles do not intersect with new line
641 bool no_line_intersects = true;
642 Vector Intersection;
643 TriangleSet triangles;
644 TriangleSet *firsttriangles = TesselStruct->GetAllTriangles(line->endpoints[0]);
645 TriangleSet *secondtriangles = TesselStruct->GetAllTriangles(line->endpoints[1]);
646 triangles.insert( firsttriangles->begin(), firsttriangles->end() );
647 triangles.insert( secondtriangles->begin(), secondtriangles->end() );
648 delete firsttriangles;
649 delete secondtriangles;
650 for (TriangleSet::const_iterator triangleiter = triangles.begin();
651 triangleiter != triangles.end(); ++triangleiter) {
652 const BoundaryTriangleSet * triangle = *triangleiter;
653 bool line_intersects = triangle->GetIntersectionInsideTriangle(
654 BPS[0]->node->getPosition(),
655 BPS[1]->node->getPosition(),
656 Intersection);
657 // switch result when coinciding with endpoint
658 bool concave_adjacent_line = false;
659 bool intersection_is_endnode = false;
660 for (int j=0;j<2;++j) {
661 if (Intersection.DistanceSquared(BPS[j]->node->getPosition()) < MYEPSILON) {
662 intersection_is_endnode = true;
663 // check whether its an adjacent triangle and if it's concavely connected
664 // only then are we in danger of cutting through it and need to check
665 // sign of normal vector and intersecting line
666 for (int i=0;i<2;++i)
667 for (int lineindex=0;lineindex < 3;++lineindex)
668 if ((triangle->lines[lineindex]->ContainsBoundaryPoint(line->endpoints[i]))
669 && (triangle->lines[lineindex]->ContainsBoundaryPoint(BPS[j]))) {
670 concave_adjacent_line = !triangle->lines[lineindex]->CheckConvexityCriterion();
671 }
672 if (concave_adjacent_line) {
673 const Vector intersector =
674 BPS[(j+1)%2]->node->getPosition() - Intersection;
675 if (triangle->NormalVector.ScalarProduct(intersector) >= -MYEPSILON) {
676 LOG(4, "ACCEPT: Intersection coincides with first endpoint "
677 << *BPS[j] << ".");
678 line_intersects = false;
679 } else {
680 LOG(4, "REJECT: Intersection ends on wrong side of triangle.");
681 }
682 } else {
683 LOG(4, "ACCEPT: Intersection coincides with first endpoint "
684 << *BPS[j] << ".");
685 line_intersects = false;
686 }
687 }
688 }
689 // if we have an intersection, check that it is within either
690 // endpoint, i.e. check that scalar product between vectors going
691 // from intersction to either endpoint has negative sign (both
692 // vectors point in opposite directions)
693 if (!intersection_is_endnode && line_intersects) {
694 const Vector firstvector = BPS[0]->node->getPosition() - Intersection;
695 const Vector secondvector = BPS[1]->node->getPosition() - Intersection;
696 if (firstvector.ScalarProduct(secondvector) >= 0)
697 line_intersects = false;
698 }
699 no_line_intersects &= !line_intersects;
700 }
701
702 if (no_line_intersects) {
703 // calculate the volume
[81ffd8]704 const double tmp = line->CalculateConvexity();
705 const double gain =
706 CalculateVolumeofGeneralTetraeder(
707 BPS[0]->node->getPosition(),
708 BPS[1]->node->getPosition(),
709 BPS[2]->node->getPosition(),
710 BPS[3]->node->getPosition());
711
712 GainMap.insert(std::make_pair(tmp, std::make_pair(line,gain) ));
[d289c6]713 LOG(2, "DEBUG: Adding concave line " << *line << " with gain of "
[81ffd8]714 << gain << ".");
[d289c6]715 } else {
716 // if 2 or 3 don't
717 LOG(2, "DEBUG: We don't added concave line " << *line
718 << " as other line intersects with adjacent triangles.");
719 }
[57066a]720 }
[1d9b7aa]721 }
722 LineRunner = LineAdvance;
723 }
[d289c6]724 // flip line with most gain
725 if (!GainMap.empty()) {
[81ffd8]726 line = GainMap.begin()->second.first;
727 const double tmp = GainMap.begin()->second.second;
[d289c6]728 volume += tmp;
[81ffd8]729
[d289c6]730// GainMap.clear();
731
732 // and flip the line
[81ffd8]733 LOG(1, "INFO: Flipping current most concave line " << *line << " with gain of "
[d289c6]734 << tmp << ".");
735 TesselStruct->FlipBaseline(line);
736 Concavity = true;
737 }
738 } while ((Concavity)); // && (run < 100)
[093645]739
[e138de]740 CalculateConcavityPerBoundaryPoint(TesselStruct);
741 StoreTrianglesinFile(mol, (const Tesselation *&)TesselStruct, filename, "");
[0077b5]742
743 // end
[d289c6]744 LOG(0, "RESULT: Added volume in convexization is " << volume << ".");
[0077b5]745 return volume;
746};
747
[6ac7ee]748
[7dea7c]749/** Stores triangles to file.
750 * \param *out output stream for debugging
751 * \param *mol molecule with atoms and bonds
[71b20e]752 * \param *TesselStruct Tesselation with boundary triangles
[7dea7c]753 * \param *filename prefix of filename
754 * \param *extraSuffix intermediate suffix
755 */
[71b20e]756void StoreTrianglesinFile(const molecule * const mol, const Tesselation * const TesselStruct, const char *filename, const char *extraSuffix)
[7dea7c]757{
[ce7bfd]758 //Info FunctionInfo(__func__);
[caa06ef]759 PointCloudAdaptor< molecule > cloud(const_cast<molecule *>(mol), mol->name);
[7dea7c]760 // 4. Store triangles in tecplot file
761 if (filename != NULL) {
762 if (DoTecplotOutput) {
763 string OutputName(filename);
764 OutputName.append(extraSuffix);
765 OutputName.append(TecplotSuffix);
766 ofstream *tecplot = new ofstream(OutputName.c_str());
[34c43a]767 WriteTecplotFile(tecplot, TesselStruct, cloud, -1);
[7dea7c]768 tecplot->close();
769 delete(tecplot);
770 }
771 if (DoRaster3DOutput) {
772 string OutputName(filename);
773 OutputName.append(extraSuffix);
774 OutputName.append(Raster3DSuffix);
775 ofstream *rasterplot = new ofstream(OutputName.c_str());
[34c43a]776 WriteRaster3dFile(rasterplot, TesselStruct, cloud);
[7dea7c]777 rasterplot->close();
778 delete(rasterplot);
779 }
780 }
781};
[03648b]782
[6ac7ee]783/** Tesselates the non convex boundary by rolling a virtual sphere along the surface of the molecule.
784 * \param *out output stream for debugging
785 * \param *mol molecule structure with Atom's and Bond's
[776b64]786 * \param *&TesselStruct Tesselation filled with points, lines and triangles on boundary on return
[6bd7e0]787 * \param *&LCList atoms in LinkedCell_deprecated list
[57066a]788 * \param RADIUS radius of the virtual sphere
[6ac7ee]789 * \param *filename filename prefix for output of vertex data
[4fc93f]790 * \return true - tesselation successful, false - tesselation failed
[6ac7ee]791 */
[6bd7e0]792bool FindNonConvexBorder(molecule* const mol, Tesselation *&TesselStruct, const LinkedCell_deprecated *&LCList, const double RADIUS, const char *filename = NULL)
[03648b]793{
[ce7bfd]794 //Info FunctionInfo(__func__);
[3d919e]795 bool freeLC = false;
[4fc93f]796 bool status = false;
[6613ec]797 CandidateForTesselation *baseline = NULL;
[1e168b]798 bool OneLoopWithoutSuccessFlag = true; // marks whether we went once through all baselines without finding any without two triangles
[a2a2f7]799// bool TesselationFailFlag = false;
[357fba]800
[a7b761b]801 mol->getAtomCount();
[357fba]802
[776b64]803 if (TesselStruct == NULL) {
[47d041]804 LOG(1, "Allocating Tesselation struct ...");
[776b64]805 TesselStruct= new Tesselation;
[ef0e6d]806 } else {
[776b64]807 delete(TesselStruct);
[47d041]808 LOG(1, "Re-Allocating Tesselation struct ...");
[776b64]809 TesselStruct = new Tesselation;
[3d919e]810 }
[ad37ab]811
[57066a]812 // initialise Linked Cell
[caa06ef]813 PointCloudAdaptor< molecule > cloud(mol, mol->name);
[3d919e]814 if (LCList == NULL) {
[6bd7e0]815 LCList = new LinkedCell_deprecated(cloud, 2.*RADIUS);
[3d919e]816 freeLC = true;
817 }
818
[57066a]819 // 1. get starting triangle
[ce70970]820 if (!TesselStruct->FindStartingTriangle(RADIUS, LCList)) {
[47d041]821 ELOG(0, "No valid starting triangle found.");
[b5c2d7]822 //performCriticalExit();
[ce70970]823 }
[f8bd7d]824 if (filename != NULL) {
825 if ((DoSingleStepOutput && ((TesselStruct->TrianglesOnBoundary.size() % SingleStepWidth == 0)))) { // if we have a new triangle and want to output each new triangle configuration
[34c43a]826 TesselStruct->Output(filename, cloud);
[f8bd7d]827 }
828 }
[3d919e]829
[57066a]830 // 2. expand from there
[1e168b]831 while ((!TesselStruct->OpenLines.empty()) && (OneLoopWithoutSuccessFlag)) {
[b40a93]832 (cerr << "There are " << TesselStruct->TrianglesOnBoundary.size() << " triangles and " << TesselStruct->OpenLines.size() << " open lines to scan for candidates." << endl);
833 // 2a. print OpenLines without candidates
[47d041]834 LOG(1, "There are the following open lines to scan for a candidates:");
[1e168b]835 for (CandidateMap::iterator Runner = TesselStruct->OpenLines.begin(); Runner != TesselStruct->OpenLines.end(); Runner++)
[b40a93]836 if (Runner->second->pointlist.empty())
[47d041]837 LOG(1, " " << *(Runner->second));
[1e168b]838
[734816]839 // 2b. find best candidate for each OpenLine
[a2a2f7]840 const bool TesselationFailFlag = TesselStruct->FindCandidatesforOpenLines(RADIUS, LCList);
841 ASSERT( TesselationFailFlag,
842 "FindNonConvexBorder() - at least one open line without candidate exists.");
[3d919e]843
[734816]844 // 2c. print OpenLines with candidates again
[47d041]845 LOG(1, "There are " << TesselStruct->OpenLines.size() << " open lines to scan for the best candidates:");
[1e168b]846 for (CandidateMap::iterator Runner = TesselStruct->OpenLines.begin(); Runner != TesselStruct->OpenLines.end(); Runner++)
[47d041]847 LOG(1, " " << *(Runner->second));
[1e168b]848
[734816]849 // 2d. search for smallest ShortestAngle among all candidates
850 double ShortestAngle = 4.*M_PI;
[1e168b]851 for (CandidateMap::iterator Runner = TesselStruct->OpenLines.begin(); Runner != TesselStruct->OpenLines.end(); Runner++) {
852 if (Runner->second->ShortestAngle < ShortestAngle) {
853 baseline = Runner->second;
854 ShortestAngle = baseline->ShortestAngle;
[47d041]855 LOG(1, "New best candidate is " << *baseline->BaseLine << " with point " << *(*baseline->pointlist.begin()) << " and angle " << baseline->ShortestAngle);
[1e168b]856 }
857 }
[734816]858 // 2e. if we found one, add candidate
[f67b6e]859 if ((ShortestAngle == 4.*M_PI) || (baseline->pointlist.empty()))
[7dea7c]860 OneLoopWithoutSuccessFlag = false;
[1e168b]861 else {
[474961]862 TesselStruct->AddCandidatePolygon(*baseline, RADIUS, LCList);
[1e168b]863 }
864
[734816]865 // 2f. write temporary envelope
[1e168b]866 if (filename != NULL) {
867 if ((DoSingleStepOutput && ((TesselStruct->TrianglesOnBoundary.size() % SingleStepWidth == 0)))) { // if we have a new triangle and want to output each new triangle configuration
[34c43a]868 TesselStruct->Output(filename, cloud);
[1e168b]869 }
[3d919e]870 }
871 }
[4fc93f]872// // check envelope for consistency
873// status = CheckListOfBaselines(TesselStruct);
874//
875// // look whether all points are inside of the convex envelope, otherwise add them via degenerated triangles
876// //->InsertStraddlingPoints(mol, LCList);
[9879f6]877// for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
[57066a]878// class TesselPoint *Runner = NULL;
[9879f6]879// Runner = *iter;
[47d041]880// LOG(1, "Checking on " << Runner->Name << " ... ");
[e138de]881// if (!->IsInnerPoint(Runner, LCList)) {
[47d041]882// LOG(2, Runner->Name << " is outside of envelope, adding via degenerated triangles.");
[e138de]883// ->AddBoundaryPointByDegeneratedTriangle(Runner, LCList);
[57066a]884// } else {
[47d041]885// LOG(2, Runner->Name << " is inside of or on envelope.");
[57066a]886// }
887// }
[357fba]888
[f67b6e]889// // Purges surplus triangles.
890// TesselStruct->RemoveDegeneratedTriangles();
[b32dbb]891//
892// // check envelope for consistency
893// status = CheckListOfBaselines(TesselStruct);
[b998c3]894
[b8d215]895// cout << "before correction" << endl;
[c72112]896
[b998c3]897 // store before correction
[c72112]898 StoreTrianglesinFile(mol, TesselStruct, filename, "");
[b998c3]899
[6613ec]900// // correct degenerated polygons
901// TesselStruct->CorrectAllDegeneratedPolygons();
902//
[e69c87]903 // check envelope for consistency
904 status = CheckListOfBaselines(TesselStruct);
[ef0e6d]905
[57066a]906 // write final envelope
[e138de]907 CalculateConcavityPerBoundaryPoint(TesselStruct);
[b8d215]908// cout << "after correction" << endl;
[c72112]909 StoreTrianglesinFile(mol, TesselStruct, filename, "");
[8c54a3]910
[3d919e]911 if (freeLC)
912 delete(LCList);
[4fc93f]913
914 return status;
[6ac7ee]915};
[03648b]916
[57066a]917
[ad37ab]918/** Finds a hole of sufficient size in \a *mols to embed \a *srcmol into it.
[ca2587]919 * \param *out output stream for debugging
[ad37ab]920 * \param *mols molecules in the domain to embed in between
921 * \param *srcmol embedding molecule
[ca2587]922 * \return *Vector new center of \a *srcmol for embedding relative to \a this
923 */
[e138de]924Vector* FindEmbeddingHole(MoleculeListClass *mols, molecule *srcmol)
[ca2587]925{
[ce7bfd]926 //Info FunctionInfo(__func__);
[ca2587]927 Vector *Center = new Vector;
928 Center->Zero();
929 // calculate volume/shape of \a *srcmol
930
931 // find embedding holes
932
933 // if more than one, let user choose
934
935 // return embedding center
936 return Center;
937};
938
Note: See TracBrowser for help on using the repository browser.