Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / gmxlib / gmx_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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <stdlib.h>
43 #include "typedefs.h"
44 #include "smalloc.h"
45 #include "vec.h"
46 #include "gmx_fatal.h"
47 #include "gmx_matrix.h"
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         ptr[i] = ptr[i-1]+m;
61     return ptr;
62 }
63
64 void free_matrix(double **a,int n)
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         fprintf(fp,"Multiplying %d x %d matrix with a %d x %d matrix\n",
80                 n,m,m,n);
81     if (fp)
82         for(i=0; (i<n); i++) 
83         {
84             for(j=0; (j<m); j++) 
85                 fprintf(fp," %7g",x[i][j]);
86             fprintf(fp,"\n");
87         }
88 #endif
89     for(i=0; (i<m); i++) 
90     {
91         for(j=0; (j<m); j++) 
92         {
93             z[i][j] = 0;
94             for(k=0; (k<n); k++)
95                 z[i][j] += x[k][i]*y[j][k];
96         }
97     }
98 }
99
100 static void dump_matrix(FILE *fp,const char *title,int n,double **a)
101 {
102     double d=1;
103     int i,j;
104
105     fprintf(fp,"%s\n",title);
106     for(i=0; (i<n); i++)  
107     {
108         d = d*a[i][i];
109         for(j=0; (j<n); j++)  
110         {
111             fprintf(fp," %8.2f",a[i][j]);
112         }
113         fprintf(fp,"\n");
114     }
115     fprintf(fp,"Prod a[i][i] = %g\n",d);
116 }
117
118 int matrix_invert(FILE *fp,int n,double **a)
119 {
120     int i,j,m,lda,*ipiv,lwork,info;
121     double **test=NULL,**id,*work;
122   
123 #ifdef DEBUG_MATRIX
124     if (fp) 
125     {
126         fprintf(fp,"Inverting %d square matrix\n",n);
127         test = alloc_matrix(n,n);  
128         for(i=0; (i<n); i++)  
129         {
130             for(j=0; (j<n); j++)  
131             {
132                 test[i][j] = a[i][j];
133             }
134         }
135         dump_matrix(fp,"before inversion",n,a);
136     }
137 #endif
138     snew(ipiv,n);
139     lwork = n*n;
140     snew(work,lwork);
141     m = lda   = n;
142     info  = 0;
143     F77_FUNC(dgetrf,DGETRF)(&n,&m,a[0],&lda,ipiv,&info);
144 #ifdef DEBUG_MATRIX
145     if (fp)
146         dump_matrix(fp,"after dgetrf",n,a);
147 #endif
148     if (info != 0)
149         return info;
150     F77_FUNC(dgetri,DGETRI)(&n,a[0],&lda,ipiv,work,&lwork,&info);
151 #ifdef DEBUG_MATRIX
152     if (fp)
153         dump_matrix(fp,"after dgetri",n,a);
154 #endif
155     if (info != 0)
156         return info;
157     
158 #ifdef DEBUG_MATRIX
159     if (fp) 
160     {
161         id = alloc_matrix(n,n);
162         matrix_multiply(fp,n,n,test,a,id);
163         dump_matrix(fp,"And here is the product of A and Ainv",n,id);
164         free_matrix(id,n);
165         free_matrix(test,n);
166     }
167 #endif
168     sfree(ipiv);
169     sfree(work);
170   
171     return 0;
172 }
173
174 double multi_regression(FILE *fp,int nrow,double *y,int ncol,
175                         double **xx,double *a0)
176 {
177   int    row,niter,i,j;
178   double ax,chi2,**a,**at,**ata,*atx;
179   
180   a   = alloc_matrix(nrow,ncol);
181   at  = alloc_matrix(ncol,nrow);
182   ata = alloc_matrix(ncol,ncol);
183   for(i=0; (i<nrow); i++)
184       for(j=0; (j<ncol); j++)
185           at[j][i] = a[i][j] = xx[j][i];
186   matrix_multiply(fp,nrow,ncol,a,at,ata);
187   if ((row = matrix_invert(fp,ncol,ata)) != 0) {
188       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.",
189                 row);
190   }
191   snew(atx,ncol);
192   
193   for(i=0; (i<ncol); i++)  
194   {
195       atx[i] = 0;
196       for(j=0; (j<nrow); j++)
197           atx[i] += at[i][j]*y[j];
198   }
199   for(i=0; (i<ncol); i++) 
200   {
201       a0[i] = 0;
202       for(j=0; (j<ncol); j++)
203           a0[i] += ata[i][j]*atx[j];
204   }
205   chi2 = 0;
206   for(j=0; (j<nrow); j++)
207   {
208       ax = 0;
209       for(i=0; (i<ncol); i++)  
210           ax += a0[i]*a[j][i];
211       chi2 += sqr(y[j]-ax);
212   }
213   
214   sfree(atx);
215   free_matrix(a,nrow);
216   free_matrix(at,ncol);
217   free_matrix(ata,ncol);
218   
219   return chi2;
220 }
221