| 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 | * vector_ops.cpp
|
|---|
| 10 | *
|
|---|
| 11 | * Created on: Apr 1, 2010
|
|---|
| 12 | * Author: crueger
|
|---|
| 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 "CodePatterns/Info.hpp"
|
|---|
| 23 | #include "CodePatterns/Log.hpp"
|
|---|
| 24 | #include "CodePatterns/Verbose.hpp"
|
|---|
| 25 | #include "LinearAlgebra/defs.hpp"
|
|---|
| 26 | #include "LinearAlgebra/fast_functions.hpp"
|
|---|
| 27 | #include "LinearAlgebra/leastsquaremin.hpp"
|
|---|
| 28 | #include "LinearAlgebra/Plane.hpp"
|
|---|
| 29 | #include "LinearAlgebra/Vector.hpp"
|
|---|
| 30 |
|
|---|
| 31 | #include <gsl/gsl_linalg.h>
|
|---|
| 32 | #include <gsl/gsl_matrix.h>
|
|---|
| 33 | #include <gsl/gsl_permutation.h>
|
|---|
| 34 | #include <gsl/gsl_vector.h>
|
|---|
| 35 | #include <gsl/gsl_multimin.h>
|
|---|
| 36 |
|
|---|
| 37 | /**
|
|---|
| 38 | * !@file
|
|---|
| 39 | * These files defines several common operation on vectors that should not
|
|---|
| 40 | * become part of the main vector class, because they are either to complex
|
|---|
| 41 | * or need methods from other subsystems that should not be moved to
|
|---|
| 42 | * the LinAlg-Subsystem
|
|---|
| 43 | */
|
|---|
| 44 |
|
|---|
| 45 | /** Creates a new vector as the one with least square distance to a given set of \a vectors.
|
|---|
| 46 | * \param *vectors set of vectors
|
|---|
| 47 | * \param num number of vectors
|
|---|
| 48 | * \return true if success, false if failed due to linear dependency
|
|---|
| 49 | */
|
|---|
| 50 | bool LSQdistance(Vector &res,const Vector **vectors, int num)
|
|---|
| 51 | {
|
|---|
| 52 | int j;
|
|---|
| 53 |
|
|---|
| 54 | for (j=0;j<num;j++) {
|
|---|
| 55 | Log() << Verbose(1) << j << "th atom's vector: " << vectors[j] << endl;
|
|---|
| 56 | }
|
|---|
| 57 |
|
|---|
| 58 | int np = 3;
|
|---|
| 59 | struct LSQ_params par;
|
|---|
| 60 |
|
|---|
| 61 | const gsl_multimin_fminimizer_type *T =
|
|---|
| 62 | gsl_multimin_fminimizer_nmsimplex;
|
|---|
| 63 | gsl_multimin_fminimizer *s = NULL;
|
|---|
| 64 | gsl_vector *ss, *y;
|
|---|
| 65 | gsl_multimin_function minex_func;
|
|---|
| 66 |
|
|---|
| 67 | size_t iter = 0, i;
|
|---|
| 68 | int status;
|
|---|
| 69 | double size;
|
|---|
| 70 |
|
|---|
| 71 | /* Initial vertex size vector */
|
|---|
| 72 | ss = gsl_vector_alloc (np);
|
|---|
| 73 | y = gsl_vector_alloc (np);
|
|---|
| 74 |
|
|---|
| 75 | /* Set all step sizes to 1 */
|
|---|
| 76 | gsl_vector_set_all (ss, 1.0);
|
|---|
| 77 |
|
|---|
| 78 | /* Starting point */
|
|---|
| 79 | par.vectors = vectors;
|
|---|
| 80 | par.num = num;
|
|---|
| 81 |
|
|---|
| 82 | for (i=NDIM;i--;)
|
|---|
| 83 | gsl_vector_set(y, i, (vectors[0]->at(i) - vectors[1]->at(i))/2.);
|
|---|
| 84 |
|
|---|
| 85 | /* Initialize method and iterate */
|
|---|
| 86 | minex_func.f = &LSQ;
|
|---|
| 87 | minex_func.n = np;
|
|---|
| 88 | minex_func.params = (void *)∥
|
|---|
| 89 |
|
|---|
| 90 | s = gsl_multimin_fminimizer_alloc (T, np);
|
|---|
| 91 | gsl_multimin_fminimizer_set (s, &minex_func, y, ss);
|
|---|
| 92 |
|
|---|
| 93 | do
|
|---|
| 94 | {
|
|---|
| 95 | iter++;
|
|---|
| 96 | status = gsl_multimin_fminimizer_iterate(s);
|
|---|
| 97 |
|
|---|
| 98 | if (status)
|
|---|
| 99 | break;
|
|---|
| 100 |
|
|---|
| 101 | size = gsl_multimin_fminimizer_size (s);
|
|---|
| 102 | status = gsl_multimin_test_size (size, 1e-2);
|
|---|
| 103 |
|
|---|
| 104 | if (status == GSL_SUCCESS)
|
|---|
| 105 | {
|
|---|
| 106 | printf ("converged to minimum at\n");
|
|---|
| 107 | }
|
|---|
| 108 |
|
|---|
| 109 | printf ("%5d ", (int)iter);
|
|---|
| 110 | for (i = 0; i < (size_t)np; i++)
|
|---|
| 111 | {
|
|---|
| 112 | printf ("%10.3e ", gsl_vector_get (s->x, i));
|
|---|
| 113 | }
|
|---|
| 114 | printf ("f() = %7.3f size = %.3f\n", s->fval, size);
|
|---|
| 115 | }
|
|---|
| 116 | while (status == GSL_CONTINUE && iter < 100);
|
|---|
| 117 |
|
|---|
| 118 | for (i=(size_t)np;i--;)
|
|---|
| 119 | res[i] = gsl_vector_get(s->x, i);
|
|---|
| 120 | gsl_vector_free(y);
|
|---|
| 121 | gsl_vector_free(ss);
|
|---|
| 122 | gsl_multimin_fminimizer_free (s);
|
|---|
| 123 |
|
|---|
| 124 | return true;
|
|---|
| 125 | };
|
|---|