Split lines with many copyright years
[gromacs.git] / src / gromacs / gmxana / eigio.cpp
blob4433fb463c248804a46d769933ef19d6c7781b01
1 /*
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.
38 #include "gmxpre.h"
40 #include "eigio.h"
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,
49 int* natoms,
50 gmx_bool* bFit,
51 rvec** xref,
52 gmx_bool* bDMR,
53 rvec** xav,
54 gmx_bool* bDMA,
55 int* nvec,
56 int** eignr,
57 rvec*** eigvec,
58 real** eigval)
60 gmx_trr_header_t head;
61 int i, snew_size;
62 struct t_fileio* status;
63 rvec* x;
64 matrix box;
65 gmx_bool bOK;
67 *bDMR = FALSE;
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;
73 snew(*xav, *natoms);
74 gmx_trr_read_frame_data(status, &head, box, *xav, nullptr, nullptr);
76 if ((head.t >= -1.1) && (head.t <= -0.9))
78 snew(*xref, *natoms);
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);
85 if (*bFit)
87 fprintf(stderr, "Read %smass weighted reference structure with %d atoms from %s\n",
88 *bDMR ? "" : "non ", *natoms, file);
90 else
92 fprintf(stderr, "Eigenvectors in %s were determined without fitting\n", file);
93 sfree(*xref);
94 *xref = nullptr;
96 gmx_trr_read_frame_header(status, &head, &bOK);
97 gmx_trr_read_frame_data(status, &head, box, *xav, nullptr, nullptr);
99 else
101 *bFit = TRUE;
102 *xref = nullptr;
104 *bDMA = (head.lambda > 0.5);
105 if ((head.t <= -0.01) || (head.t >= 0.01))
107 fprintf(stderr,
108 "WARNING: %s does not start with t=0, which should be the "
109 "average structure. This might not be a eigenvector file. "
110 "Some things might go wrong.\n",
111 file);
113 else
115 fprintf(stderr, "Read %smass weighted average/minimum structure with %d atoms from %s\n",
116 *bDMA ? "" : "non ", *natoms, file);
119 snew(x, *natoms);
120 snew_size = 10;
121 snew(*eignr, snew_size);
122 snew(*eigval, snew_size);
123 snew(*eigvec, snew_size);
125 *nvec = 0;
126 while (gmx_trr_read_frame_header(status, &head, &bOK))
128 gmx_trr_read_frame_data(status, &head, box, x, nullptr, nullptr);
129 if (*nvec >= snew_size)
131 snew_size += 10;
132 srenew(*eignr, snew_size);
133 srenew(*eigval, snew_size);
134 srenew(*eigvec, snew_size);
136 i = head.step;
137 (*eigval)[*nvec] = head.t;
138 (*eignr)[*nvec] = i - 1;
139 snew((*eigvec)[*nvec], *natoms);
140 for (i = 0; i < *natoms; i++)
142 copy_rvec(x[i], (*eigvec)[*nvec][i]);
144 (*nvec)++;
146 sfree(x);
147 gmx_trr_close(status);
148 fprintf(stderr, "Read %d eigenvectors (for %d atoms)\n\n", *nvec, *natoms);
152 void write_eigenvectors(const char* trrname,
153 int natoms,
154 const real mat[],
155 gmx_bool bReverse,
156 int begin,
157 int end,
158 int WriteXref,
159 const rvec* xref,
160 gmx_bool bDMR,
161 const rvec xav[],
162 gmx_bool bDMA,
163 const real eigval[])
165 struct t_fileio* trrout;
166 int ndim, i, j, d, vec;
167 matrix zerobox;
168 rvec* x;
170 ndim = natoms * DIM;
171 clear_mat(zerobox);
172 snew(x, natoms);
174 fprintf(stderr, "\nWriting %saverage structure & eigenvectors %d--%d to %s\n",
175 (WriteXref == eWXR_YES) ? "reference, " : "", begin, end, trrname);
177 trrout = gmx_trr_open(trrname, "w");
178 if (WriteXref == eWXR_YES)
180 /* misuse lambda: 0/1 mass weighted fit no/yes */
181 gmx_trr_write_frame(trrout, -1, -1, bDMR ? 1.0 : 0.0, zerobox, natoms, xref, nullptr, nullptr);
183 else if (WriteXref == eWXR_NOFIT)
185 /* misuse lambda: -1 no fit */
186 gmx_trr_write_frame(trrout, -1, -1, -1.0, zerobox, natoms, x, nullptr, nullptr);
189 /* misuse lambda: 0/1 mass weighted analysis no/yes */
190 gmx_trr_write_frame(trrout, 0, 0, bDMA ? 1.0 : 0.0, zerobox, natoms, xav, nullptr, nullptr);
192 for (i = 0; i <= (end - begin); i++)
195 if (!bReverse)
197 vec = i;
199 else
201 vec = ndim - i - 1;
204 for (j = 0; j < natoms; j++)
206 for (d = 0; d < DIM; d++)
208 x[j][d] = mat[vec * ndim + DIM * j + d];
212 /* Store the eigenvalue in the time field */
213 gmx_trr_write_frame(trrout, begin + i, eigval[vec], 0, zerobox, natoms, x, nullptr, nullptr);
215 gmx_trr_close(trrout);
217 sfree(x);