source: src/ellipsoid.cpp@ 68f03d

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 68f03d was 8cbb97, checked in by Tillmann Crueger <crueger@…>, 15 years ago

Merge branch 'VectorRefactoring' into StructureRefactoring

Conflicts:

molecuilder/src/Legacy/oldmenu.cpp
molecuilder/src/Makefile.am
molecuilder/src/analysis_correlation.cpp
molecuilder/src/boundary.cpp
molecuilder/src/builder.cpp
molecuilder/src/config.cpp
molecuilder/src/ellipsoid.cpp
molecuilder/src/linkedcell.cpp
molecuilder/src/molecule.cpp
molecuilder/src/molecule_fragmentation.cpp
molecuilder/src/molecule_geometry.cpp
molecuilder/src/molecule_graph.cpp
molecuilder/src/moleculelist.cpp
molecuilder/src/tesselation.cpp
molecuilder/src/tesselationhelpers.cpp
molecuilder/src/unittests/AnalysisCorrelationToSurfaceUnitTest.cpp
molecuilder/src/unittests/bondgraphunittest.cpp
molecuilder/src/vector.cpp
molecuilder/src/vector.hpp

  • Property mode set to 100644
File size: 16.9 KB
Line 
1/*
2 * ellipsoid.cpp
3 *
4 * Created on: Jan 20, 2009
5 * Author: heber
6 */
7
8#include <gsl/gsl_multimin.h>
9#include <gsl/gsl_vector.h>
10
11#include <iomanip>
12
13#include <set>
14
15#include "boundary.hpp"
16#include "ellipsoid.hpp"
17#include "linkedcell.hpp"
18#include "log.hpp"
19#include "tesselation.hpp"
20#include "vector.hpp"
21#include "verbose.hpp"
22
23/** Determines squared distance for a given point \a x to surface of ellipsoid.
24 * \param x given point
25 * \param EllipsoidCenter center of ellipsoid
26 * \param EllipsoidLength[3] three lengths of half axis of ellipsoid
27 * \param EllipsoidAngle[3] three rotation angles of ellipsoid
28 * \return squared distance from point to surface
29 */
30double SquaredDistanceToEllipsoid(Vector &x, Vector &EllipsoidCenter, double *EllipsoidLength, double *EllipsoidAngle)
31{
32 Vector helper, RefPoint;
33 double distance = -1.;
34 double Matrix[NDIM*NDIM];
35 double InverseLength[3];
36 double psi,theta,phi; // euler angles in ZX'Z'' convention
37
38 //Log() << Verbose(3) << "Begin of SquaredDistanceToEllipsoid" << endl;
39
40 for(int i=0;i<3;i++)
41 InverseLength[i] = 1./EllipsoidLength[i];
42
43 // 1. translate coordinate system so that ellipsoid center is in origin
44 RefPoint = helper = x - EllipsoidCenter;
45 //Log() << Verbose(4) << "Translated given point is at " << RefPoint << "." << endl;
46
47 // 2. transform coordinate system by inverse of rotation matrix and of diagonal matrix
48 psi = EllipsoidAngle[0];
49 theta = EllipsoidAngle[1];
50 phi = EllipsoidAngle[2];
51 Matrix[0] = cos(psi)*cos(phi) - sin(psi)*cos(theta)*sin(phi);
52 Matrix[1] = -cos(psi)*sin(phi) - sin(psi)*cos(theta)*cos(phi);
53 Matrix[2] = sin(psi)*sin(theta);
54 Matrix[3] = sin(psi)*cos(phi) + cos(psi)*cos(theta)*sin(phi);
55 Matrix[4] = cos(psi)*cos(theta)*cos(phi) - sin(psi)*sin(phi);
56 Matrix[5] = -cos(psi)*sin(theta);
57 Matrix[6] = sin(theta)*sin(phi);
58 Matrix[7] = sin(theta)*cos(phi);
59 Matrix[8] = cos(theta);
60 helper.MatrixMultiplication(Matrix);
61 helper.ScaleAll(InverseLength);
62 //Log() << Verbose(4) << "Transformed RefPoint is at " << helper << "." << endl;
63
64 // 3. construct intersection point with unit sphere and ray between origin and x
65 helper.Normalize(); // is simply normalizes vector in distance direction
66 //Log() << Verbose(4) << "Transformed intersection is at " << helper << "." << endl;
67
68 // 4. transform back the constructed intersection point
69 psi = -EllipsoidAngle[0];
70 theta = -EllipsoidAngle[1];
71 phi = -EllipsoidAngle[2];
72 helper.ScaleAll(EllipsoidLength);
73 Matrix[0] = cos(psi)*cos(phi) - sin(psi)*cos(theta)*sin(phi);
74 Matrix[1] = -cos(psi)*sin(phi) - sin(psi)*cos(theta)*cos(phi);
75 Matrix[2] = sin(psi)*sin(theta);
76 Matrix[3] = sin(psi)*cos(phi) + cos(psi)*cos(theta)*sin(phi);
77 Matrix[4] = cos(psi)*cos(theta)*cos(phi) - sin(psi)*sin(phi);
78 Matrix[5] = -cos(psi)*sin(theta);
79 Matrix[6] = sin(theta)*sin(phi);
80 Matrix[7] = sin(theta)*cos(phi);
81 Matrix[8] = cos(theta);
82 helper.MatrixMultiplication(Matrix);
83 //Log() << Verbose(4) << "Intersection is at " << helper << "." << endl;
84
85 // 5. determine distance between backtransformed point and x
86 distance = RefPoint.DistanceSquared(helper);
87 //Log() << Verbose(4) << "Squared distance between intersection and RefPoint is " << distance << "." << endl;
88
89 return distance;
90 //Log() << Verbose(3) << "End of SquaredDistanceToEllipsoid" << endl;
91};
92
93/** structure for ellipsoid minimisation containing points to fit to.
94 */
95struct EllipsoidMinimisation {
96 int N; //!< dimension of vector set
97 Vector *x; //!< array of vectors
98};
99
100/** Sum of squared distance to ellipsoid to be minimised.
101 * \param *x parameters for the ellipsoid
102 * \param *params EllipsoidMinimisation with set of data points to minimise distance to and dimension
103 * \return sum of squared distance, \sa SquaredDistanceToEllipsoid()
104 */
105double SumSquaredDistance (const gsl_vector * x, void * params)
106{
107 Vector *set= ((struct EllipsoidMinimisation *)params)->x;
108 int N = ((struct EllipsoidMinimisation *)params)->N;
109 double SumDistance = 0.;
110 double distance;
111 Vector Center;
112 double EllipsoidLength[3], EllipsoidAngle[3];
113
114 // put parameters into suitable ellipsoid form
115 for (int i=0;i<3;i++) {
116 Center[i] = gsl_vector_get(x, i+0);
117 EllipsoidLength[i] = gsl_vector_get(x, i+3);
118 EllipsoidAngle[i] = gsl_vector_get(x, i+6);
119 }
120
121 // go through all points and sum distance
122 for (int i=0;i<N;i++) {
123 distance = SquaredDistanceToEllipsoid(set[i], Center, EllipsoidLength, EllipsoidAngle);
124 if (!isnan(distance)) {
125 SumDistance += distance;
126 } else {
127 SumDistance = GSL_NAN;
128 break;
129 }
130 }
131
132 //Log() << Verbose(0) << "Current summed distance is " << SumDistance << "." << endl;
133 return SumDistance;
134};
135
136/** Finds best fitting ellipsoid parameter set in Least square sense for a given point set.
137 * \param *out output stream for debugging
138 * \param *set given point set
139 * \param N number of points in set
140 * \param EllipsoidParamter[3] three parameters in ellipsoid equation
141 * \return true - fit successful, false - fit impossible
142 */
143bool FitPointSetToEllipsoid(Vector *set, int N, Vector *EllipsoidCenter, double *EllipsoidLength, double *EllipsoidAngle)
144{
145 int status = GSL_SUCCESS;
146 DoLog(2) && (Log() << Verbose(2) << "Begin of FitPointSetToEllipsoid " << endl);
147 if (N >= 3) { // check that enough points are given (9 d.o.f.)
148 struct EllipsoidMinimisation par;
149 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
150 gsl_multimin_fminimizer *s = NULL;
151 gsl_vector *ss, *x;
152 gsl_multimin_function minex_func;
153
154 size_t iter = 0;
155 double size;
156
157 /* Starting point */
158 x = gsl_vector_alloc (9);
159 for (int i=0;i<3;i++) {
160 gsl_vector_set (x, i+0, EllipsoidCenter->at(i));
161 gsl_vector_set (x, i+3, EllipsoidLength[i]);
162 gsl_vector_set (x, i+6, EllipsoidAngle[i]);
163 }
164 par.x = set;
165 par.N = N;
166
167 /* Set initial step sizes */
168 ss = gsl_vector_alloc (9);
169 for (int i=0;i<3;i++) {
170 gsl_vector_set (ss, i+0, 0.1);
171 gsl_vector_set (ss, i+3, 1.0);
172 gsl_vector_set (ss, i+6, M_PI/20.);
173 }
174
175 /* Initialize method and iterate */
176 minex_func.n = 9;
177 minex_func.f = &SumSquaredDistance;
178 minex_func.params = (void *)&par;
179
180 s = gsl_multimin_fminimizer_alloc (T, 9);
181 gsl_multimin_fminimizer_set (s, &minex_func, x, ss);
182
183 do {
184 iter++;
185 status = gsl_multimin_fminimizer_iterate(s);
186
187 if (status)
188 break;
189
190 size = gsl_multimin_fminimizer_size (s);
191 status = gsl_multimin_test_size (size, 1e-2);
192
193 if (status == GSL_SUCCESS) {
194 for (int i=0;i<3;i++) {
195 EllipsoidCenter->at(i) = gsl_vector_get (s->x,i+0);
196 EllipsoidLength[i] = gsl_vector_get (s->x, i+3);
197 EllipsoidAngle[i] = gsl_vector_get (s->x, i+6);
198 }
199 DoLog(4) && (Log() << Verbose(4) << setprecision(3) << "Converged fit at: " << *EllipsoidCenter << ", lengths " << EllipsoidLength[0] << ", " << EllipsoidLength[1] << ", " << EllipsoidLength[2] << ", angles " << EllipsoidAngle[0] << ", " << EllipsoidAngle[1] << ", " << EllipsoidAngle[2] << " with summed distance " << s->fval << "." << endl);
200 }
201
202 } while (status == GSL_CONTINUE && iter < 1000);
203
204 gsl_vector_free(x);
205 gsl_vector_free(ss);
206 gsl_multimin_fminimizer_free (s);
207
208 } else {
209 DoLog(3) && (Log() << Verbose(3) << "Not enough points provided for fit to ellipsoid." << endl);
210 return false;
211 }
212 DoLog(2) && (Log() << Verbose(2) << "End of FitPointSetToEllipsoid" << endl);
213 if (status == GSL_SUCCESS)
214 return true;
215 else
216 return false;
217};
218
219/** Picks a number of random points from a LC neighbourhood as a fitting set.
220 * \param *out output stream for debugging
221 * \param *T Tesselation containing boundary points
222 * \param *LC linked cell list of all atoms
223 * \param *&x random point set on return (not allocated!)
224 * \param PointsToPick number of points in set to pick
225 */
226void PickRandomNeighbouredPointSet(class Tesselation *T, class LinkedCell *LC, Vector *&x, size_t PointsToPick)
227{
228 size_t PointsLeft = 0;
229 size_t PointsPicked = 0;
230 int Nlower[NDIM], Nupper[NDIM];
231 set<int> PickedAtomNrs; // ordered list of picked atoms
232 set<int>::iterator current;
233 int index;
234 TesselPoint *Candidate = NULL;
235 DoLog(2) && (Log() << Verbose(2) << "Begin of PickRandomPointSet" << endl);
236
237 // allocate array
238 if (x == NULL) {
239 x = new Vector[PointsToPick];
240 } else {
241 DoeLog(2) && (eLog()<< Verbose(2) << "Given pointer to vector array seems already allocated." << endl);
242 }
243
244 do {
245 for(int i=0;i<NDIM;i++) // pick three random indices
246 LC->n[i] = (rand() % LC->N[i]);
247 DoLog(2) && (Log() << Verbose(2) << "INFO: Center cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " ... ");
248 // get random cell
249 const LinkedCell::LinkedNodes *List = LC->GetCurrentCell();
250 if (List == NULL) { // set index to it
251 continue;
252 }
253 DoLog(2) && (Log() << Verbose(2) << "with No. " << LC->index << "." << endl);
254
255 DoLog(2) && (Log() << Verbose(2) << "LC Intervals:");
256 for (int i=0;i<NDIM;i++) {
257 Nlower[i] = ((LC->n[i]-1) >= 0) ? LC->n[i]-1 : 0;
258 Nupper[i] = ((LC->n[i]+1) < LC->N[i]) ? LC->n[i]+1 : LC->N[i]-1;
259 DoLog(0) && (Log() << Verbose(0) << " [" << Nlower[i] << "," << Nupper[i] << "] ");
260 }
261 DoLog(0) && (Log() << Verbose(0) << endl);
262
263 // count whether there are sufficient atoms in this cell+neighbors
264 PointsLeft=0;
265 for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
266 for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
267 for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
268 const LinkedCell::LinkedNodes *List = LC->GetCurrentCell();
269 PointsLeft += List->size();
270 }
271 DoLog(2) && (Log() << Verbose(2) << "There are " << PointsLeft << " atoms in this neighbourhood." << endl);
272 if (PointsLeft < PointsToPick) { // ensure that we can pick enough points in its neighbourhood at all.
273 continue;
274 }
275
276 // pre-pick a fixed number of atoms
277 PickedAtomNrs.clear();
278 do {
279 index = (rand() % PointsLeft);
280 current = PickedAtomNrs.find(index); // not present?
281 if (current == PickedAtomNrs.end()) {
282 //Log() << Verbose(2) << "Picking atom nr. " << index << "." << endl;
283 PickedAtomNrs.insert(index);
284 }
285 } while (PickedAtomNrs.size() < PointsToPick);
286
287 index = 0; // now go through all and pick those whose from PickedAtomsNr
288 PointsPicked=0;
289 current = PickedAtomNrs.begin();
290 for (LC->n[0] = Nlower[0]; LC->n[0] <= Nupper[0]; LC->n[0]++)
291 for (LC->n[1] = Nlower[1]; LC->n[1] <= Nupper[1]; LC->n[1]++)
292 for (LC->n[2] = Nlower[2]; LC->n[2] <= Nupper[2]; LC->n[2]++) {
293 const LinkedCell::LinkedNodes *List = LC->GetCurrentCell();
294// Log() << Verbose(2) << "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << " containing " << List->size() << " points." << endl;
295 if (List != NULL) {
296// if (List->begin() != List->end())
297// Log() << Verbose(2) << "Going through candidates ... " << endl;
298// else
299// Log() << Verbose(2) << "Cell is empty ... " << endl;
300 for (LinkedCell::LinkedNodes::const_iterator Runner = List->begin(); Runner != List->end(); Runner++) {
301 if ((current != PickedAtomNrs.end()) && (*current == index)) {
302 Candidate = (*Runner);
303 DoLog(2) && (Log() << Verbose(2) << "Current picked node is " << **Runner << " with index " << index << "." << endl);
304 x[PointsPicked++] = *Candidate->node; // we have one more atom picked
305 current++; // next pre-picked atom
306 }
307 index++; // next atom nr.
308 }
309// } else {
310// Log() << Verbose(2) << "List for this index not allocated!" << endl;
311 }
312 }
313 DoLog(2) && (Log() << Verbose(2) << "The following points were picked: " << endl);
314 for (size_t i=0;i<PointsPicked;i++)
315 DoLog(2) && (Log() << Verbose(2) << x[i] << endl);
316 if (PointsPicked == PointsToPick) // break out of loop if we have all
317 break;
318 } while(1);
319
320 DoLog(2) && (Log() << Verbose(2) << "End of PickRandomPointSet" << endl);
321};
322
323/** Picks a number of random points from a set of boundary points as a fitting set.
324 * \param *out output stream for debugging
325 * \param *T Tesselation containing boundary points
326 * \param *&x random point set on return (not allocated!)
327 * \param PointsToPick number of points in set to pick
328 */
329void PickRandomPointSet(class Tesselation *T, Vector *&x, size_t PointsToPick)
330{
331 size_t PointsLeft = (size_t) T->PointsOnBoundaryCount;
332 size_t PointsPicked = 0;
333 double value, threshold;
334 PointMap *List = &T->PointsOnBoundary;
335 DoLog(2) && (Log() << Verbose(2) << "Begin of PickRandomPointSet" << endl);
336
337 // allocate array
338 if (x == NULL) {
339 x = new Vector[PointsToPick];
340 } else {
341 DoeLog(2) && (eLog()<< Verbose(2) << "Given pointer to vector array seems already allocated." << endl);
342 }
343
344 if (List != NULL)
345 for (PointMap::iterator Runner = List->begin(); Runner != List->end(); Runner++) {
346 threshold = 1. - (double)(PointsToPick - PointsPicked)/(double)PointsLeft;
347 value = (double)rand()/(double)RAND_MAX;
348 //Log() << Verbose(3) << "Current node is " << *Runner->second->node << " with " << value << " ... " << threshold << ": ";
349 if (value > threshold) {
350 x[PointsPicked] = (*Runner->second->node->node);
351 PointsPicked++;
352 //Log() << Verbose(0) << "IN." << endl;
353 } else {
354 //Log() << Verbose(0) << "OUT." << endl;
355 }
356 PointsLeft--;
357 }
358 DoLog(2) && (Log() << Verbose(2) << "The following points were picked: " << endl);
359 for (size_t i=0;i<PointsPicked;i++)
360 DoLog(3) && (Log() << Verbose(3) << x[i] << endl);
361
362 DoLog(2) && (Log() << Verbose(2) << "End of PickRandomPointSet" << endl);
363};
364
365/** Finds best fitting ellipsoid parameter set in least square sense for a given point set.
366 * \param *out output stream for debugging
367 * \param *T Tesselation containing boundary points
368 * \param *LCList linked cell list of all atoms
369 * \param N number of unique points in ellipsoid fit, must be greater equal 6
370 * \param number of fits (i.e. parameter sets in output file)
371 * \param *filename name for output file
372 */
373void FindDistributionOfEllipsoids(class Tesselation *T, class LinkedCell *LCList, int N, int number, const char *filename)
374{
375 ofstream output;
376 Vector *x = NULL;
377 Vector Center;
378 Vector EllipsoidCenter;
379 double EllipsoidLength[3];
380 double EllipsoidAngle[3];
381 double distance, MaxDistance, MinDistance;
382 DoLog(0) && (Log() << Verbose(0) << "Begin of FindDistributionOfEllipsoids" << endl);
383
384 // construct center of gravity of boundary point set for initial ellipsoid center
385 Center.Zero();
386 for (PointMap::iterator Runner = T->PointsOnBoundary.begin(); Runner != T->PointsOnBoundary.end(); Runner++)
387 Center += (*Runner->second->node->node);
388 Center.Scale(1./T->PointsOnBoundaryCount);
389 DoLog(1) && (Log() << Verbose(1) << "Center is at " << Center << "." << endl);
390
391 // Output header
392 output.open(filename, ios::trunc);
393 output << "# Nr.\tCenterX\tCenterY\tCenterZ\ta\tb\tc\tpsi\ttheta\tphi" << endl;
394
395 // loop over desired number of parameter sets
396 for (;number >0;number--) {
397 DoLog(1) && (Log() << Verbose(1) << "Determining data set " << number << " ... " << endl);
398 // pick the point set
399 x = NULL;
400 //PickRandomPointSet(T, LCList, x, N);
401 PickRandomNeighbouredPointSet(T, LCList, x, N);
402
403 // calculate some sensible starting values for parameter fit
404 MaxDistance = 0.;
405 MinDistance = x[0].ScalarProduct(x[0]);
406 for (int i=0;i<N;i++) {
407 distance = x[i].ScalarProduct(x[i]);
408 if (distance > MaxDistance)
409 MaxDistance = distance;
410 if (distance < MinDistance)
411 MinDistance = distance;
412 }
413 //Log() << Verbose(2) << "MinDistance " << MinDistance << ", MaxDistance " << MaxDistance << "." << endl;
414 EllipsoidCenter = Center; // use Center of Gravity as initial center of ellipsoid
415 for (int i=0;i<3;i++)
416 EllipsoidAngle[i] = 0.;
417 EllipsoidLength[0] = sqrt(MaxDistance);
418 EllipsoidLength[1] = sqrt((MaxDistance+MinDistance)/2.);
419 EllipsoidLength[2] = sqrt(MinDistance);
420
421 // fit the parameters
422 if (FitPointSetToEllipsoid(x, N, &EllipsoidCenter, &EllipsoidLength[0], &EllipsoidAngle[0])) {
423 DoLog(1) && (Log() << Verbose(1) << "Picking succeeded!" << endl);
424 // output obtained parameter set
425 output << number << "\t";
426 for (int i=0;i<3;i++)
427 output << setprecision(9) << EllipsoidCenter[i] << "\t";
428 for (int i=0;i<3;i++)
429 output << setprecision(9) << EllipsoidLength[i] << "\t";
430 for (int i=0;i<3;i++)
431 output << setprecision(9) << EllipsoidAngle[i] << "\t";
432 output << endl;
433 } else { // increase N to pick one more
434 DoLog(1) && (Log() << Verbose(1) << "Picking failed!" << endl);
435 number++;
436 }
437 delete[](x); // free allocated memory for point set
438 }
439 // close output and finish
440 output.close();
441
442 DoLog(0) && (Log() << Verbose(0) << "End of FindDistributionOfEllipsoids" << endl);
443};
Note: See TracBrowser for help on using the repository browser.