Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / gmxana / eigio.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "smalloc.h"
40 #include "vec.h"
41 #include "eigio.h"
42 #include "trnio.h"
43 #include "tpxio.h"
44 #include "statutil.h"
45 #include "futil.h"
46
47 void read_eigenvectors(const char *file, int *natoms, gmx_bool *bFit,
48                        rvec **xref, gmx_bool *bDMR,
49                        rvec **xav, gmx_bool *bDMA,
50                        int *nvec, int **eignr,
51                        rvec ***eigvec, real **eigval)
52 {
53     t_trnheader head;
54     int         i, snew_size;
55     t_fileio   *status;
56     rvec       *x;
57     matrix      box;
58     gmx_bool    bOK;
59
60     *bDMR = FALSE;
61
62     /* read (reference (t=-1) and) average (t=0) structure */
63     status = open_trn(file, "r");
64     fread_trnheader(status, &head, &bOK);
65     *natoms = head.natoms;
66     snew(*xav, *natoms);
67     fread_htrn(status, &head, box, *xav, NULL, NULL);
68
69     if ((head.t >= -1.1) && (head.t <= -0.9))
70     {
71         snew(*xref, *natoms);
72         for (i = 0; i < *natoms; i++)
73         {
74             copy_rvec((*xav)[i], (*xref)[i]);
75         }
76         *bDMR = (head.lambda > 0.5);
77         *bFit = (head.lambda > -0.5);
78         if (*bFit)
79         {
80             fprintf(stderr, "Read %smass weighted reference structure with %d atoms from %s\n", *bDMR ? "" : "non ", *natoms, file);
81         }
82         else
83         {
84             fprintf(stderr, "Eigenvectors in %s were determined without fitting\n", file);
85             sfree(*xref);
86             *xref = NULL;
87         }
88         fread_trnheader(status, &head, &bOK);
89         fread_htrn(status, &head, box, *xav, NULL, NULL);
90     }
91     else
92     {
93         *bFit = TRUE;
94         *xref = NULL;
95     }
96     *bDMA = (head.lambda > 0.5);
97     if ((head.t <= -0.01) || (head.t >= 0.01))
98     {
99         fprintf(stderr, "WARNING: %s does not start with t=0, which should be the "
100                 "average structure. This might not be a eigenvector file. "
101                 "Some things might go wrong.\n",
102                 file);
103     }
104     else
105     {
106         fprintf(stderr,
107                 "Read %smass weighted average/minimum structure with %d atoms from %s\n",
108                 *bDMA ? "" : "non ", *natoms, file);
109     }
110
111     snew(x, *natoms);
112     snew_size = 10;
113     snew(*eignr, snew_size);
114     snew(*eigval, snew_size);
115     snew(*eigvec, snew_size);
116
117     *nvec = 0;
118     while (fread_trnheader(status, &head, &bOK))
119     {
120         fread_htrn(status, &head, box, x, NULL, NULL);
121         if (*nvec >= snew_size)
122         {
123             snew_size += 10;
124             srenew(*eignr, snew_size);
125             srenew(*eigval, snew_size);
126             srenew(*eigvec, snew_size);
127         }
128         i                = head.step;
129         (*eigval)[*nvec] = head.t;
130         (*eignr)[*nvec]  = i-1;
131         snew((*eigvec)[*nvec], *natoms);
132         for (i = 0; i < *natoms; i++)
133         {
134             copy_rvec(x[i], (*eigvec)[*nvec][i]);
135         }
136         (*nvec)++;
137     }
138     sfree(x);
139     fprintf(stderr, "Read %d eigenvectors (for %d atoms)\n\n", *nvec, *natoms);
140 }
141
142
143 void write_eigenvectors(const char *trnname, int natoms, real mat[],
144                         gmx_bool bReverse, int begin, int end,
145                         int WriteXref, rvec *xref, gmx_bool bDMR,
146                         rvec xav[], gmx_bool bDMA, real eigval[])
147 {
148     t_fileio *trnout;
149     int       ndim, i, j, d, vec;
150     matrix    zerobox;
151     rvec     *x;
152
153     ndim = natoms*DIM;
154     clear_mat(zerobox);
155     snew(x, natoms);
156
157     fprintf (stderr,
158              "\nWriting %saverage structure & eigenvectors %d--%d to %s\n",
159              (WriteXref == eWXR_YES) ? "reference, " : "",
160              begin, end, trnname);
161
162     trnout = open_tpx(trnname, "w");
163     if (WriteXref == eWXR_YES)
164     {
165         /* misuse lambda: 0/1 mass weighted fit no/yes */
166         fwrite_trn(trnout, -1, -1, bDMR ? 1.0 : 0.0, zerobox, natoms, xref, NULL, NULL);
167     }
168     else if (WriteXref == eWXR_NOFIT)
169     {
170         /* misuse lambda: -1 no fit */
171         fwrite_trn(trnout, -1, -1, -1.0, zerobox, natoms, x, NULL, NULL);
172     }
173
174     /* misuse lambda: 0/1 mass weighted analysis no/yes */
175     fwrite_trn(trnout, 0, 0, bDMA ? 1.0 : 0.0, zerobox, natoms, xav, NULL, NULL);
176
177     for (i = 0; i <= (end-begin); i++)
178     {
179
180         if (!bReverse)
181         {
182             vec = i;
183         }
184         else
185         {
186             vec = ndim-i-1;
187         }
188
189         for (j = 0; j < natoms; j++)
190         {
191             for (d = 0; d < DIM; d++)
192             {
193                 x[j][d] = mat[vec*ndim+DIM*j+d];
194             }
195         }
196
197         /* Store the eigenvalue in the time field */
198         fwrite_trn(trnout, begin+i, eigval[vec], 0, zerobox, natoms, x, NULL, NULL);
199     }
200     close_trn(trnout);
201
202     sfree(x);
203 }