4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
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
27 * Good ROcking Metal Altar for Chronical Sinners
29 static char *SRCID_constr_h = "$Id$";
36 extern int bshakef(FILE *log, /* Log file */
37 int natoms, /* Total number of atoms */
38 real invmass[], /* Atomic masses */
39 int nblocks, /* The number of shake blocks */
40 int sblock[], /* The shake blocks */
41 t_idef *idef, /* The interaction def */
42 t_inputrec *ir, /* Input record */
43 matrix box, /* The box */
44 rvec x_s[], /* Coords before update */
45 rvec xp[], /* Output coords */
46 t_nrnb *nrnb, /* Performance measure */
47 real lambda, /* FEP lambda */
48 real *dvdlambda); /* FEP force */
49 /* Shake all the atoms blockwise. It is assumed that all the constraints
50 * in the idef->shakes field are sorted, to ascending block nr. The
51 * sblock array points into the idef->shakes.iatoms field, with block 0
53 * at sblock[0] and running to ( < ) sblock[1], block n running from
54 * sblock[n] to sblock[n+1]. Array sblock should be large enough.
55 * Return 0 when OK, -1 when shake-error
57 extern void csettle(FILE *log,
58 int nshake, /* Number of water molecules */
59 int owptr[], /* pointer to Oxygen in b4 & after */
60 real b4[], /* Old coordinates */
61 real after[], /* New coords, to be settled */
62 real dOH, /* Constraint length Ox-Hyd */
63 real dHH, /* Constraint length Hyd-Hyd */
64 real mO, /* Mass of Oxygen */
65 real mH, /* Mass of Hydrogen */
68 extern void cshake(atom_id iatom[],int ncon,int *nnit,int maxnit,
69 real dist2[],real xp[],real rij[],real m2[],real omega,
70 real invmass[],real tt[],real lagr[],int *nerror);
71 /* Regular iterative shake */
73 extern void constrain(FILE *log,t_topology *top,t_inputrec *ir,int step,
74 t_mdatoms *md,int start,int homenr,
75 rvec *x,rvec *xprime,rvec *min_proj,matrix box,
76 real lambda,real *dvdlambda,t_nrnb *nrnb,
79 * When bCoordinates=TRUE constrains coordinates xprime using th
80 * directions in x, min_proj is not used.
82 * When bCoordinates=FALSE, calculates the components xprime in
83 * the constraint directions and subtracts these components from min_proj.
84 * So when min_proj=xprime, the constraint components are projected out.
86 * Init_constraints must have be called once, before calling constrain.
89 extern bool init_constraints(FILE *log,t_topology *top,t_inputrec *ir,
90 t_mdatoms *md,int start,int homenr,
92 /* Initialize constraints stuff */
94 /* C routines for LINCS algorithm */
95 extern void clincsp(rvec *x,rvec *f,rvec *fp,int ncons,
96 int *bla1,int *bla2,int *blnr,int *blbnb,
97 real *blc,real *blcc,real *blm,
98 int nrec,real *invmass,rvec *r,
99 real *vbo,real *vbn,real *vbt);
101 extern void clincs(rvec *x,rvec *xp,int ncons,
102 int *bla1,int *bla2,int *blnr,int *blbnb,real *bllen,
103 real *blc,real *blcc,real *blm,
104 int nit,int nrec,real *invmass,rvec *r,
105 real *vbo,real *vbn,real *vbt,real wangle,int *warn,
108 extern void cconerr(real *max,real *rms,int *imax,rvec *xprime,
109 int ncons,int *bla1,int *bla2,real *bllen);
111 void lincs_warning(rvec *x,rvec *xprime,
112 int ncons,int *bla1,int *bla2,real *bllen,real wangle);
116 extern void F77_FUNC(fsettle,FSETTLE)(int *nshake,int owptr[],
117 real b4[],real after[],
119 real *mO,real *mH,int *error);
120 extern void F77_FUNC(fshake,FSHAKE)(atom_id iatom[],int *ncon,
121 int *nit, int *maxnit,
122 real dist2[],real xp[],
123 real rij[],real m2[],real *omega,
124 real invmass[],real tt[],
125 real lambda[],int *error);
126 extern void F77_FUNC(flincsp,FLINCSP)(real *x,real *f,real *fp,
127 int *nc, int *bla1,int *bla2,
128 int *blnr,int *blbnb,
129 real *blc,real *blcc,real *blm,
130 int *nrec,real *invmass,
131 real *r, real *rhs1, real *rhs2,
133 extern void F77_FUNC(flincs,FLINCS)(real *x,real *xp,int *nc,
134 int *bla1,int *bla2,int *blnr,
135 int *blbnb,real *bllen,
136 real *blc,real *blcc,real *blm,
137 int *nit,int *nrec,real *invmass,
138 real *r,real *temp1,real *temp2,
139 real *temp3,real *wangle,
140 int *warn,real *lambda);