Merge branch 'release-4-5-patches' of git@git.gromacs.org:gromacs into release-4...
[gromacs.git] / include / statutil.h
blobd66f4fcc85fba7fc8d6c0622669dc5c9113eab22
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 * Gromacs Runs On Most of All Computer Systems
36 #ifndef _statutil_h
37 #define _statutil_h
39 #include <stdio.h>
40 #include "typedefs.h"
41 #include "filenm.h"
42 #include "readinp.h"
43 #include "wman.h"
44 #include "pdbio.h"
45 #include "oenv.h"
46 #include "gmxfio.h"
49 #ifdef __cplusplus
50 extern "C" {
51 #endif
52 #if 0 /* avoid screwing up indentation */
54 #endif
57 /* The code below is to facilitate controlled begin and end of
58 trajectory reading. Corresponding routines in
59 src/gmxlib/tcontrol.c
61 enum { TBEGIN, TEND, TDELTA, TNR };
63 gmx_bool bTimeSet(int tcontrol);
65 real rTimeValue(int tcontrol);
67 void setTimeValue(int tcontrol,real value);
69 /* End trajectory time control */
71 /* a dedicated status type contains fp, etc. */
72 typedef struct t_trxstatus t_trxstatus;
74 typedef int t_first_x(t_trxstatus **status,const char *fn,real *t,rvec **x,
75 matrix box);
77 typedef gmx_bool t_next_x(t_trxstatus *status,real *t,int natoms,rvec x[],
78 matrix box);
80 /* I/O function types */
83 /* LEGACY FUNCTIONS
85 The program names, command lines, etc. are now also set in the output_env
86 structure. That is now the preferred location, but the functions here
87 are still available as legacy functions. Because they all act on inherently
88 global informaion, their existence in a multi-threaded environment is not
89 a real problem. */
91 /* Return the name of the program */
92 const char *command_line(void);
93 void set_command_line(int argc, char *argv[]);
95 /* set the program name to the provided string, but note
96 * that it must be a real file - we determine the library
97 * directory from its location!
98 */
99 const char *Program(void);
100 void set_program_name(const char *argvzero);
101 /* Id. without leading directory */
102 const char *ShortProgram(void);
104 /************************************************
105 * Trajectory functions
106 ************************************************/
108 int prec2ndec(real prec);
109 /* Convert precision in 1/(nm) to number of decimal places */
111 void clear_trxframe(t_trxframe *fr,gmx_bool bFirst);
112 /* Set all content gmx_booleans to FALSE.
113 * When bFirst = TRUE, set natoms=-1, all pointers to NULL
114 * and all data to zero.
117 void set_trxframe_ePBC(t_trxframe *fr,int ePBC);
118 /* Set the type of periodic boundary conditions, ePBC=-1 is not set */
120 int nframes_read(t_trxstatus *status);
121 /* Returns the number of frames read from the trajectory */
123 int write_trxframe_indexed(t_trxstatus *status,t_trxframe *fr,int nind,
124 atom_id *ind, gmx_conect gc);
125 /* Write an indexed frame to a TRX file, see write_trxframe. gc may be NULL */
127 int write_trxframe(t_trxstatus *status,t_trxframe *fr,gmx_conect gc);
128 /* Write a frame to a TRX file.
129 * Only entries for which the gmx_boolean is TRUE will be written,
130 * except for step, time, lambda and/or box, which may not be
131 * omitted for certain trajectory formats.
132 * The precision for .xtc and .gro is fr->prec, when fr->bPrec=FALSE,
133 * the precision is set to 1000.
134 * gc is important for pdb file writing only and may be NULL.
137 int write_trx(t_trxstatus *status,int nind,atom_id *ind,t_atoms *atoms,
138 int step,real time,matrix box,rvec x[],rvec *v,
139 gmx_conect gc);
140 /* Write an indexed frame to a TRX file.
141 * v can be NULL.
142 * atoms can be NULL for file types which don't need atom names.
145 void close_trx(t_trxstatus *status);
146 /* Close trj file as opened with read_first_x, read_frist_frame
147 * or open_trx. Identical to close_trj.
150 t_trxstatus *open_trx(const char *outfile,const char *filemode);
151 /* Open a TRX file and return an allocated status pointer */
153 /* get a fileio from a trxstatus */
154 t_fileio *trx_get_fileio(t_trxstatus *status);
157 gmx_bool bRmod_fd(double a, double b, double c,gmx_bool bDouble);
158 /* Returns TRUE when (a - b) MOD c = 0, using a margin which is slightly
159 * larger than the float/double precision.
162 #ifdef GMX_DOUBLE
163 #define bRmod(a,b,c) bRmod_fd(a,b,c,TRUE)
164 #else
165 #define bRmod(a,b,c) bRmod_fd(a,b,c,FALSE)
166 #endif
168 int check_times2(real t,real t0,real tp,real tpp,gmx_bool bDouble);
169 /* This routine checkes if the read-in time is correct or not;
170 * returns -1 if t<tbegin or t MOD dt = t0,
171 * 0 if tbegin <= t <=tend+margin,
172 * 1 if t>tend
173 * where margin is 0.1*min(t-tp,tp-tpp), if this positive, 0 otherwise.
174 * tp and tpp should be the time of the previous frame and the one before.
175 * The mod is done with single or double precision accuracy depending
176 * on the value of bDouble.
179 int check_times(real t);
180 /* This routine checkes if the read-in time is correct or not;
181 * returns -1 if t<tbegin,
182 * 0 if tbegin <= t <=tend,
183 * 1 if t>tend
190 /* For trxframe.flags, used in trxframe read routines.
191 * When a READ flag is set, the field will be read when present,
192 * but a frame might be returned which does not contain the field.
193 * When a NEED flag is set, frames not containing the field will be skipped.
195 #define TRX_READ_X (1<<0)
196 #define TRX_NEED_X (1<<1)
197 #define TRX_READ_V (1<<2)
198 #define TRX_NEED_V (1<<3)
199 #define TRX_READ_F (1<<4)
200 #define TRX_NEED_F (1<<5)
201 /* Useful for reading natoms from a trajectory without skipping */
202 #define TRX_DONT_SKIP (1<<6)
204 /* For trxframe.not_ok */
205 #define HEADER_NOT_OK (1<<0)
206 #define DATA_NOT_OK (1<<1)
207 #define FRAME_NOT_OK (HEADER_NOT_OK | DATA_NOT_OK)
209 int read_first_frame(const output_env_t oenv,t_trxstatus **status,
210 const char *fn, t_trxframe *fr,int flags);
211 /* Read the first frame which is in accordance with flags, which are
212 * defined further up in this file.
213 * Returns natoms when succeeded, 0 otherwise.
214 * Memory will be allocated for flagged entries.
215 * The flags are copied to fr for subsequent calls to read_next_frame.
216 * Returns TRUE when succeeded, FALSE otherwise.
219 gmx_bool read_next_frame(const output_env_t oenv,t_trxstatus *status,
220 t_trxframe *fr);
221 /* Reads the next frame which is in accordance with fr->flags.
222 * Returns TRUE when succeeded, FALSE otherwise.
225 int read_first_x(const output_env_t oenv,t_trxstatus **status,
226 const char *fn, real *t,rvec **x,matrix box);
227 /* These routines read first coordinates and box, and allocates
228 * memory for the coordinates, for a trajectory file.
229 * The routine returns the number of atoms, or 0 when something is wrong.
230 * The integer in status should be passed to calls of read_next_x
233 gmx_bool read_next_x(const output_env_t oenv,t_trxstatus *status,real *t,
234 int natoms, rvec x[],matrix box);
235 /* Read coordinates and box from a trajectory file. Return TRUE when all well,
236 * or FALSE when end of file (or last frame requested by user).
237 * status is the integer set in read_first_x.
240 void close_trj(t_trxstatus *status);
241 /* Close trj file as opened with read_first_x, read_frist_frame
242 * or open_trx. Identical to close_trx.
245 void rewind_trj(t_trxstatus *status);
246 /* Rewind trj file as opened with read_first_x */
248 t_topology *read_top(const char *fn,int *ePBC);
249 /* Extract a topology data structure from a topology file.
250 * If ePBC!=NULL *ePBC gives the pbc type.
253 /*****************************************************
254 * Some command line parsing routines
255 *****************************************************/
257 #define PCA_CAN_VIEW (1<<5)
258 /* add option -w to view output files (must be implemented in program) */
259 #define PCA_CAN_BEGIN (1<<6)
260 #define PCA_CAN_END (1<<7)
261 #define PCA_CAN_DT (1<<14)
262 #define PCA_CAN_TIME (PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_DT)
263 /* adds options -b and -e for begin and end time for reading trajectories */
264 #define PCA_TIME_UNIT (1<<15)
265 /* set time unit for output */
266 #define PCA_KEEP_ARGS (1<<8)
267 /* keep parsed args in argv (doesn't make sense without NOEXIT_ON_ARGS) */
268 #define PCA_SILENT (1<<9)
269 /* don't print options by default */
270 #define PCA_CAN_SET_DEFFNM (1<<10)
271 /* does something for non-master mdrun nodes */
272 #define PCA_NOEXIT_ON_ARGS (1<<11)
273 /* no fatal_error when invalid options are encountered */
274 #define PCA_QUIET (1<<12)
275 /* does something for non-master mdrun nodes */
276 #define PCA_BE_NICE (1<<13)
277 /* Default to low priority, unless configured with --disable-nice */
278 #define PCA_NOT_READ_NODE (1<<16)
279 /* Is this node not reading: for parallel all nodes but the master */
281 int iscan(int argc,char *argv[],int *i);
282 /* Scan an int from the argument at *i. If the argument length
283 * is > 2, the int is assumed to be in the remainder of the arg,
284 * eg: -p32, else the int is assumed to be in the next argument
285 * eg: -p 32. If neither is the case the routine exits with an error,
286 * otherwise it returns the value found. If the value is in the next
287 * argument *i is incremented. You typically would want to pass
288 * a loop variable to this routine.
290 gmx_large_int_t istepscan(int argc,char *argv[],int *i);
291 /* Same as above, but for large integer values */
293 double dscan(int argc,char *argv[],int *i);
294 /* Routine similar to the above, but working on doubles. */
296 char *sscan(int argc,char *argv[],int *i);
297 /* Routine similar to the above, but working on strings. The pointer
298 * returned is a pointer to the argv field.
301 void vscan(int argc,char *argv[],int *i,rvec *vec);
302 /* Routine similar to the above, but working on rvecs. */
304 int nenum(const char *const enumc[]);
305 /* returns ordinal number of selected enum from args
306 * depends on enumc[0] pointing to one of the other elements
307 * array must be terminated by a NULL pointer
310 void parse_common_args(int *argc,char *argv[],unsigned long Flags,
311 int nfile,t_filenm fnm[],int npargs,t_pargs *pa,
312 int ndesc,const char **desc,
313 int nbugs,const char **bugs,
314 output_env_t *oenv);
315 /* Get arguments from the arg-list. The arguments extracted
316 * are removed from the list. If manual is NULL a default message is displayed
317 * when errors are encountered. The Flags argument, when non-0 enables
318 * some input checks. Using this routine also means that the arguments
319 * -b and -e will be used for begin and end time, whether this is
320 * appropriate or not!
323 #ifdef __cplusplus
325 #endif
327 #endif