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