2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
46 #include "gmx_fatal.h"
47 #include "gmx_matrix.h"
48 #include "gmx_lapack.h"
50 double **alloc_matrix(int n,int m)
55 /* There's always time for more pointer arithmetic! */
56 /* This is necessary in order to be able to work with LAPACK */
64 void free_matrix(double **a,int n)
73 void matrix_multiply(FILE *fp,int n,int m,double **x,double **y,double **z)
79 fprintf(fp,"Multiplying %d x %d matrix with a %d x %d matrix\n",
85 fprintf(fp," %7g",x[i][j]);
95 z[i][j] += x[k][i]*y[j][k];
100 static void dump_matrix(FILE *fp,const char *title,int n,double **a)
105 fprintf(fp,"%s\n",title);
111 fprintf(fp," %8.2f",a[i][j]);
115 fprintf(fp,"Prod a[i][i] = %g\n",d);
118 int matrix_invert(FILE *fp,int n,double **a)
120 int i,j,m,lda,*ipiv,lwork,info;
121 double **test=NULL,**id,*work;
126 fprintf(fp,"Inverting %d square matrix\n",n);
127 test = alloc_matrix(n,n);
132 test[i][j] = a[i][j];
135 dump_matrix(fp,"before inversion",n,a);
143 F77_FUNC(dgetrf,DGETRF)(&n,&m,a[0],&lda,ipiv,&info);
146 dump_matrix(fp,"after dgetrf",n,a);
150 F77_FUNC(dgetri,DGETRI)(&n,a[0],&lda,ipiv,work,&lwork,&info);
153 dump_matrix(fp,"after dgetri",n,a);
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);
174 double multi_regression(FILE *fp,int nrow,double *y,int ncol,
175 double **xx,double *a0)
178 double ax,chi2,**a,**at,**ata,*atx;
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.",
193 for(i=0; (i<ncol); i++)
196 for(j=0; (j<nrow); j++)
197 atx[i] += at[i][j]*y[j];
199 for(i=0; (i<ncol); i++)
202 for(j=0; (j<ncol); j++)
203 a0[i] += ata[i][j]*atx[j];
206 for(j=0; (j<nrow); j++)
209 for(i=0; (i<ncol); i++)
211 chi2 += sqr(y[j]-ax);
216 free_matrix(at,ncol);
217 free_matrix(ata,ncol);