Replaced all occurences of str(n)casecmp with gmx_str(n)casecmp.
[gromacs/rigid-bodies.git] / src / gmxlib / mtxio.c
blob9faacbc34bfe612c7d3cd8815ad628a27c97fbc1
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Green Red Orange Magenta Azure Cyan Skyblue
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 /* This module provides routines to read/write sparse or full storage
40 * matrices from/to files. It is normally used for the Hessian matrix
41 * in normal mode analysis.
44 #include "xdrf.h"
45 #include "smalloc.h"
46 #include "gmxfio.h"
47 #include "copyrite.h"
48 #include "gmx_fatal.h"
49 #include "mtxio.h"
50 #include "gmxfio.h"
53 /* Just a number to identify our file type */
54 #define GMX_MTXIO_MAGIC_NUMBER 0x34ce8fd2
56 #define GMX_MTXIO_FULL_MATRIX 0
57 #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
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.
82 void gmx_mtxio_write(const char * filename,
83 int nrow,
84 int ncol,
85 real * full_matrix,
86 gmx_sparsematrix_t * sparse_matrix)
88 t_fileio *fio;
89 XDR * xd;
90 int i,j,prec;
91 bool bDum = TRUE;
92 bool bRead = FALSE;
93 size_t sz;
95 if(full_matrix!=NULL && sparse_matrix!=NULL)
97 gmx_fatal(FARGS,"Both full AND sparse matrix specified to gmx_mtxio_write().\n");
100 fio = gmx_fio_open(filename,"w");
101 gmx_fio_checktype(fio);
102 xd = gmx_fio_getxdr(fio);
104 /* Write magic number */
105 i = GMX_MTXIO_MAGIC_NUMBER;
106 gmx_fio_do_int(fio, i);
108 /* Write generating Gromacs version */
109 gmx_fio_write_string(fio, GromacsVersion());
111 /* Write 1 for double, 0 for single precision */
112 if(sizeof(real)==sizeof(double))
113 prec = 1;
114 else
115 prec = 0;
116 gmx_fio_do_int(fio, prec);
118 gmx_fio_do_int(fio, nrow);
119 gmx_fio_do_int(fio, ncol);
121 if(full_matrix!=NULL)
123 /* Full matrix storage format */
124 i = GMX_MTXIO_FULL_MATRIX;
125 gmx_fio_do_int(fio, i);
126 sz = nrow*ncol;
127 bDum=gmx_fio_ndo_real(fio, full_matrix,sz);
129 else
131 /* Sparse storage */
132 i = GMX_MTXIO_SPARSE_MATRIX;
133 gmx_fio_do_int(fio, i);
135 gmx_fio_do_bool(fio, sparse_matrix->compressed_symmetric);
136 gmx_fio_do_int(fio, sparse_matrix->nrow);
137 if(sparse_matrix->nrow != nrow)
139 gmx_fatal(FARGS,"Internal inconsistency in sparse matrix.\n");
141 bDum=gmx_fio_ndo_int(fio, sparse_matrix->ndata,sparse_matrix->nrow);
142 for(i=0;i<sparse_matrix->nrow;i++)
144 for(j=0;j<sparse_matrix->ndata[i];j++)
146 gmx_fio_do_int(fio, sparse_matrix->data[i][j].col);
147 gmx_fio_do_real(fio, sparse_matrix->data[i][j].value);
151 gmx_fio_close(fio);
155 void
156 gmx_mtxio_read (const char * filename,
157 int * nrow,
158 int * ncol,
159 real ** full_matrix,
160 gmx_sparsematrix_t ** sparse_matrix)
162 t_fileio *fio;
163 XDR * xd;
164 int i,j,prec;
165 bool bDum = TRUE;
166 bool bRead = TRUE;
167 char gmxver[256];
168 size_t sz;
170 fio = gmx_fio_open(filename,"r");
171 gmx_fio_checktype(fio);
172 xd = gmx_fio_getxdr(fio);
174 /* Read and check magic number */
175 i = GMX_MTXIO_MAGIC_NUMBER;
176 gmx_fio_do_int(fio, i);
178 if(i!=GMX_MTXIO_MAGIC_NUMBER)
180 gmx_fatal(FARGS,
181 "No matrix data found in file. Note that the Hessian matrix format changed\n"
182 "in Gromacs 3.3 to enable portable files and sparse matrix storage.\n");
185 /* Read generating Gromacs version */
186 gmx_fio_do_string(fio, gmxver);
188 /* Write 1 for double, 0 for single precision */
189 if(sizeof(real)==sizeof(double))
190 prec = 1;
191 else
192 prec = 0;
193 gmx_fio_do_int(fio, prec);
195 fprintf(stderr,"Reading %s precision matrix generated by Gromacs %s\n",
196 (prec == 1) ? "double" : "single",gmxver);
198 gmx_fio_do_int(fio, i);
199 *nrow=i;
200 gmx_fio_do_int(fio, i);
201 *ncol=i;
203 gmx_fio_do_int(fio, i);
205 if(i==GMX_MTXIO_FULL_MATRIX)
207 printf("Full matrix storage format, nrow=%d, ncols=%d\n",*nrow,*ncol);
209 sz = (*nrow) * (*ncol);
210 snew((*full_matrix),sz);
211 bDum=gmx_fio_ndo_real(fio, (*full_matrix),sz);
213 else
215 /* Sparse storage */
216 printf("Sparse matrix storage format, nrow=%d, ncols=%d\n",*nrow,*ncol);
218 snew((*sparse_matrix),1);
219 gmx_fio_do_bool(fio, (*sparse_matrix)->compressed_symmetric);
220 gmx_fio_do_int(fio, (*sparse_matrix)->nrow);
221 if((*sparse_matrix)->nrow != *nrow)
223 gmx_fatal(FARGS,"Internal inconsistency in sparse matrix.\n");
225 snew((*sparse_matrix)->ndata,(*sparse_matrix)->nrow);
226 snew((*sparse_matrix)->nalloc,(*sparse_matrix)->nrow);
227 snew((*sparse_matrix)->data,(*sparse_matrix)->nrow);
228 bDum=gmx_fio_ndo_int(fio, (*sparse_matrix)->ndata,
229 (*sparse_matrix)->nrow);
231 for(i=0;i<(*sparse_matrix)->nrow;i++)
233 (*sparse_matrix)->nalloc[i] = (*sparse_matrix)->ndata[i] + 10;
234 snew(((*sparse_matrix)->data[i]),(*sparse_matrix)->nalloc[i]);
236 for(j=0;j<(*sparse_matrix)->ndata[i];j++)
238 gmx_fio_do_int(fio, (*sparse_matrix)->data[i][j].col);
239 gmx_fio_do_real(fio, (*sparse_matrix)->data[i][j].value);
243 gmx_fio_close(fio);