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-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source 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.
42 #include "gromacs/fileio/tpxio.h"
43 #include "gromacs/fileio/trrio.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/futil.h"
46 #include "gromacs/utility/smalloc.h"
48 void read_eigenvectors(const char* file,
60 gmx_trr_header_t head;
62 struct t_fileio* status;
69 /* read (reference (t=-1) and) average (t=0) structure */
70 status = gmx_trr_open(file, "r");
71 gmx_trr_read_frame_header(status, &head, &bOK);
72 *natoms = head.natoms;
74 gmx_trr_read_frame_data(status, &head, box, *xav, nullptr, nullptr);
76 if ((head.t >= -1.1) && (head.t <= -0.9))
79 for (i = 0; i < *natoms; i++)
81 copy_rvec((*xav)[i], (*xref)[i]);
83 *bDMR = (head.lambda > 0.5);
84 *bFit = (head.lambda > -0.5);
88 "Read %smass weighted reference structure with %d atoms from %s\n",
95 fprintf(stderr, "Eigenvectors in %s were determined without fitting\n", file);
99 gmx_trr_read_frame_header(status, &head, &bOK);
100 gmx_trr_read_frame_data(status, &head, box, *xav, nullptr, nullptr);
107 *bDMA = (head.lambda > 0.5);
108 if ((head.t <= -0.01) || (head.t >= 0.01))
111 "WARNING: %s does not start with t=0, which should be the "
112 "average structure. This might not be a eigenvector file. "
113 "Some things might go wrong.\n",
119 "Read %smass weighted average/minimum structure with %d atoms from %s\n",
127 snew(*eignr, snew_size);
128 snew(*eigval, snew_size);
129 snew(*eigvec, snew_size);
132 while (gmx_trr_read_frame_header(status, &head, &bOK))
134 gmx_trr_read_frame_data(status, &head, box, x, nullptr, nullptr);
135 if (*nvec >= snew_size)
138 srenew(*eignr, snew_size);
139 srenew(*eigval, snew_size);
140 srenew(*eigvec, snew_size);
143 (*eigval)[*nvec] = head.t;
144 (*eignr)[*nvec] = i - 1;
145 snew((*eigvec)[*nvec], *natoms);
146 for (i = 0; i < *natoms; i++)
148 copy_rvec(x[i], (*eigvec)[*nvec][i]);
153 gmx_trr_close(status);
154 fprintf(stderr, "Read %d eigenvectors (for %d atoms)\n\n", *nvec, *natoms);
158 void write_eigenvectors(const char* trrname,
171 struct t_fileio* trrout;
172 int ndim, i, j, d, vec;
181 "\nWriting %saverage structure & eigenvectors %d--%d to %s\n",
182 (WriteXref == eWXR_YES) ? "reference, " : "",
187 trrout = gmx_trr_open(trrname, "w");
188 if (WriteXref == eWXR_YES)
190 /* misuse lambda: 0/1 mass weighted fit no/yes */
191 gmx_trr_write_frame(trrout, -1, -1, bDMR ? 1.0 : 0.0, zerobox, natoms, xref, nullptr, nullptr);
193 else if (WriteXref == eWXR_NOFIT)
195 /* misuse lambda: -1 no fit */
196 gmx_trr_write_frame(trrout, -1, -1, -1.0, zerobox, natoms, x, nullptr, nullptr);
199 /* misuse lambda: 0/1 mass weighted analysis no/yes */
200 gmx_trr_write_frame(trrout, 0, 0, bDMA ? 1.0 : 0.0, zerobox, natoms, xav, nullptr, nullptr);
202 for (i = 0; i <= (end - begin); i++)
214 for (j = 0; j < natoms; j++)
216 for (d = 0; d < DIM; d++)
218 x[j][d] = mat[vec * ndim + DIM * j + d];
222 /* Store the eigenvalue in the time field */
223 gmx_trr_write_frame(trrout, begin + i, eigval[vec], 0, zerobox, natoms, x, nullptr, nullptr);
225 gmx_trr_close(trrout);