Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / legacyheaders / constr.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifndef _constr_h
39 #define _constr_h
40 #include "typedefs.h"
41 #include "types/commrec.h"
42
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
46
47 enum
48 {
49     econqCoord,         /* Constrain coordinates (mass weighted)           */
50     econqVeloc,         /* Constrain velocities (mass weighted)            */
51     econqDeriv,         /* Constrain a derivative (mass weighted),         *
52                          * for instance velocity or acceleration,          *
53                          * constraint virial can not be calculated.        */
54     econqDeriv_FlexCon, /* As econqDeriv, but only output flex. con.       */
55     econqForce,         /* Constrain forces (non mass-weighted)            */
56     econqForceDispl     /* Constrain forces (mass-weighted 1/0 for freeze) */
57 };
58
59 int n_flexible_constraints(struct gmx_constr *constr);
60 /* Returns the total number of flexible constraints in the system */
61
62 void too_many_constraint_warnings(int eConstrAlg, int warncount);
63 /* Generate a fatal error because of too many LINCS/SETTLE warnings */
64
65 gmx_shakedata_t shake_init();
66 /* Initializes and return the SHAKE data structure */
67
68 gmx_bool bshakef(FILE           *log,          /* Log file                      */
69                  gmx_shakedata_t shaked,       /* Total number of atoms */
70                  real            invmass[],    /* Atomic masses         */
71                  int             nblocks,      /* The number of shake blocks    */
72                  int             sblock[],     /* The shake blocks             */
73                  t_idef         *idef,         /* The interaction def           */
74                  t_inputrec     *ir,           /* Input record                  */
75                  rvec            x_s[],        /* Coords before update          */
76                  rvec            prime[],      /* Output coords         */
77                  t_nrnb         *nrnb,         /* Performance measure          */
78                  real           *lagr,         /* The Lagrange multipliers     */
79                  real            lambda,       /* FEP lambda                   */
80                  real           *dvdlambda,    /* FEP force                    */
81                  real            invdt,        /* 1/delta_t                    */
82                  rvec           *v,            /* Also constrain v if v!=NULL  */
83                  gmx_bool        bCalcVir,     /* Calculate r x m delta_r      */
84                  tensor          vir_r_m_dr,   /* sum r x m delta_r            */
85                  gmx_bool        bDumpOnError, /* Dump debugging stuff on error*/
86                  int             econq,        /* which type of constrainint is occurring */
87                  t_vetavars     *vetavar);     /* veta for pressure control */
88 /* Shake all the atoms blockwise. It is assumed that all the constraints
89  * in the idef->shakes field are sorted, to ascending block nr. The
90  * sblock array points into the idef->shakes.iatoms field, with block 0
91  * starting
92  * at sblock[0] and running to ( < ) sblock[1], block n running from
93  * sblock[n] to sblock[n+1]. Array sblock should be large enough.
94  * Return TRUE when OK, FALSE when shake-error
95  */
96
97 gmx_settledata_t settle_init(real mO, real mH, real invmO, real invmH,
98                              real dOH, real dHH);
99 /* Initializes and returns a structure with SETTLE parameters */
100
101 void csettle(gmx_settledata_t settled,
102              int              nsettle,          /* Number of settles            */
103              t_iatom          iatoms[],         /* The settle iatom list        */
104              const t_pbc     *pbc,              /* PBC data pointer, can be NULL  */
105              real             b4[],             /* Old coordinates              */
106              real             after[],          /* New coords, to be settled    */
107              real             invdt,            /* 1/delta_t                    */
108              real            *v,                /* Also constrain v if v!=NULL  */
109              int              calcvir_atom_end, /* Calculate r x m delta_r up to this atom */
110              tensor           vir_r_m_dr,       /* sum r x m delta_r            */
111              int             *xerror,
112              t_vetavars      *vetavar           /* variables for pressure control */
113              );
114
115 void settle_proj(gmx_settledata_t settled, int econq,
116                  int nsettle, t_iatom iatoms[],
117                  const t_pbc *pbc,   /* PBC data pointer, can be NULL  */
118                  rvec x[],
119                  rvec *der, rvec *derp,
120                  int CalcVirAtomEnd, tensor vir_r_m_dder,
121                  t_vetavars *vetavar);
122 /* Analytical algorithm to subtract the components of derivatives
123  * of coordinates working on settle type constraint.
124  */
125
126 void cshake(atom_id iatom[], int ncon, int *nnit, int maxnit,
127             real dist2[], real xp[], real rij[], real m2[], real omega,
128             real invmass[], real tt[], real lagr[], int *nerror);
129 /* Regular iterative shake */
130
131 void crattle(atom_id iatom[], int ncon, int *nnit, int maxnit,
132              real dist2[], real vp[], real rij[], real m2[], real omega,
133              real invmass[], real tt[], real lagr[], int *nerror, real invdt, t_vetavars *vetavar);
134
135 gmx_bool constrain(FILE *log, gmx_bool bLog, gmx_bool bEner,
136                    gmx_constr_t constr,
137                    t_idef *idef,
138                    t_inputrec *ir,
139                    gmx_ekindata_t *ekind,
140                    t_commrec *cr,
141                    gmx_int64_t step, int delta_step,
142                    t_mdatoms *md,
143                    rvec *x, rvec *xprime, rvec *min_proj,
144                    gmx_bool bMolPBC, matrix box,
145                    real lambda, real *dvdlambda,
146                    rvec *v, tensor *vir,
147                    t_nrnb *nrnb, int econq, gmx_bool bPscal, real veta, real vetanew);
148 /*
149  * When econq=econqCoord constrains coordinates xprime using th
150  * directions in x, min_proj is not used.
151  *
152  * When econq=econqDeriv, calculates the components xprime in
153  * the constraint directions and subtracts these components from min_proj.
154  * So when min_proj=xprime, the constraint components are projected out.
155  *
156  * When econq=econqDeriv_FlexCon, the same is done as with econqDeriv,
157  * but only the components of the flexible constraints are stored.
158  *
159  * When bMolPBC=TRUE, assume that molecules might be broken: correct PBC.
160  *
161  * delta_step is used for determining the constraint reference lengths
162  * when lenA != lenB or will the pull code with a pulling rate.
163  * step + delta_step is the step at which the final configuration
164  * is meant to be; for update delta_step = 1.
165  *
166  * If v!=NULL also constrain v by adding the constraint corrections / dt.
167  *
168  * If vir!=NULL calculate the constraint virial.
169  *
170  * if veta != NULL, constraints are done assuming isotropic pressure control
171  * (i.e. constraining rdot.r = (v + veta*r).r = 0 instead of v
172  *
173  * Init_constraints must have be called once, before calling constrain.
174  *
175  * Return TRUE if OK, FALSE in case of shake error
176  *
177  */
178
179 gmx_constr_t init_constraints(FILE *log,
180                               gmx_mtop_t *mtop, t_inputrec *ir,
181                               gmx_edsam_t ed, t_state *state,
182                               t_commrec *cr);
183 /* Initialize constraints stuff */
184
185 void set_constraints(gmx_constr_t    constr,
186                      gmx_localtop_t *top,
187                      t_inputrec     *ir,
188                      t_mdatoms      *md,
189                      t_commrec      *cr);
190 /* Set up all the local constraints for the node */
191
192 /* The at2con t_blocka struct returned by the routines below
193  * contains a list of constraints per atom.
194  * The F_CONSTRNC constraints in this structure number consecutively
195  * after the F_CONSTR constraints.
196  */
197
198 t_blocka make_at2con(int start, int natoms,
199                      t_ilist *ilist, t_iparams *iparams,
200                      gmx_bool bDynamics, int *nflexiblecons);
201 /* Returns a block struct to go from atoms to constraints */
202
203 const t_blocka *atom2constraints_moltype(gmx_constr_t constr);
204 /* Returns the an array of atom to constraints lists for the moltypes */
205
206 const int **atom2settle_moltype(gmx_constr_t constr);
207 /* Returns the an array of atom to settle for the moltypes */
208
209 #define constr_iatomptr(nconstr, iatom_constr, iatom_constrnc, con) ((con) < (nconstr) ? (iatom_constr)+(con)*3 : (iatom_constrnc)+(con-nconstr)*3)
210 /* Macro for getting the constraint iatoms for a constraint number con
211  * which comes from a list where F_CONSTR and F_CONSTRNC constraints
212  * are concatenated.
213  */
214
215 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop);
216 /* Returns if there are inter charge group constraints */
217
218 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop);
219 /* Returns if there are inter charge group settles */
220
221 real *constr_rmsd_data(gmx_constr_t constr);
222 /* Return the data for determining constraint RMS relative deviations.
223  * Returns NULL when LINCS is not used.
224  */
225
226 real constr_rmsd(gmx_constr_t constr, gmx_bool bSD2);
227 /* Return the RMSD of the constraint, bSD2 selects the second SD step */
228
229 real *lincs_rmsd_data(gmx_lincsdata_t lincsd);
230 /* Return the data for determining constraint RMS relative deviations */
231
232 real lincs_rmsd(gmx_lincsdata_t lincsd, gmx_bool bSD2);
233 /* Return the RMSD of the constraint, bSD2 selects the second SD step */
234
235 gmx_lincsdata_t init_lincs(FILE *fplog, gmx_mtop_t *mtop,
236                            int nflexcon_global, t_blocka *at2con,
237                            gmx_bool bPLINCS, int nIter, int nProjOrder);
238 /* Initializes and returns the lincs data struct */
239
240 void set_lincs(t_idef *idef, t_mdatoms *md,
241                gmx_bool bDynamics, t_commrec *cr,
242                gmx_lincsdata_t li);
243 /* Initialize lincs stuff */
244
245 void set_lincs_matrix(gmx_lincsdata_t li, real *invmass, real lambda);
246 /* Sets the elements of the LINCS constraint coupling matrix */
247
248 real constr_r_max(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir);
249 /* Returns an estimate of the maximum distance between atoms
250  * required for LINCS.
251  */
252
253 gmx_bool
254 constrain_lincs(FILE *log, gmx_bool bLog, gmx_bool bEner,
255                 t_inputrec *ir,
256                 gmx_int64_t step,
257                 gmx_lincsdata_t lincsd, t_mdatoms *md,
258                 t_commrec *cr,
259                 rvec *x, rvec *xprime, rvec *min_proj,
260                 matrix box, t_pbc *pbc,
261                 real lambda, real *dvdlambda,
262                 real invdt, rvec *v,
263                 gmx_bool bCalcVir, tensor vir_r_m_dr,
264                 int econ,
265                 t_nrnb *nrnb,
266                 int maxwarn, int *warncount);
267 /* Returns if the constraining succeeded */
268
269
270 /* helper functions for andersen temperature control, because the
271  * gmx_constr construct is only defined in constr.c. Return the list
272  * of blocks (get_sblock) and the number of blocks (get_nblocks).  */
273
274 int *get_sblock(struct gmx_constr *constr);
275
276 int get_nblocks(struct gmx_constr *constr);
277
278 #ifdef __cplusplus
279 }
280 #endif
281 #endif