Merge branch release-5-0 into release-5-1
[alexxy/gromacs.git] / src / gromacs / fileio / mtxio.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) 2012,2013,2014,2015, 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 "mtxio.h"
40
41 /* This module provides routines to read/write sparse or full storage
42  * matrices from/to files. It is normally used for the Hessian matrix
43  * in normal mode analysis.
44  */
45
46 #include "gromacs/fileio/gmxfio.h"
47 #include "gromacs/fileio/gmxfio-xdr.h"
48 #include "gromacs/linearalgebra/sparsematrix.h"
49 #include "gromacs/utility/baseversion.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
52
53 /* Just a number to identify our file type */
54 #define GMX_MTXIO_MAGIC_NUMBER  0x34ce8fd2
55
56 #define GMX_MTXIO_FULL_MATRIX     0
57 #define GMX_MTXIO_SPARSE_MATRIX   1
58
59
60
61 /* Matrix file format definition:
62  *
63  * All entries are stored in XDR format.
64  *
65  * 1. Magic number integer, should be GMX_MTXIO_MAGIC_NUMBER
66  * 2. An XDR string specifying the Gromacs version used to generate the file.
67  * 3. Integer to denote precision. 1 if double, 0 if single precision.
68  * 4. Two integers specifying number of rows and columns.
69  * 5. Integer to denote storage type:
70  *    GMX_MTXIO_FULL_MATRIX or GMX_MTXIO_SPARSE_MATRIX
71  *
72  * 6. Matrix data.
73  *    a) In case of full matrix, this is nrow*ncol floating-point values.
74  *    b) In case of sparse matrix the data is:
75  *       - Integer specifying compressed_symmetric format (1=yes, 0=no)
76  *       - Integer specifying number of rows (again)
77  *       - nrow integers specifying the number of data entries on each row ("ndata")
78  *       - All the actual entries for each row, stored contiguous.
79  *         Each entry consists of an integer column index and floating-point data value.
80  */
81
82 void gmx_mtxio_write(const char *             filename,
83                      int                      nrow,
84                      int                      ncol,
85                      real *                   full_matrix,
86                      gmx_sparsematrix_t *     sparse_matrix)
87 {
88     t_fileio   *fio;
89     int         i, j, prec;
90     size_t      sz;
91
92     if (full_matrix != NULL && sparse_matrix != NULL)
93     {
94         gmx_fatal(FARGS, "Both full AND sparse matrix specified to gmx_mtxio_write().\n");
95     }
96
97     fio = gmx_fio_open(filename, "w");
98
99     /* Write magic number */
100     i = GMX_MTXIO_MAGIC_NUMBER;
101     gmx_fio_do_int(fio, i);
102
103     /* Write generating Gromacs version */
104     gmx_fio_write_string(fio, gmx_version());
105
106     /* Write 1 for double, 0 for single precision */
107     if (sizeof(real) == sizeof(double))
108     {
109         prec = 1;
110     }
111     else
112     {
113         prec = 0;
114     }
115     gmx_fio_do_int(fio, prec);
116
117     gmx_fio_do_int(fio, nrow);
118     gmx_fio_do_int(fio, ncol);
119
120     if (full_matrix != NULL)
121     {
122         /* Full matrix storage format */
123         i = GMX_MTXIO_FULL_MATRIX;
124         gmx_fio_do_int(fio, i);
125         sz   = nrow*ncol;
126         gmx_fio_ndo_real(fio, full_matrix, sz);
127     }
128     else
129     {
130         /* Sparse storage */
131         i = GMX_MTXIO_SPARSE_MATRIX;
132         gmx_fio_do_int(fio, i);
133
134         gmx_fio_do_gmx_bool(fio, sparse_matrix->compressed_symmetric);
135         gmx_fio_do_int(fio, sparse_matrix->nrow);
136         if (sparse_matrix->nrow != nrow)
137         {
138             gmx_fatal(FARGS, "Internal inconsistency in sparse matrix.\n");
139         }
140         gmx_fio_ndo_int(fio, sparse_matrix->ndata, sparse_matrix->nrow);
141         for (i = 0; i < sparse_matrix->nrow; i++)
142         {
143             for (j = 0; j < sparse_matrix->ndata[i]; j++)
144             {
145                 gmx_fio_do_int(fio, sparse_matrix->data[i][j].col);
146                 gmx_fio_do_real(fio, sparse_matrix->data[i][j].value);
147             }
148         }
149     }
150     gmx_fio_close(fio);
151 }
152
153
154 void
155 gmx_mtxio_read (const char *            filename,
156                 int *                   nrow,
157                 int *                   ncol,
158                 real **                 full_matrix,
159                 gmx_sparsematrix_t **   sparse_matrix)
160 {
161     t_fileio   *fio;
162     int         i, j, prec;
163     char        gmxver[256];
164     size_t      sz;
165
166     fio = gmx_fio_open(filename, "r");
167
168     /* Read and check magic number */
169     i = GMX_MTXIO_MAGIC_NUMBER;
170     gmx_fio_do_int(fio, i);
171
172     if (i != GMX_MTXIO_MAGIC_NUMBER)
173     {
174         gmx_fatal(FARGS,
175                   "No matrix data found in file. Note that the Hessian matrix format changed\n"
176                   "in GROMACS 3.3 to enable portable files and sparse matrix storage.\n");
177     }
178
179     /* Read generating Gromacs version */
180     gmx_fio_do_string(fio, gmxver);
181
182     /* Write 1 for double, 0 for single precision */
183     if (sizeof(real) == sizeof(double))
184     {
185         prec = 1;
186     }
187     else
188     {
189         prec = 0;
190     }
191     gmx_fio_do_int(fio, prec);
192
193     fprintf(stderr, "Reading %s precision matrix generated by GROMACS %s\n",
194             (prec == 1) ? "double" : "single", gmxver);
195
196     gmx_fio_do_int(fio, i);
197     *nrow = i;
198     gmx_fio_do_int(fio, i);
199     *ncol = i;
200
201     gmx_fio_do_int(fio, i);
202
203     if (i == GMX_MTXIO_FULL_MATRIX && NULL != full_matrix)
204     {
205         printf("Full matrix storage format, nrow=%d, ncols=%d\n", *nrow, *ncol);
206
207         sz = (*nrow) * (*ncol);
208         snew((*full_matrix), sz);
209         gmx_fio_ndo_real(fio, (*full_matrix), sz);
210     }
211     else if (NULL != sparse_matrix)
212     {
213         /* Sparse storage */
214         printf("Sparse matrix storage format, nrow=%d, ncols=%d\n", *nrow, *ncol);
215
216         snew((*sparse_matrix), 1);
217         gmx_fio_do_gmx_bool(fio, (*sparse_matrix)->compressed_symmetric);
218         gmx_fio_do_int(fio, (*sparse_matrix)->nrow);
219         if ((*sparse_matrix)->nrow != *nrow)
220         {
221             gmx_fatal(FARGS, "Internal inconsistency in sparse matrix.\n");
222         }
223         snew((*sparse_matrix)->ndata, (*sparse_matrix)->nrow);
224         snew((*sparse_matrix)->nalloc, (*sparse_matrix)->nrow);
225         snew((*sparse_matrix)->data, (*sparse_matrix)->nrow);
226         gmx_fio_ndo_int(fio, (*sparse_matrix)->ndata,
227                         (*sparse_matrix)->nrow);
228
229         for (i = 0; i < (*sparse_matrix)->nrow; i++)
230         {
231             (*sparse_matrix)->nalloc[i] = (*sparse_matrix)->ndata[i] + 10;
232             snew(((*sparse_matrix)->data[i]), (*sparse_matrix)->nalloc[i]);
233
234             for (j = 0; j < (*sparse_matrix)->ndata[i]; j++)
235             {
236                 gmx_fio_do_int(fio, (*sparse_matrix)->data[i][j].col);
237                 gmx_fio_do_real(fio, (*sparse_matrix)->data[i][j].value);
238             }
239         }
240     }
241     gmx_fio_close(fio);
242 }