Apply re-formatting to C++ in src/ tree.
[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 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "eigio.h"
41
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"
47
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)
59 {
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;
66
67     *bDMR = FALSE;
68
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);
75
76     if ((head.t >= -1.1) && (head.t <= -0.9))
77     {
78         snew(*xref, *natoms);
79         for (i = 0; i < *natoms; i++)
80         {
81             copy_rvec((*xav)[i], (*xref)[i]);
82         }
83         *bDMR = (head.lambda > 0.5);
84         *bFit = (head.lambda > -0.5);
85         if (*bFit)
86         {
87             fprintf(stderr,
88                     "Read %smass weighted reference structure with %d atoms from %s\n",
89                     *bDMR ? "" : "non ",
90                     *natoms,
91                     file);
92         }
93         else
94         {
95             fprintf(stderr, "Eigenvectors in %s were determined without fitting\n", file);
96             sfree(*xref);
97             *xref = nullptr;
98         }
99         gmx_trr_read_frame_header(status, &head, &bOK);
100         gmx_trr_read_frame_data(status, &head, box, *xav, nullptr, nullptr);
101     }
102     else
103     {
104         *bFit = TRUE;
105         *xref = nullptr;
106     }
107     *bDMA = (head.lambda > 0.5);
108     if ((head.t <= -0.01) || (head.t >= 0.01))
109     {
110         fprintf(stderr,
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",
114                 file);
115     }
116     else
117     {
118         fprintf(stderr,
119                 "Read %smass weighted average/minimum structure with %d atoms from %s\n",
120                 *bDMA ? "" : "non ",
121                 *natoms,
122                 file);
123     }
124
125     snew(x, *natoms);
126     snew_size = 10;
127     snew(*eignr, snew_size);
128     snew(*eigval, snew_size);
129     snew(*eigvec, snew_size);
130
131     *nvec = 0;
132     while (gmx_trr_read_frame_header(status, &head, &bOK))
133     {
134         gmx_trr_read_frame_data(status, &head, box, x, nullptr, nullptr);
135         if (*nvec >= snew_size)
136         {
137             snew_size += 10;
138             srenew(*eignr, snew_size);
139             srenew(*eigval, snew_size);
140             srenew(*eigvec, snew_size);
141         }
142         i                = head.step;
143         (*eigval)[*nvec] = head.t;
144         (*eignr)[*nvec]  = i - 1;
145         snew((*eigvec)[*nvec], *natoms);
146         for (i = 0; i < *natoms; i++)
147         {
148             copy_rvec(x[i], (*eigvec)[*nvec][i]);
149         }
150         (*nvec)++;
151     }
152     sfree(x);
153     gmx_trr_close(status);
154     fprintf(stderr, "Read %d eigenvectors (for %d atoms)\n\n", *nvec, *natoms);
155 }
156
157
158 void write_eigenvectors(const char* trrname,
159                         int         natoms,
160                         const real  mat[],
161                         gmx_bool    bReverse,
162                         int         begin,
163                         int         end,
164                         int         WriteXref,
165                         const rvec* xref,
166                         gmx_bool    bDMR,
167                         const rvec  xav[],
168                         gmx_bool    bDMA,
169                         const real  eigval[])
170 {
171     struct t_fileio* trrout;
172     int              ndim, i, j, d, vec;
173     matrix           zerobox;
174     rvec*            x;
175
176     ndim = natoms * DIM;
177     clear_mat(zerobox);
178     snew(x, natoms);
179
180     fprintf(stderr,
181             "\nWriting %saverage structure & eigenvectors %d--%d to %s\n",
182             (WriteXref == eWXR_YES) ? "reference, " : "",
183             begin,
184             end,
185             trrname);
186
187     trrout = gmx_trr_open(trrname, "w");
188     if (WriteXref == eWXR_YES)
189     {
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);
192     }
193     else if (WriteXref == eWXR_NOFIT)
194     {
195         /* misuse lambda: -1 no fit */
196         gmx_trr_write_frame(trrout, -1, -1, -1.0, zerobox, natoms, x, nullptr, nullptr);
197     }
198
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);
201
202     for (i = 0; i <= (end - begin); i++)
203     {
204
205         if (!bReverse)
206         {
207             vec = i;
208         }
209         else
210         {
211             vec = ndim - i - 1;
212         }
213
214         for (j = 0; j < natoms; j++)
215         {
216             for (d = 0; d < DIM; d++)
217             {
218                 x[j][d] = mat[vec * ndim + DIM * j + d];
219             }
220         }
221
222         /* Store the eigenvalue in the time field */
223         gmx_trr_write_frame(trrout, begin + i, eigval[vec], 0, zerobox, natoms, x, nullptr, nullptr);
224     }
225     gmx_trr_close(trrout);
226
227     sfree(x);
228 }