Fixed the config.h inclusion
[gromacs.git] / include / constr.h
blob95b106add21933342d9864cdc29fffd5e4b1fecb
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.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Gromacs Runs On Most of All Computer Systems
36 #ifdef HAVE_CONFIG_H
37 #include<config.h>
38 #endif
39 #include<callf77.h>
41 extern bool bshakef(FILE *log, /* Log file */
42 int natoms, /* Total number of atoms */
43 real invmass[], /* Atomic masses */
44 int nblocks, /* The number of shake blocks */
45 int sblock[], /* The shake blocks */
46 t_idef *idef, /* The interaction def */
47 t_inputrec *ir, /* Input record */
48 matrix box, /* The box */
49 rvec x_s[], /* Coords before update */
50 rvec xp[], /* Output coords */
51 t_nrnb *nrnb, /* Performance measure */
52 real lambda, /* FEP lambda */
53 real *dvdlambda, /* FEP force */
54 bool bDumpOnError); /* Dump debugging stuff on error*/
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 TRUE when OK, FALSE when shake-error
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);
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 */
79 extern bool 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);
85 * When bCoordinates=TRUE constrains coordinates xprime using th
86 * directions in x, min_proj is not used.
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.
92 * Init_constraints must have be called once, before calling constrain.
94 * Return TRUE if OK, FALSE in case of shake error
98 extern int count_constraints(t_topology *top,t_commrec *cr);
99 /* Returns the total number of constraints in the system */
101 extern int init_constraints(FILE *log,t_topology *top,t_inputrec *ir,
102 t_mdatoms *md,int start,int homenr,
103 bool bOnlyCoords,t_commrec *cr);
104 /* Initialize constraints stuff */
106 /* C routines for LINCS algorithm */
107 extern void clincsp(rvec *x,rvec *f,rvec *fp,int ncons,
108 int *bla1,int *bla2,int *blnr,int *blbnb,
109 real *blc,real *blcc,real *blm,
110 int nrec,real *invmass,rvec *r,
111 real *vbo,real *vbn,real *vbt);
113 extern void clincs(rvec *x,rvec *xp,int ncons,
114 int *bla1,int *bla2,int *blnr,int *blbnb,real *bllen,
115 real *blc,real *blcc,real *blm,
116 int nit,int nrec,real *invmass,rvec *r,
117 real *vbo,real *vbn,real *vbt,real wangle,int *warn,
118 real *lambda);
120 extern void cconerr(real *max,real *rms,int *imax,rvec *xprime,
121 int ncons,int *bla1,int *bla2,real *bllen);
123 void lincs_warning(rvec *x,rvec *xprime,
124 int ncons,int *bla1,int *bla2,real *bllen,real wangle);
127 #ifdef USE_FORTRAN
128 extern void F77_FUNC(fsettle,FSETTLE)(int *nshake,int owptr[],
129 real b4[],real after[],
130 real *dOH,real *dHH,
131 real *mO,real *mH,int *error);
132 extern void F77_FUNC(fshake,FSHAKE)(atom_id iatom[],int *ncon,
133 int *nit, int *maxnit,
134 real dist2[],real xp[],
135 real rij[],real m2[],real *omega,
136 real invmass[],real tt[],
137 real lambda[],int *error);
138 extern void F77_FUNC(flincsp,FLINCSP)(real *x,real *f,real *fp,
139 int *nc, int *bla1,int *bla2,
140 int *blnr,int *blbnb,
141 real *blc,real *blcc,real *blm,
142 int *nrec,real *invmass,
143 real *r, real *rhs1, real *rhs2,
144 real *sol);
145 extern void F77_FUNC(flincs,FLINCS)(real *x,real *xp,int *nc,
146 int *bla1,int *bla2,int *blnr,
147 int *blbnb,real *bllen,
148 real *blc,real *blcc,real *blm,
149 int *nit,int *nrec,real *invmass,
150 real *r,real *temp1,real *temp2,
151 real *temp3,real *wangle,
152 int *warn,real *lambda);
153 #endif