Default SD integrator is now sd1 by Nicu Goga
[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,2014, 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
42 #ifdef __cplusplus
43 extern "C" {
44 #endif
45
46 enum
47 {
48     econqCoord,         /* Constrain coordinates (mass weighted)           */
49     econqVeloc,         /* Constrain velocities (mass weighted)            */
50     econqDeriv,         /* Constrain a derivative (mass weighted),         *
51                          * for instance velocity or acceleration,          *
52                          * constraint virial can not be calculated.        */
53     econqDeriv_FlexCon, /* As econqDeriv, but only output flex. con.       */
54     econqForce,         /* Constrain forces (non mass-weighted)            */
55     econqForceDispl     /* Constrain forces (mass-weighted 1/0 for freeze) */
56 };
57
58 int n_flexible_constraints(struct gmx_constr *constr);
59 /* Returns the total number of flexible constraints in the system */
60
61 void too_many_constraint_warnings(int eConstrAlg, int warncount);
62 /* Generate a fatal error because of too many LINCS/SETTLE warnings */
63
64 gmx_shakedata_t shake_init();
65 /* Initializes and return the SHAKE data structure */
66
67 gmx_bool bshakef(FILE           *log,          /* Log file                      */
68                  gmx_shakedata_t shaked,       /* Total number of atoms */
69                  real            invmass[],    /* Atomic masses         */
70                  int             nblocks,      /* The number of shake blocks    */
71                  int             sblock[],     /* The shake blocks             */
72                  t_idef         *idef,         /* The interaction def           */
73                  t_inputrec     *ir,           /* Input record                  */
74                  rvec            x_s[],        /* Coords before update          */
75                  rvec            prime[],      /* Output coords         */
76                  t_nrnb         *nrnb,         /* Performance measure          */
77                  real           *lagr,         /* The Lagrange multipliers     */
78                  real            lambda,       /* FEP lambda                   */
79                  real           *dvdlambda,    /* FEP force                    */
80                  real            invdt,        /* 1/delta_t                    */
81                  rvec           *v,            /* Also constrain v if v!=NULL  */
82                  gmx_bool        bCalcVir,     /* Calculate r x m delta_r      */
83                  tensor          vir_r_m_dr,   /* sum r x m delta_r            */
84                  gmx_bool        bDumpOnError, /* Dump debugging stuff on error*/
85                  int             econq,        /* which type of constrainint is occurring */
86                  t_vetavars     *vetavar);     /* veta for pressure control */
87 /* Shake all the atoms blockwise. It is assumed that all the constraints
88  * in the idef->shakes field are sorted, to ascending block nr. The
89  * sblock array points into the idef->shakes.iatoms field, with block 0
90  * starting
91  * at sblock[0] and running to ( < ) sblock[1], block n running from
92  * sblock[n] to sblock[n+1]. Array sblock should be large enough.
93  * Return TRUE when OK, FALSE when shake-error
94  */
95
96 gmx_settledata_t settle_init(real mO, real mH, real invmO, real invmH,
97                              real dOH, real dHH);
98 /* Initializes and returns a structure with SETTLE parameters */
99
100 void csettle(gmx_settledata_t settled,
101              int              nsettle,          /* Number of settles            */
102              t_iatom          iatoms[],         /* The settle iatom list        */
103              const t_pbc     *pbc,              /* PBC data pointer, can be NULL  */
104              real             b4[],             /* Old coordinates              */
105              real             after[],          /* New coords, to be settled    */
106              real             invdt,            /* 1/delta_t                    */
107              real            *v,                /* Also constrain v if v!=NULL  */
108              int              calcvir_atom_end, /* Calculate r x m delta_r up to this atom */
109              tensor           vir_r_m_dr,       /* sum r x m delta_r            */
110              int             *xerror,
111              t_vetavars      *vetavar           /* variables for pressure control */
112              );
113
114 void settle_proj(gmx_settledata_t settled, int econq,
115                  int nsettle, t_iatom iatoms[],
116                  const t_pbc *pbc,   /* PBC data pointer, can be NULL  */
117                  rvec x[],
118                  rvec *der, rvec *derp,
119                  int CalcVirAtomEnd, tensor vir_r_m_dder,
120                  t_vetavars *vetavar);
121 /* Analytical algorithm to subtract the components of derivatives
122  * of coordinates working on settle type constraint.
123  */
124
125 void cshake(atom_id iatom[], int ncon, int *nnit, int maxnit,
126             real dist2[], real xp[], real rij[], real m2[], real omega,
127             real invmass[], real tt[], real lagr[], int *nerror);
128 /* Regular iterative shake */
129
130 void crattle(atom_id iatom[], int ncon, int *nnit, int maxnit,
131              real dist2[], real vp[], real rij[], real m2[], real omega,
132              real invmass[], real tt[], real lagr[], int *nerror, real invdt, t_vetavars *vetavar);
133
134 gmx_bool constrain(FILE *log, gmx_bool bLog, gmx_bool bEner,
135                    gmx_constr_t constr,
136                    t_idef *idef,
137                    t_inputrec *ir,
138                    gmx_ekindata_t *ekind,
139                    t_commrec *cr,
140                    gmx_int64_t step, int delta_step,
141                    real step_scaling,
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  * step_scaling can be used to update coordinates based on the time
167  * step multiplied by this factor. Thus, normally 1.0 is passed. The
168  * SD1 integrator uses 0.5 in one of its calls, to correct positions
169  * for half a step of changed velocities.
170  *
171  * If v!=NULL also constrain v by adding the constraint corrections / dt.
172  *
173  * If vir!=NULL calculate the constraint virial.
174  *
175  * if veta != NULL, constraints are done assuming isotropic pressure control
176  * (i.e. constraining rdot.r = (v + veta*r).r = 0 instead of v
177  *
178  * Init_constraints must have be called once, before calling constrain.
179  *
180  * Return TRUE if OK, FALSE in case of shake error
181  *
182  */
183
184 gmx_constr_t init_constraints(FILE *log,
185                               gmx_mtop_t *mtop, t_inputrec *ir,
186                               gmx_edsam_t ed, t_state *state,
187                               t_commrec *cr);
188 /* Initialize constraints stuff */
189
190 void set_constraints(gmx_constr_t    constr,
191                      gmx_localtop_t *top,
192                      t_inputrec     *ir,
193                      t_mdatoms      *md,
194                      t_commrec      *cr);
195 /* Set up all the local constraints for the node */
196
197 /* The at2con t_blocka struct returned by the routines below
198  * contains a list of constraints per atom.
199  * The F_CONSTRNC constraints in this structure number consecutively
200  * after the F_CONSTR constraints.
201  */
202
203 t_blocka make_at2con(int start, int natoms,
204                      t_ilist *ilist, t_iparams *iparams,
205                      gmx_bool bDynamics, int *nflexiblecons);
206 /* Returns a block struct to go from atoms to constraints */
207
208 const t_blocka *atom2constraints_moltype(gmx_constr_t constr);
209 /* Returns the an array of atom to constraints lists for the moltypes */
210
211 const int **atom2settle_moltype(gmx_constr_t constr);
212 /* Returns the an array of atom to settle for the moltypes */
213
214 #define constr_iatomptr(nconstr, iatom_constr, iatom_constrnc, con) ((con) < (nconstr) ? (iatom_constr)+(con)*3 : (iatom_constrnc)+(con-nconstr)*3)
215 /* Macro for getting the constraint iatoms for a constraint number con
216  * which comes from a list where F_CONSTR and F_CONSTRNC constraints
217  * are concatenated.
218  */
219
220 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop);
221 /* Returns if there are inter charge group constraints */
222
223 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop);
224 /* Returns if there are inter charge group settles */
225
226 real *constr_rmsd_data(gmx_constr_t constr);
227 /* Return the data for determining constraint RMS relative deviations.
228  * Returns NULL when LINCS is not used.
229  */
230
231 real constr_rmsd(gmx_constr_t constr, gmx_bool bSD2);
232 /* Return the RMSD of the constraint, bSD2 selects the second SD step */
233
234 real *lincs_rmsd_data(gmx_lincsdata_t lincsd);
235 /* Return the data for determining constraint RMS relative deviations */
236
237 real lincs_rmsd(gmx_lincsdata_t lincsd, gmx_bool bSD2);
238 /* Return the RMSD of the constraint, bSD2 selects the second SD step */
239
240 gmx_lincsdata_t init_lincs(FILE *fplog, gmx_mtop_t *mtop,
241                            int nflexcon_global, t_blocka *at2con,
242                            gmx_bool bPLINCS, int nIter, int nProjOrder);
243 /* Initializes and returns the lincs data struct */
244
245 void set_lincs(t_idef *idef, t_mdatoms *md,
246                gmx_bool bDynamics, t_commrec *cr,
247                gmx_lincsdata_t li);
248 /* Initialize lincs stuff */
249
250 void set_lincs_matrix(gmx_lincsdata_t li, real *invmass, real lambda);
251 /* Sets the elements of the LINCS constraint coupling matrix */
252
253 real constr_r_max(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir);
254 /* Returns an estimate of the maximum distance between atoms
255  * required for LINCS.
256  */
257
258 gmx_bool
259 constrain_lincs(FILE *log, gmx_bool bLog, gmx_bool bEner,
260                 t_inputrec *ir,
261                 gmx_int64_t step,
262                 gmx_lincsdata_t lincsd, t_mdatoms *md,
263                 t_commrec *cr,
264                 rvec *x, rvec *xprime, rvec *min_proj,
265                 matrix box, t_pbc *pbc,
266                 real lambda, real *dvdlambda,
267                 real invdt, rvec *v,
268                 gmx_bool bCalcVir, tensor vir_r_m_dr,
269                 int econ,
270                 t_nrnb *nrnb,
271                 int maxwarn, int *warncount);
272 /* Returns if the constraining succeeded */
273
274
275 /* helper functions for andersen temperature control, because the
276  * gmx_constr construct is only defined in constr.c. Return the list
277  * of blocks (get_sblock) and the number of blocks (get_nblocks).  */
278
279 int *get_sblock(struct gmx_constr *constr);
280
281 int get_nblocks(struct gmx_constr *constr);
282
283 #ifdef __cplusplus
284 }
285 #endif
286 #endif