changed reading hint
[gromacs/adressmacs.git] / src / gmxlib / xtcio.c
blobd526d42e944f3518e7091f8b2d49d4fb47ab0fb0
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 2.0
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
16 * Please refer to:
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
23 * or e-mail to:
24 * gromacs@chem.rug.nl
26 * And Hey:
27 * Green Red Orange Magenta Azure Cyan Skyblue
29 static char *SRCID_xtcio_c = "$Id$";
31 #include <string.h>
32 #include "typedefs.h"
33 #include "xdrf.h"
34 #include "gmxfio.h"
35 #include "xtcio.h"
36 #include "smalloc.h"
37 #include "vec.h"
38 #include "futil.h"
39 #include "fatal.h"
41 #define XTC_MAGIC 1995
43 int open_xtc(char *fn,char *mode)
45 return fio_open(fn,mode);
48 void close_xtc(int fp)
50 fio_close(fp);
53 static void check_xtc_magic(int magic)
55 if (magic != XTC_MAGIC)
56 fatal_error(0,"Magic Number Error in XTC file (read %d, should be %d)",
57 magic,XTC_MAGIC);
60 int xtc_check(char *str,bool bResult,char *file,int line)
62 if (!bResult) {
63 if (debug)
64 fprintf(debug,"\nXTC error: read/write of %s failed, "
65 "source file %s, line %d\n",str,file,line);
66 return 0;
68 return 1;
71 void xtc_check_fat_err(char *str,bool bResult,char *file,int line)
73 if (!bResult) {
74 fatal_error(0,"XTC read/write of %s failed, "
75 "source file %s, line %d\n",str,file,line);
79 static int xtc_header(XDR *xd,int *magic,int *natoms,int *step,real *time,
80 bool *bOK)
82 int result;
84 if (xdr_int(xd,magic) == 0)
85 return 0;
86 result=XTC_CHECK("natoms", xdr_int(xd,natoms)); /* number of atoms */
87 if (result)
88 result=XTC_CHECK("step", xdr_int(xd,step)); /* frame number */
89 if (result)
90 result=XTC_CHECK("time", xdr_real(xd,time)); /* time */
91 *bOK=(result!=0);
93 return result;
96 static int xtc_coord(XDR *xd,int *natoms,matrix box,rvec *x,real *prec)
98 int i,j,result;
100 /* box */
101 result=1;
102 for(i=0; ((i<DIM) && result); i++)
103 for(j=0; ((j<DIM) && result); j++)
104 result=XTC_CHECK("box",xdr_real(xd,&(box[i][j])));
106 if (result)
107 /* coordinates */
108 result=XTC_CHECK("x",xdr3drcoord(xd,x[0],natoms,prec));
110 return result;
113 static int xtc_io(XDR *xd,int *magic,
114 int *natoms,int *step,real *time,
115 matrix box,rvec *x,real *prec,bool *bOK)
117 if (!xtc_header(xd,magic,natoms,step,time,bOK))
118 return 0;
119 return xtc_coord(xd,natoms,box,x,prec);
122 int write_xtc(int fp,
123 int natoms,int step,real time,
124 matrix box,rvec *x,real prec)
126 int magic_number = XTC_MAGIC;
127 XDR *xd;
128 bool bDum;
130 xd = fio_getxdr(fp);
131 /* write magic number and xtc identidier */
132 if (!xtc_header(xd,&magic_number,&natoms,&step,&time,&bDum))
133 return 0;
135 /* write data */
136 return xtc_coord(xd,&natoms,box,x,&prec);
139 int read_first_xtc(int fp,int *natoms,int *step,real *time,
140 matrix box,rvec **x,real *prec,bool *bOK)
142 int magic;
143 XDR *xd;
145 *bOK=TRUE;
146 xd = fio_getxdr(fp);
148 /* read header and malloc x */
149 if ( !xtc_header(xd,&magic,natoms,step,time,bOK))
150 return 0;
152 /* Check magic number */
153 check_xtc_magic(magic);
155 snew(*x,*natoms);
157 *bOK=xtc_coord(xd,natoms,box,*x,prec);
159 return *bOK;
162 int read_next_xtc(int fp,
163 int *natoms,int *step,real *time,
164 matrix box,rvec *x,real *prec,bool *bOK)
166 int magic;
167 XDR *xd;
169 *bOK=TRUE;
170 xd = fio_getxdr(fp);
172 /* read header */
173 if (!xtc_header(xd,&magic,natoms,step,time,bOK))
174 return 0;
176 /* Check magic number */
177 check_xtc_magic(magic);
179 *bOK=xtc_coord(xd,natoms,box,x,prec);
181 return *bOK;