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