source: src/boundary.cpp@ 89c8b2

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 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 89c8b2 was 34e0592, checked in by Frederik Heber <heber@…>, 15 years ago

Merge branch 'ConcaveHull' of ssh://stud64d-02/home/metzler/workspace/espack into Ticket14

Conflicts:

molecuilder/src/boundary.cpp
molecuilder/src/tesselation.cpp

  • Property mode set to 100755
File size: 58.4 KB
Line 
1/** \file boundary.hpp
2 *
3 * Implementations and super-function for envelopes
4 */
5
6
7#include "boundary.hpp"
8
9#include<gsl/gsl_poly.h>
10
11// ========================================== F U N C T I O N S =================================
12
13
14/** Determines greatest diameters of a cluster defined by its convex envelope.
15 * Looks at lines parallel to one axis and where they intersect on the projected planes
16 * \param *out output stream for debugging
17 * \param *BoundaryPoints NDIM set of boundary points defining the convex envelope on each projected plane
18 * \param *mol molecule structure representing the cluster
19 * \param IsAngstroem whether we have angstroem or atomic units
20 * \return NDIM array of the diameters
21 */
22double *
23GetDiametersOfCluster(ofstream *out, Boundaries *BoundaryPtr, molecule *mol,
24 bool IsAngstroem)
25{
26 // get points on boundary of NULL was given as parameter
27 bool BoundaryFreeFlag = false;
28 Boundaries *BoundaryPoints = BoundaryPtr;
29 if (BoundaryPoints == NULL) {
30 BoundaryFreeFlag = true;
31 BoundaryPoints = GetBoundaryPoints(out, mol);
32 } else {
33 *out << Verbose(1) << "Using given boundary points set." << endl;
34 }
35 // determine biggest "diameter" of cluster for each axis
36 Boundaries::iterator Neighbour, OtherNeighbour;
37 double *GreatestDiameter = new double[NDIM];
38 for (int i = 0; i < NDIM; i++)
39 GreatestDiameter[i] = 0.;
40 double OldComponent, tmp, w1, w2;
41 Vector DistanceVector, OtherVector;
42 int component, Othercomponent;
43 for (int axis = 0; axis < NDIM; axis++)
44 { // regard each projected plane
45 //*out << Verbose(1) << "Current axis is " << axis << "." << endl;
46 for (int j = 0; j < 2; j++)
47 { // and for both axis on the current plane
48 component = (axis + j + 1) % NDIM;
49 Othercomponent = (axis + 1 + ((j + 1) & 1)) % NDIM;
50 //*out << Verbose(1) << "Current component is " << component << ", Othercomponent is " << Othercomponent << "." << endl;
51 for (Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner
52 != BoundaryPoints[axis].end(); runner++)
53 {
54 //*out << Verbose(2) << "Current runner is " << *(runner->second.second) << "." << endl;
55 // seek for the neighbours pair where the Othercomponent sign flips
56 Neighbour = runner;
57 Neighbour++;
58 if (Neighbour == BoundaryPoints[axis].end()) // make it wrap around
59 Neighbour = BoundaryPoints[axis].begin();
60 DistanceVector.CopyVector(&runner->second.second->x);
61 DistanceVector.SubtractVector(&Neighbour->second.second->x);
62 do
63 { // seek for neighbour pair where it flips
64 OldComponent = DistanceVector.x[Othercomponent];
65 Neighbour++;
66 if (Neighbour == BoundaryPoints[axis].end()) // make it wrap around
67 Neighbour = BoundaryPoints[axis].begin();
68 DistanceVector.CopyVector(&runner->second.second->x);
69 DistanceVector.SubtractVector(&Neighbour->second.second->x);
70 //*out << Verbose(3) << "OldComponent is " << OldComponent << ", new one is " << DistanceVector.x[Othercomponent] << "." << endl;
71 }
72 while ((runner != Neighbour) && (fabs(OldComponent / fabs(
73 OldComponent) - DistanceVector.x[Othercomponent] / fabs(
74 DistanceVector.x[Othercomponent])) < MYEPSILON)); // as long as sign does not flip
75 if (runner != Neighbour)
76 {
77 OtherNeighbour = Neighbour;
78 if (OtherNeighbour == BoundaryPoints[axis].begin()) // make it wrap around
79 OtherNeighbour = BoundaryPoints[axis].end();
80 OtherNeighbour--;
81 //*out << Verbose(2) << "The pair, where the sign of OtherComponent flips, is: " << *(Neighbour->second.second) << " and " << *(OtherNeighbour->second.second) << "." << endl;
82 // now we have found the pair: Neighbour and OtherNeighbour
83 OtherVector.CopyVector(&runner->second.second->x);
84 OtherVector.SubtractVector(&OtherNeighbour->second.second->x);
85 //*out << Verbose(2) << "Distances to Neighbour and OtherNeighbour are " << DistanceVector.x[component] << " and " << OtherVector.x[component] << "." << endl;
86 //*out << Verbose(2) << "OtherComponents to Neighbour and OtherNeighbour are " << DistanceVector.x[Othercomponent] << " and " << OtherVector.x[Othercomponent] << "." << endl;
87 // do linear interpolation between points (is exact) to extract exact intersection between Neighbour and OtherNeighbour
88 w1 = fabs(OtherVector.x[Othercomponent]);
89 w2 = fabs(DistanceVector.x[Othercomponent]);
90 tmp = fabs((w1 * DistanceVector.x[component] + w2
91 * OtherVector.x[component]) / (w1 + w2));
92 // mark if it has greater diameter
93 //*out << Verbose(2) << "Comparing current greatest " << GreatestDiameter[component] << " to new " << tmp << "." << endl;
94 GreatestDiameter[component] = (GreatestDiameter[component]
95 > tmp) ? GreatestDiameter[component] : tmp;
96 } //else
97 //*out << Verbose(2) << "Saw no sign flip, probably top or bottom node." << endl;
98 }
99 }
100 }
101 *out << Verbose(0) << "RESULT: The biggest diameters are "
102 << GreatestDiameter[0] << " and " << GreatestDiameter[1] << " and "
103 << GreatestDiameter[2] << " " << (IsAngstroem ? "angstrom"
104 : "atomiclength") << "." << endl;
105
106 // free reference lists
107 if (BoundaryFreeFlag)
108 delete[] (BoundaryPoints);
109
110 return GreatestDiameter;
111}
112;
113
114/** Creates the objects in a VRML file.
115 * \param *out output stream for debugging
116 * \param *vrmlfile output stream for tecplot data
117 * \param *Tess Tesselation structure with constructed triangles
118 * \param *mol molecule structure with atom positions
119 */
120void WriteVrmlFile(ofstream *out, ofstream *vrmlfile, class Tesselation *Tess, class molecule *mol)
121{
122 atom *Walker = mol->start;
123 bond *Binder = mol->first;
124 int i;
125 Vector *center = mol->DetermineCenterOfAll(out);
126 if (vrmlfile != NULL) {
127 //cout << Verbose(1) << "Writing Raster3D file ... ";
128 *vrmlfile << "#VRML V2.0 utf8" << endl;
129 *vrmlfile << "#Created by molecuilder" << endl;
130 *vrmlfile << "#All atoms as spheres" << endl;
131 while (Walker->next != mol->end) {
132 Walker = Walker->next;
133 *vrmlfile << "Sphere {" << endl << " "; // 2 is sphere type
134 for (i=0;i<NDIM;i++)
135 *vrmlfile << Walker->x.x[i]-center->x[i] << " ";
136 *vrmlfile << "\t0.1\t1. 1. 1." << endl; // radius 0.05 and white as colour
137 }
138
139 *vrmlfile << "# All bonds as vertices" << endl;
140 while (Binder->next != mol->last) {
141 Binder = Binder->next;
142 *vrmlfile << "3" << endl << " "; // 2 is round-ended cylinder type
143 for (i=0;i<NDIM;i++)
144 *vrmlfile << Binder->leftatom->x.x[i]-center->x[i] << " ";
145 *vrmlfile << "\t0.03\t";
146 for (i=0;i<NDIM;i++)
147 *vrmlfile << Binder->rightatom->x.x[i]-center->x[i] << " ";
148 *vrmlfile << "\t0.03\t0. 0. 1." << endl; // radius 0.05 and blue as colour
149 }
150
151 *vrmlfile << "# All tesselation triangles" << endl;
152 for (TriangleMap::iterator TriangleRunner = Tess->TrianglesOnBoundary.begin(); TriangleRunner != Tess->TrianglesOnBoundary.end(); TriangleRunner++) {
153 *vrmlfile << "1" << endl << " "; // 1 is triangle type
154 for (i=0;i<3;i++) { // print each node
155 for (int j=0;j<NDIM;j++) // and for each node all NDIM coordinates
156 *vrmlfile << TriangleRunner->second->endpoints[i]->node->node->x[j]-center->x[j] << " ";
157 *vrmlfile << "\t";
158 }
159 *vrmlfile << "1. 0. 0." << endl; // red as colour
160 *vrmlfile << "18" << endl << " 0.5 0.5 0.5" << endl; // 18 is transparency type for previous object
161 }
162 } else {
163 cerr << "ERROR: Given vrmlfile is " << vrmlfile << "." << endl;
164 }
165 delete(center);
166};
167
168/** Creates the objects in a raster3d file (renderable with a header.r3d).
169 * \param *out output stream for debugging
170 * \param *rasterfile output stream for tecplot data
171 * \param *Tess Tesselation structure with constructed triangles
172 * \param *mol molecule structure with atom positions
173 */
174void WriteRaster3dFile(ofstream *out, ofstream *rasterfile, class Tesselation *Tess, class molecule *mol)
175{
176 atom *Walker = mol->start;
177 bond *Binder = mol->first;
178 int i;
179 Vector *center = mol->DetermineCenterOfAll(out);
180 if (rasterfile != NULL) {
181 //cout << Verbose(1) << "Writing Raster3D file ... ";
182 *rasterfile << "# Raster3D object description, created by MoleCuilder" << endl;
183 *rasterfile << "@header.r3d" << endl;
184 *rasterfile << "# All atoms as spheres" << endl;
185 while (Walker->next != mol->end) {
186 Walker = Walker->next;
187 *rasterfile << "2" << endl << " "; // 2 is sphere type
188 for (i=0;i<NDIM;i++)
189 *rasterfile << Walker->x.x[i]-center->x[i] << " ";
190 *rasterfile << "\t0.1\t1. 1. 1." << endl; // radius 0.05 and white as colour
191 }
192
193 *rasterfile << "# All bonds as vertices" << endl;
194 while (Binder->next != mol->last) {
195 Binder = Binder->next;
196 *rasterfile << "3" << endl << " "; // 2 is round-ended cylinder type
197 for (i=0;i<NDIM;i++)
198 *rasterfile << Binder->leftatom->x.x[i]-center->x[i] << " ";
199 *rasterfile << "\t0.03\t";
200 for (i=0;i<NDIM;i++)
201 *rasterfile << Binder->rightatom->x.x[i]-center->x[i] << " ";
202 *rasterfile << "\t0.03\t0. 0. 1." << endl; // radius 0.05 and blue as colour
203 }
204
205 *rasterfile << "# All tesselation triangles" << endl;
206 *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";
207 for (TriangleMap::iterator TriangleRunner = Tess->TrianglesOnBoundary.begin(); TriangleRunner != Tess->TrianglesOnBoundary.end(); TriangleRunner++) {
208 *rasterfile << "1" << endl << " "; // 1 is triangle type
209 for (i=0;i<3;i++) { // print each node
210 for (int j=0;j<NDIM;j++) // and for each node all NDIM coordinates
211 *rasterfile << TriangleRunner->second->endpoints[i]->node->node->x[j]-center->x[j] << " ";
212 *rasterfile << "\t";
213 }
214 *rasterfile << "1. 0. 0." << endl; // red as colour
215 //*rasterfile << "18" << endl << " 0.5 0.5 0.5" << endl; // 18 is transparency type for previous object
216 }
217 *rasterfile << "9\n# terminating special property\n";
218 } else {
219 cerr << "ERROR: Given rasterfile is " << rasterfile << "." << endl;
220 }
221 delete(center);
222};
223
224/** This function creates the tecplot file, displaying the tesselation of the hull.
225 * \param *out output stream for debugging
226 * \param *tecplot output stream for tecplot data
227 * \param N arbitrary number to differentiate various zones in the tecplot format
228 */
229void WriteTecplotFile(ofstream *out, ofstream *tecplot, class Tesselation *TesselStruct, class molecule *mol, int N)
230{
231 if ((tecplot != NULL) && (TesselStruct != NULL)) {
232 // write header
233 *tecplot << "TITLE = \"3D CONVEX SHELL\"" << endl;
234 *tecplot << "VARIABLES = \"X\" \"Y\" \"Z\" \"U\"" << endl;
235 *tecplot << "ZONE T=\"TRIANGLES" << N << "\", N=" << TesselStruct->PointsOnBoundary.size() << ", E=" << TesselStruct->TrianglesOnBoundary.size() << ", DATAPACKING=POINT, ZONETYPE=FETRIANGLE" << endl;
236 int *LookupList = new int[mol->AtomCount];
237 for (int i = 0; i < mol->AtomCount; i++)
238 LookupList[i] = -1;
239
240 // print atom coordinates
241 *out << Verbose(2) << "The following triangles were created:";
242 int Counter = 1;
243 TesselPoint *Walker = NULL;
244 for (PointMap::iterator target = TesselStruct->PointsOnBoundary.begin(); target != TesselStruct->PointsOnBoundary.end(); target++) {
245 Walker = target->second->node;
246 LookupList[Walker->nr] = Counter++;
247 *tecplot << Walker->node->x[0] << " " << Walker->node->x[1] << " " << Walker->node->x[2] << " " << target->second->value << endl;
248 }
249 *tecplot << endl;
250 // print connectivity
251 for (TriangleMap::iterator runner = TesselStruct->TrianglesOnBoundary.begin(); runner != TesselStruct->TrianglesOnBoundary.end(); runner++) {
252 *out << " " << runner->second->endpoints[0]->node->Name << "<->" << runner->second->endpoints[1]->node->Name << "<->" << runner->second->endpoints[2]->node->Name;
253 *tecplot << LookupList[runner->second->endpoints[0]->node->nr] << " " << LookupList[runner->second->endpoints[1]->node->nr] << " " << LookupList[runner->second->endpoints[2]->node->nr] << endl;
254 }
255 delete[] (LookupList);
256 *out << endl;
257 }
258}
259
260
261/** Determines the boundary points of a cluster.
262 * Does a projection per axis onto the orthogonal plane, transforms into spherical coordinates, sorts them by the angle
263 * and looks at triples: if the middle has less a distance than the allowed maximum height of the triangle formed by the plane's
264 * center and first and last point in the triple, it is thrown out.
265 * \param *out output stream for debugging
266 * \param *mol molecule structure representing the cluster
267 */
268Boundaries *GetBoundaryPoints(ofstream *out, molecule *mol)
269{
270 atom *Walker = NULL;
271 PointMap PointsOnBoundary;
272 LineMap LinesOnBoundary;
273 TriangleMap TrianglesOnBoundary;
274 Vector *MolCenter = mol->DetermineCenterOfAll(out);
275 Vector helper;
276
277 *out << Verbose(1) << "Finding all boundary points." << endl;
278 Boundaries *BoundaryPoints = new Boundaries[NDIM]; // first is alpha, second is (r, nr)
279 BoundariesTestPair BoundaryTestPair;
280 Vector AxisVector, AngleReferenceVector, AngleReferenceNormalVector;
281 double radius, angle;
282 // 3a. Go through every axis
283 for (int axis = 0; axis < NDIM; axis++) {
284 AxisVector.Zero();
285 AngleReferenceVector.Zero();
286 AngleReferenceNormalVector.Zero();
287 AxisVector.x[axis] = 1.;
288 AngleReferenceVector.x[(axis + 1) % NDIM] = 1.;
289 AngleReferenceNormalVector.x[(axis + 2) % NDIM] = 1.;
290
291 *out << Verbose(1) << "Axisvector is " << AxisVector << " and AngleReferenceVector is " << AngleReferenceVector << ", and AngleReferenceNormalVector is " << AngleReferenceNormalVector << "." << endl;
292
293 // 3b. construct set of all points, transformed into cylindrical system and with left and right neighbours
294 Walker = mol->start;
295 while (Walker->next != mol->end) {
296 Walker = Walker->next;
297 Vector ProjectedVector;
298 ProjectedVector.CopyVector(&Walker->x);
299 ProjectedVector.SubtractVector(MolCenter);
300 ProjectedVector.ProjectOntoPlane(&AxisVector);
301
302 // correct for negative side
303 radius = ProjectedVector.NormSquared();
304 if (fabs(radius) > MYEPSILON)
305 angle = ProjectedVector.Angle(&AngleReferenceVector);
306 else
307 angle = 0.; // otherwise it's a vector in Axis Direction and unimportant for boundary issues
308
309 //*out << "Checking sign in quadrant : " << ProjectedVector.Projection(&AngleReferenceNormalVector) << "." << endl;
310 if (ProjectedVector.ScalarProduct(&AngleReferenceNormalVector) > 0) {
311 angle = 2. * M_PI - angle;
312 }
313 *out << Verbose(2) << "Inserting " << *Walker << ": (r, alpha) = (" << radius << "," << angle << "): " << ProjectedVector << endl;
314 BoundaryTestPair = BoundaryPoints[axis].insert(BoundariesPair(angle, DistancePair (radius, Walker)));
315 if (!BoundaryTestPair.second) { // same point exists, check first r, then distance of original vectors to center of gravity
316 *out << Verbose(2) << "Encountered two vectors whose projection onto axis " << axis << " is equal: " << endl;
317 *out << Verbose(2) << "Present vector: " << *BoundaryTestPair.first->second.second << endl;
318 *out << Verbose(2) << "New vector: " << *Walker << endl;
319 double tmp = ProjectedVector.NormSquared();
320 if ((tmp - BoundaryTestPair.first->second.first) > MYEPSILON) {
321 BoundaryTestPair.first->second.first = tmp;
322 BoundaryTestPair.first->second.second = Walker;
323 *out << Verbose(2) << "Keeping new vector due to larger projected distance " << tmp << "." << endl;
324 } else if (fabs(tmp - BoundaryTestPair.first->second.first) < MYEPSILON) {
325 helper.CopyVector(&Walker->x);
326 helper.SubtractVector(MolCenter);
327 tmp = helper.NormSquared();
328 helper.CopyVector(&BoundaryTestPair.first->second.second->x);
329 helper.SubtractVector(MolCenter);
330 if (helper.NormSquared() < tmp) {
331 BoundaryTestPair.first->second.second = Walker;
332 *out << Verbose(2) << "Keeping new vector due to larger distance to molecule center " << helper.NormSquared() << "." << endl;
333 } else {
334 *out << Verbose(2) << "Keeping present vector due to larger distance to molecule center " << tmp << "." << endl;
335 }
336 } else {
337 *out << Verbose(2) << "Keeping present vector due to larger projected distance " << tmp << "." << endl;
338 }
339 }
340 }
341 // printing all inserted for debugging
342 // {
343 // *out << Verbose(2) << "Printing list of candidates for axis " << axis << " which we have inserted so far." << endl;
344 // int i=0;
345 // for(Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
346 // if (runner != BoundaryPoints[axis].begin())
347 // *out << ", " << i << ": " << *runner->second.second;
348 // else
349 // *out << i << ": " << *runner->second.second;
350 // i++;
351 // }
352 // *out << endl;
353 // }
354 // 3c. throw out points whose distance is less than the mean of left and right neighbours
355 bool flag = false;
356 *out << Verbose(1) << "Looking for candidates to kick out by convex condition ... " << endl;
357 do { // do as long as we still throw one out per round
358 flag = false;
359 Boundaries::iterator left = BoundaryPoints[axis].end();
360 Boundaries::iterator right = BoundaryPoints[axis].end();
361 for (Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
362 // set neighbours correctly
363 if (runner == BoundaryPoints[axis].begin()) {
364 left = BoundaryPoints[axis].end();
365 } else {
366 left = runner;
367 }
368 left--;
369 right = runner;
370 right++;
371 if (right == BoundaryPoints[axis].end()) {
372 right = BoundaryPoints[axis].begin();
373 }
374 // check distance
375
376 // construct the vector of each side of the triangle on the projected plane (defined by normal vector AxisVector)
377 {
378 Vector SideA, SideB, SideC, SideH;
379 SideA.CopyVector(&left->second.second->x);
380 SideA.SubtractVector(MolCenter);
381 SideA.ProjectOntoPlane(&AxisVector);
382 // *out << "SideA: ";
383 // SideA.Output(out);
384 // *out << endl;
385
386 SideB.CopyVector(&right->second.second->x);
387 SideB.SubtractVector(MolCenter);
388 SideB.ProjectOntoPlane(&AxisVector);
389 // *out << "SideB: ";
390 // SideB.Output(out);
391 // *out << endl;
392
393 SideC.CopyVector(&left->second.second->x);
394 SideC.SubtractVector(&right->second.second->x);
395 SideC.ProjectOntoPlane(&AxisVector);
396 // *out << "SideC: ";
397 // SideC.Output(out);
398 // *out << endl;
399
400 SideH.CopyVector(&runner->second.second->x);
401 SideH.SubtractVector(MolCenter);
402 SideH.ProjectOntoPlane(&AxisVector);
403 // *out << "SideH: ";
404 // SideH.Output(out);
405 // *out << endl;
406
407 // calculate each length
408 double a = SideA.Norm();
409 //double b = SideB.Norm();
410 //double c = SideC.Norm();
411 double h = SideH.Norm();
412 // calculate the angles
413 double alpha = SideA.Angle(&SideH);
414 double beta = SideA.Angle(&SideC);
415 double gamma = SideB.Angle(&SideH);
416 double delta = SideC.Angle(&SideH);
417 double MinDistance = a * sin(beta) / (sin(delta)) * (((alpha < M_PI / 2.) || (gamma < M_PI / 2.)) ? 1. : -1.);
418 //*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;
419 *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;
420 if ((fabs(h / fabs(h) - MinDistance / fabs(MinDistance)) < MYEPSILON) && ((h - MinDistance)) < -MYEPSILON) {
421 // throw out point
422 *out << Verbose(1) << "Throwing out " << *runner->second.second << "." << endl;
423 BoundaryPoints[axis].erase(runner);
424 flag = true;
425 }
426 }
427 }
428 } while (flag);
429 }
430 delete(MolCenter);
431 return BoundaryPoints;
432};
433
434/** Tesselates the convex boundary by finding all boundary points.
435 * \param *out output stream for debugging
436 * \param *mol molecule structure with Atom's and Bond's
437 * \param *TesselStruct Tesselation filled with points, lines and triangles on boundary on return
438 * \param *LCList atoms in LinkedCell list
439 * \param *filename filename prefix for output of vertex data
440 * \return *TesselStruct is filled with convex boundary and tesselation is stored under \a *filename.
441 */
442void FindConvexBorder(ofstream *out, molecule* mol, class LinkedCell *LCList, const char *filename)
443{
444 bool BoundaryFreeFlag = false;
445 Boundaries *BoundaryPoints = NULL;
446
447 cout << Verbose(1) << "Begin of FindConvexBorder" << endl;
448
449 if (mol->TesselStruct != NULL) // free if allocated
450 delete(mol->TesselStruct);
451 mol->TesselStruct = new class Tesselation;
452
453 // 1. Find all points on the boundary
454 if (BoundaryPoints == NULL) {
455 BoundaryFreeFlag = true;
456 BoundaryPoints = GetBoundaryPoints(out, mol);
457 } else {
458 *out << Verbose(1) << "Using given boundary points set." << endl;
459 }
460
461// printing all inserted for debugging
462 for (int axis=0; axis < NDIM; axis++)
463 {
464 *out << Verbose(2) << "Printing list of candidates for axis " << axis << " which we have inserted so far." << endl;
465 int i=0;
466 for(Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++) {
467 if (runner != BoundaryPoints[axis].begin())
468 *out << ", " << i << ": " << *runner->second.second;
469 else
470 *out << i << ": " << *runner->second.second;
471 i++;
472 }
473 *out << endl;
474 }
475
476 // 2. fill the boundary point list
477 for (int axis = 0; axis < NDIM; axis++)
478 for (Boundaries::iterator runner = BoundaryPoints[axis].begin(); runner != BoundaryPoints[axis].end(); runner++)
479 if (!mol->TesselStruct->AddBoundaryPoint(runner->second.second, 0))
480 *out << Verbose(3) << "WARNING: Point " << *(runner->second.second) << " is already present!" << endl;
481
482 *out << Verbose(2) << "I found " << mol->TesselStruct->PointsOnBoundaryCount << " points on the convex boundary." << endl;
483 // now we have the whole set of edge points in the BoundaryList
484
485 // listing for debugging
486 // *out << Verbose(1) << "Listing PointsOnBoundary:";
487 // for(PointMap::iterator runner = PointsOnBoundary.begin(); runner != PointsOnBoundary.end(); runner++) {
488 // *out << " " << *runner->second;
489 // }
490 // *out << endl;
491
492 // 3a. guess starting triangle
493 mol->TesselStruct->GuessStartingTriangle(out);
494
495 // 3b. go through all lines, that are not yet part of two triangles (only of one so far)
496 mol->TesselStruct->TesselateOnBoundary(out, mol);
497
498 // 3c. check whether all atoms lay inside the boundary, if not, add to boundary points, segment triangle into three with the new point
499 if (!mol->TesselStruct->InsertStraddlingPoints(out, mol, LCList))
500 *out << Verbose(1) << "Insertion of straddling points failed!" << endl;
501
502 *out << Verbose(2) << "I created " << mol->TesselStruct->TrianglesOnBoundary.size() << " intermediate triangles with " << mol->TesselStruct->LinesOnBoundary.size() << " lines and " << mol->TesselStruct->PointsOnBoundary.size() << " points." << endl;
503
504 // 4. Store triangles in tecplot file
505 if (filename != NULL) {
506 if (DoTecplotOutput) {
507 string OutputName(filename);
508 OutputName.append("_intermed");
509 OutputName.append(TecplotSuffix);
510 ofstream *tecplot = new ofstream(OutputName.c_str());
511 WriteTecplotFile(out, tecplot, mol->TesselStruct, mol, 0);
512 tecplot->close();
513 delete(tecplot);
514 }
515 if (DoRaster3DOutput) {
516 string OutputName(filename);
517 OutputName.append("_intermed");
518 OutputName.append(Raster3DSuffix);
519 ofstream *rasterplot = new ofstream(OutputName.c_str());
520 WriteRaster3dFile(out, rasterplot, mol->TesselStruct, mol);
521 rasterplot->close();
522 delete(rasterplot);
523 }
524 }
525
526 // 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
527 bool AllConvex;
528 class BoundaryLineSet *line = NULL;
529 do {
530 AllConvex = true;
531 for (LineMap::iterator LineRunner = mol->TesselStruct->LinesOnBoundary.begin(); LineRunner != mol->TesselStruct->LinesOnBoundary.end(); LineRunner++) {
532 line = LineRunner->second;
533 *out << Verbose(1) << "INFO: Current line is " << *line << "." << endl;
534 if (!line->CheckConvexityCriterion(out)) {
535 *out << Verbose(1) << "... line " << *line << " is concave, flipping it." << endl;
536
537 // flip the line
538 if (!mol->TesselStruct->PickFarthestofTwoBaselines(out, line))
539 *out << Verbose(1) << "ERROR: Correction of concave baselines failed!" << endl;
540 else
541 *out << Verbose(1) << "INFO: Correction of concave baselines worked." << endl;
542 }
543 }
544 } while (!AllConvex);
545
546 // 3e. we need another correction here, for TesselPoints that are below the surface (i.e. have an odd number of concave triangles surrounding it)
547// if (!mol->TesselStruct->CorrectConcaveTesselPoints(out))
548// *out << Verbose(1) << "Correction of concave tesselpoints failed!" << endl;
549
550 *out << Verbose(2) << "I created " << mol->TesselStruct->TrianglesOnBoundary.size() << " triangles with " << mol->TesselStruct->LinesOnBoundary.size() << " lines and " << mol->TesselStruct->PointsOnBoundary.size() << " points." << endl;
551
552 // 4. Store triangles in tecplot file
553 if (filename != NULL) {
554 if (DoTecplotOutput) {
555 string OutputName(filename);
556 OutputName.append(TecplotSuffix);
557 ofstream *tecplot = new ofstream(OutputName.c_str());
558 WriteTecplotFile(out, tecplot, mol->TesselStruct, mol, 0);
559 tecplot->close();
560 delete(tecplot);
561 }
562 if (DoRaster3DOutput) {
563 string OutputName(filename);
564 OutputName.append(Raster3DSuffix);
565 ofstream *rasterplot = new ofstream(OutputName.c_str());
566 WriteRaster3dFile(out, rasterplot, mol->TesselStruct, mol);
567 rasterplot->close();
568 delete(rasterplot);
569 }
570 }
571
572
573 // free reference lists
574 if (BoundaryFreeFlag)
575 delete[] (BoundaryPoints);
576
577 cout << Verbose(1) << "End of FindConvexBorder" << endl;
578};
579
580/** Creates a convex envelope from a given non-convex one.
581 * -# First step, remove concave spots, i.e. singular "dents"
582 * -# We go through all PointsOnBoundary.
583 * -# We CheckConvexityCriterion() for all its lines.
584 * -# If all its lines are concave, it cannot be on the convex envelope.
585 * -# Hence, we remove it and re-create all its triangles from its getCircleOfConnectedPoints()
586 * -# We calculate the additional volume.
587 * -# We go over all lines until none yields a concavity anymore.
588 * -# Second step, remove concave lines, i.e. line-shape "dents"
589 * -# We go through all LinesOnBoundary
590 * -# We CheckConvexityCriterion()
591 * -# If it returns concave, we flip the line in this quadrupel of points (abusing the degeneracy of the tesselation)
592 * -# We CheckConvexityCriterion(),
593 * -# if it's concave, we continue
594 * -# if not, we mark an error and stop
595 * Note: This routine - for free - calculates the difference in volume between convex and
596 * non-convex envelope, as the former is easy to calculate - VolumeOfConvexEnvelope() - it
597 * can be used to compute volumes of arbitrary shapes.
598 * \param *out output stream for debugging
599 * \param *TesselStruct non-convex envelope, is changed in return!
600 * \param *mol molecule
601 * \param *filename name of file
602 * \return volume difference between the non- and the created convex envelope
603 */
604double ConvexizeNonconvexEnvelope(ofstream *out, class Tesselation *TesselStruct, molecule *mol, char *filename)
605{
606 double volume = 0;
607 class BoundaryPointSet *point = NULL;
608 class BoundaryLineSet *line = NULL;
609 bool Concavity;
610 PointMap::iterator PointRunner, PointAdvance;
611 LineMap::iterator LineRunner, LineAdvance;
612 TriangleMap::iterator TriangleRunner, TriangleAdvance;
613
614 *out << Verbose(0) << "Begin of ConvexizeNonconvexEnvelope" << endl;
615
616 // check whether there is something to work on
617 if (TesselStruct == NULL) {
618 *out << Verbose(1) << "ERROR: TesselStruct is empty!" << endl;
619 return volume;
620 }
621
622 //CalculateConcavityPerBoundaryPoint(out, TesselStruct);
623 StoreTrianglesinFile(out, mol, filename, "-first");
624
625 // First step: RemovePointFromTesselatedSurface
626 do {
627 Concavity = false;
628 PointRunner = TesselStruct->PointsOnBoundary.begin();
629 PointAdvance = PointRunner; // we need an advanced point, as the PointRunner might get removed
630 while (PointRunner != TesselStruct->PointsOnBoundary.end()) {
631 PointAdvance++;
632 point = PointRunner->second;
633 *out << Verbose(1) << "INFO: Current point is " << *point << "." << endl;
634 for (LineMap::iterator LineRunner = point->lines.begin(); LineRunner != point->lines.end(); LineRunner++) {
635 line = LineRunner->second;
636 *out << Verbose(2) << "INFO: Current line of point " << *point << " is " << *line << "." << endl;
637 }
638 if (!line->CheckConvexityCriterion(out)) {
639 *out << Verbose(1) << "... point " << *point << " cannot be on convex envelope." << endl;
640 // remove the point
641 Concavity = true;
642 TesselStruct->RemovePointFromTesselatedSurface(out, point);
643 }
644 PointRunner = PointAdvance;
645 }
646
647 //CalculateConcavityPerBoundaryPoint(out, TesselStruct);
648 //StoreTrianglesinFile(out, mol, filename, "-second");
649
650 // second step: PickFarthestofTwoBaselines
651 LineRunner = TesselStruct->LinesOnBoundary.begin();
652 LineAdvance = LineRunner; // we need an advanced line, as the LineRunner might get removed
653 while (LineRunner != TesselStruct->LinesOnBoundary.end()) {
654 LineAdvance++;
655 line = LineRunner->second;
656 *out << Verbose(1) << "INFO: Picking farthest baseline for line is " << *line << "." << endl;
657 // take highest of both lines
658 if (TesselStruct->IsConvexRectangle(out, line) == NULL) {
659 TesselStruct->PickFarthestofTwoBaselines(out, line);
660 Concavity = true;
661 }
662 LineRunner = LineAdvance;
663 }
664 } while (Concavity);
665 CalculateConcavityPerBoundaryPoint(out, TesselStruct);
666 StoreTrianglesinFile(out, mol, filename, "-third");
667
668 // third step: IsConvexRectangle
669 LineRunner = TesselStruct->LinesOnBoundary.begin();
670 LineAdvance = LineRunner; // we need an advanced line, as the LineRunner might get removed
671 while (LineRunner != TesselStruct->LinesOnBoundary.end()) {
672 LineAdvance++;
673 line = LineRunner->second;
674 *out << Verbose(1) << "INFO: Current line is " << *line << "." << endl;
675 //if (LineAdvance != TesselStruct->LinesOnBoundary.end())
676 //*out << Verbose(1) << "INFO: Next line will be " << *(LineAdvance->second) << "." << endl;
677 if (!line->CheckConvexityCriterion(out)) {
678 *out << Verbose(1) << "... line " << *line << " is concave, flipping it." << endl;
679
680 // take highest of both lines
681 point = TesselStruct->IsConvexRectangle(out, line);
682 if (point != NULL)
683 TesselStruct->RemovePointFromTesselatedSurface(out, point);
684 }
685 LineRunner = LineAdvance;
686 }
687
688 CalculateConcavityPerBoundaryPoint(out, TesselStruct);
689 StoreTrianglesinFile(out, mol, filename, "-fourth");
690
691 // end
692 *out << Verbose(0) << "End of ConvexizeNonconvexEnvelope" << endl;
693 return volume;
694};
695
696/** Calculates the concavity for each of the BoundaryPointSet's in a Tesselation.
697 * Sets BoundaryPointSet::value equal to the number of connected lines that are not convex.
698 * \param *out output stream for debugging
699 * \param *TesselStruct pointer to Tesselation structure
700 */
701void CalculateConcavityPerBoundaryPoint(ofstream *out, class Tesselation *TesselStruct)
702{
703 class BoundaryPointSet *point = NULL;
704 class BoundaryLineSet *line = NULL;
705 // calculate remaining concavity
706 for (PointMap::iterator PointRunner = TesselStruct->PointsOnBoundary.begin(); PointRunner != TesselStruct->PointsOnBoundary.end(); PointRunner++) {
707 point = PointRunner->second;
708 *out << Verbose(1) << "INFO: Current point is " << *point << "." << endl;
709 point->value = 0;
710 for (LineMap::iterator LineRunner = point->lines.begin(); LineRunner != point->lines.end(); LineRunner++) {
711 line = LineRunner->second;
712 *out << Verbose(2) << "INFO: Current line of point " << *point << " is " << *line << "." << endl;
713 if (!line->CheckConvexityCriterion(out))
714 point->value += 1;
715 }
716 }
717};
718
719/** Stores triangles to file.
720 * \param *out output stream for debugging
721 * \param *mol molecule with atoms and bonds
722 * \param *filename prefix of filename
723 * \param *extraSuffix intermediate suffix
724 */
725void StoreTrianglesinFile(ofstream *out, molecule *mol, const char *filename, const char *extraSuffix)
726{
727 // 4. Store triangles in tecplot file
728 if (filename != NULL) {
729 if (DoTecplotOutput) {
730 string OutputName(filename);
731 OutputName.append(extraSuffix);
732 OutputName.append(TecplotSuffix);
733 ofstream *tecplot = new ofstream(OutputName.c_str());
734 WriteTecplotFile(out, tecplot, mol->TesselStruct, mol, 0);
735 tecplot->close();
736 delete(tecplot);
737 }
738 if (DoRaster3DOutput) {
739 string OutputName(filename);
740 OutputName.append(extraSuffix);
741 OutputName.append(Raster3DSuffix);
742 ofstream *rasterplot = new ofstream(OutputName.c_str());
743 WriteRaster3dFile(out, rasterplot, mol->TesselStruct, mol);
744 rasterplot->close();
745 delete(rasterplot);
746 }
747 }
748};
749
750/** Determines the volume of a cluster.
751 * Determines first the convex envelope, then tesselates it and calculates its volume.
752 * \param *out output stream for debugging
753 * \param *TesselStruct Tesselation filled with points, lines and triangles on boundary on return
754 * \param *configuration needed for path to store convex envelope file
755 * \return determined volume of the cluster in cubed config:GetIsAngstroem()
756 */
757double VolumeOfConvexEnvelope(ofstream *out, class Tesselation *TesselStruct, class config *configuration)
758{
759 bool IsAngstroem = configuration->GetIsAngstroem();
760 double volume = 0.;
761 double PyramidVolume = 0.;
762 double G, h;
763 Vector x, y;
764 double a, b, c;
765
766 // 6a. Every triangle forms a pyramid with the center of gravity as its peak, sum up the volumes
767 *out << Verbose(1)
768 << "Calculating the volume of the pyramids formed out of triangles and center of gravity."
769 << endl;
770 for (TriangleMap::iterator runner = TesselStruct->TrianglesOnBoundary.begin(); runner != TesselStruct->TrianglesOnBoundary.end(); runner++)
771 { // go through every triangle, calculate volume of its pyramid with CoG as peak
772 x.CopyVector(runner->second->endpoints[0]->node->node);
773 x.SubtractVector(runner->second->endpoints[1]->node->node);
774 y.CopyVector(runner->second->endpoints[0]->node->node);
775 y.SubtractVector(runner->second->endpoints[2]->node->node);
776 a = sqrt(runner->second->endpoints[0]->node->node->DistanceSquared(runner->second->endpoints[1]->node->node));
777 b = sqrt(runner->second->endpoints[0]->node->node->DistanceSquared(runner->second->endpoints[2]->node->node));
778 c = sqrt(runner->second->endpoints[2]->node->node->DistanceSquared(runner->second->endpoints[1]->node->node));
779 G = sqrt(((a + b + c) * (a + b + c) - 2 * (a * a + b * b + c * c)) / 16.); // area of tesselated triangle
780 x.MakeNormalVector(runner->second->endpoints[0]->node->node, runner->second->endpoints[1]->node->node, runner->second->endpoints[2]->node->node);
781 x.Scale(runner->second->endpoints[1]->node->node->ScalarProduct(&x));
782 h = x.Norm(); // distance of CoG to triangle
783 PyramidVolume = (1. / 3.) * G * h; // this formula holds for _all_ pyramids (independent of n-edge base or (not) centered peak)
784 *out << Verbose(2) << "Area of triangle is " << G << " "
785 << (IsAngstroem ? "angstrom" : "atomiclength") << "^2, height is "
786 << h << " and the volume is " << PyramidVolume << " "
787 << (IsAngstroem ? "angstrom" : "atomiclength") << "^3." << endl;
788 volume += PyramidVolume;
789 }
790 *out << Verbose(0) << "RESULT: The summed volume is " << setprecision(10)
791 << volume << " " << (IsAngstroem ? "angstrom" : "atomiclength") << "^3."
792 << endl;
793
794 return volume;
795}
796;
797
798/** Creates multiples of the by \a *mol given cluster and suspends them in water with a given final density.
799 * We get cluster volume by VolumeOfConvexEnvelope() and its diameters by GetDiametersOfCluster()
800 * \param *out output stream for debugging
801 * \param *configuration needed for path to store convex envelope file
802 * \param *mol molecule structure representing the cluster
803 * \param ClusterVolume guesstimated cluster volume, if equal 0 we used VolumeOfConvexEnvelope() instead.
804 * \param celldensity desired average density in final cell
805 */
806void
807PrepareClustersinWater(ofstream *out, config *configuration, molecule *mol,
808 double ClusterVolume, double celldensity)
809{
810 // transform to PAS
811 mol->PrincipalAxisSystem(out, true);
812
813 // some preparations beforehand
814 bool IsAngstroem = configuration->GetIsAngstroem();
815 Boundaries *BoundaryPoints = GetBoundaryPoints(out, mol);
816 class Tesselation *TesselStruct = NULL;
817 LinkedCell LCList(mol, 10.);
818 FindConvexBorder(out, mol, &LCList, NULL);
819 double clustervolume;
820 if (ClusterVolume == 0)
821 clustervolume = VolumeOfConvexEnvelope(out, TesselStruct, configuration);
822 else
823 clustervolume = ClusterVolume;
824 double *GreatestDiameter = GetDiametersOfCluster(out, BoundaryPoints, mol, IsAngstroem);
825 Vector BoxLengths;
826 int repetition[NDIM] =
827 { 1, 1, 1 };
828 int TotalNoClusters = 1;
829 for (int i = 0; i < NDIM; i++)
830 TotalNoClusters *= repetition[i];
831
832 // sum up the atomic masses
833 double totalmass = 0.;
834 atom *Walker = mol->start;
835 while (Walker->next != mol->end)
836 {
837 Walker = Walker->next;
838 totalmass += Walker->type->mass;
839 }
840 *out << Verbose(0) << "RESULT: The summed mass is " << setprecision(10)
841 << totalmass << " atomicmassunit." << endl;
842
843 *out << Verbose(0) << "RESULT: The average density is " << setprecision(10)
844 << totalmass / clustervolume << " atomicmassunit/"
845 << (IsAngstroem ? "angstrom" : "atomiclength") << "^3." << endl;
846
847 // solve cubic polynomial
848 *out << Verbose(1) << "Solving equidistant suspension in water problem ..."
849 << endl;
850 double cellvolume;
851 if (IsAngstroem)
852 cellvolume = (TotalNoClusters * totalmass / SOLVENTDENSITY_A - (totalmass
853 / clustervolume)) / (celldensity - 1);
854 else
855 cellvolume = (TotalNoClusters * totalmass / SOLVENTDENSITY_a0 - (totalmass
856 / clustervolume)) / (celldensity - 1);
857 *out << Verbose(1) << "Cellvolume needed for a density of " << celldensity
858 << " g/cm^3 is " << cellvolume << " " << (IsAngstroem ? "angstrom"
859 : "atomiclength") << "^3." << endl;
860
861 double minimumvolume = TotalNoClusters * (GreatestDiameter[0]
862 * GreatestDiameter[1] * GreatestDiameter[2]);
863 *out << Verbose(1)
864 << "Minimum volume of the convex envelope contained in a rectangular box is "
865 << minimumvolume << " atomicmassunit/" << (IsAngstroem ? "angstrom"
866 : "atomiclength") << "^3." << endl;
867 if (minimumvolume > cellvolume)
868 {
869 cerr << Verbose(0)
870 << "ERROR: the containing box already has a greater volume than the envisaged cell volume!"
871 << endl;
872 cout << Verbose(0)
873 << "Setting Box dimensions to minimum possible, the greatest diameters."
874 << endl;
875 for (int i = 0; i < NDIM; i++)
876 BoxLengths.x[i] = GreatestDiameter[i];
877 mol->CenterEdge(out, &BoxLengths);
878 }
879 else
880 {
881 BoxLengths.x[0] = (repetition[0] * GreatestDiameter[0] + repetition[1]
882 * GreatestDiameter[1] + repetition[2] * GreatestDiameter[2]);
883 BoxLengths.x[1] = (repetition[0] * repetition[1] * GreatestDiameter[0]
884 * GreatestDiameter[1] + repetition[0] * repetition[2]
885 * GreatestDiameter[0] * GreatestDiameter[2] + repetition[1]
886 * repetition[2] * GreatestDiameter[1] * GreatestDiameter[2]);
887 BoxLengths.x[2] = minimumvolume - cellvolume;
888 double x0 = 0., x1 = 0., x2 = 0.;
889 if (gsl_poly_solve_cubic(BoxLengths.x[0], BoxLengths.x[1],
890 BoxLengths.x[2], &x0, &x1, &x2) == 1) // either 1 or 3 on return
891 *out << Verbose(0) << "RESULT: The resulting spacing is: " << x0
892 << " ." << endl;
893 else
894 {
895 *out << Verbose(0) << "RESULT: The resulting spacings are: " << x0
896 << " and " << x1 << " and " << x2 << " ." << endl;
897 x0 = x2; // sorted in ascending order
898 }
899
900 cellvolume = 1;
901 for (int i = 0; i < NDIM; i++)
902 {
903 BoxLengths.x[i] = repetition[i] * (x0 + GreatestDiameter[i]);
904 cellvolume *= BoxLengths.x[i];
905 }
906
907 // set new box dimensions
908 *out << Verbose(0) << "Translating to box with these boundaries." << endl;
909 mol->SetBoxDimension(&BoxLengths);
910 mol->CenterInBox((ofstream *) &cout);
911 }
912 // update Box of atoms by boundary
913 mol->SetBoxDimension(&BoxLengths);
914 *out << Verbose(0) << "RESULT: The resulting cell dimensions are: "
915 << BoxLengths.x[0] << " and " << BoxLengths.x[1] << " and "
916 << BoxLengths.x[2] << " with total volume of " << cellvolume << " "
917 << (IsAngstroem ? "angstrom" : "atomiclength") << "^3." << endl;
918}
919;
920
921
922/** Fills the empty space of the simulation box with water/
923 * \param *out output stream for debugging
924 * \param *List list of molecules already present in box
925 * \param *TesselStruct contains tesselated surface
926 * \param *filler molecule which the box is to be filled with
927 * \param configuration contains box dimensions
928 * \param distance[NDIM] distance between filling molecules in each direction
929 * \param RandAtomDisplacement maximum distance for random displacement per atom
930 * \param RandMolDisplacement maximum distance for random displacement per filler molecule
931 * \param DoRandomRotation true - do random rotiations, false - don't
932 * \return *mol pointer to new molecule with filled atoms
933 */
934molecule * FillBoxWithMolecule(ofstream *out, MoleculeListClass *List, molecule *filler, config &configuration, double distance[NDIM], double RandomAtomDisplacement, double RandomMolDisplacement, bool DoRandomRotation)
935{
936 molecule *Filling = new molecule(filler->elemente);
937 Vector CurrentPosition;
938 int N[NDIM];
939 int n[NDIM];
940 double *M = filler->ReturnFullMatrixforSymmetric(filler->cell_size);
941 double Rotations[NDIM*NDIM];
942 Vector AtomTranslations;
943 Vector FillerTranslations;
944 Vector FillerDistance;
945 double FillIt = false;
946 atom *Walker = NULL;
947 bond *Binder = NULL;
948 int i;
949 LinkedCell *LCList[List->ListOfMolecules.size()];
950
951 *out << Verbose(0) << "Begin of FillBoxWithMolecule" << endl;
952
953 i=0;
954 for (MoleculeList::iterator ListRunner = List->ListOfMolecules.begin(); ListRunner != List->ListOfMolecules.end(); ListRunner++) {
955 *out << Verbose(1) << "Pre-creating linked cell lists for molecule " << *ListRunner << "." << endl;
956 LCList[i] = new LinkedCell((*ListRunner), 5.); // get linked cell list
957 if ((*ListRunner)->TesselStruct == NULL) {
958 *out << Verbose(1) << "Pre-creating tesselation for molecule " << *ListRunner << "." << endl;
959 FindNonConvexBorder((ofstream *)&cout, (*ListRunner), LCList[i], NULL, 5.);
960 }
961 i++;
962 }
963
964 // Center filler at origin
965 filler->CenterOrigin(out);
966 filler->Center.Zero();
967
968 filler->CountAtoms(out);
969 atom * CopyAtoms[filler->AtomCount];
970 int nr = 0;
971
972 // calculate filler grid in [0,1]^3
973 FillerDistance.Init(distance[0], distance[1], distance[2]);
974 FillerDistance.InverseMatrixMultiplication(M);
975 *out << Verbose(1) << "INFO: Grid steps are ";
976 for(int i=0;i<NDIM;i++) {
977 N[i] = (int) ceil(1./FillerDistance.x[i]);
978 *out << N[i];
979 if (i != NDIM-1)
980 *out<< ", ";
981 else
982 *out << "." << endl;
983 }
984
985 // go over [0,1]^3 filler grid
986 for (n[0] = 0; n[0] < N[0]; n[0]++)
987 for (n[1] = 0; n[1] < N[1]; n[1]++)
988 for (n[2] = 0; n[2] < N[2]; n[2]++) {
989 // calculate position of current grid vector in untransformed box
990 CurrentPosition.Init((double)n[0]/(double)N[0], (double)n[1]/(double)N[1], (double)n[2]/(double)N[2]);
991 CurrentPosition.MatrixMultiplication(M);
992 *out << Verbose(2) << "INFO: Current Position is " << CurrentPosition << "." << endl;
993 // Check whether point is in- or outside
994 FillIt = true;
995 i=0;
996 for (MoleculeList::iterator ListRunner = List->ListOfMolecules.begin(); ListRunner != List->ListOfMolecules.end(); ListRunner++) {
997 // get linked cell list
998 if ((*ListRunner)->TesselStruct == NULL) {
999 *out << Verbose(1) << "ERROR: TesselStruct of " << (*ListRunner) << " is NULL. Didn't we pre-create it?" << endl;
1000 FillIt = false;
1001 } else
1002 FillIt = FillIt && (!(*ListRunner)->TesselStruct->IsInnerPoint(out, CurrentPosition, LCList[i++]));
1003 }
1004
1005 if (FillIt) {
1006 // fill in Filler
1007 *out << Verbose(2) << "Space at " << CurrentPosition << " is unoccupied by any molecule, filling in." << endl;
1008
1009 // create molecule random translation vector ...
1010 for (int i=0;i<NDIM;i++)
1011 FillerTranslations.x[i] = RandomMolDisplacement*(rand()/(RAND_MAX/2.) - 1.);
1012 *out << Verbose(3) << "INFO: Translating this filler by " << FillerTranslations << "." << endl;
1013
1014 // go through all atoms
1015 nr=0;
1016 Walker = filler->start;
1017 while (Walker->next != filler->end) {
1018 Walker = Walker->next;
1019 // copy atom ...
1020 CopyAtoms[Walker->nr] = new atom(Walker);
1021
1022 // create atomic random translation vector ...
1023 for (int i=0;i<NDIM;i++)
1024 AtomTranslations.x[i] = RandomAtomDisplacement*(rand()/(RAND_MAX/2.) - 1.);
1025
1026 // ... and rotation matrix
1027 if (DoRandomRotation) {
1028 double phi[NDIM];
1029 for (int i=0;i<NDIM;i++) {
1030 phi[i] = rand()/(RAND_MAX/(2.*M_PI));
1031 }
1032
1033 Rotations[0] = cos(phi[0]) *cos(phi[2]) + (sin(phi[0])*sin(phi[1])*sin(phi[2]));
1034 Rotations[3] = sin(phi[0]) *cos(phi[2]) - (cos(phi[0])*sin(phi[1])*sin(phi[2]));
1035 Rotations[6] = cos(phi[1])*sin(phi[2]) ;
1036 Rotations[1] = - sin(phi[0])*cos(phi[1]) ;
1037 Rotations[4] = cos(phi[0])*cos(phi[1]) ;
1038 Rotations[7] = sin(phi[1]) ;
1039 Rotations[3] = - cos(phi[0]) *sin(phi[2]) + (sin(phi[0])*sin(phi[1])*cos(phi[2]));
1040 Rotations[5] = - sin(phi[0]) *sin(phi[2]) - (cos(phi[0])*sin(phi[1])*cos(phi[2]));
1041 Rotations[8] = cos(phi[1])*cos(phi[2]) ;
1042 }
1043
1044 // ... and put at new position
1045 if (DoRandomRotation)
1046 CopyAtoms[Walker->nr]->x.MatrixMultiplication(Rotations);
1047 CopyAtoms[Walker->nr]->x.AddVector(&AtomTranslations);
1048 CopyAtoms[Walker->nr]->x.AddVector(&FillerTranslations);
1049 CopyAtoms[Walker->nr]->x.AddVector(&CurrentPosition);
1050
1051 // insert into Filling
1052
1053 // FIXME: gives completely different results if CopyAtoms[..] used instead of Walker, why???
1054 *out << Verbose(4) << "Filling atom " << *Walker << ", translated to " << AtomTranslations << ", at final position is " << (CopyAtoms[Walker->nr]->x) << "." << endl;
1055 Filling->AddAtom(CopyAtoms[Walker->nr]);
1056 }
1057
1058 // go through all bonds and add as well
1059 Binder = filler->first;
1060 while(Binder->next != filler->last) {
1061 Binder = Binder->next;
1062 *out << Verbose(3) << "Adding Bond between " << *CopyAtoms[Binder->leftatom->nr] << " and " << *CopyAtoms[Binder->rightatom->nr]<< "." << endl;
1063 Filling->AddBond(CopyAtoms[Binder->leftatom->nr], CopyAtoms[Binder->rightatom->nr], Binder->BondDegree);
1064 }
1065 } else {
1066 // leave empty
1067 *out << Verbose(2) << "Space at " << CurrentPosition << " is occupied." << endl;
1068 }
1069 }
1070 *out << Verbose(0) << "End of FillBoxWithMolecule" << endl;
1071
1072 return Filling;
1073};
1074
1075
1076/** Tesselates the non convex boundary by rolling a virtual sphere along the surface of the molecule.
1077 * \param *out output stream for debugging
1078 * \param *mol molecule structure with Atom's and Bond's
1079 * \param *Tess Tesselation filled with points, lines and triangles on boundary on return
1080 * \param *LCList atoms in LinkedCell list
1081 * \param *filename filename prefix for output of vertex data
1082 * \para RADIUS radius of the virtual sphere
1083 */
1084void FindNonConvexBorder(ofstream *out, molecule* mol, class LinkedCell *LCList, const char *filename, const double RADIUS)
1085{
1086 int N = 0;
1087 bool freeLC = false;
1088 ofstream *tempstream = NULL;
1089 char NumberName[255];
1090 int TriangleFilesWritten = 0;
1091
1092 *out << Verbose(1) << "Entering search for non convex hull. " << endl;
1093 if (mol->TesselStruct == NULL) {
1094 *out << Verbose(1) << "Allocating Tesselation struct ..." << endl;
1095 mol->TesselStruct = new Tesselation;
1096 } else {
1097 delete(mol->TesselStruct);
1098 *out << Verbose(1) << "Re-Allocating Tesselation struct ..." << endl;
1099 mol->TesselStruct = new Tesselation;
1100 }
1101 LineMap::iterator baseline;
1102 LineMap::iterator testline;
1103 *out << Verbose(0) << "Begin of FindNonConvexBorder\n";
1104 bool flag = false; // marks whether we went once through all baselines without finding any without two triangles
1105 bool failflag = false;
1106
1107 if (LCList == NULL) {
1108 LCList = new LinkedCell(mol, 2.*RADIUS);
1109 freeLC = true;
1110 }
1111
1112 mol->TesselStruct->FindStartingTriangle(out, RADIUS, LCList);
1113
1114 baseline = mol->TesselStruct->LinesOnBoundary.begin();
1115 // the outward most line is dangerous, as we may end up with wrapping up the starting triangle, hence
1116 // terminating the algorithm too early.
1117 if (baseline != mol->TesselStruct->LinesOnBoundary.end()) // skip first line as it its the outwardmost!
1118 baseline++;
1119 while ((baseline != mol->TesselStruct->LinesOnBoundary.end()) || (flag)) {
1120 if (baseline->second->triangles.size() == 1) {
1121 failflag = mol->TesselStruct->FindNextSuitableTriangle(out, *(baseline->second), *(((baseline->second->triangles.begin()))->second), RADIUS, N, LCList); //the line is there, so there is a triangle, but only one.
1122 flag = flag || failflag;
1123 if (!failflag)
1124 cerr << "WARNING: FindNextSuitableTriangle failed." << endl;
1125 // write temporary envelope
1126 if ((DoSingleStepOutput && (mol->TesselStruct->TrianglesOnBoundaryCount % SingleStepWidth == 0))) { // if we have a new triangle and want to output each new triangle configuration
1127 TriangleMap::iterator runner = mol->TesselStruct->TrianglesOnBoundary.end();
1128 runner--;
1129 class BoundaryTriangleSet *triangle = runner->second;
1130 if (triangle != NULL) {
1131 sprintf(NumberName, "-%04d-%s_%s_%s", TriangleFilesWritten, triangle->endpoints[0]->node->Name, triangle->endpoints[1]->node->Name, triangle->endpoints[2]->node->Name);
1132 if (DoTecplotOutput) {
1133 string NameofTempFile(filename);
1134 NameofTempFile.append(NumberName);
1135 for(size_t npos = NameofTempFile.find_first_of(' '); npos != string::npos; npos = NameofTempFile.find(' ', npos))
1136 NameofTempFile.erase(npos, 1);
1137 NameofTempFile.append(TecplotSuffix);
1138 *out << Verbose(1) << "Writing temporary non convex hull to file " << NameofTempFile << ".\n";
1139 tempstream = new ofstream(NameofTempFile.c_str(), ios::trunc);
1140 WriteTecplotFile(out, tempstream, mol->TesselStruct, mol, TriangleFilesWritten);
1141 tempstream->close();
1142 tempstream->flush();
1143 delete(tempstream);
1144 }
1145
1146 if (DoRaster3DOutput) {
1147 string NameofTempFile(filename);
1148 NameofTempFile.append(NumberName);
1149 for(size_t npos = NameofTempFile.find_first_of(' '); npos != string::npos; npos = NameofTempFile.find(' ', npos))
1150 NameofTempFile.erase(npos, 1);
1151 NameofTempFile.append(Raster3DSuffix);
1152 *out << Verbose(1) << "Writing temporary non convex hull to file " << NameofTempFile << ".\n";
1153 tempstream = new ofstream(NameofTempFile.c_str(), ios::trunc);
1154 WriteRaster3dFile(out, tempstream, mol->TesselStruct, mol);
1155 // // include the current position of the virtual sphere in the temporary raster3d file
1156 // // make the circumsphere's center absolute again
1157 // helper.CopyVector(BaseRay->endpoints[0]->node->node);
1158 // helper.AddVector(BaseRay->endpoints[1]->node->node);
1159 // helper.Scale(0.5);
1160 // (*it)->OptCenter.AddVector(&helper);
1161 // Vector *center = mol->DetermineCenterOfAll(out);
1162 // (*it)->OptCenter.SubtractVector(center);
1163 // delete(center);
1164 // // and add to file plus translucency object
1165 // *tempstream << "# current virtual sphere\n";
1166 // *tempstream << "8\n 25.0 0.6 -1.0 -1.0 -1.0 0.2 0 0 0 0\n";
1167 // *tempstream << "2\n " << (*it)->OptCenter.x[0] << " "
1168 // << (*it)->OptCenter.x[1] << " " << (*it)->OptCenter.x[2]
1169 // << "\t" << RADIUS << "\t1 0 0\n";
1170 // *tempstream << "9\n terminating special property\n";
1171 tempstream->close();
1172 tempstream->flush();
1173 delete(tempstream);
1174 }
1175 }
1176 if (DoTecplotOutput || DoRaster3DOutput)
1177 TriangleFilesWritten++;
1178 }
1179 } else {
1180 //cout << Verbose(1) << "Line " << *baseline->second << " has " << baseline->second->triangles.size() << " triangles adjacent" << endl;
1181 if (baseline->second->triangles.size() != 2)
1182 *out << Verbose(1) << "ERROR: TESSELATION FINISHED WITH INVALID TRIANGLE COUNT!" << endl;
1183 }
1184
1185 N++;
1186 baseline++;
1187 if ((baseline == mol->TesselStruct->LinesOnBoundary.end()) && (flag)) {
1188 baseline = mol->TesselStruct->LinesOnBoundary.begin(); // restart if we reach end due to newly inserted lines
1189 flag = false;
1190 }
1191 }
1192
1193 // Purges surplus triangles.
1194 mol->TesselStruct->RemoveDegeneratedTriangles();
1195
1196 // write final envelope
1197 if (filename != 0) {
1198 *out << Verbose(1) << "Writing final tecplot file\n";
1199 if (DoTecplotOutput) {
1200 string OutputName(filename);
1201 OutputName.append(TecplotSuffix);
1202 ofstream *tecplot = new ofstream(OutputName.c_str());
1203 WriteTecplotFile(out, tecplot, mol->TesselStruct, mol, -1);
1204 tecplot->close();
1205 delete(tecplot);
1206 }
1207 if (DoRaster3DOutput) {
1208 string OutputName(filename);
1209 OutputName.append(Raster3DSuffix);
1210 ofstream *raster = new ofstream(OutputName.c_str());
1211 WriteRaster3dFile(out, raster, mol->TesselStruct, mol);
1212 raster->close();
1213 delete(raster);
1214 }
1215 } else {
1216 cerr << "ERROR: Could definitively not find all necessary triangles!" << endl;
1217 }
1218
1219 cout << Verbose(2) << "Check: List of Baselines with not two connected triangles:" << endl;
1220 int counter = 0;
1221 for (testline = mol->TesselStruct->LinesOnBoundary.begin(); testline != mol->TesselStruct->LinesOnBoundary.end(); testline++) {
1222 if (testline->second->triangles.size() != 2) {
1223 cout << Verbose(2) << *testline->second << "\t" << testline->second->triangles.size() << endl;
1224 counter++;
1225 }
1226 }
1227 if (counter == 0)
1228 *out << Verbose(2) << "None." << endl;
1229
1230// // Tests the IsInnerAtom() function.
1231// Vector x (0, 0, 0);
1232// *out << Verbose(0) << "Point to check: " << x << endl;
1233// *out << Verbose(0) << "Check: IsInnerPoint() returns " << mol->TesselStruct->IsInnerPoint(out, x, LCList)
1234// << "for vector " << x << "." << endl;
1235// TesselPoint* a = mol->TesselStruct->PointsOnBoundary.begin()->second->node;
1236// *out << Verbose(0) << "Point to check: " << *a << " (on boundary)." << endl;
1237// *out << Verbose(0) << "Check: IsInnerAtom() returns " << mol->TesselStruct->IsInnerPoint(out, a, LCList)
1238// << "for atom " << a << " (on boundary)." << endl;
1239// LinkedNodes *List = NULL;
1240// for (int i=0;i<NDIM;i++) { // each axis
1241// LCList->n[i] = LCList->N[i]-1; // current axis is topmost cell
1242// for (LCList->n[(i+1)%NDIM]=0;LCList->n[(i+1)%NDIM]<LCList->N[(i+1)%NDIM];LCList->n[(i+1)%NDIM]++)
1243// for (LCList->n[(i+2)%NDIM]=0;LCList->n[(i+2)%NDIM]<LCList->N[(i+2)%NDIM];LCList->n[(i+2)%NDIM]++) {
1244// List = LCList->GetCurrentCell();
1245// //cout << Verbose(2) << "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << "." << endl;
1246// if (List != NULL) {
1247// for (LinkedNodes::iterator Runner = List->begin();Runner != List->end();Runner++) {
1248// if (mol->TesselStruct->PointsOnBoundary.find((*Runner)->nr) == mol->TesselStruct->PointsOnBoundary.end()) {
1249// a = *Runner;
1250// i=3;
1251// for (int j=0;j<NDIM;j++)
1252// LCList->n[j] = LCList->N[j];
1253// break;
1254// }
1255// }
1256// }
1257// }
1258// }
1259// *out << Verbose(0) << "Check: IsInnerPoint() returns " << mol->TesselStruct->IsInnerPoint(out, a, LCList)
1260// << "for atom " << a << " (inside)." << endl;
1261
1262 if (freeLC)
1263 delete(LCList);
1264 *out << Verbose(0) << "End of FindNonConvexBorder\n";
1265};
1266
1267/** Finds a hole of sufficient size in \a this molecule to embed \a *srcmol into it.
1268 * \param *out output stream for debugging
1269 * \param *srcmol molecule to embed into
1270 * \return *Vector new center of \a *srcmol for embedding relative to \a this
1271 */
1272Vector* molecule::FindEmbeddingHole(ofstream *out, molecule *srcmol)
1273{
1274 Vector *Center = new Vector;
1275 Center->Zero();
1276 // calculate volume/shape of \a *srcmol
1277
1278 // find embedding holes
1279
1280 // if more than one, let user choose
1281
1282 // return embedding center
1283 return Center;
1284};
1285
Note: See TracBrowser for help on using the repository browser.