| 1 | /* | 
|---|
| 2 | * Project: MoleCuilder | 
|---|
| 3 | * Description: creates and alters molecular systems | 
|---|
| 4 | * Copyright (C)  2010-2012 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 "CodePatterns/Log.hpp" | 
|---|
| 30 | #include "ellipsoid.hpp" | 
|---|
| 31 | #include "LinearAlgebra/Vector.hpp" | 
|---|
| 32 | #include "LinearAlgebra/RealSpaceMatrix.hpp" | 
|---|
| 33 | #include "LinkedCell/linkedcell.hpp" | 
|---|
| 34 | #include "Tesselation/BoundaryPointSet.hpp" | 
|---|
| 35 | #include "Tesselation/boundary.hpp" | 
|---|
| 36 | #include "Tesselation/tesselation.hpp" | 
|---|
| 37 |  | 
|---|
| 38 | #include "RandomNumbers/RandomNumberGeneratorFactory.hpp" | 
|---|
| 39 | #include "RandomNumbers/RandomNumberGenerator.hpp" | 
|---|
| 40 |  | 
|---|
| 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 | */ | 
|---|
| 48 | double SquaredDistanceToEllipsoid(Vector &x, Vector &EllipsoidCenter, double *EllipsoidLength, double *EllipsoidAngle) | 
|---|
| 49 | { | 
|---|
| 50 | Vector helper, RefPoint; | 
|---|
| 51 | double distance = -1.; | 
|---|
| 52 | RealSpaceMatrix Matrix; | 
|---|
| 53 | double InverseLength[3]; | 
|---|
| 54 | double psi,theta,phi; // euler angles in ZX'Z'' convention | 
|---|
| 55 |  | 
|---|
| 56 | //LOG(3, "Begin of SquaredDistanceToEllipsoid"); | 
|---|
| 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 | 
|---|
| 62 | RefPoint = helper = x - EllipsoidCenter; | 
|---|
| 63 | //LOG(4, "Translated given point is at " << RefPoint << "."); | 
|---|
| 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]; | 
|---|
| 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)); | 
|---|
| 78 | helper *= Matrix; | 
|---|
| 79 | helper.ScaleAll(InverseLength); | 
|---|
| 80 | //LOG(4, "Transformed RefPoint is at " << helper << "."); | 
|---|
| 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 | 
|---|
| 84 | //LOG(4, "Transformed intersection is at " << helper << "."); | 
|---|
| 85 |  | 
|---|
| 86 | // 4. transform back the constructed intersection point | 
|---|
| 87 | psi = -EllipsoidAngle[0]; | 
|---|
| 88 | theta = -EllipsoidAngle[1]; | 
|---|
| 89 | phi = -EllipsoidAngle[2]; | 
|---|
| 90 | helper.ScaleAll(EllipsoidLength); | 
|---|
| 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)); | 
|---|
| 100 | helper *= Matrix; | 
|---|
| 101 | //LOG(4, "Intersection is at " << helper << "."); | 
|---|
| 102 |  | 
|---|
| 103 | // 5. determine distance between backtransformed point and x | 
|---|
| 104 | distance = RefPoint.DistanceSquared(helper); | 
|---|
| 105 | //LOG(4, "Squared distance between intersection and RefPoint is " << distance << "."); | 
|---|
| 106 |  | 
|---|
| 107 | return distance; | 
|---|
| 108 | //LOG(3, "End of SquaredDistanceToEllipsoid"); | 
|---|
| 109 | }; | 
|---|
| 110 |  | 
|---|
| 111 | /** structure for ellipsoid minimisation containing points to fit to. | 
|---|
| 112 | */ | 
|---|
| 113 | struct EllipsoidMinimisation { | 
|---|
| 114 | int N;      //!< dimension of vector set | 
|---|
| 115 | Vector *x;  //!< array of vectors | 
|---|
| 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 | */ | 
|---|
| 123 | double SumSquaredDistance (const gsl_vector * x, void * params) | 
|---|
| 124 | { | 
|---|
| 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++) { | 
|---|
| 134 | Center[i] = gsl_vector_get(x, i+0); | 
|---|
| 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 |  | 
|---|
| 150 | //LOG(0, "Current summed distance is " << SumDistance << "."); | 
|---|
| 151 | return SumDistance; | 
|---|
| 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 | */ | 
|---|
| 161 | bool FitPointSetToEllipsoid(Vector *set, int N, Vector *EllipsoidCenter, double *EllipsoidLength, double *EllipsoidAngle) | 
|---|
| 162 | { | 
|---|
| 163 | int status = GSL_SUCCESS; | 
|---|
| 164 | LOG(2, "Begin of FitPointSetToEllipsoid "); | 
|---|
| 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++) { | 
|---|
| 178 | gsl_vector_set (x, i+0, EllipsoidCenter->at(i)); | 
|---|
| 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 *)∥ | 
|---|
| 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++) { | 
|---|
| 213 | EllipsoidCenter->at(i) = gsl_vector_get (s->x,i+0); | 
|---|
| 214 | EllipsoidLength[i] = gsl_vector_get (s->x, i+3); | 
|---|
| 215 | EllipsoidAngle[i] = gsl_vector_get (s->x, i+6); | 
|---|
| 216 | } | 
|---|
| 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 << "."); | 
|---|
| 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 { | 
|---|
| 227 | LOG(3, "Not enough points provided for fit to ellipsoid."); | 
|---|
| 228 | return false; | 
|---|
| 229 | } | 
|---|
| 230 | LOG(2, "End of FitPointSetToEllipsoid"); | 
|---|
| 231 | if (status == GSL_SUCCESS) | 
|---|
| 232 | return true; | 
|---|
| 233 | else | 
|---|
| 234 | return false; | 
|---|
| 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 | */ | 
|---|
| 244 | void PickRandomNeighbouredPointSet(class Tesselation *T, class LinkedCell_deprecated *LC, Vector *&x, size_t PointsToPick) | 
|---|
| 245 | { | 
|---|
| 246 | size_t PointsLeft = 0; | 
|---|
| 247 | size_t PointsPicked = 0; | 
|---|
| 248 | int Nlower[NDIM], Nupper[NDIM]; | 
|---|
| 249 | set<int> PickedAtomNrs;   // ordered list of picked atoms | 
|---|
| 250 | set<int>::iterator current; | 
|---|
| 251 | int index; | 
|---|
| 252 | TesselPoint *Candidate = NULL; | 
|---|
| 253 | LOG(2, "Begin of PickRandomPointSet"); | 
|---|
| 254 |  | 
|---|
| 255 | // allocate array | 
|---|
| 256 | if (x == NULL) { | 
|---|
| 257 | x = new Vector[PointsToPick]; | 
|---|
| 258 | } else { | 
|---|
| 259 | ELOG(2, "Given pointer to vector array seems already allocated."); | 
|---|
| 260 | } | 
|---|
| 261 |  | 
|---|
| 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 |  | 
|---|
| 280 | do { | 
|---|
| 281 | for(int i=0;i<NDIM;i++) // pick three random indices | 
|---|
| 282 | LC->n[i] = ((int)random() % LC->N[i]); | 
|---|
| 283 | LOG(2, "INFO: Center cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << "."); | 
|---|
| 284 | // get random cell | 
|---|
| 285 | const TesselPointSTLList *List = LC->GetCurrentCell(); | 
|---|
| 286 | if (List == NULL) {  // set index to it | 
|---|
| 287 | continue; | 
|---|
| 288 | } | 
|---|
| 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 | } | 
|---|
| 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]++) { | 
|---|
| 309 | const TesselPointSTLList *List = LC->GetCurrentCell(); | 
|---|
| 310 | PointsLeft += List->size(); | 
|---|
| 311 | } | 
|---|
| 312 | LOG(2, "There are " << PointsLeft << " atoms in this neighbourhood."); | 
|---|
| 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 { | 
|---|
| 320 | index = (((int)random()) % PointsLeft); | 
|---|
| 321 | current = PickedAtomNrs.find(index);  // not present? | 
|---|
| 322 | if (current == PickedAtomNrs.end()) { | 
|---|
| 323 | //LOG(2, "Picking atom Nr. " << index << "."); | 
|---|
| 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]++) { | 
|---|
| 334 | const TesselPointSTLList *List = LC->GetCurrentCell(); | 
|---|
| 335 | //          LOG(2, "Current cell is " << LC->n[0] << ", " << LC->n[1] << ", " << LC->n[2] << " with No. " << LC->index << " containing " << List->size() << " points."); | 
|---|
| 336 | if (List != NULL) { | 
|---|
| 337 | //            if (List->begin() != List->end()) | 
|---|
| 338 | //              LOG(2, "Going through candidates ... "); | 
|---|
| 339 | //            else | 
|---|
| 340 | //              LOG(2, "Cell is empty ... "); | 
|---|
| 341 | for (TesselPointSTLList::const_iterator Runner = List->begin(); Runner != List->end(); Runner++) { | 
|---|
| 342 | if ((current != PickedAtomNrs.end()) && (*current == index)) { | 
|---|
| 343 | Candidate = (*Runner); | 
|---|
| 344 | LOG(2, "Current picked node is " << (*Runner)->getName() << " with index " << index << "."); | 
|---|
| 345 | x[PointsPicked++] = Candidate->getPosition();    // we have one more atom picked | 
|---|
| 346 | current++;    // next pre-picked atom | 
|---|
| 347 | } | 
|---|
| 348 | index++;  // next atom Nr. | 
|---|
| 349 | } | 
|---|
| 350 | //          } else { | 
|---|
| 351 | //            LOG(2, "List for this index not allocated!"); | 
|---|
| 352 | } | 
|---|
| 353 | } | 
|---|
| 354 | LOG(2, "The following points were picked: "); | 
|---|
| 355 | for (size_t i=0;i<PointsPicked;i++) | 
|---|
| 356 | LOG(2, x[i]); | 
|---|
| 357 | if (PointsPicked == PointsToPick)  // break out of loop if we have all | 
|---|
| 358 | break; | 
|---|
| 359 | } while(1); | 
|---|
| 360 |  | 
|---|
| 361 | LOG(2, "End of PickRandomPointSet"); | 
|---|
| 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 | */ | 
|---|
| 370 | void PickRandomPointSet(class Tesselation *T, Vector *&x, size_t PointsToPick) | 
|---|
| 371 | { | 
|---|
| 372 | size_t PointsLeft = (size_t) T->PointsOnBoundaryCount; | 
|---|
| 373 | size_t PointsPicked = 0; | 
|---|
| 374 | double value, threshold; | 
|---|
| 375 | PointMap *List = &T->PointsOnBoundary; | 
|---|
| 376 | LOG(2, "Begin of PickRandomPointSet"); | 
|---|
| 377 |  | 
|---|
| 378 | // allocate array | 
|---|
| 379 | if (x == NULL) { | 
|---|
| 380 | x = new Vector[PointsToPick]; | 
|---|
| 381 | } else { | 
|---|
| 382 | ELOG(2, "Given pointer to vector array seems already allocated."); | 
|---|
| 383 | } | 
|---|
| 384 |  | 
|---|
| 385 | RandomNumberGenerator &random = RandomNumberGeneratorFactory::getInstance().makeRandomNumberGenerator("mt19937", "uniform_int"); | 
|---|
| 386 | const double rng_min = random.min(); | 
|---|
| 387 | const double rng_max = random.max(); | 
|---|
| 388 | if (List != NULL) | 
|---|
| 389 | for (PointMap::iterator Runner = List->begin(); Runner != List->end(); Runner++) { | 
|---|
| 390 | threshold = 1. - (double)(PointsToPick - PointsPicked)/(double)PointsLeft; | 
|---|
| 391 | value = (double)random()/(double)(rng_max-rng_min); | 
|---|
| 392 | if (value > threshold) { | 
|---|
| 393 | x[PointsPicked] = (Runner->second->node->getPosition()); | 
|---|
| 394 | PointsPicked++; | 
|---|
| 395 | //LOG(3, "Current node is " << *Runner->second->node << " with " << value << " ... " << threshold << ": IN."); | 
|---|
| 396 | } else { | 
|---|
| 397 | //LOG(3, "Current node is " << *Runner->second->node << " with " << value << " ... " << threshold << ": OUT."); | 
|---|
| 398 | } | 
|---|
| 399 | PointsLeft--; | 
|---|
| 400 | } | 
|---|
| 401 | LOG(2, "The following points were picked: "); | 
|---|
| 402 | for (size_t i=0;i<PointsPicked;i++) | 
|---|
| 403 | LOG(3, x[i]); | 
|---|
| 404 |  | 
|---|
| 405 | LOG(2, "End of PickRandomPointSet"); | 
|---|
| 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 | */ | 
|---|
| 416 | void FindDistributionOfEllipsoids(class Tesselation *T, class LinkedCell_deprecated *LCList, int N, int number, const char *filename) | 
|---|
| 417 | { | 
|---|
| 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; | 
|---|
| 425 | LOG(0, "Begin of FindDistributionOfEllipsoids"); | 
|---|
| 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++) | 
|---|
| 430 | Center += (Runner->second->node->getPosition()); | 
|---|
| 431 | Center.Scale(1./T->PointsOnBoundaryCount); | 
|---|
| 432 | LOG(4, "DEBUG: Center of PointsOnBoundary is at " << Center << "."); | 
|---|
| 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--) { | 
|---|
| 440 | LOG(1, "Determining data set " << number << " ... "); | 
|---|
| 441 | // pick the point set | 
|---|
| 442 | x = NULL; | 
|---|
| 443 | //PickRandomPointSet(T, LCList, x, N); | 
|---|
| 444 | PickRandomNeighbouredPointSet(T, LCList, x, N); | 
|---|
| 445 |  | 
|---|
| 446 | // calculate some sensible starting values for parameter fit | 
|---|
| 447 | MaxDistance = 0.; | 
|---|
| 448 | MinDistance = x[0].ScalarProduct(x[0]); | 
|---|
| 449 | for (int i=0;i<N;i++) { | 
|---|
| 450 | distance = x[i].ScalarProduct(x[i]); | 
|---|
| 451 | if (distance > MaxDistance) | 
|---|
| 452 | MaxDistance = distance; | 
|---|
| 453 | if (distance < MinDistance) | 
|---|
| 454 | MinDistance = distance; | 
|---|
| 455 | } | 
|---|
| 456 | //LOG(2, "MinDistance " << MinDistance << ", MaxDistance " << MaxDistance << "."); | 
|---|
| 457 | EllipsoidCenter = Center;  // use Center of Gravity as initial center of ellipsoid | 
|---|
| 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 | 
|---|
| 465 | if (FitPointSetToEllipsoid(x, N, &EllipsoidCenter, &EllipsoidLength[0], &EllipsoidAngle[0])) { | 
|---|
| 466 | LOG(1, "Picking succeeded!"); | 
|---|
| 467 | // output obtained parameter set | 
|---|
| 468 | output << number << "\t"; | 
|---|
| 469 | for (int i=0;i<3;i++) | 
|---|
| 470 | output << setprecision(9) << EllipsoidCenter[i] << "\t"; | 
|---|
| 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 | 
|---|
| 477 | LOG(1, "Picking failed!"); | 
|---|
| 478 | number++; | 
|---|
| 479 | } | 
|---|
| 480 | delete[](x);  // free allocated memory for point set | 
|---|
| 481 | } | 
|---|
| 482 | // close output and finish | 
|---|
| 483 | output.close(); | 
|---|
| 484 |  | 
|---|
| 485 | LOG(0, "End of FindDistributionOfEllipsoids"); | 
|---|
| 486 | }; | 
|---|