Another quick bug fix.
[gromacs/rigid-bodies.git] / include / statutil.h
blobf24aed8c6ed6be191bfb2ae1f27a6c89996e55ee
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 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
46 #if 0 /* avoid screwing up indentation */
48 #endif
50 #include <stdio.h>
51 #include "typedefs.h"
52 #include "filenm.h"
53 #include "readinp.h"
54 #include "wman.h"
55 #include "pdbio.h"
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 extern bool bTimeSet(int tcontrol);
65 extern real rTimeValue(int tcontrol);
67 extern void setTimeValue(int tcontrol,real value);
69 /* End trajectory time control */
71 typedef int t_first_x(int *status,const char *fn,real *t,rvec **x,matrix box);
73 typedef bool t_next_x(int status,real *t,int natoms,rvec x[],matrix box);
75 /* I/O function types */
78 /* This is for re-entrant versions of the functions below. The env_info
79 structure holds information about program name, cmd line, default times,
80 etc. Right now, there is a global variable of this type in statutil.c to
81 enable the old functions, but this should change. All the functions below
82 now have re-entrant versions listed with them. */
83 extern void set_output_env(output_env_t oenv, bool view, bool xvgr_codes,
84 const char *timenm);
85 extern void init_output_env(output_env_t oenv, int argc, char *argv[],
86 bool view, bool xvgr_codes, const char *timenm);
90 /* Return the name of the program */
91 extern const char *command_line(void);
92 extern void set_command_line(int argc, char *argv[]);
94 /* set the program name to the provided string, but note
95 * that it must be a real file - we determine the library
96 * directory from its location!
97 */
98 extern const char *Program(void);
99 extern void set_program_name(const char *argvzero);
100 /* Id. without leading directory */
101 extern const char *ShortProgram(void);
103 /************************************************
104 * Trajectory functions
105 ************************************************/
107 extern int prec2ndec(real prec);
108 /* Convert precision in 1/(nm) to number of decimal places */
110 extern void clear_trxframe(t_trxframe *fr,bool bFirst);
111 /* Set all content booleans to FALSE.
112 * When bFirst = TRUE, set natoms=-1, all pointers to NULL
113 * and all data to zero.
116 extern void set_trxframe_ePBC(t_trxframe *fr,int ePBC);
117 /* Set the type of periodic boundary conditions, ePBC=-1 is not set */
119 extern int nframes_read(void);
120 /* Returns the number of frames read from the trajectory */
122 int write_trxframe_indexed(int status,t_trxframe *fr,int nind,atom_id *ind,
123 gmx_conect gc);
124 /* Write an indexed frame to a TRX file, see write_trxframe. gc may be NULL */
126 int write_trxframe(int status,t_trxframe *fr,gmx_conect gc);
127 /* Write a frame to a TRX file.
128 * Only entries for which the boolean is TRUE will be written,
129 * except for step, time, lambda and/or box, which may not be
130 * omitted for certain trajectory formats.
131 * The precision for .xtc and .gro is fr->prec, when fr->bPrec=FALSE,
132 * the precision is set to 1000.
133 * gc is important for pdb file writing only and may be NULL.
136 int write_trx(int status,int nind,atom_id *ind,t_atoms *atoms,
137 int step,real time,matrix box,rvec x[],rvec *v,
138 gmx_conect gc);
139 /* Write an indexed frame to a TRX file.
140 * v can be NULL.
141 * atoms can be NULL for file types which don't need atom names.
144 void close_trx(int status);
145 /* Close trj file as opened with read_first_x, read_frist_frame
146 * or open_trx. Identical to close_trj.
149 int open_trx(const char *outfile,const char *filemode);
150 /* Open a TRX file and return the file number */
152 extern bool bRmod_fd(double a, double b, double c,bool bDouble);
153 /* Returns TRUE when (a - b) MOD c = 0, using a margin which is slightly
154 * larger than the float/double precision.
157 #ifdef GMX_DOUBLE
158 #define bRmod(a,b,c) bRmod_fd(a,b,c,TRUE)
159 #else
160 #define bRmod(a,b,c) bRmod_fd(a,b,c,FALSE)
161 #endif
163 extern int check_times2(real t,real t0,real tp,real tpp,bool bDouble);
164 /* This routine checkes if the read-in time is correct or not;
165 * returns -1 if t<tbegin or t MOD dt = t0,
166 * 0 if tbegin <= t <=tend+margin,
167 * 1 if t>tend
168 * where margin is 0.1*min(t-tp,tp-tpp), if this positive, 0 otherwise.
169 * tp and tpp should be the time of the previous frame and the one before.
170 * The mod is done with single or double precision accuracy depending
171 * on the value of bDouble.
174 extern int check_times(real t);
175 /* This routine checkes if the read-in time is correct or not;
176 * returns -1 if t<tbegin,
177 * 0 if tbegin <= t <=tend,
178 * 1 if t>tend
181 extern const char *get_time_unit(const output_env_t oenv);
182 /* return time unit (e.g. ps or ns) */
184 extern const char *get_time_label(const output_env_t oenv);
185 /* return time unit label (e.g. "Time (ps)") */
187 extern const char *get_xvgr_tlabel(const output_env_t oenv);
188 /* retrun x-axis time label for xmgr */
190 extern real get_time_factor(const output_env_t oenv);
191 /* return time conversion factor from ps (i.e. 1e-3 for ps->ns) */
193 extern real get_time_invfactor(const output_env_t oenv);
194 /* return inverse time conversion factor from ps (i.e. 1e3 for ps->ns) */
196 extern real conv_time(const output_env_t oenv, real time);
197 /* return converted time */
199 extern void conv_times(const output_env_t oenv, int n, real *time);
200 /* convert array of times */
202 extern bool get_view(const output_env_t oenv);
203 /* Return TRUE when user requested viewing of the file */
205 extern bool get_print_xvgr_codes(const output_env_t oenv);
206 /* Return TRUE when user wants printing of legends etc. in the file. */
211 /* For trxframe.flags, used in trxframe read routines.
212 * When a READ flag is set, the field will be read when present,
213 * but a frame might be returned which does not contain the field.
214 * When a NEED flag is set, frames not containing the field will be skipped.
216 #define TRX_READ_X (1<<0)
217 #define TRX_NEED_X (1<<1)
218 #define TRX_READ_V (1<<2)
219 #define TRX_NEED_V (1<<3)
220 #define TRX_READ_F (1<<4)
221 #define TRX_NEED_F (1<<5)
222 /* Useful for reading natoms from a trajectory without skipping */
223 #define TRX_DONT_SKIP (1<<6)
225 /* For trxframe.not_ok */
226 #define HEADER_NOT_OK (1<<0)
227 #define DATA_NOT_OK (1<<1)
228 #define FRAME_NOT_OK (HEADER_NOT_OK | DATA_NOT_OK)
230 extern int read_first_frame(const output_env_t oenv,int *status,const char *fn,
231 t_trxframe *fr,int flags);
232 /* Read the first frame which is in accordance with flags, which are
233 * defined further up in this file.
234 * Returns natoms when succeeded, 0 otherwise.
235 * Memory will be allocated for flagged entries.
236 * The flags are copied to fr for subsequent calls to read_next_frame.
237 * Returns TRUE when succeeded, FALSE otherwise.
240 extern bool read_next_frame(const output_env_t oenv,int status,t_trxframe *fr);
241 /* Reads the next frame which is in accordance with fr->flags.
242 * Returns TRUE when succeeded, FALSE otherwise.
245 extern int read_first_x(const output_env_t oenv,int *status,const char *fn,
246 real *t,rvec **x,matrix box);
247 /* These routines read first coordinates and box, and allocates
248 * memory for the coordinates, for a trajectory file.
249 * The routine returns the number of atoms, or 0 when something is wrong.
250 * The integer in status should be passed to calls of read_next_x
253 extern bool read_next_x(const output_env_t oenv,int status,real *t,int natoms,
254 rvec x[],matrix box);
255 /* Read coordinates and box from a trajectory file. Return TRUE when all well,
256 * or FALSE when end of file (or last frame requested by user).
257 * status is the integer set in read_first_x.
260 extern void close_trj(int status);
261 /* Close trj file as opened with read_first_x, read_frist_frame
262 * or open_trx. Identical to close_trx.
265 extern void rewind_trj(int status);
266 /* Rewind trj file as opened with read_first_x */
268 extern t_topology *read_top(const char *fn,int *ePBC);
269 /* Extract a topology data structure from a topology file.
270 * If ePBC!=NULL *ePBC gives the pbc type.
273 /*****************************************************
274 * Some command line parsing routines
275 *****************************************************/
277 #define PCA_CAN_VIEW (1<<5)
278 /* add option -w to view output files (must be implemented in program) */
279 #define PCA_CAN_BEGIN (1<<6)
280 #define PCA_CAN_END (1<<7)
281 #define PCA_CAN_DT (1<<14)
282 #define PCA_CAN_TIME (PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_DT)
283 /* adds options -b and -e for begin and end time for reading trajectories */
284 #define PCA_TIME_UNIT (1<<15)
285 /* set time unit for output */
286 #define PCA_KEEP_ARGS (1<<8)
287 /* keep parsed args in argv (doesn't make sense without NOEXIT_ON_ARGS) */
288 #define PCA_SILENT (1<<9)
289 /* don't print options by default */
290 #define PCA_CAN_SET_DEFFNM (1<<10)
291 /* does something for non-master mdrun nodes */
292 #define PCA_NOEXIT_ON_ARGS (1<<11)
293 /* no fatal_error when invalid options are encountered */
294 #define PCA_QUIET (1<<12)
295 /* does something for non-master mdrun nodes */
296 #define PCA_BE_NICE (1<<13)
297 /* Default to low priority, unless configured with --disable-nice */
298 #define PCA_NOT_READ_NODE (1<<16)
299 /* Is this node not reading: for parallel all nodes but the master */
301 extern int iscan(int argc,char *argv[],int *i);
302 /* Scan an int from the argument at *i. If the argument length
303 * is > 2, the int is assumed to be in the remainder of the arg,
304 * eg: -p32, else the int is assumed to be in the next argument
305 * eg: -p 32. If neither is the case the routine exits with an error,
306 * otherwise it returns the value found. If the value is in the next
307 * argument *i is incremented. You typically would want to pass
308 * a loop variable to this routine.
310 extern gmx_large_int_t istepscan(int argc,char *argv[],int *i);
311 /* Same as above, but for large integer values */
313 extern double dscan(int argc,char *argv[],int *i);
314 /* Routine similar to the above, but working on doubles. */
316 extern char *sscan(int argc,char *argv[],int *i);
317 /* Routine similar to the above, but working on strings. The pointer
318 * returned is a pointer to the argv field.
321 extern void vscan(int argc,char *argv[],int *i,rvec *vec);
322 /* Routine similar to the above, but working on rvecs. */
324 extern int nenum(const char *const enumc[]);
325 /* returns ordinal number of selected enum from args
326 * depends on enumc[0] pointing to one of the other elements
327 * array must be terminated by a NULL pointer
330 #ifdef HAVE_MOTIF
331 extern void gmx_gui(int *argc,char *argv[],
332 int nfile,t_filenm fnm[],int npargs,t_pargs pa[],
333 int ndesc,const char *desc[],
334 int nbugs,const char *bugs[]);
335 /* This function plops up a Motif dialog box in which the command-line options
336 * can be changed.
338 #endif
340 extern void parse_common_args(int *argc,char *argv[],unsigned long Flags,
341 int nfile,t_filenm fnm[],int npargs,t_pargs *pa,
342 int ndesc,const char **desc,
343 int nbugs,const char **bugs,
344 output_env_t *oenv);
345 /* Get arguments from the arg-list. The arguments extracted
346 * are removed from the list. If manual is NULL a default message is displayed
347 * when errors are encountered. The Flags argument, when non-0 enables
348 * some input checks. Using this routine also means that the arguments
349 * -b and -e will be used for begin and end time, whether this is
350 * appropriate or not!
353 #ifdef __cplusplus
355 #endif
357 #endif