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