source: src/Tesselation/ellipsoid.cpp@ 2a0271

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 2a0271 was 0aa122, checked in by Frederik Heber <heber@…>, 13 years ago

Updated all source files's copyright note to current year 2012.

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