3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
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)
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;
67 fread_htrn(status, &head, box, *xav, NULL, NULL);
69 if ((head.t >= -1.1) && (head.t <= -0.9))
72 for (i = 0; i < *natoms; i++)
74 copy_rvec((*xav)[i], (*xref)[i]);
76 *bDMR = (head.lambda > 0.5);
77 *bFit = (head.lambda > -0.5);
80 fprintf(stderr, "Read %smass weighted reference structure with %d atoms from %s\n", *bDMR ? "" : "non ", *natoms, file);
84 fprintf(stderr, "Eigenvectors in %s were determined without fitting\n", file);
88 fread_trnheader(status, &head, &bOK);
89 fread_htrn(status, &head, box, *xav, NULL, NULL);
96 *bDMA = (head.lambda > 0.5);
97 if ((head.t <= -0.01) || (head.t >= 0.01))
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",
107 "Read %smass weighted average/minimum structure with %d atoms from %s\n",
108 *bDMA ? "" : "non ", *natoms, file);
113 snew(*eignr, snew_size);
114 snew(*eigval, snew_size);
115 snew(*eigvec, snew_size);
118 while (fread_trnheader(status, &head, &bOK))
120 fread_htrn(status, &head, box, x, NULL, NULL);
121 if (*nvec >= snew_size)
124 srenew(*eignr, snew_size);
125 srenew(*eigval, snew_size);
126 srenew(*eigvec, snew_size);
129 (*eigval)[*nvec] = head.t;
130 (*eignr)[*nvec] = i-1;
131 snew((*eigvec)[*nvec], *natoms);
132 for (i = 0; i < *natoms; i++)
134 copy_rvec(x[i], (*eigvec)[*nvec][i]);
139 fprintf(stderr, "Read %d eigenvectors (for %d atoms)\n\n", *nvec, *natoms);
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[])
149 int ndim, i, j, d, vec;
158 "\nWriting %saverage structure & eigenvectors %d--%d to %s\n",
159 (WriteXref == eWXR_YES) ? "reference, " : "",
160 begin, end, trnname);
162 trnout = open_tpx(trnname, "w");
163 if (WriteXref == eWXR_YES)
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);
168 else if (WriteXref == eWXR_NOFIT)
170 /* misuse lambda: -1 no fit */
171 fwrite_trn(trnout, -1, -1, -1.0, zerobox, natoms, x, NULL, NULL);
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);
177 for (i = 0; i <= (end-begin); i++)
189 for (j = 0; j < natoms; j++)
191 for (d = 0; d < DIM; d++)
193 x[j][d] = mat[vec*ndim+DIM*j+d];
197 /* Store the eigenvalue in the time field */
198 fwrite_trn(trnout, begin+i, eigval[vec], 0, zerobox, natoms, x, NULL, NULL);