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,2015,2017 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
42 /* This module provides routines to read/write sparse or full storage
43 * matrices from/to files. It is normally used for the Hessian matrix
44 * in normal mode analysis.
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/gmxfio_xdr.h"
49 #include "gromacs/linearalgebra/sparsematrix.h"
50 #include "gromacs/utility/baseversion.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
54 /* Just a number to identify our file type */
55 #define GMX_MTXIO_MAGIC_NUMBER 0x34ce8fd2
57 #define GMX_MTXIO_FULL_MATRIX 0
58 #define GMX_MTXIO_SPARSE_MATRIX 1
61 /* Matrix file format definition:
63 * All entries are stored in XDR format.
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
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.
82 void gmx_mtxio_write(const char* filename, int nrow, int ncol, real* full_matrix, gmx_sparsematrix_t* sparse_matrix)
88 if (full_matrix != nullptr && sparse_matrix != nullptr)
90 gmx_fatal(FARGS, "Both full AND sparse matrix specified to gmx_mtxio_write().\n");
93 fio = gmx_fio_open(filename, "w");
95 /* Write magic number */
96 i = GMX_MTXIO_MAGIC_NUMBER;
97 gmx_fio_do_int(fio, i);
99 /* Write generating Gromacs version */
100 gmx_fio_write_string(fio, gmx_version());
102 /* Write 1 for double, 0 for single precision */
103 if (sizeof(real) == sizeof(double))
111 gmx_fio_do_int(fio, prec);
113 gmx_fio_do_int(fio, nrow);
114 gmx_fio_do_int(fio, ncol);
116 if (full_matrix != nullptr)
118 /* Full matrix storage format */
119 i = GMX_MTXIO_FULL_MATRIX;
120 gmx_fio_do_int(fio, i);
122 gmx_fio_ndo_real(fio, full_matrix, sz);
127 i = GMX_MTXIO_SPARSE_MATRIX;
128 gmx_fio_do_int(fio, i);
130 gmx_fio_do_gmx_bool(fio, sparse_matrix->compressed_symmetric);
131 gmx_fio_do_int(fio, sparse_matrix->nrow);
132 if (sparse_matrix->nrow != nrow)
134 gmx_fatal(FARGS, "Internal inconsistency in sparse matrix.\n");
136 gmx_fio_ndo_int(fio, sparse_matrix->ndata, sparse_matrix->nrow);
137 for (i = 0; i < sparse_matrix->nrow; i++)
139 for (j = 0; j < sparse_matrix->ndata[i]; j++)
141 gmx_fio_do_int(fio, sparse_matrix->data[i][j].col);
142 gmx_fio_do_real(fio, sparse_matrix->data[i][j].value);
150 void gmx_mtxio_read(const char* filename, int* nrow, int* ncol, real** full_matrix, gmx_sparsematrix_t** sparse_matrix)
157 fio = gmx_fio_open(filename, "r");
159 /* Read and check magic number */
160 i = GMX_MTXIO_MAGIC_NUMBER;
161 gmx_fio_do_int(fio, i);
163 if (i != GMX_MTXIO_MAGIC_NUMBER)
166 "No matrix data found in file. Note that the Hessian matrix format changed\n"
167 "in GROMACS 3.3 to enable portable files and sparse matrix storage.\n");
170 /* Read generating Gromacs version */
171 gmx_fio_do_string(fio, gmxver);
173 /* Write 1 for double, 0 for single precision */
174 if (sizeof(real) == sizeof(double))
182 gmx_fio_do_int(fio, prec);
185 "Reading %s precision matrix generated by GROMACS %s\n",
186 (prec == 1) ? "double" : "single",
189 gmx_fio_do_int(fio, i);
191 gmx_fio_do_int(fio, i);
194 gmx_fio_do_int(fio, i);
196 if (i == GMX_MTXIO_FULL_MATRIX && nullptr != full_matrix)
198 printf("Full matrix storage format, nrow=%d, ncols=%d\n", *nrow, *ncol);
200 sz = (*nrow) * (*ncol);
201 snew((*full_matrix), sz);
202 gmx_fio_ndo_real(fio, (*full_matrix), sz);
204 else if (nullptr != sparse_matrix)
207 printf("Sparse matrix storage format, nrow=%d, ncols=%d\n", *nrow, *ncol);
209 snew((*sparse_matrix), 1);
210 gmx_fio_do_gmx_bool(fio, (*sparse_matrix)->compressed_symmetric);
211 gmx_fio_do_int(fio, (*sparse_matrix)->nrow);
212 if ((*sparse_matrix)->nrow != *nrow)
214 gmx_fatal(FARGS, "Internal inconsistency in sparse matrix.\n");
216 snew((*sparse_matrix)->ndata, (*sparse_matrix)->nrow);
217 snew((*sparse_matrix)->nalloc, (*sparse_matrix)->nrow);
218 snew((*sparse_matrix)->data, (*sparse_matrix)->nrow);
219 gmx_fio_ndo_int(fio, (*sparse_matrix)->ndata, (*sparse_matrix)->nrow);
221 for (i = 0; i < (*sparse_matrix)->nrow; i++)
223 (*sparse_matrix)->nalloc[i] = (*sparse_matrix)->ndata[i] + 10;
224 snew(((*sparse_matrix)->data[i]), (*sparse_matrix)->nalloc[i]);
226 for (j = 0; j < (*sparse_matrix)->ndata[i]; j++)
228 gmx_fio_do_int(fio, (*sparse_matrix)->data[i][j].col);
229 gmx_fio_do_real(fio, (*sparse_matrix)->data[i][j].value);