changed reading hint
[gromacs/adressmacs.git] / include / constr.h
blob68877b2a0e0b89eceaff6af7d495f9befd69ecd6
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 * Good ROcking Metal Altar for Chronical Sinners
29 static char *SRCID_constr_h = "$Id$";
31 extern int bshakef(FILE *log, /* Log file */
32 int natoms, /* Total number of atoms */
33 real invmass[], /* Atomic masses */
34 int nblocks, /* The number of shake blocks */
35 int sblock[], /* The shake blocks */
36 t_idef *idef, /* The interaction def */
37 t_inputrec *ir, /* Input record */
38 matrix box, /* The box */
39 rvec x_s[], /* Coords before update */
40 rvec xp[], /* Output coords */
41 t_nrnb *nrnb, /* Performance measure */
42 real *dvdlambda); /* FEP force */
43 /* Shake all the atoms blockwise. It is assumed that all the constraints
44 * in the idef->shakes field are sorted, to ascending block nr. The
45 * sblock array points into the idef->shakes.iatoms field, with block 0
46 * starting
47 * at sblock[0] and running to ( < ) sblock[1], block n running from
48 * sblock[n] to sblock[n+1]. Array sblock should be large enough.
49 * Return 0 when OK, -1 when shake-error
51 extern void csettle(FILE *log,
52 int nshake, /* Number of water molecules */
53 int owptr[], /* pointer to Oxygen in b4 & after */
54 real b4[], /* Old coordinates */
55 real after[], /* New coords, to be settled */
56 real dOH, /* Constraint length Ox-Hyd */
57 real dHH, /* Constraint length Hyd-Hyd */
58 real mO, /* Mass of Oxygen */
59 real mH, /* Mass of Hydrogen */
60 int *xerror);
62 extern void cshake(atom_id iatom[],int ncon,int *nnit,int maxnit,
63 real dist2[],real xp[],real rij[],real m2[],
64 real invmass[],real tt[],int *nerror);
65 /* Regular iterative shake */
67 extern void constrain(FILE *log,t_topology *top,t_inputrec *ir,int step,
68 t_mdatoms *md,int start,int homenr,
69 rvec *x,rvec *xprime,matrix box,
70 real lambda,real *dvdlambda,t_nrnb *nrnb);
71 /* Constrain coordinates xprime using the directions in x,
72 * init_constraints must have be called once, before calling constrain
75 extern bool init_constraints(FILE *log,t_topology *top,t_inputrec *ir,
76 t_mdatoms *md,int start,int homenr);
77 /* Initialize constraints stuff */
79 /* C routines for LINCS algorithm */
80 extern void clincs(rvec *x,rvec *xp,int ncons,
81 int *bla1,int *bla2,int *blnr,int *blbnb,real *bllen,
82 real *blc,real *blcc,real *blm,
83 int nit,int nrec,real *invmass,rvec * r,
84 real *vbo,real *vbn,real *vbt,real wangle,int *warn,
85 real *lambda);
87 extern void cconerr(real *max,real *rms,int *imax,rvec *xprime,
88 int ncons,int *bla1,int *bla2,real *bllen);
90 void lincs_warning(rvec *x,rvec *xprime,
91 int ncons,int *bla1,int *bla2,real *bllen,real wangle);