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"
54 #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"
65 typedef struct gmx_constr {
66 int ncon_tot; /* The total number of constraints */
67 int nflexcon; /* The number of flexible constraints */
68 int n_at2con_mt; /* The size of at2con = #moltypes */
69 t_blocka *at2con_mt; /* A list of atoms to constraints */
70 int n_at2settle_mt; /* The size of at2settle = #moltypes */
71 int **at2settle_mt; /* A list of atoms to settles */
72 gmx_bool bInterCGsettles;
73 gmx_lincsdata_t lincsd; /* LINCS data */
74 gmx_shakedata_t shaked; /* SHAKE data */
75 gmx_settledata_t settled; /* SETTLE data */
76 int nblocks; /* The number of SHAKE blocks */
77 int *sblock; /* The SHAKE blocks */
78 int sblock_nalloc; /* The allocation size of sblock */
79 real *lagr; /* Lagrange multipliers for SHAKE */
80 int lagr_nalloc; /* The allocation size of lagr */
81 int maxwarn; /* The maximum number of warnings */
84 gmx_edsam_t ed; /* The essential dynamics data */
86 tensor *vir_r_m_dr_th; /* Thread local working data */
87 int *settle_error; /* Thread local working data */
89 gmx_mtop_t *warn_mtop; /* Only used for printing warnings */
97 static void *init_vetavars(t_vetavars *vars,
98 gmx_bool constr_deriv,
99 real veta, real vetanew, t_inputrec *ir, gmx_ekindata_t *ekind, gmx_bool bPscal)
104 /* first, set the alpha integrator variable */
105 if ((ir->opts.nrdf[0] > 0) && bPscal)
107 vars->alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
113 g = 0.5*veta*ir->delta_t;
114 vars->rscale = exp(g)*series_sinhx(g);
115 g = -0.25*vars->alpha*veta*ir->delta_t;
116 vars->vscale = exp(g)*series_sinhx(g);
117 vars->rvscale = vars->vscale*vars->rscale;
118 vars->veta = vetanew;
122 snew(vars->vscale_nhc, ir->opts.ngtc);
123 if ((ekind == NULL) || (!bPscal))
125 for (i = 0; i < ir->opts.ngtc; i++)
127 vars->vscale_nhc[i] = 1;
132 for (i = 0; i < ir->opts.ngtc; i++)
134 vars->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
140 vars->vscale_nhc = NULL;
146 static void free_vetavars(t_vetavars *vars)
148 if (vars->vscale_nhc != NULL)
150 sfree(vars->vscale_nhc);
154 static int pcomp(const void *p1, const void *p2)
157 atom_id min1, min2, max1, max2;
158 t_sortblock *a1 = (t_sortblock *)p1;
159 t_sortblock *a2 = (t_sortblock *)p2;
161 db = a1->blocknr-a2->blocknr;
168 min1 = min(a1->iatom[1], a1->iatom[2]);
169 max1 = max(a1->iatom[1], a1->iatom[2]);
170 min2 = min(a2->iatom[1], a2->iatom[2]);
171 max2 = max(a2->iatom[1], a2->iatom[2]);
183 int n_flexible_constraints(struct gmx_constr *constr)
189 nflexcon = constr->nflexcon;
199 void too_many_constraint_warnings(int eConstrAlg, int warncount)
201 const char *abort = "- aborting to avoid logfile runaway.\n"
202 "This normally happens when your system is not sufficiently equilibrated,"
203 "or if you are changing lambda too fast in free energy simulations.\n";
206 "Too many %s warnings (%d)\n"
207 "If you know what you are doing you can %s"
208 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
209 "but normally it is better to fix the problem",
210 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
211 (eConstrAlg == econtLINCS) ?
212 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
215 static void write_constr_pdb(const char *fn, const char *title,
217 int start, int homenr, t_commrec *cr,
218 rvec x[], matrix box)
220 char fname[STRLEN], format[STRLEN];
222 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
229 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
230 if (DOMAINDECOMP(cr))
233 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
240 sprintf(fname, "%s.pdb", fn);
242 sprintf(format, "%s\n", get_pdbformat());
244 out = gmx_fio_fopen(fname, "w");
246 fprintf(out, "TITLE %s\n", title);
247 gmx_write_pdb_box(out, -1, box);
248 for (i = start; i < start+homenr; i++)
252 if (i >= dd->nat_home && i < dd_ac0)
256 ii = dd->gatindex[i];
262 gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
263 fprintf(out, format, "ATOM", (ii+1)%100000,
264 anm, resnm, ' ', resnr%10000, ' ',
265 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ]);
267 fprintf(out, "TER\n");
272 static void dump_confs(FILE *fplog, gmx_int64_t step, gmx_mtop_t *mtop,
273 int start, int homenr, t_commrec *cr,
274 rvec x[], rvec xprime[], matrix box)
276 char buf[256], buf2[22];
278 char *env = getenv("GMX_SUPPRESS_DUMP");
284 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
285 write_constr_pdb(buf, "initial coordinates",
286 mtop, start, homenr, cr, x, box);
287 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
288 write_constr_pdb(buf, "coordinates after constraining",
289 mtop, start, homenr, cr, xprime, box);
292 fprintf(fplog, "Wrote pdb files with previous and current coordinates\n");
294 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
297 static void pr_sortblock(FILE *fp, const char *title, int nsb, t_sortblock sb[])
301 fprintf(fp, "%s\n", title);
302 for (i = 0; (i < nsb); i++)
304 fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
305 i, sb[i].iatom[0], sb[i].iatom[1], sb[i].iatom[2],
310 gmx_bool constrain(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
311 struct gmx_constr *constr,
312 t_idef *idef, t_inputrec *ir, gmx_ekindata_t *ekind,
314 gmx_int64_t step, int delta_step,
316 rvec *x, rvec *xprime, rvec *min_proj,
317 gmx_bool bMolPBC, matrix box,
318 real lambda, real *dvdlambda,
319 rvec *v, tensor *vir,
320 t_nrnb *nrnb, int econq, gmx_bool bPscal,
321 real veta, real vetanew)
324 int start, homenr, nrend;
326 int ncons, settle_error;
329 real invdt, vir_fac, t;
332 t_pbc pbc, *pbc_null;
337 if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
339 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");
347 nrend = start+homenr;
349 /* set constants for pressure control integration */
350 init_vetavars(&vetavar, econq != econqCoord,
351 veta, vetanew, ir, ekind, bPscal);
353 if (ir->delta_t == 0)
359 invdt = 1/ir->delta_t;
362 if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
364 /* Set the constraint lengths for the step at which this configuration
365 * is meant to be. The invmasses should not be changed.
367 lambda += delta_step*ir->fepvals->delta_lambda;
372 clear_mat(vir_r_m_dr);
377 settle = &idef->il[F_SETTLE];
378 nsettle = settle->nr/(1+NRAL(F_SETTLE));
382 nth = gmx_omp_nthreads_get(emntSETTLE);
389 if (nth > 1 && constr->vir_r_m_dr_th == NULL)
391 snew(constr->vir_r_m_dr_th, nth);
392 snew(constr->settle_error, nth);
397 /* We do not need full pbc when constraints do not cross charge groups,
398 * i.e. when dd->constraint_comm==NULL.
399 * Note that PBC for constraints is different from PBC for bondeds.
400 * For constraints there is both forward and backward communication.
402 if (ir->ePBC != epbcNONE &&
403 (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm == NULL))
405 /* With pbc=screw the screw has been changed to a shift
406 * by the constraint coordinate communication routine,
407 * so that here we can use normal pbc.
409 pbc_null = set_pbc_dd(&pbc, ir->ePBC, cr->dd, FALSE, box);
416 /* Communicate the coordinates required for the non-local constraints
417 * for LINCS and/or SETTLE.
421 dd_move_x_constraints(cr->dd, box, x, xprime);
423 else if (PARTDECOMP(cr))
425 pd_move_x_constraints(cr, x, xprime);
428 if (constr->lincsd != NULL)
430 bOK = constrain_lincs(fplog, bLog, bEner, ir, step, constr->lincsd, md, cr,
432 box, pbc_null, lambda, dvdlambda,
433 invdt, v, vir != NULL, vir_r_m_dr,
435 constr->maxwarn, &constr->warncount_lincs);
436 if (!bOK && constr->maxwarn >= 0)
440 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
441 econstr_names[econtLINCS], gmx_step_str(step, buf));
447 if (constr->nblocks > 0)
452 bOK = bshakef(fplog, constr->shaked,
453 md->invmass, constr->nblocks, constr->sblock,
454 idef, ir, x, xprime, nrnb,
455 constr->lagr, lambda, dvdlambda,
456 invdt, v, vir != NULL, vir_r_m_dr,
457 constr->maxwarn >= 0, econq, &vetavar);
460 bOK = bshakef(fplog, constr->shaked,
461 md->invmass, constr->nblocks, constr->sblock,
462 idef, ir, x, min_proj, nrnb,
463 constr->lagr, lambda, dvdlambda,
464 invdt, NULL, vir != NULL, vir_r_m_dr,
465 constr->maxwarn >= 0, econq, &vetavar);
468 gmx_fatal(FARGS, "Internal error, SHAKE called for constraining something else than coordinates");
472 if (!bOK && constr->maxwarn >= 0)
476 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
477 econstr_names[econtSHAKE], gmx_step_str(step, buf));
485 int calcvir_atom_end;
489 calcvir_atom_end = 0;
493 calcvir_atom_end = md->start + md->homenr;
499 #pragma omp parallel for num_threads(nth) schedule(static)
500 for (th = 0; th < nth; th++)
502 int start_th, end_th;
506 clear_mat(constr->vir_r_m_dr_th[th]);
509 start_th = (nsettle* th )/nth;
510 end_th = (nsettle*(th+1))/nth;
511 if (start_th >= 0 && end_th - start_th > 0)
513 csettle(constr->settled,
515 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
518 invdt, v ? v[0] : NULL, calcvir_atom_end,
519 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
520 th == 0 ? &settle_error : &constr->settle_error[th],
524 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
527 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
531 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
537 case econqForceDispl:
538 #pragma omp parallel for num_threads(nth) schedule(static)
539 for (th = 0; th < nth; th++)
541 int start_th, end_th;
545 clear_mat(constr->vir_r_m_dr_th[th]);
548 start_th = (nsettle* th )/nth;
549 end_th = (nsettle*(th+1))/nth;
551 if (start_th >= 0 && end_th - start_th > 0)
553 settle_proj(constr->settled, econq,
555 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
558 xprime, min_proj, calcvir_atom_end,
559 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
563 /* This is an overestimate */
564 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
566 case econqDeriv_FlexCon:
567 /* Nothing to do, since the are no flexible constraints in settles */
570 gmx_incons("Unknown constraint quantity for settle");
576 /* Combine virial and error info of the other threads */
577 for (i = 1; i < nth; i++)
579 m_add(vir_r_m_dr, constr->vir_r_m_dr_th[i], vir_r_m_dr);
580 settle_error = constr->settle_error[i];
583 if (econq == econqCoord && settle_error >= 0)
586 if (constr->maxwarn >= 0)
590 "\nstep " "%"GMX_PRId64 ": Water molecule starting at atom %d can not be "
591 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
592 step, ddglatnr(cr->dd, settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
595 fprintf(fplog, "%s", buf);
597 fprintf(stderr, "%s", buf);
598 constr->warncount_settle++;
599 if (constr->warncount_settle > constr->maxwarn)
601 too_many_constraint_warnings(-1, constr->warncount_settle);
608 free_vetavars(&vetavar);
615 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
618 vir_fac = 0.5/ir->delta_t;
621 case econqForceDispl:
626 gmx_incons("Unsupported constraint quantity for virial");
631 vir_fac *= 2; /* only constraining over half the distance here */
633 for (i = 0; i < DIM; i++)
635 for (j = 0; j < DIM; j++)
637 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
644 dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
647 if (econq == econqCoord)
649 if (ir->ePull == epullCONSTRAINT)
651 if (EI_DYNAMICS(ir->eI))
653 t = ir->init_t + (step + delta_step)*ir->delta_t;
659 set_pbc(&pbc, ir->ePBC, box);
660 pull_constraint(ir->pull, md, &pbc, cr, ir->delta_t, t, x, xprime, v, *vir);
662 if (constr->ed && delta_step > 0)
664 /* apply the essential dynamcs constraints here */
665 do_edsam(ir, step, cr, xprime, v, box, constr->ed);
672 real *constr_rmsd_data(struct gmx_constr *constr)
676 return lincs_rmsd_data(constr->lincsd);
684 real constr_rmsd(struct gmx_constr *constr, gmx_bool bSD2)
688 return lincs_rmsd(constr->lincsd, bSD2);
696 static void make_shake_sblock_pd(struct gmx_constr *constr,
697 t_idef *idef, t_mdatoms *md)
706 /* Since we are processing the local topology,
707 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
709 ncons = idef->il[F_CONSTR].nr/3;
711 init_blocka(&sblocks);
712 gen_sblocks(NULL, md->start, md->start+md->homenr, idef, &sblocks, FALSE);
715 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
716 nblocks=blocks->multinr[idef->nodeid] - bstart;
719 constr->nblocks = sblocks.nr;
722 fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n",
723 ncons, bstart, constr->nblocks);
726 /* Calculate block number for each atom */
727 inv_sblock = make_invblocka(&sblocks, md->nr);
729 done_blocka(&sblocks);
731 /* Store the block number in temp array and
732 * sort the constraints in order of the sblock number
733 * and the atom numbers, really sorting a segment of the array!
736 pr_idef(fplog, 0, "Before Sort", idef);
738 iatom = idef->il[F_CONSTR].iatoms;
740 for (i = 0; (i < ncons); i++, iatom += 3)
742 for (m = 0; (m < 3); m++)
744 sb[i].iatom[m] = iatom[m];
746 sb[i].blocknr = inv_sblock[iatom[1]];
749 /* Now sort the blocks */
752 pr_sortblock(debug, "Before sorting", ncons, sb);
753 fprintf(debug, "Going to sort constraints\n");
756 qsort(sb, ncons, (size_t)sizeof(*sb), pcomp);
760 pr_sortblock(debug, "After sorting", ncons, sb);
763 iatom = idef->il[F_CONSTR].iatoms;
764 for (i = 0; (i < ncons); i++, iatom += 3)
766 for (m = 0; (m < 3); m++)
768 iatom[m] = sb[i].iatom[m];
772 pr_idef(fplog, 0, "After Sort", idef);
776 snew(constr->sblock, constr->nblocks+1);
778 for (i = 0; (i < ncons); i++)
780 if (sb[i].blocknr != bnr)
783 constr->sblock[j++] = 3*i;
787 constr->sblock[j++] = 3*ncons;
789 if (j != (constr->nblocks+1))
791 fprintf(stderr, "bstart: %d\n", bstart);
792 fprintf(stderr, "j: %d, nblocks: %d, ncons: %d\n",
793 j, constr->nblocks, ncons);
794 for (i = 0; (i < ncons); i++)
796 fprintf(stderr, "i: %5d sb[i].blocknr: %5u\n", i, sb[i].blocknr);
798 for (j = 0; (j <= constr->nblocks); j++)
800 fprintf(stderr, "sblock[%3d]=%5d\n", j, (int)constr->sblock[j]);
802 gmx_fatal(FARGS, "DEATH HORROR: "
803 "sblocks does not match idef->il[F_CONSTR]");
809 static void make_shake_sblock_dd(struct gmx_constr *constr,
810 t_ilist *ilcon, t_block *cgs,
816 if (dd->ncg_home+1 > constr->sblock_nalloc)
818 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
819 srenew(constr->sblock, constr->sblock_nalloc);
823 iatom = ilcon->iatoms;
826 for (c = 0; c < ncons; c++)
828 if (c == 0 || iatom[1] >= cgs->index[cg+1])
830 constr->sblock[constr->nblocks++] = 3*c;
831 while (iatom[1] >= cgs->index[cg+1])
838 constr->sblock[constr->nblocks] = 3*ncons;
841 t_blocka make_at2con(int start, int natoms,
842 t_ilist *ilist, t_iparams *iparams,
843 gmx_bool bDynamics, int *nflexiblecons)
845 int *count, ncon, con, con_tot, nflexcon, ftype, i, a;
852 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
854 ncon = ilist[ftype].nr/3;
855 ia = ilist[ftype].iatoms;
856 for (con = 0; con < ncon; con++)
858 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
859 iparams[ia[0]].constr.dB == 0);
864 if (bDynamics || !bFlexCon)
866 for (i = 1; i < 3; i++)
875 *nflexiblecons = nflexcon;
878 at2con.nalloc_index = at2con.nr+1;
879 snew(at2con.index, at2con.nalloc_index);
881 for (a = 0; a < natoms; a++)
883 at2con.index[a+1] = at2con.index[a] + count[a];
886 at2con.nra = at2con.index[natoms];
887 at2con.nalloc_a = at2con.nra;
888 snew(at2con.a, at2con.nalloc_a);
890 /* The F_CONSTRNC constraints have constraint numbers
891 * that continue after the last F_CONSTR constraint.
894 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
896 ncon = ilist[ftype].nr/3;
897 ia = ilist[ftype].iatoms;
898 for (con = 0; con < ncon; con++)
900 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
901 iparams[ia[0]].constr.dB == 0);
902 if (bDynamics || !bFlexCon)
904 for (i = 1; i < 3; i++)
907 at2con.a[at2con.index[a]+count[a]++] = con_tot;
920 static int *make_at2settle(int natoms, const t_ilist *ilist)
926 /* Set all to no settle */
927 for (a = 0; a < natoms; a++)
932 stride = 1 + NRAL(F_SETTLE);
934 for (s = 0; s < ilist->nr; s += stride)
936 at2s[ilist->iatoms[s+1]] = s/stride;
937 at2s[ilist->iatoms[s+2]] = s/stride;
938 at2s[ilist->iatoms[s+3]] = s/stride;
944 void set_constraints(struct gmx_constr *constr,
945 gmx_localtop_t *top, t_inputrec *ir,
946 t_mdatoms *md, t_commrec *cr)
955 if (constr->ncon_tot > 0)
957 /* We are using the local topology,
958 * so there are only F_CONSTR constraints.
960 ncons = idef->il[F_CONSTR].nr/3;
962 /* With DD we might also need to call LINCS with ncons=0 for
963 * communicating coordinates to other nodes that do have constraints.
965 if (ir->eConstrAlg == econtLINCS)
967 set_lincs(idef, md, EI_DYNAMICS(ir->eI), cr, constr->lincsd);
969 if (ir->eConstrAlg == econtSHAKE)
973 make_shake_sblock_dd(constr, &idef->il[F_CONSTR], &top->cgs, cr->dd);
977 make_shake_sblock_pd(constr, idef, md);
979 if (ncons > constr->lagr_nalloc)
981 constr->lagr_nalloc = over_alloc_dd(ncons);
982 srenew(constr->lagr, constr->lagr_nalloc);
987 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
989 settle = &idef->il[F_SETTLE];
990 iO = settle->iatoms[1];
991 iH = settle->iatoms[2];
993 settle_init(md->massT[iO], md->massT[iH],
994 md->invmass[iO], md->invmass[iH],
995 idef->iparams[settle->iatoms[0]].settle.doh,
996 idef->iparams[settle->iatoms[0]].settle.dhh);
999 /* Make a selection of the local atoms for essential dynamics */
1000 if (constr->ed && cr->dd)
1002 dd_make_local_ed_indices(cr->dd, constr->ed);
1006 static void constr_recur(t_blocka *at2con,
1007 t_ilist *ilist, t_iparams *iparams, gmx_bool bTopB,
1008 int at, int depth, int nc, int *path,
1009 real r0, real r1, real *r2max,
1021 ncon1 = ilist[F_CONSTR].nr/3;
1022 ia1 = ilist[F_CONSTR].iatoms;
1023 ia2 = ilist[F_CONSTRNC].iatoms;
1025 /* Loop over all constraints connected to this atom */
1026 for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
1029 /* Do not walk over already used constraints */
1031 for (a1 = 0; a1 < depth; a1++)
1033 if (con == path[a1])
1040 ia = constr_iatomptr(ncon1, ia1, ia2, con);
1041 /* Flexible constraints currently have length 0, which is incorrect */
1044 len = iparams[ia[0]].constr.dA;
1048 len = iparams[ia[0]].constr.dB;
1050 /* In the worst case the bond directions alternate */
1061 /* Assume angles of 120 degrees between all bonds */
1062 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
1064 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
1067 fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
1068 for (a1 = 0; a1 < depth; a1++)
1070 fprintf(debug, " %d %5.3f",
1072 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
1074 fprintf(debug, " %d %5.3f\n", con, len);
1077 /* Limit the number of recursions to 1000*nc,
1078 * so a call does not take more than a second,
1079 * even for highly connected systems.
1081 if (depth + 1 < nc && *count < 1000*nc)
1093 constr_recur(at2con, ilist, iparams,
1094 bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
1101 static real constr_r_max_moltype(gmx_moltype_t *molt, t_iparams *iparams,
1104 int natoms, nflexcon, *path, at, count;
1107 real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
1109 if (molt->ilist[F_CONSTR].nr == 0 &&
1110 molt->ilist[F_CONSTRNC].nr == 0)
1115 natoms = molt->atoms.nr;
1117 at2con = make_at2con(0, natoms, molt->ilist, iparams,
1118 EI_DYNAMICS(ir->eI), &nflexcon);
1119 snew(path, 1+ir->nProjOrder);
1120 for (at = 0; at < 1+ir->nProjOrder; at++)
1126 for (at = 0; at < natoms; at++)
1132 constr_recur(&at2con, molt->ilist, iparams,
1133 FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
1135 if (ir->efep == efepNO)
1137 rmax = sqrt(r2maxA);
1142 for (at = 0; at < natoms; at++)
1147 constr_recur(&at2con, molt->ilist, iparams,
1148 TRUE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxB, &count);
1150 lam0 = ir->fepvals->init_lambda;
1151 if (EI_DYNAMICS(ir->eI))
1153 lam0 += ir->init_step*ir->fepvals->delta_lambda;
1155 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1156 if (EI_DYNAMICS(ir->eI))
1158 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1159 rmax = max(rmax, (1 - lam1)*sqrt(r2maxA) + lam1*sqrt(r2maxB));
1163 done_blocka(&at2con);
1169 real constr_r_max(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir)
1175 for (mt = 0; mt < mtop->nmoltype; mt++)
1178 constr_r_max_moltype(&mtop->moltype[mt],
1179 mtop->ffparams.iparams, ir));
1184 fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
1190 gmx_constr_t init_constraints(FILE *fplog,
1191 gmx_mtop_t *mtop, t_inputrec *ir,
1192 gmx_edsam_t ed, t_state *state,
1195 int ncon, nset, nmol, settle_type, i, natoms, mt, nflexcon;
1196 struct gmx_constr *constr;
1199 gmx_mtop_ilistloop_t iloop;
1202 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1203 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1204 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
1206 if (ncon+nset == 0 && ir->ePull != epullCONSTRAINT && ed == NULL)
1213 constr->ncon_tot = ncon;
1214 constr->nflexcon = 0;
1217 constr->n_at2con_mt = mtop->nmoltype;
1218 snew(constr->at2con_mt, constr->n_at2con_mt);
1219 for (mt = 0; mt < mtop->nmoltype; mt++)
1221 constr->at2con_mt[mt] = make_at2con(0, mtop->moltype[mt].atoms.nr,
1222 mtop->moltype[mt].ilist,
1223 mtop->ffparams.iparams,
1224 EI_DYNAMICS(ir->eI), &nflexcon);
1225 for (i = 0; i < mtop->nmolblock; i++)
1227 if (mtop->molblock[i].type == mt)
1229 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1234 if (constr->nflexcon > 0)
1238 fprintf(fplog, "There are %d flexible constraints\n",
1240 if (ir->fc_stepsize == 0)
1243 "WARNING: step size for flexible constraining = 0\n"
1244 " All flexible constraints will be rigid.\n"
1245 " Will try to keep all flexible constraints at their original length,\n"
1246 " but the lengths may exhibit some drift.\n\n");
1247 constr->nflexcon = 0;
1250 if (constr->nflexcon > 0)
1252 please_cite(fplog, "Hess2002");
1256 if (ir->eConstrAlg == econtLINCS)
1258 constr->lincsd = init_lincs(fplog, mtop,
1259 constr->nflexcon, constr->at2con_mt,
1260 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1261 ir->nLincsIter, ir->nProjOrder);
1264 if (ir->eConstrAlg == econtSHAKE)
1266 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1268 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1270 if (constr->nflexcon)
1272 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");
1274 please_cite(fplog, "Ryckaert77a");
1277 please_cite(fplog, "Barth95a");
1280 constr->shaked = shake_init();
1286 please_cite(fplog, "Miyamoto92a");
1288 constr->bInterCGsettles = inter_charge_group_settles(mtop);
1290 /* Check that we have only one settle type */
1292 iloop = gmx_mtop_ilistloop_init(mtop);
1293 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
1295 for (i = 0; i < ilist[F_SETTLE].nr; i += 4)
1297 if (settle_type == -1)
1299 settle_type = ilist[F_SETTLE].iatoms[i];
1301 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1304 "The [molecules] section of your topology specifies more than one block of\n"
1305 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1306 "are trying to partition your solvent into different *groups* (e.g. for\n"
1307 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1308 "files specify groups. Otherwise, you may wish to change the least-used\n"
1309 "block of molecules with SETTLE constraints into 3 normal constraints.");
1314 constr->n_at2settle_mt = mtop->nmoltype;
1315 snew(constr->at2settle_mt, constr->n_at2settle_mt);
1316 for (mt = 0; mt < mtop->nmoltype; mt++)
1318 constr->at2settle_mt[mt] =
1319 make_at2settle(mtop->moltype[mt].atoms.nr,
1320 &mtop->moltype[mt].ilist[F_SETTLE]);
1324 constr->maxwarn = 999;
1325 env = getenv("GMX_MAXCONSTRWARN");
1328 constr->maxwarn = 0;
1329 sscanf(env, "%d", &constr->maxwarn);
1333 "Setting the maximum number of constraint warnings to %d\n",
1339 "Setting the maximum number of constraint warnings to %d\n",
1343 if (constr->maxwarn < 0 && fplog)
1345 fprintf(fplog, "maxwarn < 0, will not stop on constraint errors\n");
1347 constr->warncount_lincs = 0;
1348 constr->warncount_settle = 0;
1350 /* Initialize the essential dynamics sampling.
1351 * Put the pointer to the ED struct in constr */
1353 if (ed != NULL || state->edsamstate.nED > 0)
1355 init_edsam(mtop, ir, cr, ed, state->x, state->box, &state->edsamstate);
1358 constr->warn_mtop = mtop;
1363 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1365 return constr->at2con_mt;
1368 const int **atom2settle_moltype(gmx_constr_t constr)
1370 return (const int **)constr->at2settle_mt;
1374 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1376 const gmx_moltype_t *molt;
1380 int nat, *at2cg, cg, a, ftype, i;
1384 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1386 molt = &mtop->moltype[mtop->molblock[mb].type];
1388 if (molt->ilist[F_CONSTR].nr > 0 ||
1389 molt->ilist[F_CONSTRNC].nr > 0 ||
1390 molt->ilist[F_SETTLE].nr > 0)
1393 snew(at2cg, molt->atoms.nr);
1394 for (cg = 0; cg < cgs->nr; cg++)
1396 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1402 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1404 il = &molt->ilist[ftype];
1405 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(ftype))
1407 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1421 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1423 const gmx_moltype_t *molt;
1427 int nat, *at2cg, cg, a, ftype, i;
1431 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1433 molt = &mtop->moltype[mtop->molblock[mb].type];
1435 if (molt->ilist[F_SETTLE].nr > 0)
1438 snew(at2cg, molt->atoms.nr);
1439 for (cg = 0; cg < cgs->nr; cg++)
1441 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1447 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1449 il = &molt->ilist[ftype];
1450 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(F_SETTLE))
1452 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1453 at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1467 /* helper functions for andersen temperature control, because the
1468 * gmx_constr construct is only defined in constr.c. Return the list
1469 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1471 extern int *get_sblock(struct gmx_constr *constr)
1473 return constr->sblock;
1476 extern int get_nblocks(struct gmx_constr *constr)
1478 return constr->nblocks;