3cea7e812a1227272e822d51b101a0cd751d89fa
[alexxy/gromacs.git] / src / gromacs / gmxana / eigio.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "eigio.h"
40
41 #include "gromacs/fileio/tpxio.h"
42 #include "gromacs/fileio/trrio.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/utility/futil.h"
45 #include "gromacs/utility/smalloc.h"
46
47 void read_eigenvectors(const char* file,
48                        int*        natoms,
49                        gmx_bool*   bFit,
50                        rvec**      xref,
51                        gmx_bool*   bDMR,
52                        rvec**      xav,
53                        gmx_bool*   bDMA,
54                        int*        nvec,
55                        int**       eignr,
56                        rvec***     eigvec,
57                        real**      eigval)
58 {
59     gmx_trr_header_t head;
60     int              i, snew_size;
61     struct t_fileio* status;
62     rvec*            x;
63     matrix           box;
64     gmx_bool         bOK;
65
66     *bDMR = FALSE;
67
68     /* read (reference (t=-1) and) average (t=0) structure */
69     status = gmx_trr_open(file, "r");
70     gmx_trr_read_frame_header(status, &head, &bOK);
71     *natoms = head.natoms;
72     snew(*xav, *natoms);
73     gmx_trr_read_frame_data(status, &head, box, *xav, nullptr, nullptr);
74
75     if ((head.t >= -1.1) && (head.t <= -0.9))
76     {
77         snew(*xref, *natoms);
78         for (i = 0; i < *natoms; i++)
79         {
80             copy_rvec((*xav)[i], (*xref)[i]);
81         }
82         *bDMR = (head.lambda > 0.5);
83         *bFit = (head.lambda > -0.5);
84         if (*bFit)
85         {
86             fprintf(stderr, "Read %smass weighted reference structure with %d atoms from %s\n",
87                     *bDMR ? "" : "non ", *natoms, file);
88         }
89         else
90         {
91             fprintf(stderr, "Eigenvectors in %s were determined without fitting\n", file);
92             sfree(*xref);
93             *xref = nullptr;
94         }
95         gmx_trr_read_frame_header(status, &head, &bOK);
96         gmx_trr_read_frame_data(status, &head, box, *xav, nullptr, nullptr);
97     }
98     else
99     {
100         *bFit = TRUE;
101         *xref = nullptr;
102     }
103     *bDMA = (head.lambda > 0.5);
104     if ((head.t <= -0.01) || (head.t >= 0.01))
105     {
106         fprintf(stderr,
107                 "WARNING: %s does not start with t=0, which should be the "
108                 "average structure. This might not be a eigenvector file. "
109                 "Some things might go wrong.\n",
110                 file);
111     }
112     else
113     {
114         fprintf(stderr, "Read %smass weighted average/minimum structure with %d atoms from %s\n",
115                 *bDMA ? "" : "non ", *natoms, file);
116     }
117
118     snew(x, *natoms);
119     snew_size = 10;
120     snew(*eignr, snew_size);
121     snew(*eigval, snew_size);
122     snew(*eigvec, snew_size);
123
124     *nvec = 0;
125     while (gmx_trr_read_frame_header(status, &head, &bOK))
126     {
127         gmx_trr_read_frame_data(status, &head, box, x, nullptr, nullptr);
128         if (*nvec >= snew_size)
129         {
130             snew_size += 10;
131             srenew(*eignr, snew_size);
132             srenew(*eigval, snew_size);
133             srenew(*eigvec, snew_size);
134         }
135         i                = head.step;
136         (*eigval)[*nvec] = head.t;
137         (*eignr)[*nvec]  = i - 1;
138         snew((*eigvec)[*nvec], *natoms);
139         for (i = 0; i < *natoms; i++)
140         {
141             copy_rvec(x[i], (*eigvec)[*nvec][i]);
142         }
143         (*nvec)++;
144     }
145     sfree(x);
146     gmx_trr_close(status);
147     fprintf(stderr, "Read %d eigenvectors (for %d atoms)\n\n", *nvec, *natoms);
148 }
149
150
151 void write_eigenvectors(const char* trrname,
152                         int         natoms,
153                         const real  mat[],
154                         gmx_bool    bReverse,
155                         int         begin,
156                         int         end,
157                         int         WriteXref,
158                         const rvec* xref,
159                         gmx_bool    bDMR,
160                         const rvec  xav[],
161                         gmx_bool    bDMA,
162                         const real  eigval[])
163 {
164     struct t_fileio* trrout;
165     int              ndim, i, j, d, vec;
166     matrix           zerobox;
167     rvec*            x;
168
169     ndim = natoms * DIM;
170     clear_mat(zerobox);
171     snew(x, natoms);
172
173     fprintf(stderr, "\nWriting %saverage structure & eigenvectors %d--%d to %s\n",
174             (WriteXref == eWXR_YES) ? "reference, " : "", begin, end, trrname);
175
176     trrout = gmx_trr_open(trrname, "w");
177     if (WriteXref == eWXR_YES)
178     {
179         /* misuse lambda: 0/1 mass weighted fit no/yes */
180         gmx_trr_write_frame(trrout, -1, -1, bDMR ? 1.0 : 0.0, zerobox, natoms, xref, nullptr, nullptr);
181     }
182     else if (WriteXref == eWXR_NOFIT)
183     {
184         /* misuse lambda: -1 no fit */
185         gmx_trr_write_frame(trrout, -1, -1, -1.0, zerobox, natoms, x, nullptr, nullptr);
186     }
187
188     /* misuse lambda: 0/1 mass weighted analysis no/yes */
189     gmx_trr_write_frame(trrout, 0, 0, bDMA ? 1.0 : 0.0, zerobox, natoms, xav, nullptr, nullptr);
190
191     for (i = 0; i <= (end - begin); i++)
192     {
193
194         if (!bReverse)
195         {
196             vec = i;
197         }
198         else
199         {
200             vec = ndim - i - 1;
201         }
202
203         for (j = 0; j < natoms; j++)
204         {
205             for (d = 0; d < DIM; d++)
206             {
207                 x[j][d] = mat[vec * ndim + DIM * j + d];
208             }
209         }
210
211         /* Store the eigenvalue in the time field */
212         gmx_trr_write_frame(trrout, begin + i, eigval[vec], 0, zerobox, natoms, x, nullptr, nullptr);
213     }
214     gmx_trr_close(trrout);
215
216     sfree(x);
217 }