82aa7bcc83f7cdf149daa3e3ddfa73509a26b450
[alexxy/gromacs.git] / src / gromacs / linearalgebra / matrix.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2008, The GROMACS development team.
6  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "matrix.h"
38
39 #include "config.h"
40
41 #include <stdio.h>
42
43 #include "gromacs/utility/fatalerror.h"
44 #include "gromacs/utility/smalloc.h"
45
46 #include "gmx_lapack.h"
47
48 double **alloc_matrix(int n, int m)
49 {
50     double **ptr;
51     int      i;
52
53     /* There's always time for more pointer arithmetic! */
54     /* This is necessary in order to be able to work with LAPACK */
55     snew(ptr, n);
56     snew(ptr[0], n*m);
57     for (i = 1; (i < n); i++)
58     {
59         ptr[i] = ptr[i-1]+m;
60     }
61     return ptr;
62 }
63
64 void free_matrix(double **a)
65 {
66     int i;
67
68     sfree(a[0]);
69     sfree(a);
70 }
71
72 #define DEBUG_MATRIX
73 void matrix_multiply(FILE *fp, int n, int m, double **x, double **y, double **z)
74 {
75     int i, j, k;
76
77 #ifdef DEBUG_MATRIX
78     if (fp)
79     {
80         fprintf(fp, "Multiplying %d x %d matrix with a %d x %d matrix\n",
81                 n, m, m, n);
82     }
83     if (fp)
84     {
85         for (i = 0; (i < n); i++)
86         {
87             for (j = 0; (j < m); j++)
88             {
89                 fprintf(fp, " %7g", x[i][j]);
90             }
91             fprintf(fp, "\n");
92         }
93     }
94 #endif
95     for (i = 0; (i < m); i++)
96     {
97         for (j = 0; (j < m); j++)
98         {
99             z[i][j] = 0;
100             for (k = 0; (k < n); k++)
101             {
102                 z[i][j] += x[k][i]*y[j][k];
103             }
104         }
105     }
106 }
107
108 static void dump_matrix(FILE *fp, const char *title, int n, double **a)
109 {
110     double d = 1;
111     int    i, j;
112
113     fprintf(fp, "%s\n", title);
114     for (i = 0; (i < n); i++)
115     {
116         d = d*a[i][i];
117         for (j = 0; (j < n); j++)
118         {
119             fprintf(fp, " %8.2f", a[i][j]);
120         }
121         fprintf(fp, "\n");
122     }
123     fprintf(fp, "Prod a[i][i] = %g\n", d);
124 }
125
126 int matrix_invert(FILE *fp, int n, double **a)
127 {
128     int      i, j, m, lda, *ipiv, lwork, info;
129     double **test = NULL, **id, *work;
130
131 #ifdef DEBUG_MATRIX
132     if (fp)
133     {
134         fprintf(fp, "Inverting %d square matrix\n", n);
135         test = alloc_matrix(n, n);
136         for (i = 0; (i < n); i++)
137         {
138             for (j = 0; (j < n); j++)
139             {
140                 test[i][j] = a[i][j];
141             }
142         }
143         dump_matrix(fp, "before inversion", n, a);
144     }
145 #endif
146     snew(ipiv, n);
147     lwork = n*n;
148     snew(work, lwork);
149     m     = lda   = n;
150     info  = 0;
151     F77_FUNC(dgetrf, DGETRF) (&n, &m, a[0], &lda, ipiv, &info);
152 #ifdef DEBUG_MATRIX
153     if (fp)
154     {
155         dump_matrix(fp, "after dgetrf", n, a);
156     }
157 #endif
158     if (info != 0)
159     {
160         return info;
161     }
162     F77_FUNC(dgetri, DGETRI) (&n, a[0], &lda, ipiv, work, &lwork, &info);
163 #ifdef DEBUG_MATRIX
164     if (fp)
165     {
166         dump_matrix(fp, "after dgetri", n, a);
167     }
168 #endif
169     if (info != 0)
170     {
171         return info;
172     }
173
174 #ifdef DEBUG_MATRIX
175     if (fp)
176     {
177         id = alloc_matrix(n, n);
178         matrix_multiply(fp, n, n, test, a, id);
179         dump_matrix(fp, "And here is the product of A and Ainv", n, id);
180         free_matrix(id);
181         free_matrix(test);
182     }
183 #endif
184     sfree(ipiv);
185     sfree(work);
186
187     return 0;
188 }
189
190 double multi_regression(FILE *fp, int nrow, double *y, int ncol,
191                         double **xx, double *a0)
192 {
193     int    row, niter, i, j;
194     double ax, chi2, **a, **at, **ata, *atx;
195
196     a   = alloc_matrix(nrow, ncol);
197     at  = alloc_matrix(ncol, nrow);
198     ata = alloc_matrix(ncol, ncol);
199     for (i = 0; (i < nrow); i++)
200     {
201         for (j = 0; (j < ncol); j++)
202         {
203             at[j][i] = a[i][j] = xx[j][i];
204         }
205     }
206     matrix_multiply(fp, nrow, ncol, a, at, ata);
207     if ((row = matrix_invert(fp, ncol, ata)) != 0)
208     {
209         gmx_fatal(FARGS, "Matrix inversion failed. Incorrect row = %d.\nThis probably indicates that you do not have sufficient data points, or that some parameters are linearly dependent.",
210                   row);
211     }
212     snew(atx, ncol);
213
214     for (i = 0; (i < ncol); i++)
215     {
216         atx[i] = 0;
217         for (j = 0; (j < nrow); j++)
218         {
219             atx[i] += at[i][j]*y[j];
220         }
221     }
222     for (i = 0; (i < ncol); i++)
223     {
224         a0[i] = 0;
225         for (j = 0; (j < ncol); j++)
226         {
227             a0[i] += ata[i][j]*atx[j];
228         }
229     }
230     chi2 = 0;
231     for (j = 0; (j < nrow); j++)
232     {
233         ax = 0;
234         for (i = 0; (i < ncol); i++)
235         {
236             ax += a0[i]*a[j][i];
237         }
238         chi2 += (y[j] - ax) * (y[j] - ax);
239     }
240
241     sfree(atx);
242     free_matrix(a);
243     free_matrix(at);
244     free_matrix(ata);
245
246     return chi2;
247 }