Partial commit of the project to remove all static variables.
[gromacs.git] / src / gmxlib / xtcio.c
blob0008cc41adc83396d230cf37f6b202dd1be37642
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.1
11 * Copyright (c) 1991-2001, University of Groningen, The Netherlands
12 * This program is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU General Public License
14 * as published by the Free Software Foundation; either version 2
15 * of the License, or (at your option) any later version.
17 * If you want to redistribute modifications, please consider that
18 * scientific software is very special. Version control is crucial -
19 * bugs must be traceable. We will be happy to consider code for
20 * inclusion in the official distribution, but derived work must not
21 * be called official GROMACS. Details are found in the README & COPYING
22 * files - if they are missing, get the official version at www.gromacs.org.
24 * To help us fund GROMACS development, we humbly ask that you cite
25 * the papers on the package - you can find them in the top README file.
27 * For more info, check our website at http://www.gromacs.org
29 * And Hey:
30 * Gnomes, ROck Monsters And Chili Sauce
33 #include <string.h>
34 #include "typedefs.h"
35 #include "xdrf.h"
36 #include "gmxfio.h"
37 #include "xtcio.h"
38 #include "smalloc.h"
39 #include "vec.h"
40 #include "futil.h"
41 #include "fatal.h"
43 #define XTC_MAGIC 1995
45 int open_xtc(char *fn,char *mode)
47 return fio_open(fn,mode);
50 void close_xtc(int fp)
52 fio_close(fp);
55 static void check_xtc_magic(int magic)
57 if (magic != XTC_MAGIC)
58 fatal_error(0,"Magic Number Error in XTC file (read %d, should be %d)",
59 magic,XTC_MAGIC);
62 int xtc_check(char *str,bool bResult,char *file,int line)
64 if (!bResult) {
65 if (debug)
66 fprintf(debug,"\nXTC error: read/write of %s failed, "
67 "source file %s, line %d\n",str,file,line);
68 return 0;
70 return 1;
73 void xtc_check_fat_err(char *str,bool bResult,char *file,int line)
75 if (!bResult) {
76 fatal_error(0,"XTC read/write of %s failed, "
77 "source file %s, line %d\n",str,file,line);
81 static int xtc_header(XDR *xd,int *magic,int *natoms,int *step,real *time,
82 bool *bOK)
84 int result;
86 if (xdr_int(xd,magic) == 0)
87 return 0;
88 result=XTC_CHECK("natoms", xdr_int(xd,natoms)); /* number of atoms */
89 if (result)
90 result=XTC_CHECK("step", xdr_int(xd,step)); /* frame number */
91 if (result)
92 result=XTC_CHECK("time", xdr_real(xd,time)); /* time */
93 *bOK=(result!=0);
95 return result;
98 static int xtc_coord(XDR *xd,int *natoms,matrix box,rvec *x,real *prec)
100 int i,j,result;
102 /* box */
103 result=1;
104 for(i=0; ((i<DIM) && result); i++)
105 for(j=0; ((j<DIM) && result); j++)
106 result=XTC_CHECK("box",xdr_real(xd,&(box[i][j])));
108 if (result)
109 /* coordinates */
110 result=XTC_CHECK("x",xdr3drcoord(xd,x[0],natoms,prec));
112 return result;
115 static int xtc_io(XDR *xd,int *magic,
116 int *natoms,int *step,real *time,
117 matrix box,rvec *x,real *prec,bool *bOK)
119 if (!xtc_header(xd,magic,natoms,step,time,bOK))
120 return 0;
121 return xtc_coord(xd,natoms,box,x,prec);
124 int write_xtc(int fp,
125 int natoms,int step,real time,
126 matrix box,rvec *x,real prec)
128 int magic_number = XTC_MAGIC;
129 XDR *xd;
130 bool bDum;
132 xd = fio_getxdr(fp);
133 /* write magic number and xtc identidier */
134 if (!xtc_header(xd,&magic_number,&natoms,&step,&time,&bDum))
135 return 0;
137 /* write data */
138 return xtc_coord(xd,&natoms,box,x,&prec);
141 int read_first_xtc(int fp,int *natoms,int *step,real *time,
142 matrix box,rvec **x,real *prec,bool *bOK)
144 int magic;
145 XDR *xd;
147 *bOK=TRUE;
148 xd = fio_getxdr(fp);
150 /* read header and malloc x */
151 if ( !xtc_header(xd,&magic,natoms,step,time,bOK))
152 return 0;
154 /* Check magic number */
155 check_xtc_magic(magic);
157 snew(*x,*natoms);
159 *bOK=xtc_coord(xd,natoms,box,*x,prec);
161 return *bOK;
164 int read_next_xtc(int fp,
165 int natoms,int *step,real *time,
166 matrix box,rvec *x,real *prec,bool *bOK)
168 int magic;
169 int n;
170 XDR *xd;
172 *bOK=TRUE;
173 xd = fio_getxdr(fp);
175 /* read header */
176 if (!xtc_header(xd,&magic,&n,step,time,bOK))
177 return 0;
178 if (n>natoms)
179 fatal_error(0, "Frame contains more atoms (%d) than expected (%d)",
180 n, natoms);
182 /* Check magic number */
183 check_xtc_magic(magic);
185 *bOK=xtc_coord(xd,&natoms,box,x,prec);
187 return *bOK;