0d78627b552bd662fa5855f4e9df792b2e7bc7b3
[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, 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 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <stdio.h>
44
45 #include "gromacs/legacyheaders/gmx_fatal.h"
46 #include "gromacs/legacyheaders/smalloc.h"
47 #include "gromacs/legacyheaders/vec.h"
48
49 #include "gmx_lapack.h"
50
51 double **alloc_matrix(int n, int m)
52 {
53     double **ptr;
54     int      i;
55
56     /* There's always time for more pointer arithmetic! */
57     /* This is necessary in order to be able to work with LAPACK */
58     snew(ptr, n);
59     snew(ptr[0], n*m);
60     for (i = 1; (i < n); i++)
61     {
62         ptr[i] = ptr[i-1]+m;
63     }
64     return ptr;
65 }
66
67 void free_matrix(double **a)
68 {
69     int i;
70
71     sfree(a[0]);
72     sfree(a);
73 }
74
75 #define DEBUG_MATRIX
76 void matrix_multiply(FILE *fp, int n, int m, double **x, double **y, double **z)
77 {
78     int i, j, k;
79
80 #ifdef DEBUG_MATRIX
81     if (fp)
82     {
83         fprintf(fp, "Multiplying %d x %d matrix with a %d x %d matrix\n",
84                 n, m, m, n);
85     }
86     if (fp)
87     {
88         for (i = 0; (i < n); i++)
89         {
90             for (j = 0; (j < m); j++)
91             {
92                 fprintf(fp, " %7g", x[i][j]);
93             }
94             fprintf(fp, "\n");
95         }
96     }
97 #endif
98     for (i = 0; (i < m); i++)
99     {
100         for (j = 0; (j < m); j++)
101         {
102             z[i][j] = 0;
103             for (k = 0; (k < n); k++)
104             {
105                 z[i][j] += x[k][i]*y[j][k];
106             }
107         }
108     }
109 }
110
111 static void dump_matrix(FILE *fp, const char *title, int n, double **a)
112 {
113     double d = 1;
114     int    i, j;
115
116     fprintf(fp, "%s\n", title);
117     for (i = 0; (i < n); i++)
118     {
119         d = d*a[i][i];
120         for (j = 0; (j < n); j++)
121         {
122             fprintf(fp, " %8.2f", a[i][j]);
123         }
124         fprintf(fp, "\n");
125     }
126     fprintf(fp, "Prod a[i][i] = %g\n", d);
127 }
128
129 int matrix_invert(FILE *fp, int n, double **a)
130 {
131     int      i, j, m, lda, *ipiv, lwork, info;
132     double **test = NULL, **id, *work;
133
134 #ifdef DEBUG_MATRIX
135     if (fp)
136     {
137         fprintf(fp, "Inverting %d square matrix\n", n);
138         test = alloc_matrix(n, n);
139         for (i = 0; (i < n); i++)
140         {
141             for (j = 0; (j < n); j++)
142             {
143                 test[i][j] = a[i][j];
144             }
145         }
146         dump_matrix(fp, "before inversion", n, a);
147     }
148 #endif
149     snew(ipiv, n);
150     lwork = n*n;
151     snew(work, lwork);
152     m     = lda   = n;
153     info  = 0;
154     F77_FUNC(dgetrf, DGETRF) (&n, &m, a[0], &lda, ipiv, &info);
155 #ifdef DEBUG_MATRIX
156     if (fp)
157     {
158         dump_matrix(fp, "after dgetrf", n, a);
159     }
160 #endif
161     if (info != 0)
162     {
163         return info;
164     }
165     F77_FUNC(dgetri, DGETRI) (&n, a[0], &lda, ipiv, work, &lwork, &info);
166 #ifdef DEBUG_MATRIX
167     if (fp)
168     {
169         dump_matrix(fp, "after dgetri", n, a);
170     }
171 #endif
172     if (info != 0)
173     {
174         return info;
175     }
176
177 #ifdef DEBUG_MATRIX
178     if (fp)
179     {
180         id = alloc_matrix(n, n);
181         matrix_multiply(fp, n, n, test, a, id);
182         dump_matrix(fp, "And here is the product of A and Ainv", n, id);
183         free_matrix(id);
184         free_matrix(test);
185     }
186 #endif
187     sfree(ipiv);
188     sfree(work);
189
190     return 0;
191 }
192
193 double multi_regression(FILE *fp, int nrow, double *y, int ncol,
194                         double **xx, double *a0)
195 {
196     int    row, niter, i, j;
197     double ax, chi2, **a, **at, **ata, *atx;
198
199     a   = alloc_matrix(nrow, ncol);
200     at  = alloc_matrix(ncol, nrow);
201     ata = alloc_matrix(ncol, ncol);
202     for (i = 0; (i < nrow); i++)
203     {
204         for (j = 0; (j < ncol); j++)
205         {
206             at[j][i] = a[i][j] = xx[j][i];
207         }
208     }
209     matrix_multiply(fp, nrow, ncol, a, at, ata);
210     if ((row = matrix_invert(fp, ncol, ata)) != 0)
211     {
212         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.",
213                   row);
214     }
215     snew(atx, ncol);
216
217     for (i = 0; (i < ncol); i++)
218     {
219         atx[i] = 0;
220         for (j = 0; (j < nrow); j++)
221         {
222             atx[i] += at[i][j]*y[j];
223         }
224     }
225     for (i = 0; (i < ncol); i++)
226     {
227         a0[i] = 0;
228         for (j = 0; (j < ncol); j++)
229         {
230             a0[i] += ata[i][j]*atx[j];
231         }
232     }
233     chi2 = 0;
234     for (j = 0; (j < nrow); j++)
235     {
236         ax = 0;
237         for (i = 0; (i < ncol); i++)
238         {
239             ax += a0[i]*a[j][i];
240         }
241         chi2 += sqr(y[j]-ax);
242     }
243
244     sfree(atx);
245     free_matrix(a);
246     free_matrix(at);
247     free_matrix(ata);
248
249     return chi2;
250 }