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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * 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.
51 #include "gmx_fatal.h"
56 /* Just a number to identify our file type */
57 #define GMX_MTXIO_MAGIC_NUMBER 0x34ce8fd2
59 #define GMX_MTXIO_FULL_MATRIX 0
60 #define GMX_MTXIO_SPARSE_MATRIX 1
64 /* Matrix file format definition:
66 * All entries are stored in XDR format.
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
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.
85 void gmx_mtxio_write(const char * filename,
89 gmx_sparsematrix_t * sparse_matrix)
95 gmx_bool bRead = FALSE;
98 if(full_matrix!=NULL && sparse_matrix!=NULL)
100 gmx_fatal(FARGS,"Both full AND sparse matrix specified to gmx_mtxio_write().\n");
103 fio = gmx_fio_open(filename,"w");
104 gmx_fio_checktype(fio);
105 xd = gmx_fio_getxdr(fio);
107 /* Write magic number */
108 i = GMX_MTXIO_MAGIC_NUMBER;
109 gmx_fio_do_int(fio, i);
111 /* Write generating Gromacs version */
112 gmx_fio_write_string(fio, GromacsVersion());
114 /* Write 1 for double, 0 for single precision */
115 if(sizeof(real)==sizeof(double))
119 gmx_fio_do_int(fio, prec);
121 gmx_fio_do_int(fio, nrow);
122 gmx_fio_do_int(fio, ncol);
124 if(full_matrix!=NULL)
126 /* Full matrix storage format */
127 i = GMX_MTXIO_FULL_MATRIX;
128 gmx_fio_do_int(fio, i);
130 bDum=gmx_fio_ndo_real(fio, full_matrix,sz);
135 i = GMX_MTXIO_SPARSE_MATRIX;
136 gmx_fio_do_int(fio, i);
138 gmx_fio_do_gmx_bool(fio, sparse_matrix->compressed_symmetric);
139 gmx_fio_do_int(fio, sparse_matrix->nrow);
140 if(sparse_matrix->nrow != nrow)
142 gmx_fatal(FARGS,"Internal inconsistency in sparse matrix.\n");
144 bDum=gmx_fio_ndo_int(fio, sparse_matrix->ndata,sparse_matrix->nrow);
145 for(i=0;i<sparse_matrix->nrow;i++)
147 for(j=0;j<sparse_matrix->ndata[i];j++)
149 gmx_fio_do_int(fio, sparse_matrix->data[i][j].col);
150 gmx_fio_do_real(fio, sparse_matrix->data[i][j].value);
159 gmx_mtxio_read (const char * filename,
163 gmx_sparsematrix_t ** sparse_matrix)
168 gmx_bool bDum = TRUE;
169 gmx_bool bRead = TRUE;
173 fio = gmx_fio_open(filename,"r");
174 gmx_fio_checktype(fio);
175 xd = gmx_fio_getxdr(fio);
177 /* Read and check magic number */
178 i = GMX_MTXIO_MAGIC_NUMBER;
179 gmx_fio_do_int(fio, i);
181 if(i!=GMX_MTXIO_MAGIC_NUMBER)
184 "No matrix data found in file. Note that the Hessian matrix format changed\n"
185 "in Gromacs 3.3 to enable portable files and sparse matrix storage.\n");
188 /* Read generating Gromacs version */
189 gmx_fio_do_string(fio, gmxver);
191 /* Write 1 for double, 0 for single precision */
192 if(sizeof(real)==sizeof(double))
196 gmx_fio_do_int(fio, prec);
198 fprintf(stderr,"Reading %s precision matrix generated by Gromacs %s\n",
199 (prec == 1) ? "double" : "single",gmxver);
201 gmx_fio_do_int(fio, i);
203 gmx_fio_do_int(fio, i);
206 gmx_fio_do_int(fio, i);
208 if(i==GMX_MTXIO_FULL_MATRIX && NULL != full_matrix)
210 printf("Full matrix storage format, nrow=%d, ncols=%d\n",*nrow,*ncol);
212 sz = (*nrow) * (*ncol);
213 snew((*full_matrix),sz);
214 bDum=gmx_fio_ndo_real(fio, (*full_matrix),sz);
216 else if (NULL != sparse_matrix)
219 printf("Sparse matrix storage format, nrow=%d, ncols=%d\n",*nrow,*ncol);
221 snew((*sparse_matrix),1);
222 gmx_fio_do_gmx_bool(fio, (*sparse_matrix)->compressed_symmetric);
223 gmx_fio_do_int(fio, (*sparse_matrix)->nrow);
224 if((*sparse_matrix)->nrow != *nrow)
226 gmx_fatal(FARGS,"Internal inconsistency in sparse matrix.\n");
228 snew((*sparse_matrix)->ndata,(*sparse_matrix)->nrow);
229 snew((*sparse_matrix)->nalloc,(*sparse_matrix)->nrow);
230 snew((*sparse_matrix)->data,(*sparse_matrix)->nrow);
231 bDum=gmx_fio_ndo_int(fio, (*sparse_matrix)->ndata,
232 (*sparse_matrix)->nrow);
234 for(i=0;i<(*sparse_matrix)->nrow;i++)
236 (*sparse_matrix)->nalloc[i] = (*sparse_matrix)->ndata[i] + 10;
237 snew(((*sparse_matrix)->data[i]),(*sparse_matrix)->nalloc[i]);
239 for(j=0;j<(*sparse_matrix)->ndata[i];j++)
241 gmx_fio_do_int(fio, (*sparse_matrix)->data[i][j].col);
242 gmx_fio_do_real(fio, (*sparse_matrix)->data[i][j].value);