source: src/Tesselation/ellipsoid.cpp@ 41a467

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 41a467 was 47d041, checked in by Frederik Heber <heber@…>, 13 years ago

HUGE: Removed all calls to Log(), eLog(), replaced by LOG() and ELOG().

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