ae9b15d9e8a79bda91737282740a5e521ccc6d8a
[alexxy/gromacs.git] / src / gromacs / fileio / mtxio.c
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     gmx_bool    bDum  = TRUE;
91     gmx_bool    bRead = FALSE;
92     size_t      sz;
93
94     if (full_matrix != NULL && sparse_matrix != NULL)
95     {
96         gmx_fatal(FARGS, "Both full AND sparse matrix specified to gmx_mtxio_write().\n");
97     }
98
99     fio = gmx_fio_open(filename, "w");
100
101     /* Write magic number */
102     i = GMX_MTXIO_MAGIC_NUMBER;
103     gmx_fio_do_int(fio, i);
104
105     /* Write generating Gromacs version */
106     gmx_fio_write_string(fio, gmx_version());
107
108     /* Write 1 for double, 0 for single precision */
109     if (sizeof(real) == sizeof(double))
110     {
111         prec = 1;
112     }
113     else
114     {
115         prec = 0;
116     }
117     gmx_fio_do_int(fio, prec);
118
119     gmx_fio_do_int(fio, nrow);
120     gmx_fio_do_int(fio, ncol);
121
122     if (full_matrix != NULL)
123     {
124         /* Full matrix storage format */
125         i = GMX_MTXIO_FULL_MATRIX;
126         gmx_fio_do_int(fio, i);
127         sz   = nrow*ncol;
128         bDum = gmx_fio_ndo_real(fio, full_matrix, sz);
129     }
130     else
131     {
132         /* Sparse storage */
133         i = GMX_MTXIO_SPARSE_MATRIX;
134         gmx_fio_do_int(fio, i);
135
136         gmx_fio_do_gmx_bool(fio, sparse_matrix->compressed_symmetric);
137         gmx_fio_do_int(fio, sparse_matrix->nrow);
138         if (sparse_matrix->nrow != nrow)
139         {
140             gmx_fatal(FARGS, "Internal inconsistency in sparse matrix.\n");
141         }
142         bDum = gmx_fio_ndo_int(fio, sparse_matrix->ndata, sparse_matrix->nrow);
143         for (i = 0; i < sparse_matrix->nrow; i++)
144         {
145             for (j = 0; j < sparse_matrix->ndata[i]; j++)
146             {
147                 gmx_fio_do_int(fio, sparse_matrix->data[i][j].col);
148                 gmx_fio_do_real(fio, sparse_matrix->data[i][j].value);
149             }
150         }
151     }
152     gmx_fio_close(fio);
153 }
154
155
156 void
157 gmx_mtxio_read (const char *            filename,
158                 int *                   nrow,
159                 int *                   ncol,
160                 real **                 full_matrix,
161                 gmx_sparsematrix_t **   sparse_matrix)
162 {
163     t_fileio   *fio;
164     int         i, j, prec;
165     gmx_bool    bDum  = TRUE;
166     gmx_bool    bRead = TRUE;
167     char        gmxver[256];
168     size_t      sz;
169
170     fio = gmx_fio_open(filename, "r");
171
172     /* Read and check magic number */
173     i = GMX_MTXIO_MAGIC_NUMBER;
174     gmx_fio_do_int(fio, i);
175
176     if (i != GMX_MTXIO_MAGIC_NUMBER)
177     {
178         gmx_fatal(FARGS,
179                   "No matrix data found in file. Note that the Hessian matrix format changed\n"
180                   "in GROMACS 3.3 to enable portable files and sparse matrix storage.\n");
181     }
182
183     /* Read generating Gromacs version */
184     gmx_fio_do_string(fio, gmxver);
185
186     /* Write 1 for double, 0 for single precision */
187     if (sizeof(real) == sizeof(double))
188     {
189         prec = 1;
190     }
191     else
192     {
193         prec = 0;
194     }
195     gmx_fio_do_int(fio, prec);
196
197     fprintf(stderr, "Reading %s precision matrix generated by GROMACS %s\n",
198             (prec == 1) ? "double" : "single", gmxver);
199
200     gmx_fio_do_int(fio, i);
201     *nrow = i;
202     gmx_fio_do_int(fio, i);
203     *ncol = i;
204
205     gmx_fio_do_int(fio, i);
206
207     if (i == GMX_MTXIO_FULL_MATRIX && NULL != full_matrix)
208     {
209         printf("Full matrix storage format, nrow=%d, ncols=%d\n", *nrow, *ncol);
210
211         sz = (*nrow) * (*ncol);
212         snew((*full_matrix), sz);
213         bDum = gmx_fio_ndo_real(fio, (*full_matrix), sz);
214     }
215     else if (NULL != sparse_matrix)
216     {
217         /* Sparse storage */
218         printf("Sparse matrix storage format, nrow=%d, ncols=%d\n", *nrow, *ncol);
219
220         snew((*sparse_matrix), 1);
221         gmx_fio_do_gmx_bool(fio, (*sparse_matrix)->compressed_symmetric);
222         gmx_fio_do_int(fio, (*sparse_matrix)->nrow);
223         if ((*sparse_matrix)->nrow != *nrow)
224         {
225             gmx_fatal(FARGS, "Internal inconsistency in sparse matrix.\n");
226         }
227         snew((*sparse_matrix)->ndata, (*sparse_matrix)->nrow);
228         snew((*sparse_matrix)->nalloc, (*sparse_matrix)->nrow);
229         snew((*sparse_matrix)->data, (*sparse_matrix)->nrow);
230         bDum = gmx_fio_ndo_int(fio, (*sparse_matrix)->ndata,
231                                (*sparse_matrix)->nrow);
232
233         for (i = 0; i < (*sparse_matrix)->nrow; i++)
234         {
235             (*sparse_matrix)->nalloc[i] = (*sparse_matrix)->ndata[i] + 10;
236             snew(((*sparse_matrix)->data[i]), (*sparse_matrix)->nalloc[i]);
237
238             for (j = 0; j < (*sparse_matrix)->ndata[i]; j++)
239             {
240                 gmx_fio_do_int(fio, (*sparse_matrix)->data[i][j].col);
241                 gmx_fio_do_real(fio, (*sparse_matrix)->data[i][j].value);
242             }
243         }
244     }
245     gmx_fio_close(fio);
246 }