1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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
34 * Gromacs Runs On Most of All Computer Systems
40 #include "types/commrec.h"
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) */
58 int n_flexible_constraints(struct gmx_constr *constr);
59 /* Returns the total number of flexible constraints in the system */
61 void too_many_constraint_warnings(int eConstrAlg, int warncount);
62 /* Generate a fatal error because of too many LINCS/SETTLE warnings */
64 gmx_shakedata_t shake_init();
65 /* Initializes and return the SHAKE data structure */
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
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
96 gmx_settledata_t settle_init(real mO, real mH, real invmO, real invmH,
98 /* Initializes and returns a structure with SETTLE parameters */
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 */
111 t_vetavars *vetavar /* variables for pressure control */
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 */
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.
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 */
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);
134 gmx_bool constrain(FILE *log, gmx_bool bLog, gmx_bool bEner,
138 gmx_ekindata_t *ekind,
140 gmx_large_int_t step, int delta_step,
142 rvec *x, rvec *xprime, rvec *min_proj,
143 gmx_bool bMolPBC, matrix box,
144 real lambda, real *dvdlambda,
145 rvec *v, tensor *vir,
146 t_nrnb *nrnb, int econq, gmx_bool bPscal, real veta, real vetanew);
148 * When econq=econqCoord constrains coordinates xprime using th
149 * directions in x, min_proj is not used.
151 * When econq=econqDeriv, calculates the components xprime in
152 * the constraint directions and subtracts these components from min_proj.
153 * So when min_proj=xprime, the constraint components are projected out.
155 * When econq=econqDeriv_FlexCon, the same is done as with econqDeriv,
156 * but only the components of the flexible constraints are stored.
158 * When bMolPBC=TRUE, assume that molecules might be broken: correct PBC.
160 * delta_step is used for determining the constraint reference lengths
161 * when lenA != lenB or will the pull code with a pulling rate.
162 * step + delta_step is the step at which the final configuration
163 * is meant to be; for update delta_step = 1.
165 * If v!=NULL also constrain v by adding the constraint corrections / dt.
167 * If vir!=NULL calculate the constraint virial.
169 * if veta != NULL, constraints are done assuming isotropic pressure control
170 * (i.e. constraining rdot.r = (v + veta*r).r = 0 instead of v
172 * Init_constraints must have be called once, before calling constrain.
174 * Return TRUE if OK, FALSE in case of shake error
178 gmx_constr_t init_constraints(FILE *log,
179 gmx_mtop_t *mtop, t_inputrec *ir,
180 gmx_edsam_t ed, t_state *state,
182 /* Initialize constraints stuff */
184 void set_constraints(gmx_constr_t constr,
189 /* Set up all the local constraints for the node */
191 /* The at2con t_blocka struct returned by the routines below
192 * contains a list of constraints per atom.
193 * The F_CONSTRNC constraints in this structure number consecutively
194 * after the F_CONSTR constraints.
197 t_blocka make_at2con(int start, int natoms,
198 t_ilist *ilist, t_iparams *iparams,
199 gmx_bool bDynamics, int *nflexiblecons);
200 /* Returns a block struct to go from atoms to constraints */
202 const t_blocka *atom2constraints_moltype(gmx_constr_t constr);
203 /* Returns the an array of atom to constraints lists for the moltypes */
205 const int **atom2settle_moltype(gmx_constr_t constr);
206 /* Returns the an array of atom to settle for the moltypes */
208 #define constr_iatomptr(nconstr, iatom_constr, iatom_constrnc, con) ((con) < (nconstr) ? (iatom_constr)+(con)*3 : (iatom_constrnc)+(con-nconstr)*3)
209 /* Macro for getting the constraint iatoms for a constraint number con
210 * which comes from a list where F_CONSTR and F_CONSTRNC constraints
214 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop);
215 /* Returns if there are inter charge group constraints */
217 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop);
218 /* Returns if there are inter charge group settles */
220 real *constr_rmsd_data(gmx_constr_t constr);
221 /* Return the data for determining constraint RMS relative deviations.
222 * Returns NULL when LINCS is not used.
225 real constr_rmsd(gmx_constr_t constr, gmx_bool bSD2);
226 /* Return the RMSD of the constraint, bSD2 selects the second SD step */
228 real *lincs_rmsd_data(gmx_lincsdata_t lincsd);
229 /* Return the data for determining constraint RMS relative deviations */
231 real lincs_rmsd(gmx_lincsdata_t lincsd, gmx_bool bSD2);
232 /* Return the RMSD of the constraint, bSD2 selects the second SD step */
234 gmx_lincsdata_t init_lincs(FILE *fplog, gmx_mtop_t *mtop,
235 int nflexcon_global, t_blocka *at2con,
236 gmx_bool bPLINCS, int nIter, int nProjOrder);
237 /* Initializes and returns the lincs data struct */
239 void set_lincs(t_idef *idef, t_mdatoms *md,
240 gmx_bool bDynamics, t_commrec *cr,
242 /* Initialize lincs stuff */
244 void set_lincs_matrix(gmx_lincsdata_t li, real *invmass, real lambda);
245 /* Sets the elements of the LINCS constraint coupling matrix */
247 real constr_r_max(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir);
248 /* Returns an estimate of the maximum distance between atoms
249 * required for LINCS.
253 constrain_lincs(FILE *log, gmx_bool bLog, gmx_bool bEner,
255 gmx_large_int_t step,
256 gmx_lincsdata_t lincsd, t_mdatoms *md,
258 rvec *x, rvec *xprime, rvec *min_proj,
259 matrix box, t_pbc *pbc,
260 real lambda, real *dvdlambda,
262 gmx_bool bCalcVir, tensor vir_r_m_dr,
265 int maxwarn, int *warncount);
266 /* Returns if the constraining succeeded */
269 /* helper functions for andersen temperature control, because the
270 * gmx_constr construct is only defined in constr.c. Return the list
271 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
273 int *get_sblock(struct gmx_constr *constr);
275 int get_nblocks(struct gmx_constr *constr);