Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / tools / eigio.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-2004, 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 "smalloc.h"
43 #include "vec.h"
44 #include "eigio.h"
45 #include "trnio.h"
46 #include "tpxio.h"
47 #include "statutil.h"
48 #include "futil.h"
49
50 void read_eigenvectors(const char *file,int *natoms,gmx_bool *bFit,
51                        rvec **xref,gmx_bool *bDMR,
52                        rvec **xav,gmx_bool *bDMA,
53                        int *nvec, int **eignr, 
54                        rvec ***eigvec,real **eigval)
55 {
56   t_trnheader head;
57   int    i,snew_size;
58   t_fileio *status;
59   rvec   *x;
60   matrix box;
61   gmx_bool   bOK;
62
63   *bDMR=FALSE;
64
65   /* read (reference (t=-1) and) average (t=0) structure */
66   status=open_trn(file,"r");
67   fread_trnheader(status,&head,&bOK);
68   *natoms=head.natoms;
69   snew(*xav,*natoms);
70   fread_htrn(status,&head,box,*xav,NULL,NULL);
71
72   if ((head.t>=-1.1) && (head.t<=-0.9)) 
73   {
74       snew(*xref,*natoms);
75       for(i=0; i<*natoms; i++)
76           copy_rvec((*xav)[i],(*xref)[i]);
77       *bDMR = (head.lambda > 0.5);
78       *bFit = (head.lambda > -0.5);
79       if (*bFit)
80       {
81           fprintf(stderr,"Read %smass weighted reference structure with %d atoms from %s\n", *bDMR ? "" : "non ",*natoms,file);
82       }
83       else 
84       {
85           fprintf(stderr,"Eigenvectors in %s were determined without fitting\n",file);
86           sfree(*xref);
87           *xref=NULL;
88       }
89       fread_trnheader(status,&head,&bOK);
90       fread_htrn(status,&head,box,*xav,NULL,NULL);
91   }
92   else
93   {
94       *bFit=TRUE;
95       *xref=NULL;
96   }
97   *bDMA = (head.lambda > 0.5);
98   if ((head.t<=-0.01) || (head.t>=0.01))
99   {
100       fprintf(stderr,"WARNING: %s does not start with t=0, which should be the "
101               "average structure. This might not be a eigenvector file. "
102               "Some things might go wrong.\n",
103               file);
104   }
105   else
106   {
107       fprintf(stderr,
108               "Read %smass weighted average/minimum structure with %d atoms from %s\n",
109               *bDMA ? "" : "non ",*natoms,file);
110   }
111   
112   snew(x,*natoms);
113   snew_size=10;
114   snew(*eignr,snew_size);
115   snew(*eigval,snew_size);
116   snew(*eigvec,snew_size);
117         
118   *nvec=0;
119   while (fread_trnheader(status,&head,&bOK)) 
120   {
121       fread_htrn(status,&head,box,x,NULL,NULL);
122       if (*nvec >= snew_size) 
123       {
124           snew_size+=10;
125           srenew(*eignr,snew_size);
126           srenew(*eigval,snew_size);
127           srenew(*eigvec,snew_size);
128       }
129       i=head.step;
130       (*eigval)[*nvec]=head.t;
131       (*eignr)[*nvec]=i-1;
132       snew((*eigvec)[*nvec],*natoms);
133       for(i=0; i<*natoms; i++)
134       {
135           copy_rvec(x[i],(*eigvec)[*nvec][i]);
136       }
137       (*nvec)++;
138   }
139   sfree(x);
140   fprintf(stderr,"Read %d eigenvectors (for %d atoms)\n\n",*nvec,*natoms);
141 }
142
143
144 void write_eigenvectors(const char *trnname,int natoms,real mat[],
145                         gmx_bool bReverse,int begin,int end,
146                         int WriteXref,rvec *xref,gmx_bool bDMR,
147                         rvec xav[], gmx_bool bDMA,real eigval[])
148 {
149     t_fileio *trnout;
150     int    ndim,i,j,d,vec;
151     matrix zerobox;
152     rvec   *x;
153     
154     ndim = natoms*DIM;
155     clear_mat(zerobox);
156     snew(x,natoms);
157     
158     fprintf (stderr,
159              "\nWriting %saverage structure & eigenvectors %d--%d to %s\n",
160              (WriteXref==eWXR_YES) ? "reference, " : "",
161              begin,end,trnname);
162     
163     trnout = open_tpx(trnname,"w");
164     if (WriteXref == eWXR_YES)
165     {
166         /* misuse lambda: 0/1 mass weighted fit no/yes */  
167         fwrite_trn(trnout,-1,-1,bDMR ? 1.0 : 0.0,zerobox,natoms,xref,NULL,NULL);
168     }
169     else if (WriteXref == eWXR_NOFIT)
170     {
171         /* misuse lambda: -1 no fit */  
172         fwrite_trn(trnout,-1,-1,-1.0,zerobox,natoms,x,NULL,NULL);
173     }
174     
175     /* misuse lambda: 0/1 mass weighted analysis no/yes */ 
176     fwrite_trn(trnout,0,0,bDMA ? 1.0 : 0.0,zerobox,natoms,xav,NULL,NULL);
177
178     for(i=0; i<=(end-begin); i++) 
179     {
180         
181         if (!bReverse)
182             vec = i;
183         else
184             vec = ndim-i-1;
185         
186         for (j=0; j<natoms; j++)
187             for(d=0; d<DIM; d++)
188                 x[j][d]=mat[vec*ndim+DIM*j+d];
189         
190         /* Store the eigenvalue in the time field */
191         fwrite_trn(trnout,begin+i,eigval[vec],0,zerobox,natoms,x,NULL,NULL);
192     }
193     close_trn(trnout);
194     
195     sfree(x);
196 }
197
198