2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
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.
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/fileio/xdrf.h"
50 #include "gromacs/linearalgebra/sparsematrix.h"
51 #include "gromacs/utility/baseversion.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/smalloc.h"
55 /* Just a number to identify our file type */
56 #define GMX_MTXIO_MAGIC_NUMBER 0x34ce8fd2
58 #define GMX_MTXIO_FULL_MATRIX 0
59 #define GMX_MTXIO_SPARSE_MATRIX 1
63 /* Matrix file format definition:
65 * All entries are stored in XDR format.
67 * 1. Magic number integer, should be GMX_MTXIO_MAGIC_NUMBER
68 * 2. An XDR string specifying the Gromacs version used to generate the file.
69 * 3. Integer to denote precision. 1 if double, 0 if single precision.
70 * 4. Two integers specifying number of rows and columns.
71 * 5. Integer to denote storage type:
72 * GMX_MTXIO_FULL_MATRIX or GMX_MTXIO_SPARSE_MATRIX
75 * a) In case of full matrix, this is nrow*ncol floating-point values.
76 * b) In case of sparse matrix the data is:
77 * - Integer specifying compressed_symmetric format (1=yes, 0=no)
78 * - Integer specifying number of rows (again)
79 * - nrow integers specifying the number of data entries on each row ("ndata")
80 * - All the actual entries for each row, stored contiguous.
81 * Each entry consists of an integer column index and floating-point data value.
84 void gmx_mtxio_write(const char * filename,
88 gmx_sparsematrix_t * sparse_matrix)
94 gmx_bool bRead = FALSE;
97 if (full_matrix != NULL && sparse_matrix != NULL)
99 gmx_fatal(FARGS, "Both full AND sparse matrix specified to gmx_mtxio_write().\n");
102 fio = gmx_fio_open(filename, "w");
103 gmx_fio_checktype(fio);
104 xd = gmx_fio_getxdr(fio);
106 /* Write magic number */
107 i = GMX_MTXIO_MAGIC_NUMBER;
108 gmx_fio_do_int(fio, i);
110 /* Write generating Gromacs version */
111 gmx_fio_write_string(fio, gmx_version());
113 /* Write 1 for double, 0 for single precision */
114 if (sizeof(real) == sizeof(double))
122 gmx_fio_do_int(fio, prec);
124 gmx_fio_do_int(fio, nrow);
125 gmx_fio_do_int(fio, ncol);
127 if (full_matrix != NULL)
129 /* Full matrix storage format */
130 i = GMX_MTXIO_FULL_MATRIX;
131 gmx_fio_do_int(fio, i);
133 bDum = gmx_fio_ndo_real(fio, full_matrix, sz);
138 i = GMX_MTXIO_SPARSE_MATRIX;
139 gmx_fio_do_int(fio, i);
141 gmx_fio_do_gmx_bool(fio, sparse_matrix->compressed_symmetric);
142 gmx_fio_do_int(fio, sparse_matrix->nrow);
143 if (sparse_matrix->nrow != nrow)
145 gmx_fatal(FARGS, "Internal inconsistency in sparse matrix.\n");
147 bDum = gmx_fio_ndo_int(fio, sparse_matrix->ndata, sparse_matrix->nrow);
148 for (i = 0; i < sparse_matrix->nrow; i++)
150 for (j = 0; j < sparse_matrix->ndata[i]; j++)
152 gmx_fio_do_int(fio, sparse_matrix->data[i][j].col);
153 gmx_fio_do_real(fio, sparse_matrix->data[i][j].value);
162 gmx_mtxio_read (const char * filename,
166 gmx_sparsematrix_t ** sparse_matrix)
171 gmx_bool bDum = TRUE;
172 gmx_bool bRead = TRUE;
176 fio = gmx_fio_open(filename, "r");
177 gmx_fio_checktype(fio);
178 xd = gmx_fio_getxdr(fio);
180 /* Read and check magic number */
181 i = GMX_MTXIO_MAGIC_NUMBER;
182 gmx_fio_do_int(fio, i);
184 if (i != GMX_MTXIO_MAGIC_NUMBER)
187 "No matrix data found in file. Note that the Hessian matrix format changed\n"
188 "in Gromacs 3.3 to enable portable files and sparse matrix storage.\n");
191 /* Read generating Gromacs version */
192 gmx_fio_do_string(fio, gmxver);
194 /* Write 1 for double, 0 for single precision */
195 if (sizeof(real) == sizeof(double))
203 gmx_fio_do_int(fio, prec);
205 fprintf(stderr, "Reading %s precision matrix generated by Gromacs %s\n",
206 (prec == 1) ? "double" : "single", gmxver);
208 gmx_fio_do_int(fio, i);
210 gmx_fio_do_int(fio, i);
213 gmx_fio_do_int(fio, i);
215 if (i == GMX_MTXIO_FULL_MATRIX && NULL != full_matrix)
217 printf("Full matrix storage format, nrow=%d, ncols=%d\n", *nrow, *ncol);
219 sz = (*nrow) * (*ncol);
220 snew((*full_matrix), sz);
221 bDum = gmx_fio_ndo_real(fio, (*full_matrix), sz);
223 else if (NULL != sparse_matrix)
226 printf("Sparse matrix storage format, nrow=%d, ncols=%d\n", *nrow, *ncol);
228 snew((*sparse_matrix), 1);
229 gmx_fio_do_gmx_bool(fio, (*sparse_matrix)->compressed_symmetric);
230 gmx_fio_do_int(fio, (*sparse_matrix)->nrow);
231 if ((*sparse_matrix)->nrow != *nrow)
233 gmx_fatal(FARGS, "Internal inconsistency in sparse matrix.\n");
235 snew((*sparse_matrix)->ndata, (*sparse_matrix)->nrow);
236 snew((*sparse_matrix)->nalloc, (*sparse_matrix)->nrow);
237 snew((*sparse_matrix)->data, (*sparse_matrix)->nrow);
238 bDum = gmx_fio_ndo_int(fio, (*sparse_matrix)->ndata,
239 (*sparse_matrix)->nrow);
241 for (i = 0; i < (*sparse_matrix)->nrow; i++)
243 (*sparse_matrix)->nalloc[i] = (*sparse_matrix)->ndata[i] + 10;
244 snew(((*sparse_matrix)->data[i]), (*sparse_matrix)->nalloc[i]);
246 for (j = 0; j < (*sparse_matrix)->ndata[i]; j++)
248 gmx_fio_do_int(fio, (*sparse_matrix)->data[i][j].col);
249 gmx_fio_do_real(fio, (*sparse_matrix)->data[i][j].value);