Tagged files with gromacs 3.0 header and added some license info
[alexxy/gromacs.git] / include / constr.h
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 3.0
11  * 
12  * Copyright (c) 1991-2001
13  * BIOSON Research Institute, Dept. of Biophysical Chemistry
14  * University of Groningen, The Netherlands
15  * 
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * Do check out http://www.gromacs.org , or mail us at gromacs@gromacs.org .
32  * 
33  * And Hey:
34  * GRoups of Organic Molecules in ACtion for Science
35  */
36 static char *SRCID_constr_h = "$Id$";
37 #ifdef HAVE_CONFIG_H
38 #include<config.h>
39 #endif
40 #include<callf77.h>
41
42 extern int bshakef(FILE *log,           /* Log file                     */
43                    int natoms,          /* Total number of atoms        */
44                    real invmass[],      /* Atomic masses                */
45                    int nblocks,         /* The number of shake blocks   */
46                    int sblock[],        /* The shake blocks             */
47                    t_idef *idef,        /* The interaction def          */
48                    t_inputrec *ir,      /* Input record                 */
49                    matrix box,          /* The box                      */
50                    rvec x_s[],          /* Coords before update         */
51                    rvec xp[],           /* Output coords                */
52                    t_nrnb *nrnb,        /* Performance measure          */
53                    real lambda,         /* FEP lambda                   */
54                    real *dvdlambda);    /* FEP force                    */
55 /* Shake all the atoms blockwise. It is assumed that all the constraints
56  * in the idef->shakes field are sorted, to ascending block nr. The
57  * sblock array points into the idef->shakes.iatoms field, with block 0 
58  * starting
59  * at sblock[0] and running to ( < ) sblock[1], block n running from 
60  * sblock[n] to sblock[n+1]. Array sblock should be large enough.
61  * Return 0 when OK, -1 when shake-error
62  */
63 extern void csettle(FILE *log,
64                     int nshake,         /* Number of water molecules    */
65                     int owptr[],        /* pointer to Oxygen in b4 & after */
66                     real b4[],          /* Old coordinates              */
67                     real after[],       /* New coords, to be settled    */
68                     real dOH,           /* Constraint length Ox-Hyd     */
69                     real dHH,           /* Constraint length Hyd-Hyd    */
70                     real mO,            /* Mass of Oxygen               */
71                     real mH,            /* Mass of Hydrogen             */
72                     int *xerror);
73
74 extern void cshake(atom_id iatom[],int ncon,int *nnit,int maxnit,
75                    real dist2[],real xp[],real rij[],real m2[],real omega,
76                    real invmass[],real tt[],real lagr[],int *nerror);
77 /* Regular iterative shake */
78
79 extern void constrain(FILE *log,t_topology *top,t_inputrec *ir,int step,
80                       t_mdatoms *md,int start,int homenr,
81                       rvec *x,rvec *xprime,rvec *min_proj,matrix box,
82                       real lambda,real *dvdlambda,t_nrnb *nrnb,
83                       bool bCoordinates);
84 /*
85  * When bCoordinates=TRUE constrains coordinates xprime using th
86  * directions in x, min_proj is not used.
87  *
88  * When bCoordinates=FALSE, calculates the components xprime in
89  * the constraint directions and subtracts these components from min_proj.
90  * So when min_proj=xprime, the constraint components are projected out.
91  *
92  * Init_constraints must have be called once, before calling constrain.
93  */
94
95 extern bool init_constraints(FILE *log,t_topology *top,t_inputrec *ir,
96                              t_mdatoms *md,int start,int homenr,
97                              bool bOnlyCoords);
98 /* Initialize constraints stuff */
99
100 /* C routines for LINCS algorithm */ 
101 extern void clincsp(rvec *x,rvec *f,rvec *fp,int ncons,
102                     int *bla1,int *bla2,int *blnr,int *blbnb,
103                     real *blc,real *blcc,real *blm,
104                     int nrec,real *invmass,rvec *r,
105                     real *vbo,real *vbn,real *vbt);
106
107 extern void clincs(rvec *x,rvec *xp,int ncons,
108                    int *bla1,int *bla2,int *blnr,int *blbnb,real *bllen,
109                    real *blc,real *blcc,real *blm,
110                    int nit,int nrec,real *invmass,rvec *r,
111                    real *vbo,real *vbn,real *vbt,real wangle,int *warn,
112                    real *lambda);
113
114 extern void cconerr(real *max,real *rms,int *imax,rvec *xprime,
115                     int ncons,int *bla1,int *bla2,real *bllen);
116
117 void lincs_warning(rvec *x,rvec *xprime,
118                    int ncons,int *bla1,int *bla2,real *bllen,real wangle);
119
120
121 #ifdef USE_FORTRAN
122 extern void F77_FUNC(fsettle,FSETTLE)(int *nshake,int owptr[],
123                                       real b4[],real after[],
124                                       real *dOH,real *dHH,
125                                       real *mO,real *mH,int *error);
126 extern void F77_FUNC(fshake,FSHAKE)(atom_id iatom[],int *ncon,
127                                     int *nit, int *maxnit,
128                                     real dist2[],real xp[],
129                                     real rij[],real m2[],real *omega,
130                                     real invmass[],real tt[],
131                                     real lambda[],int *error);
132 extern void F77_FUNC(flincsp,FLINCSP)(real *x,real *f,real *fp,
133                                       int *nc, int *bla1,int *bla2,
134                                       int *blnr,int *blbnb,
135                                       real *blc,real *blcc,real *blm,
136                                       int *nrec,real *invmass,
137                                       real *r, real *rhs1, real *rhs2,
138                                       real *sol);
139 extern void F77_FUNC(flincs,FLINCS)(real *x,real *xp,int *nc,
140                                     int *bla1,int *bla2,int *blnr,
141                                     int *blbnb,real *bllen,
142                                     real *blc,real *blcc,real *blm,
143                                     int *nit,int *nrec,real *invmass,
144                                     real *r,real *temp1,real *temp2,
145                                     real *temp3,real *wangle,
146                                     int *warn,real *lambda);
147 #endif