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, 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"
64 typedef struct gmx_constr {
65 int ncon_tot; /* The total number of constraints */
66 int nflexcon; /* The number of flexible constraints */
67 int n_at2con_mt; /* The size of at2con = #moltypes */
68 t_blocka *at2con_mt; /* A list of atoms to constraints */
69 int n_at2settle_mt; /* The size of at2settle = #moltypes */
70 int **at2settle_mt; /* A list of atoms to settles */
71 gmx_bool bInterCGsettles;
72 gmx_lincsdata_t lincsd; /* LINCS data */
73 gmx_shakedata_t shaked; /* SHAKE data */
74 gmx_settledata_t settled; /* SETTLE data */
75 int nblocks; /* The number of SHAKE blocks */
76 int *sblock; /* The SHAKE blocks */
77 int sblock_nalloc; /* The allocation size of sblock */
78 real *lagr; /* Lagrange multipliers for SHAKE */
79 int lagr_nalloc; /* The allocation size of lagr */
80 int maxwarn; /* The maximum number of warnings */
83 gmx_edsam_t ed; /* The essential dynamics data */
85 tensor *vir_r_m_dr_th; /* Thread local working data */
86 int *settle_error; /* Thread local working data */
88 gmx_mtop_t *warn_mtop; /* Only used for printing warnings */
96 static void *init_vetavars(t_vetavars *vars,
97 gmx_bool constr_deriv,
98 real veta, real vetanew, t_inputrec *ir, gmx_ekindata_t *ekind, gmx_bool bPscal)
103 /* first, set the alpha integrator variable */
104 if ((ir->opts.nrdf[0] > 0) && bPscal)
106 vars->alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
112 g = 0.5*veta*ir->delta_t;
113 vars->rscale = exp(g)*series_sinhx(g);
114 g = -0.25*vars->alpha*veta*ir->delta_t;
115 vars->vscale = exp(g)*series_sinhx(g);
116 vars->rvscale = vars->vscale*vars->rscale;
117 vars->veta = vetanew;
121 snew(vars->vscale_nhc, ir->opts.ngtc);
122 if ((ekind == NULL) || (!bPscal))
124 for (i = 0; i < ir->opts.ngtc; i++)
126 vars->vscale_nhc[i] = 1;
131 for (i = 0; i < ir->opts.ngtc; i++)
133 vars->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
139 vars->vscale_nhc = NULL;
145 static void free_vetavars(t_vetavars *vars)
147 if (vars->vscale_nhc != NULL)
149 sfree(vars->vscale_nhc);
153 static int pcomp(const void *p1, const void *p2)
156 atom_id min1, min2, max1, max2;
157 t_sortblock *a1 = (t_sortblock *)p1;
158 t_sortblock *a2 = (t_sortblock *)p2;
160 db = a1->blocknr-a2->blocknr;
167 min1 = min(a1->iatom[1], a1->iatom[2]);
168 max1 = max(a1->iatom[1], a1->iatom[2]);
169 min2 = min(a2->iatom[1], a2->iatom[2]);
170 max2 = max(a2->iatom[1], a2->iatom[2]);
182 int n_flexible_constraints(struct gmx_constr *constr)
188 nflexcon = constr->nflexcon;
198 void too_many_constraint_warnings(int eConstrAlg, int warncount)
200 const char *abort = "- aborting to avoid logfile runaway.\n"
201 "This normally happens when your system is not sufficiently equilibrated,"
202 "or if you are changing lambda too fast in free energy simulations.\n";
205 "Too many %s warnings (%d)\n"
206 "If you know what you are doing you can %s"
207 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
208 "but normally it is better to fix the problem",
209 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
210 (eConstrAlg == econtLINCS) ?
211 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
214 static void write_constr_pdb(const char *fn, const char *title,
216 int start, int homenr, t_commrec *cr,
217 rvec x[], matrix box)
219 char fname[STRLEN], format[STRLEN];
221 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
228 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
229 if (DOMAINDECOMP(cr))
232 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
239 sprintf(fname, "%s.pdb", fn);
241 sprintf(format, "%s\n", get_pdbformat());
243 out = gmx_fio_fopen(fname, "w");
245 fprintf(out, "TITLE %s\n", title);
246 gmx_write_pdb_box(out, -1, box);
247 for (i = start; i < start+homenr; i++)
251 if (i >= dd->nat_home && i < dd_ac0)
255 ii = dd->gatindex[i];
261 gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
262 fprintf(out, format, "ATOM", (ii+1)%100000,
263 anm, resnm, ' ', resnr%10000, ' ',
264 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ]);
266 fprintf(out, "TER\n");
271 static void dump_confs(FILE *fplog, gmx_int64_t step, gmx_mtop_t *mtop,
272 int start, int homenr, t_commrec *cr,
273 rvec x[], rvec xprime[], matrix box)
275 char buf[256], buf2[22];
277 char *env = getenv("GMX_SUPPRESS_DUMP");
283 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
284 write_constr_pdb(buf, "initial coordinates",
285 mtop, start, homenr, cr, x, box);
286 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
287 write_constr_pdb(buf, "coordinates after constraining",
288 mtop, start, homenr, cr, xprime, box);
291 fprintf(fplog, "Wrote pdb files with previous and current coordinates\n");
293 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
296 static void pr_sortblock(FILE *fp, const char *title, int nsb, t_sortblock sb[])
300 fprintf(fp, "%s\n", title);
301 for (i = 0; (i < nsb); i++)
303 fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
304 i, sb[i].iatom[0], sb[i].iatom[1], sb[i].iatom[2],
309 gmx_bool constrain(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
310 struct gmx_constr *constr,
311 t_idef *idef, t_inputrec *ir, gmx_ekindata_t *ekind,
313 gmx_int64_t step, int delta_step,
315 rvec *x, rvec *xprime, rvec *min_proj,
316 gmx_bool bMolPBC, matrix box,
317 real lambda, real *dvdlambda,
318 rvec *v, tensor *vir,
319 t_nrnb *nrnb, int econq, gmx_bool bPscal,
320 real veta, real vetanew)
323 int start, homenr, nrend;
325 int ncons, settle_error;
328 real invdt, vir_fac, t;
331 t_pbc pbc, *pbc_null;
336 if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
338 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");
346 nrend = start+homenr;
348 /* set constants for pressure control integration */
349 init_vetavars(&vetavar, econq != econqCoord,
350 veta, vetanew, ir, ekind, bPscal);
352 if (ir->delta_t == 0)
358 invdt = 1/ir->delta_t;
361 if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
363 /* Set the constraint lengths for the step at which this configuration
364 * is meant to be. The invmasses should not be changed.
366 lambda += delta_step*ir->fepvals->delta_lambda;
371 clear_mat(vir_r_m_dr);
376 settle = &idef->il[F_SETTLE];
377 nsettle = settle->nr/(1+NRAL(F_SETTLE));
381 nth = gmx_omp_nthreads_get(emntSETTLE);
388 if (nth > 1 && constr->vir_r_m_dr_th == NULL)
390 snew(constr->vir_r_m_dr_th, nth);
391 snew(constr->settle_error, nth);
396 /* We do not need full pbc when constraints do not cross charge groups,
397 * i.e. when dd->constraint_comm==NULL.
398 * Note that PBC for constraints is different from PBC for bondeds.
399 * For constraints there is both forward and backward communication.
401 if (ir->ePBC != epbcNONE &&
402 (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm == NULL))
404 /* With pbc=screw the screw has been changed to a shift
405 * by the constraint coordinate communication routine,
406 * so that here we can use normal pbc.
408 pbc_null = set_pbc_dd(&pbc, ir->ePBC, cr->dd, FALSE, box);
415 /* Communicate the coordinates required for the non-local constraints
416 * for LINCS and/or SETTLE.
420 dd_move_x_constraints(cr->dd, box, x, xprime);
422 else if (PARTDECOMP(cr))
424 pd_move_x_constraints(cr, x, xprime);
427 if (constr->lincsd != NULL)
429 bOK = constrain_lincs(fplog, bLog, bEner, ir, step, constr->lincsd, md, cr,
431 box, pbc_null, lambda, dvdlambda,
432 invdt, v, vir != NULL, vir_r_m_dr,
434 constr->maxwarn, &constr->warncount_lincs);
435 if (!bOK && constr->maxwarn >= 0)
439 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
440 econstr_names[econtLINCS], gmx_step_str(step, buf));
446 if (constr->nblocks > 0)
451 bOK = bshakef(fplog, constr->shaked,
452 md->invmass, constr->nblocks, constr->sblock,
453 idef, ir, x, xprime, nrnb,
454 constr->lagr, lambda, dvdlambda,
455 invdt, v, vir != NULL, vir_r_m_dr,
456 constr->maxwarn >= 0, econq, &vetavar);
459 bOK = bshakef(fplog, constr->shaked,
460 md->invmass, constr->nblocks, constr->sblock,
461 idef, ir, x, min_proj, nrnb,
462 constr->lagr, lambda, dvdlambda,
463 invdt, NULL, vir != NULL, vir_r_m_dr,
464 constr->maxwarn >= 0, econq, &vetavar);
467 gmx_fatal(FARGS, "Internal error, SHAKE called for constraining something else than coordinates");
471 if (!bOK && constr->maxwarn >= 0)
475 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
476 econstr_names[econtSHAKE], gmx_step_str(step, buf));
484 int calcvir_atom_end;
488 calcvir_atom_end = 0;
492 calcvir_atom_end = md->start + md->homenr;
498 #pragma omp parallel for num_threads(nth) schedule(static)
499 for (th = 0; th < nth; th++)
501 int start_th, end_th;
505 clear_mat(constr->vir_r_m_dr_th[th]);
508 start_th = (nsettle* th )/nth;
509 end_th = (nsettle*(th+1))/nth;
510 if (start_th >= 0 && end_th - start_th > 0)
512 csettle(constr->settled,
514 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
517 invdt, v ? v[0] : NULL, calcvir_atom_end,
518 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
519 th == 0 ? &settle_error : &constr->settle_error[th],
523 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
526 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
530 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
536 case econqForceDispl:
537 #pragma omp parallel for num_threads(nth) schedule(static)
538 for (th = 0; th < nth; th++)
540 int start_th, end_th;
544 clear_mat(constr->vir_r_m_dr_th[th]);
547 start_th = (nsettle* th )/nth;
548 end_th = (nsettle*(th+1))/nth;
550 if (start_th >= 0 && end_th - start_th > 0)
552 settle_proj(constr->settled, econq,
554 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
557 xprime, min_proj, calcvir_atom_end,
558 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
562 /* This is an overestimate */
563 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
565 case econqDeriv_FlexCon:
566 /* Nothing to do, since the are no flexible constraints in settles */
569 gmx_incons("Unknown constraint quantity for settle");
575 /* Combine virial and error info of the other threads */
576 for (i = 1; i < nth; i++)
578 m_add(vir_r_m_dr, constr->vir_r_m_dr_th[i], vir_r_m_dr);
579 settle_error = constr->settle_error[i];
582 if (econq == econqCoord && settle_error >= 0)
585 if (constr->maxwarn >= 0)
589 "\nstep " "%"GMX_PRId64 ": Water molecule starting at atom %d can not be "
590 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
591 step, ddglatnr(cr->dd, settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
594 fprintf(fplog, "%s", buf);
596 fprintf(stderr, "%s", buf);
597 constr->warncount_settle++;
598 if (constr->warncount_settle > constr->maxwarn)
600 too_many_constraint_warnings(-1, constr->warncount_settle);
607 free_vetavars(&vetavar);
614 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
617 vir_fac = 0.5/ir->delta_t;
620 case econqForceDispl:
625 gmx_incons("Unsupported constraint quantity for virial");
630 vir_fac *= 2; /* only constraining over half the distance here */
632 for (i = 0; i < DIM; i++)
634 for (j = 0; j < DIM; j++)
636 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
643 dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
646 if (econq == econqCoord)
648 if (ir->ePull == epullCONSTRAINT)
650 if (EI_DYNAMICS(ir->eI))
652 t = ir->init_t + (step + delta_step)*ir->delta_t;
658 set_pbc(&pbc, ir->ePBC, box);
659 pull_constraint(ir->pull, md, &pbc, cr, ir->delta_t, t, x, xprime, v, *vir);
661 if (constr->ed && delta_step > 0)
663 /* apply the essential dynamcs constraints here */
664 do_edsam(ir, step, cr, xprime, v, box, constr->ed);
671 real *constr_rmsd_data(struct gmx_constr *constr)
675 return lincs_rmsd_data(constr->lincsd);
683 real constr_rmsd(struct gmx_constr *constr, gmx_bool bSD2)
687 return lincs_rmsd(constr->lincsd, bSD2);
695 static void make_shake_sblock_pd(struct gmx_constr *constr,
696 t_idef *idef, t_mdatoms *md)
705 /* Since we are processing the local topology,
706 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
708 ncons = idef->il[F_CONSTR].nr/3;
710 init_blocka(&sblocks);
711 gen_sblocks(NULL, md->start, md->start+md->homenr, idef, &sblocks, FALSE);
714 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
715 nblocks=blocks->multinr[idef->nodeid] - bstart;
718 constr->nblocks = sblocks.nr;
721 fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n",
722 ncons, bstart, constr->nblocks);
725 /* Calculate block number for each atom */
726 inv_sblock = make_invblocka(&sblocks, md->nr);
728 done_blocka(&sblocks);
730 /* Store the block number in temp array and
731 * sort the constraints in order of the sblock number
732 * and the atom numbers, really sorting a segment of the array!
735 pr_idef(fplog, 0, "Before Sort", idef);
737 iatom = idef->il[F_CONSTR].iatoms;
739 for (i = 0; (i < ncons); i++, iatom += 3)
741 for (m = 0; (m < 3); m++)
743 sb[i].iatom[m] = iatom[m];
745 sb[i].blocknr = inv_sblock[iatom[1]];
748 /* Now sort the blocks */
751 pr_sortblock(debug, "Before sorting", ncons, sb);
752 fprintf(debug, "Going to sort constraints\n");
755 qsort(sb, ncons, (size_t)sizeof(*sb), pcomp);
759 pr_sortblock(debug, "After sorting", ncons, sb);
762 iatom = idef->il[F_CONSTR].iatoms;
763 for (i = 0; (i < ncons); i++, iatom += 3)
765 for (m = 0; (m < 3); m++)
767 iatom[m] = sb[i].iatom[m];
771 pr_idef(fplog, 0, "After Sort", idef);
775 snew(constr->sblock, constr->nblocks+1);
777 for (i = 0; (i < ncons); i++)
779 if (sb[i].blocknr != bnr)
782 constr->sblock[j++] = 3*i;
786 constr->sblock[j++] = 3*ncons;
788 if (j != (constr->nblocks+1))
790 fprintf(stderr, "bstart: %d\n", bstart);
791 fprintf(stderr, "j: %d, nblocks: %d, ncons: %d\n",
792 j, constr->nblocks, ncons);
793 for (i = 0; (i < ncons); i++)
795 fprintf(stderr, "i: %5d sb[i].blocknr: %5u\n", i, sb[i].blocknr);
797 for (j = 0; (j <= constr->nblocks); j++)
799 fprintf(stderr, "sblock[%3d]=%5d\n", j, (int)constr->sblock[j]);
801 gmx_fatal(FARGS, "DEATH HORROR: "
802 "sblocks does not match idef->il[F_CONSTR]");
808 static void make_shake_sblock_dd(struct gmx_constr *constr,
809 t_ilist *ilcon, t_block *cgs,
815 if (dd->ncg_home+1 > constr->sblock_nalloc)
817 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
818 srenew(constr->sblock, constr->sblock_nalloc);
822 iatom = ilcon->iatoms;
825 for (c = 0; c < ncons; c++)
827 if (c == 0 || iatom[1] >= cgs->index[cg+1])
829 constr->sblock[constr->nblocks++] = 3*c;
830 while (iatom[1] >= cgs->index[cg+1])
837 constr->sblock[constr->nblocks] = 3*ncons;
840 t_blocka make_at2con(int start, int natoms,
841 t_ilist *ilist, t_iparams *iparams,
842 gmx_bool bDynamics, int *nflexiblecons)
844 int *count, ncon, con, con_tot, nflexcon, ftype, i, a;
851 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
853 ncon = ilist[ftype].nr/3;
854 ia = ilist[ftype].iatoms;
855 for (con = 0; con < ncon; con++)
857 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
858 iparams[ia[0]].constr.dB == 0);
863 if (bDynamics || !bFlexCon)
865 for (i = 1; i < 3; i++)
874 *nflexiblecons = nflexcon;
877 at2con.nalloc_index = at2con.nr+1;
878 snew(at2con.index, at2con.nalloc_index);
880 for (a = 0; a < natoms; a++)
882 at2con.index[a+1] = at2con.index[a] + count[a];
885 at2con.nra = at2con.index[natoms];
886 at2con.nalloc_a = at2con.nra;
887 snew(at2con.a, at2con.nalloc_a);
889 /* The F_CONSTRNC constraints have constraint numbers
890 * that continue after the last F_CONSTR constraint.
893 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
895 ncon = ilist[ftype].nr/3;
896 ia = ilist[ftype].iatoms;
897 for (con = 0; con < ncon; con++)
899 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
900 iparams[ia[0]].constr.dB == 0);
901 if (bDynamics || !bFlexCon)
903 for (i = 1; i < 3; i++)
906 at2con.a[at2con.index[a]+count[a]++] = con_tot;
919 static int *make_at2settle(int natoms, const t_ilist *ilist)
925 /* Set all to no settle */
926 for (a = 0; a < natoms; a++)
931 stride = 1 + NRAL(F_SETTLE);
933 for (s = 0; s < ilist->nr; s += stride)
935 at2s[ilist->iatoms[s+1]] = s/stride;
936 at2s[ilist->iatoms[s+2]] = s/stride;
937 at2s[ilist->iatoms[s+3]] = s/stride;
943 void set_constraints(struct gmx_constr *constr,
944 gmx_localtop_t *top, t_inputrec *ir,
945 t_mdatoms *md, t_commrec *cr)
954 if (constr->ncon_tot > 0)
956 /* We are using the local topology,
957 * so there are only F_CONSTR constraints.
959 ncons = idef->il[F_CONSTR].nr/3;
961 /* With DD we might also need to call LINCS with ncons=0 for
962 * communicating coordinates to other nodes that do have constraints.
964 if (ir->eConstrAlg == econtLINCS)
966 set_lincs(idef, md, EI_DYNAMICS(ir->eI), cr, constr->lincsd);
968 if (ir->eConstrAlg == econtSHAKE)
972 make_shake_sblock_dd(constr, &idef->il[F_CONSTR], &top->cgs, cr->dd);
976 make_shake_sblock_pd(constr, idef, md);
978 if (ncons > constr->lagr_nalloc)
980 constr->lagr_nalloc = over_alloc_dd(ncons);
981 srenew(constr->lagr, constr->lagr_nalloc);
986 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
988 settle = &idef->il[F_SETTLE];
989 iO = settle->iatoms[1];
990 iH = settle->iatoms[2];
992 settle_init(md->massT[iO], md->massT[iH],
993 md->invmass[iO], md->invmass[iH],
994 idef->iparams[settle->iatoms[0]].settle.doh,
995 idef->iparams[settle->iatoms[0]].settle.dhh);
998 /* Make a selection of the local atoms for essential dynamics */
999 if (constr->ed && cr->dd)
1001 dd_make_local_ed_indices(cr->dd, constr->ed);
1005 static void constr_recur(t_blocka *at2con,
1006 t_ilist *ilist, t_iparams *iparams, gmx_bool bTopB,
1007 int at, int depth, int nc, int *path,
1008 real r0, real r1, real *r2max,
1020 ncon1 = ilist[F_CONSTR].nr/3;
1021 ia1 = ilist[F_CONSTR].iatoms;
1022 ia2 = ilist[F_CONSTRNC].iatoms;
1024 /* Loop over all constraints connected to this atom */
1025 for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
1028 /* Do not walk over already used constraints */
1030 for (a1 = 0; a1 < depth; a1++)
1032 if (con == path[a1])
1039 ia = constr_iatomptr(ncon1, ia1, ia2, con);
1040 /* Flexible constraints currently have length 0, which is incorrect */
1043 len = iparams[ia[0]].constr.dA;
1047 len = iparams[ia[0]].constr.dB;
1049 /* In the worst case the bond directions alternate */
1060 /* Assume angles of 120 degrees between all bonds */
1061 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
1063 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
1066 fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
1067 for (a1 = 0; a1 < depth; a1++)
1069 fprintf(debug, " %d %5.3f",
1071 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
1073 fprintf(debug, " %d %5.3f\n", con, len);
1076 /* Limit the number of recursions to 1000*nc,
1077 * so a call does not take more than a second,
1078 * even for highly connected systems.
1080 if (depth + 1 < nc && *count < 1000*nc)
1092 constr_recur(at2con, ilist, iparams,
1093 bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
1100 static real constr_r_max_moltype(gmx_moltype_t *molt, t_iparams *iparams,
1103 int natoms, nflexcon, *path, at, count;
1106 real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
1108 if (molt->ilist[F_CONSTR].nr == 0 &&
1109 molt->ilist[F_CONSTRNC].nr == 0)
1114 natoms = molt->atoms.nr;
1116 at2con = make_at2con(0, natoms, molt->ilist, iparams,
1117 EI_DYNAMICS(ir->eI), &nflexcon);
1118 snew(path, 1+ir->nProjOrder);
1119 for (at = 0; at < 1+ir->nProjOrder; at++)
1125 for (at = 0; at < natoms; at++)
1131 constr_recur(&at2con, molt->ilist, iparams,
1132 FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
1134 if (ir->efep == efepNO)
1136 rmax = sqrt(r2maxA);
1141 for (at = 0; at < natoms; at++)
1146 constr_recur(&at2con, molt->ilist, iparams,
1147 TRUE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxB, &count);
1149 lam0 = ir->fepvals->init_lambda;
1150 if (EI_DYNAMICS(ir->eI))
1152 lam0 += ir->init_step*ir->fepvals->delta_lambda;
1154 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1155 if (EI_DYNAMICS(ir->eI))
1157 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1158 rmax = max(rmax, (1 - lam1)*sqrt(r2maxA) + lam1*sqrt(r2maxB));
1162 done_blocka(&at2con);
1168 real constr_r_max(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir)
1174 for (mt = 0; mt < mtop->nmoltype; mt++)
1177 constr_r_max_moltype(&mtop->moltype[mt],
1178 mtop->ffparams.iparams, ir));
1183 fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
1189 gmx_constr_t init_constraints(FILE *fplog,
1190 gmx_mtop_t *mtop, t_inputrec *ir,
1191 gmx_edsam_t ed, t_state *state,
1194 int ncon, nset, nmol, settle_type, i, natoms, mt, nflexcon;
1195 struct gmx_constr *constr;
1198 gmx_mtop_ilistloop_t iloop;
1201 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1202 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1203 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
1205 if (ncon+nset == 0 && ir->ePull != epullCONSTRAINT && ed == NULL)
1212 constr->ncon_tot = ncon;
1213 constr->nflexcon = 0;
1216 constr->n_at2con_mt = mtop->nmoltype;
1217 snew(constr->at2con_mt, constr->n_at2con_mt);
1218 for (mt = 0; mt < mtop->nmoltype; mt++)
1220 constr->at2con_mt[mt] = make_at2con(0, mtop->moltype[mt].atoms.nr,
1221 mtop->moltype[mt].ilist,
1222 mtop->ffparams.iparams,
1223 EI_DYNAMICS(ir->eI), &nflexcon);
1224 for (i = 0; i < mtop->nmolblock; i++)
1226 if (mtop->molblock[i].type == mt)
1228 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1233 if (constr->nflexcon > 0)
1237 fprintf(fplog, "There are %d flexible constraints\n",
1239 if (ir->fc_stepsize == 0)
1242 "WARNING: step size for flexible constraining = 0\n"
1243 " All flexible constraints will be rigid.\n"
1244 " Will try to keep all flexible constraints at their original length,\n"
1245 " but the lengths may exhibit some drift.\n\n");
1246 constr->nflexcon = 0;
1249 if (constr->nflexcon > 0)
1251 please_cite(fplog, "Hess2002");
1255 if (ir->eConstrAlg == econtLINCS)
1257 constr->lincsd = init_lincs(fplog, mtop,
1258 constr->nflexcon, constr->at2con_mt,
1259 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1260 ir->nLincsIter, ir->nProjOrder);
1263 if (ir->eConstrAlg == econtSHAKE)
1265 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1267 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1269 if (constr->nflexcon)
1271 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");
1273 please_cite(fplog, "Ryckaert77a");
1276 please_cite(fplog, "Barth95a");
1279 constr->shaked = shake_init();
1285 please_cite(fplog, "Miyamoto92a");
1287 constr->bInterCGsettles = inter_charge_group_settles(mtop);
1289 /* Check that we have only one settle type */
1291 iloop = gmx_mtop_ilistloop_init(mtop);
1292 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
1294 for (i = 0; i < ilist[F_SETTLE].nr; i += 4)
1296 if (settle_type == -1)
1298 settle_type = ilist[F_SETTLE].iatoms[i];
1300 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1303 "The [molecules] section of your topology specifies more than one block of\n"
1304 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1305 "are trying to partition your solvent into different *groups* (e.g. for\n"
1306 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1307 "files specify groups. Otherwise, you may wish to change the least-used\n"
1308 "block of molecules with SETTLE constraints into 3 normal constraints.");
1313 constr->n_at2settle_mt = mtop->nmoltype;
1314 snew(constr->at2settle_mt, constr->n_at2settle_mt);
1315 for (mt = 0; mt < mtop->nmoltype; mt++)
1317 constr->at2settle_mt[mt] =
1318 make_at2settle(mtop->moltype[mt].atoms.nr,
1319 &mtop->moltype[mt].ilist[F_SETTLE]);
1323 constr->maxwarn = 999;
1324 env = getenv("GMX_MAXCONSTRWARN");
1327 constr->maxwarn = 0;
1328 sscanf(env, "%d", &constr->maxwarn);
1332 "Setting the maximum number of constraint warnings to %d\n",
1338 "Setting the maximum number of constraint warnings to %d\n",
1342 if (constr->maxwarn < 0 && fplog)
1344 fprintf(fplog, "maxwarn < 0, will not stop on constraint errors\n");
1346 constr->warncount_lincs = 0;
1347 constr->warncount_settle = 0;
1349 /* Initialize the essential dynamics sampling.
1350 * Put the pointer to the ED struct in constr */
1352 if (ed != NULL || state->edsamstate.nED > 0)
1354 init_edsam(mtop, ir, cr, ed, state->x, state->box, &state->edsamstate);
1357 constr->warn_mtop = mtop;
1362 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1364 return constr->at2con_mt;
1367 const int **atom2settle_moltype(gmx_constr_t constr)
1369 return (const int **)constr->at2settle_mt;
1373 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1375 const gmx_moltype_t *molt;
1379 int nat, *at2cg, cg, a, ftype, i;
1383 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1385 molt = &mtop->moltype[mtop->molblock[mb].type];
1387 if (molt->ilist[F_CONSTR].nr > 0 ||
1388 molt->ilist[F_CONSTRNC].nr > 0 ||
1389 molt->ilist[F_SETTLE].nr > 0)
1392 snew(at2cg, molt->atoms.nr);
1393 for (cg = 0; cg < cgs->nr; cg++)
1395 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1401 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1403 il = &molt->ilist[ftype];
1404 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(ftype))
1406 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1420 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1422 const gmx_moltype_t *molt;
1426 int nat, *at2cg, cg, a, ftype, i;
1430 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1432 molt = &mtop->moltype[mtop->molblock[mb].type];
1434 if (molt->ilist[F_SETTLE].nr > 0)
1437 snew(at2cg, molt->atoms.nr);
1438 for (cg = 0; cg < cgs->nr; cg++)
1440 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1446 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1448 il = &molt->ilist[ftype];
1449 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(F_SETTLE))
1451 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1452 at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1466 /* helper functions for andersen temperature control, because the
1467 * gmx_constr construct is only defined in constr.c. Return the list
1468 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1470 extern int *get_sblock(struct gmx_constr *constr)
1472 return constr->sblock;
1475 extern int get_nblocks(struct gmx_constr *constr)
1477 return constr->nblocks;