2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
41 #include "gromacs/fileio/confio.h"
42 #include "types/commrec.h"
49 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/fileio/pdbio.h"
57 #include "mtop_util.h"
58 #include "gromacs/fileio/gmxfio.h"
60 #include "gmx_omp_nthreads.h"
61 #include "gromacs/essentialdynamics/edsam.h"
62 #include "gromacs/pulling/pull.h"
64 #include "gmx_fatal.h"
66 typedef struct gmx_constr {
67 int ncon_tot; /* The total number of constraints */
68 int nflexcon; /* The number of flexible constraints */
69 int n_at2con_mt; /* The size of at2con = #moltypes */
70 t_blocka *at2con_mt; /* A list of atoms to constraints */
71 int n_at2settle_mt; /* The size of at2settle = #moltypes */
72 int **at2settle_mt; /* A list of atoms to settles */
73 gmx_bool bInterCGsettles;
74 gmx_lincsdata_t lincsd; /* LINCS data */
75 gmx_shakedata_t shaked; /* SHAKE data */
76 gmx_settledata_t settled; /* SETTLE data */
77 int nblocks; /* The number of SHAKE blocks */
78 int *sblock; /* The SHAKE blocks */
79 int sblock_nalloc; /* The allocation size of sblock */
80 real *lagr; /* Lagrange multipliers for SHAKE */
81 int lagr_nalloc; /* The allocation size of lagr */
82 int maxwarn; /* The maximum number of warnings */
85 gmx_edsam_t ed; /* The essential dynamics data */
87 tensor *vir_r_m_dr_th; /* Thread local working data */
88 int *settle_error; /* Thread local working data */
90 gmx_mtop_t *warn_mtop; /* Only used for printing warnings */
98 static void *init_vetavars(t_vetavars *vars,
99 gmx_bool constr_deriv,
100 real veta, real vetanew, t_inputrec *ir, gmx_ekindata_t *ekind, gmx_bool bPscal)
105 /* first, set the alpha integrator variable */
106 if ((ir->opts.nrdf[0] > 0) && bPscal)
108 vars->alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
114 g = 0.5*veta*ir->delta_t;
115 vars->rscale = exp(g)*series_sinhx(g);
116 g = -0.25*vars->alpha*veta*ir->delta_t;
117 vars->vscale = exp(g)*series_sinhx(g);
118 vars->rvscale = vars->vscale*vars->rscale;
119 vars->veta = vetanew;
123 snew(vars->vscale_nhc, ir->opts.ngtc);
124 if ((ekind == NULL) || (!bPscal))
126 for (i = 0; i < ir->opts.ngtc; i++)
128 vars->vscale_nhc[i] = 1;
133 for (i = 0; i < ir->opts.ngtc; i++)
135 vars->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
141 vars->vscale_nhc = NULL;
147 static void free_vetavars(t_vetavars *vars)
149 if (vars->vscale_nhc != NULL)
151 sfree(vars->vscale_nhc);
155 static int pcomp(const void *p1, const void *p2)
158 atom_id min1, min2, max1, max2;
159 t_sortblock *a1 = (t_sortblock *)p1;
160 t_sortblock *a2 = (t_sortblock *)p2;
162 db = a1->blocknr-a2->blocknr;
169 min1 = min(a1->iatom[1], a1->iatom[2]);
170 max1 = max(a1->iatom[1], a1->iatom[2]);
171 min2 = min(a2->iatom[1], a2->iatom[2]);
172 max2 = max(a2->iatom[1], a2->iatom[2]);
184 int n_flexible_constraints(struct gmx_constr *constr)
190 nflexcon = constr->nflexcon;
200 void too_many_constraint_warnings(int eConstrAlg, int warncount)
202 const char *abort = "- aborting to avoid logfile runaway.\n"
203 "This normally happens when your system is not sufficiently equilibrated,"
204 "or if you are changing lambda too fast in free energy simulations.\n";
207 "Too many %s warnings (%d)\n"
208 "If you know what you are doing you can %s"
209 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
210 "but normally it is better to fix the problem",
211 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
212 (eConstrAlg == econtLINCS) ?
213 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
216 static void write_constr_pdb(const char *fn, const char *title,
218 int start, int homenr, t_commrec *cr,
219 rvec x[], matrix box)
221 char fname[STRLEN], format[STRLEN];
223 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
228 if (DOMAINDECOMP(cr))
231 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
238 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
242 sprintf(fname, "%s.pdb", fn);
244 sprintf(format, "%s\n", get_pdbformat());
246 out = gmx_fio_fopen(fname, "w");
248 fprintf(out, "TITLE %s\n", title);
249 gmx_write_pdb_box(out, -1, box);
250 for (i = start; i < start+homenr; i++)
254 if (i >= dd->nat_home && i < dd_ac0)
258 ii = dd->gatindex[i];
264 gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
265 fprintf(out, format, "ATOM", (ii+1)%100000,
266 anm, resnm, ' ', resnr%10000, ' ',
267 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ]);
269 fprintf(out, "TER\n");
274 static void dump_confs(FILE *fplog, gmx_int64_t step, gmx_mtop_t *mtop,
275 int start, int homenr, t_commrec *cr,
276 rvec x[], rvec xprime[], matrix box)
278 char buf[256], buf2[22];
280 char *env = getenv("GMX_SUPPRESS_DUMP");
286 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
287 write_constr_pdb(buf, "initial coordinates",
288 mtop, start, homenr, cr, x, box);
289 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
290 write_constr_pdb(buf, "coordinates after constraining",
291 mtop, start, homenr, cr, xprime, box);
294 fprintf(fplog, "Wrote pdb files with previous and current coordinates\n");
296 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
299 static void pr_sortblock(FILE *fp, const char *title, int nsb, t_sortblock sb[])
303 fprintf(fp, "%s\n", title);
304 for (i = 0; (i < nsb); i++)
306 fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
307 i, sb[i].iatom[0], sb[i].iatom[1], sb[i].iatom[2],
312 gmx_bool constrain(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
313 struct gmx_constr *constr,
314 t_idef *idef, t_inputrec *ir, gmx_ekindata_t *ekind,
316 gmx_int64_t step, int delta_step,
318 rvec *x, rvec *xprime, rvec *min_proj,
319 gmx_bool bMolPBC, matrix box,
320 real lambda, real *dvdlambda,
321 rvec *v, tensor *vir,
322 t_nrnb *nrnb, int econq, gmx_bool bPscal,
323 real veta, real vetanew)
326 int start, homenr, nrend;
328 int ncons, settle_error;
331 real invdt, vir_fac, t;
334 t_pbc pbc, *pbc_null;
339 if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
341 gmx_incons("constrain called for forces displacements while not doing energy minimization, can not do this while the LINCS and SETTLE constraint connection matrices are mass weighted");
349 nrend = start+homenr;
351 /* set constants for pressure control integration */
352 init_vetavars(&vetavar, econq != econqCoord,
353 veta, vetanew, ir, ekind, bPscal);
355 if (ir->delta_t == 0)
361 invdt = 1/ir->delta_t;
364 if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
366 /* Set the constraint lengths for the step at which this configuration
367 * is meant to be. The invmasses should not be changed.
369 lambda += delta_step*ir->fepvals->delta_lambda;
374 clear_mat(vir_r_m_dr);
379 settle = &idef->il[F_SETTLE];
380 nsettle = settle->nr/(1+NRAL(F_SETTLE));
384 nth = gmx_omp_nthreads_get(emntSETTLE);
391 if (nth > 1 && constr->vir_r_m_dr_th == NULL)
393 snew(constr->vir_r_m_dr_th, nth);
394 snew(constr->settle_error, nth);
399 /* We do not need full pbc when constraints do not cross charge groups,
400 * i.e. when dd->constraint_comm==NULL.
401 * Note that PBC for constraints is different from PBC for bondeds.
402 * For constraints there is both forward and backward communication.
404 if (ir->ePBC != epbcNONE &&
405 (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm == NULL))
407 /* With pbc=screw the screw has been changed to a shift
408 * by the constraint coordinate communication routine,
409 * so that here we can use normal pbc.
411 pbc_null = set_pbc_dd(&pbc, ir->ePBC, cr->dd, FALSE, box);
418 /* Communicate the coordinates required for the non-local constraints
419 * for LINCS and/or SETTLE.
423 dd_move_x_constraints(cr->dd, box, x, xprime, econq == econqCoord);
426 if (constr->lincsd != NULL)
428 bOK = constrain_lincs(fplog, bLog, bEner, ir, step, constr->lincsd, md, cr,
430 box, pbc_null, lambda, dvdlambda,
431 invdt, v, vir != NULL, vir_r_m_dr,
433 constr->maxwarn, &constr->warncount_lincs);
434 if (!bOK && constr->maxwarn >= 0)
438 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
439 econstr_names[econtLINCS], gmx_step_str(step, buf));
445 if (constr->nblocks > 0)
450 bOK = bshakef(fplog, constr->shaked,
451 md->invmass, constr->nblocks, constr->sblock,
452 idef, ir, x, xprime, nrnb,
453 constr->lagr, lambda, dvdlambda,
454 invdt, v, vir != NULL, vir_r_m_dr,
455 constr->maxwarn >= 0, econq, &vetavar);
458 bOK = bshakef(fplog, constr->shaked,
459 md->invmass, constr->nblocks, constr->sblock,
460 idef, ir, x, min_proj, nrnb,
461 constr->lagr, lambda, dvdlambda,
462 invdt, NULL, vir != NULL, vir_r_m_dr,
463 constr->maxwarn >= 0, econq, &vetavar);
466 gmx_fatal(FARGS, "Internal error, SHAKE called for constraining something else than coordinates");
470 if (!bOK && constr->maxwarn >= 0)
474 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
475 econstr_names[econtSHAKE], gmx_step_str(step, buf));
483 int calcvir_atom_end;
487 calcvir_atom_end = 0;
491 calcvir_atom_end = md->homenr;
497 #pragma omp parallel for num_threads(nth) schedule(static)
498 for (th = 0; th < nth; th++)
500 int start_th, end_th;
504 clear_mat(constr->vir_r_m_dr_th[th]);
507 start_th = (nsettle* th )/nth;
508 end_th = (nsettle*(th+1))/nth;
509 if (start_th >= 0 && end_th - start_th > 0)
511 csettle(constr->settled,
513 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
516 invdt, v ? v[0] : NULL, calcvir_atom_end,
517 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
518 th == 0 ? &settle_error : &constr->settle_error[th],
522 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
525 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
529 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
535 case econqForceDispl:
536 #pragma omp parallel for num_threads(nth) schedule(static)
537 for (th = 0; th < nth; th++)
539 int start_th, end_th;
543 clear_mat(constr->vir_r_m_dr_th[th]);
546 start_th = (nsettle* th )/nth;
547 end_th = (nsettle*(th+1))/nth;
549 if (start_th >= 0 && end_th - start_th > 0)
551 settle_proj(constr->settled, econq,
553 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
556 xprime, min_proj, calcvir_atom_end,
557 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
561 /* This is an overestimate */
562 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
564 case econqDeriv_FlexCon:
565 /* Nothing to do, since the are no flexible constraints in settles */
568 gmx_incons("Unknown constraint quantity for settle");
574 /* Combine virial and error info of the other threads */
575 for (i = 1; i < nth; i++)
577 m_add(vir_r_m_dr, constr->vir_r_m_dr_th[i], vir_r_m_dr);
578 settle_error = constr->settle_error[i];
581 if (econq == econqCoord && settle_error >= 0)
584 if (constr->maxwarn >= 0)
588 "\nstep " "%"GMX_PRId64 ": Water molecule starting at atom %d can not be "
589 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
590 step, ddglatnr(cr->dd, settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
593 fprintf(fplog, "%s", buf);
595 fprintf(stderr, "%s", buf);
596 constr->warncount_settle++;
597 if (constr->warncount_settle > constr->maxwarn)
599 too_many_constraint_warnings(-1, constr->warncount_settle);
606 free_vetavars(&vetavar);
613 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
616 vir_fac = 0.5/ir->delta_t;
619 case econqForceDispl:
624 gmx_incons("Unsupported constraint quantity for virial");
629 vir_fac *= 2; /* only constraining over half the distance here */
631 for (i = 0; i < DIM; i++)
633 for (j = 0; j < DIM; j++)
635 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
642 dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
645 if (econq == econqCoord)
647 if (ir->ePull == epullCONSTRAINT)
649 if (EI_DYNAMICS(ir->eI))
651 t = ir->init_t + (step + delta_step)*ir->delta_t;
657 set_pbc(&pbc, ir->ePBC, box);
658 pull_constraint(ir->pull, md, &pbc, cr, ir->delta_t, t, x, xprime, v, *vir);
660 if (constr->ed && delta_step > 0)
662 /* apply the essential dynamcs constraints here */
663 do_edsam(ir, step, cr, xprime, v, box, constr->ed);
670 real *constr_rmsd_data(struct gmx_constr *constr)
674 return lincs_rmsd_data(constr->lincsd);
682 real constr_rmsd(struct gmx_constr *constr, gmx_bool bSD2)
686 return lincs_rmsd(constr->lincsd, bSD2);
694 static void make_shake_sblock_serial(struct gmx_constr *constr,
695 t_idef *idef, t_mdatoms *md)
704 /* Since we are processing the local topology,
705 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
707 ncons = idef->il[F_CONSTR].nr/3;
709 init_blocka(&sblocks);
710 gen_sblocks(NULL, 0, md->homenr, idef, &sblocks, FALSE);
713 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
714 nblocks=blocks->multinr[idef->nodeid] - bstart;
717 constr->nblocks = sblocks.nr;
720 fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n",
721 ncons, bstart, constr->nblocks);
724 /* Calculate block number for each atom */
725 inv_sblock = make_invblocka(&sblocks, md->nr);
727 done_blocka(&sblocks);
729 /* Store the block number in temp array and
730 * sort the constraints in order of the sblock number
731 * and the atom numbers, really sorting a segment of the array!
734 pr_idef(fplog, 0, "Before Sort", idef);
736 iatom = idef->il[F_CONSTR].iatoms;
738 for (i = 0; (i < ncons); i++, iatom += 3)
740 for (m = 0; (m < 3); m++)
742 sb[i].iatom[m] = iatom[m];
744 sb[i].blocknr = inv_sblock[iatom[1]];
747 /* Now sort the blocks */
750 pr_sortblock(debug, "Before sorting", ncons, sb);
751 fprintf(debug, "Going to sort constraints\n");
754 qsort(sb, ncons, (size_t)sizeof(*sb), pcomp);
758 pr_sortblock(debug, "After sorting", ncons, sb);
761 iatom = idef->il[F_CONSTR].iatoms;
762 for (i = 0; (i < ncons); i++, iatom += 3)
764 for (m = 0; (m < 3); m++)
766 iatom[m] = sb[i].iatom[m];
770 pr_idef(fplog, 0, "After Sort", idef);
774 snew(constr->sblock, constr->nblocks+1);
776 for (i = 0; (i < ncons); i++)
778 if (sb[i].blocknr != bnr)
781 constr->sblock[j++] = 3*i;
785 constr->sblock[j++] = 3*ncons;
787 if (j != (constr->nblocks+1))
789 fprintf(stderr, "bstart: %d\n", bstart);
790 fprintf(stderr, "j: %d, nblocks: %d, ncons: %d\n",
791 j, constr->nblocks, ncons);
792 for (i = 0; (i < ncons); i++)
794 fprintf(stderr, "i: %5d sb[i].blocknr: %5u\n", i, sb[i].blocknr);
796 for (j = 0; (j <= constr->nblocks); j++)
798 fprintf(stderr, "sblock[%3d]=%5d\n", j, (int)constr->sblock[j]);
800 gmx_fatal(FARGS, "DEATH HORROR: "
801 "sblocks does not match idef->il[F_CONSTR]");
807 static void make_shake_sblock_dd(struct gmx_constr *constr,
808 t_ilist *ilcon, t_block *cgs,
814 if (dd->ncg_home+1 > constr->sblock_nalloc)
816 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
817 srenew(constr->sblock, constr->sblock_nalloc);
821 iatom = ilcon->iatoms;
824 for (c = 0; c < ncons; c++)
826 if (c == 0 || iatom[1] >= cgs->index[cg+1])
828 constr->sblock[constr->nblocks++] = 3*c;
829 while (iatom[1] >= cgs->index[cg+1])
836 constr->sblock[constr->nblocks] = 3*ncons;
839 t_blocka make_at2con(int start, int natoms,
840 t_ilist *ilist, t_iparams *iparams,
841 gmx_bool bDynamics, int *nflexiblecons)
843 int *count, ncon, con, con_tot, nflexcon, ftype, i, a;
850 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
852 ncon = ilist[ftype].nr/3;
853 ia = ilist[ftype].iatoms;
854 for (con = 0; con < ncon; con++)
856 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
857 iparams[ia[0]].constr.dB == 0);
862 if (bDynamics || !bFlexCon)
864 for (i = 1; i < 3; i++)
873 *nflexiblecons = nflexcon;
876 at2con.nalloc_index = at2con.nr+1;
877 snew(at2con.index, at2con.nalloc_index);
879 for (a = 0; a < natoms; a++)
881 at2con.index[a+1] = at2con.index[a] + count[a];
884 at2con.nra = at2con.index[natoms];
885 at2con.nalloc_a = at2con.nra;
886 snew(at2con.a, at2con.nalloc_a);
888 /* The F_CONSTRNC constraints have constraint numbers
889 * that continue after the last F_CONSTR constraint.
892 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
894 ncon = ilist[ftype].nr/3;
895 ia = ilist[ftype].iatoms;
896 for (con = 0; con < ncon; con++)
898 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
899 iparams[ia[0]].constr.dB == 0);
900 if (bDynamics || !bFlexCon)
902 for (i = 1; i < 3; i++)
905 at2con.a[at2con.index[a]+count[a]++] = con_tot;
918 static int *make_at2settle(int natoms, const t_ilist *ilist)
924 /* Set all to no settle */
925 for (a = 0; a < natoms; a++)
930 stride = 1 + NRAL(F_SETTLE);
932 for (s = 0; s < ilist->nr; s += stride)
934 at2s[ilist->iatoms[s+1]] = s/stride;
935 at2s[ilist->iatoms[s+2]] = s/stride;
936 at2s[ilist->iatoms[s+3]] = s/stride;
942 void set_constraints(struct gmx_constr *constr,
943 gmx_localtop_t *top, t_inputrec *ir,
944 t_mdatoms *md, t_commrec *cr)
953 if (constr->ncon_tot > 0)
955 /* We are using the local topology,
956 * so there are only F_CONSTR constraints.
958 ncons = idef->il[F_CONSTR].nr/3;
960 /* With DD we might also need to call LINCS with ncons=0 for
961 * communicating coordinates to other nodes that do have constraints.
963 if (ir->eConstrAlg == econtLINCS)
965 set_lincs(idef, md, EI_DYNAMICS(ir->eI), cr, constr->lincsd);
967 if (ir->eConstrAlg == econtSHAKE)
971 make_shake_sblock_dd(constr, &idef->il[F_CONSTR], &top->cgs, cr->dd);
975 make_shake_sblock_serial(constr, idef, md);
977 if (ncons > constr->lagr_nalloc)
979 constr->lagr_nalloc = over_alloc_dd(ncons);
980 srenew(constr->lagr, constr->lagr_nalloc);
985 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
987 settle = &idef->il[F_SETTLE];
988 iO = settle->iatoms[1];
989 iH = settle->iatoms[2];
991 settle_init(md->massT[iO], md->massT[iH],
992 md->invmass[iO], md->invmass[iH],
993 idef->iparams[settle->iatoms[0]].settle.doh,
994 idef->iparams[settle->iatoms[0]].settle.dhh);
997 /* Make a selection of the local atoms for essential dynamics */
998 if (constr->ed && cr->dd)
1000 dd_make_local_ed_indices(cr->dd, constr->ed);
1004 static void constr_recur(t_blocka *at2con,
1005 t_ilist *ilist, t_iparams *iparams, gmx_bool bTopB,
1006 int at, int depth, int nc, int *path,
1007 real r0, real r1, real *r2max,
1019 ncon1 = ilist[F_CONSTR].nr/3;
1020 ia1 = ilist[F_CONSTR].iatoms;
1021 ia2 = ilist[F_CONSTRNC].iatoms;
1023 /* Loop over all constraints connected to this atom */
1024 for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
1027 /* Do not walk over already used constraints */
1029 for (a1 = 0; a1 < depth; a1++)
1031 if (con == path[a1])
1038 ia = constr_iatomptr(ncon1, ia1, ia2, con);
1039 /* Flexible constraints currently have length 0, which is incorrect */
1042 len = iparams[ia[0]].constr.dA;
1046 len = iparams[ia[0]].constr.dB;
1048 /* In the worst case the bond directions alternate */
1059 /* Assume angles of 120 degrees between all bonds */
1060 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
1062 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
1065 fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
1066 for (a1 = 0; a1 < depth; a1++)
1068 fprintf(debug, " %d %5.3f",
1070 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
1072 fprintf(debug, " %d %5.3f\n", con, len);
1075 /* Limit the number of recursions to 1000*nc,
1076 * so a call does not take more than a second,
1077 * even for highly connected systems.
1079 if (depth + 1 < nc && *count < 1000*nc)
1091 constr_recur(at2con, ilist, iparams,
1092 bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
1099 static real constr_r_max_moltype(gmx_moltype_t *molt, t_iparams *iparams,
1102 int natoms, nflexcon, *path, at, count;
1105 real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
1107 if (molt->ilist[F_CONSTR].nr == 0 &&
1108 molt->ilist[F_CONSTRNC].nr == 0)
1113 natoms = molt->atoms.nr;
1115 at2con = make_at2con(0, natoms, molt->ilist, iparams,
1116 EI_DYNAMICS(ir->eI), &nflexcon);
1117 snew(path, 1+ir->nProjOrder);
1118 for (at = 0; at < 1+ir->nProjOrder; at++)
1124 for (at = 0; at < natoms; at++)
1130 constr_recur(&at2con, molt->ilist, iparams,
1131 FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
1133 if (ir->efep == efepNO)
1135 rmax = sqrt(r2maxA);
1140 for (at = 0; at < natoms; at++)
1145 constr_recur(&at2con, molt->ilist, iparams,
1146 TRUE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxB, &count);
1148 lam0 = ir->fepvals->init_lambda;
1149 if (EI_DYNAMICS(ir->eI))
1151 lam0 += ir->init_step*ir->fepvals->delta_lambda;
1153 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1154 if (EI_DYNAMICS(ir->eI))
1156 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1157 rmax = max(rmax, (1 - lam1)*sqrt(r2maxA) + lam1*sqrt(r2maxB));
1161 done_blocka(&at2con);
1167 real constr_r_max(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir)
1173 for (mt = 0; mt < mtop->nmoltype; mt++)
1176 constr_r_max_moltype(&mtop->moltype[mt],
1177 mtop->ffparams.iparams, ir));
1182 fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
1188 gmx_constr_t init_constraints(FILE *fplog,
1189 gmx_mtop_t *mtop, t_inputrec *ir,
1190 gmx_edsam_t ed, t_state *state,
1193 int ncon, nset, nmol, settle_type, i, natoms, mt, nflexcon;
1194 struct gmx_constr *constr;
1197 gmx_mtop_ilistloop_t iloop;
1200 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1201 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1202 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
1204 if (ncon+nset == 0 && ir->ePull != epullCONSTRAINT && ed == NULL)
1211 constr->ncon_tot = ncon;
1212 constr->nflexcon = 0;
1215 constr->n_at2con_mt = mtop->nmoltype;
1216 snew(constr->at2con_mt, constr->n_at2con_mt);
1217 for (mt = 0; mt < mtop->nmoltype; mt++)
1219 constr->at2con_mt[mt] = make_at2con(0, mtop->moltype[mt].atoms.nr,
1220 mtop->moltype[mt].ilist,
1221 mtop->ffparams.iparams,
1222 EI_DYNAMICS(ir->eI), &nflexcon);
1223 for (i = 0; i < mtop->nmolblock; i++)
1225 if (mtop->molblock[i].type == mt)
1227 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1232 if (constr->nflexcon > 0)
1236 fprintf(fplog, "There are %d flexible constraints\n",
1238 if (ir->fc_stepsize == 0)
1241 "WARNING: step size for flexible constraining = 0\n"
1242 " All flexible constraints will be rigid.\n"
1243 " Will try to keep all flexible constraints at their original length,\n"
1244 " but the lengths may exhibit some drift.\n\n");
1245 constr->nflexcon = 0;
1248 if (constr->nflexcon > 0)
1250 please_cite(fplog, "Hess2002");
1254 if (ir->eConstrAlg == econtLINCS)
1256 constr->lincsd = init_lincs(fplog, mtop,
1257 constr->nflexcon, constr->at2con_mt,
1258 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1259 ir->nLincsIter, ir->nProjOrder);
1262 if (ir->eConstrAlg == econtSHAKE)
1264 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1266 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1268 if (constr->nflexcon)
1270 gmx_fatal(FARGS, "For this system also velocities and/or forces need to be constrained, this can not be done with SHAKE, you should select LINCS");
1272 please_cite(fplog, "Ryckaert77a");
1275 please_cite(fplog, "Barth95a");
1278 constr->shaked = shake_init();
1284 please_cite(fplog, "Miyamoto92a");
1286 constr->bInterCGsettles = inter_charge_group_settles(mtop);
1288 /* Check that we have only one settle type */
1290 iloop = gmx_mtop_ilistloop_init(mtop);
1291 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
1293 for (i = 0; i < ilist[F_SETTLE].nr; i += 4)
1295 if (settle_type == -1)
1297 settle_type = ilist[F_SETTLE].iatoms[i];
1299 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1302 "The [molecules] section of your topology specifies more than one block of\n"
1303 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1304 "are trying to partition your solvent into different *groups* (e.g. for\n"
1305 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1306 "files specify groups. Otherwise, you may wish to change the least-used\n"
1307 "block of molecules with SETTLE constraints into 3 normal constraints.");
1312 constr->n_at2settle_mt = mtop->nmoltype;
1313 snew(constr->at2settle_mt, constr->n_at2settle_mt);
1314 for (mt = 0; mt < mtop->nmoltype; mt++)
1316 constr->at2settle_mt[mt] =
1317 make_at2settle(mtop->moltype[mt].atoms.nr,
1318 &mtop->moltype[mt].ilist[F_SETTLE]);
1322 constr->maxwarn = 999;
1323 env = getenv("GMX_MAXCONSTRWARN");
1326 constr->maxwarn = 0;
1327 sscanf(env, "%d", &constr->maxwarn);
1331 "Setting the maximum number of constraint warnings to %d\n",
1337 "Setting the maximum number of constraint warnings to %d\n",
1341 if (constr->maxwarn < 0 && fplog)
1343 fprintf(fplog, "maxwarn < 0, will not stop on constraint errors\n");
1345 constr->warncount_lincs = 0;
1346 constr->warncount_settle = 0;
1348 /* Initialize the essential dynamics sampling.
1349 * Put the pointer to the ED struct in constr */
1351 if (ed != NULL || state->edsamstate.nED > 0)
1353 init_edsam(mtop, ir, cr, ed, state->x, state->box, &state->edsamstate);
1356 constr->warn_mtop = mtop;
1361 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1363 return constr->at2con_mt;
1366 const int **atom2settle_moltype(gmx_constr_t constr)
1368 return (const int **)constr->at2settle_mt;
1372 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1374 const gmx_moltype_t *molt;
1378 int nat, *at2cg, cg, a, ftype, i;
1382 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1384 molt = &mtop->moltype[mtop->molblock[mb].type];
1386 if (molt->ilist[F_CONSTR].nr > 0 ||
1387 molt->ilist[F_CONSTRNC].nr > 0 ||
1388 molt->ilist[F_SETTLE].nr > 0)
1391 snew(at2cg, molt->atoms.nr);
1392 for (cg = 0; cg < cgs->nr; cg++)
1394 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1400 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1402 il = &molt->ilist[ftype];
1403 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(ftype))
1405 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1419 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1421 const gmx_moltype_t *molt;
1425 int nat, *at2cg, cg, a, ftype, i;
1429 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1431 molt = &mtop->moltype[mtop->molblock[mb].type];
1433 if (molt->ilist[F_SETTLE].nr > 0)
1436 snew(at2cg, molt->atoms.nr);
1437 for (cg = 0; cg < cgs->nr; cg++)
1439 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1445 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1447 il = &molt->ilist[ftype];
1448 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(F_SETTLE))
1450 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1451 at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1465 /* helper functions for andersen temperature control, because the
1466 * gmx_constr construct is only defined in constr.c. Return the list
1467 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1469 extern int *get_sblock(struct gmx_constr *constr)
1471 return constr->sblock;
1474 extern int get_nblocks(struct gmx_constr *constr)
1476 return constr->nblocks;