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