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,2015,2016, 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.
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/essentialdynamics/edsam.h"
51 #include "gromacs/fileio/confio.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/pdbio.h"
54 #include "gromacs/gmxlib/nrnb.h"
55 #include "gromacs/math/utilities.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdlib/gmx_omp_nthreads.h"
58 #include "gromacs/mdlib/mdrun.h"
59 #include "gromacs/mdlib/splitter.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/pulling/pull.h"
65 #include "gromacs/topology/block.h"
66 #include "gromacs/topology/invblock.h"
67 #include "gromacs/topology/mtop_util.h"
68 #include "gromacs/utility/exceptions.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/pleasecite.h"
71 #include "gromacs/utility/smalloc.h"
72 #include "gromacs/utility/txtdump.h"
74 typedef struct gmx_constr {
75 int ncon_tot; /* The total number of constraints */
76 int nflexcon; /* The number of flexible constraints */
77 int n_at2con_mt; /* The size of at2con = #moltypes */
78 t_blocka *at2con_mt; /* A list of atoms to constraints */
79 int n_at2settle_mt; /* The size of at2settle = #moltypes */
80 int **at2settle_mt; /* A list of atoms to settles */
81 gmx_bool bInterCGsettles;
82 gmx_lincsdata_t lincsd; /* LINCS data */
83 gmx_shakedata_t shaked; /* SHAKE data */
84 gmx_settledata_t settled; /* SETTLE data */
85 int nblocks; /* The number of SHAKE blocks */
86 int *sblock; /* The SHAKE blocks */
87 int sblock_nalloc; /* The allocation size of sblock */
88 real *lagr; /* -2 times the Lagrange multipliers for SHAKE */
89 int lagr_nalloc; /* The allocation size of lagr */
90 int maxwarn; /* The maximum number of warnings */
93 gmx_edsam_t ed; /* The essential dynamics data */
95 tensor *vir_r_m_dr_th; /* Thread local working data */
96 int *settle_error; /* Thread local working data */
98 const gmx_mtop_t *warn_mtop; /* Only used for printing warnings */
106 static int pcomp(const void *p1, const void *p2)
109 int min1, min2, max1, max2;
110 t_sortblock *a1 = (t_sortblock *)p1;
111 t_sortblock *a2 = (t_sortblock *)p2;
113 db = a1->blocknr-a2->blocknr;
120 min1 = std::min(a1->iatom[1], a1->iatom[2]);
121 max1 = std::max(a1->iatom[1], a1->iatom[2]);
122 min2 = std::min(a2->iatom[1], a2->iatom[2]);
123 max2 = std::max(a2->iatom[1], a2->iatom[2]);
135 int n_flexible_constraints(struct gmx_constr *constr)
141 nflexcon = constr->nflexcon;
151 static void clear_constraint_quantity_nonlocal(gmx_domdec_t *dd, rvec *q)
153 int nonlocal_at_start, nonlocal_at_end, at;
155 dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
157 for (at = nonlocal_at_start; at < nonlocal_at_end; at++)
163 void too_many_constraint_warnings(int eConstrAlg, int warncount)
166 "Too many %s warnings (%d)\n"
167 "If you know what you are doing you can %s"
168 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
169 "but normally it is better to fix the problem",
170 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
171 (eConstrAlg == econtLINCS) ?
172 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
175 static void write_constr_pdb(const char *fn, const char *title,
176 const gmx_mtop_t *mtop,
177 int start, int homenr, t_commrec *cr,
178 rvec x[], matrix box)
182 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
187 if (DOMAINDECOMP(cr))
190 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
197 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
201 sprintf(fname, "%s.pdb", fn);
204 out = gmx_fio_fopen(fname, "w");
206 fprintf(out, "TITLE %s\n", title);
207 gmx_write_pdb_box(out, -1, box);
208 for (i = start; i < start+homenr; i++)
212 if (i >= dd->nat_home && i < dd_ac0)
216 ii = dd->gatindex[i];
222 gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
223 gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, anm, ' ', resnm, ' ', resnr, ' ',
224 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 0.0, "");
226 fprintf(out, "TER\n");
231 static void dump_confs(FILE *fplog, gmx_int64_t step, const gmx_mtop_t *mtop,
232 int start, int homenr, t_commrec *cr,
233 rvec x[], rvec xprime[], matrix box)
235 char buf[256], buf2[22];
237 char *env = getenv("GMX_SUPPRESS_DUMP");
243 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
244 write_constr_pdb(buf, "initial coordinates",
245 mtop, start, homenr, cr, x, box);
246 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
247 write_constr_pdb(buf, "coordinates after constraining",
248 mtop, start, homenr, cr, xprime, box);
251 fprintf(fplog, "Wrote pdb files with previous and current coordinates\n");
253 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
256 static void pr_sortblock(FILE *fp, const char *title, int nsb, t_sortblock sb[])
260 fprintf(fp, "%s\n", title);
261 for (i = 0; (i < nsb); i++)
263 fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
264 i, sb[i].iatom[0], sb[i].iatom[1], sb[i].iatom[2],
269 gmx_bool constrain(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
270 struct gmx_constr *constr,
271 t_idef *idef, t_inputrec *ir,
273 gmx_int64_t step, int delta_step,
276 rvec *x, rvec *xprime, rvec *min_proj,
277 gmx_bool bMolPBC, matrix box,
278 real lambda, real *dvdlambda,
279 rvec *v, tensor *vir,
280 t_nrnb *nrnb, int econq)
288 real invdt, vir_fac = 0, t;
291 t_pbc pbc, *pbc_null;
295 if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
297 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");
306 scaled_delta_t = step_scaling * ir->delta_t;
308 /* Prepare time step for use in constraint implementations, and
309 avoid generating inf when ir->delta_t = 0. */
310 if (ir->delta_t == 0)
316 invdt = 1.0/scaled_delta_t;
319 if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
321 /* Set the constraint lengths for the step at which this configuration
322 * is meant to be. The invmasses should not be changed.
324 lambda += delta_step*ir->fepvals->delta_lambda;
329 clear_mat(vir_r_m_dr);
334 settle = &idef->il[F_SETTLE];
335 nsettle = settle->nr/(1+NRAL(F_SETTLE));
339 nth = gmx_omp_nthreads_get(emntSETTLE);
346 if (nth > 1 && constr->vir_r_m_dr_th == NULL)
348 snew(constr->vir_r_m_dr_th, nth);
349 snew(constr->settle_error, nth);
354 /* We do not need full pbc when constraints do not cross charge groups,
355 * i.e. when dd->constraint_comm==NULL.
356 * Note that PBC for constraints is different from PBC for bondeds.
357 * For constraints there is both forward and backward communication.
359 if (ir->ePBC != epbcNONE &&
360 (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm == NULL))
362 /* With pbc=screw the screw has been changed to a shift
363 * by the constraint coordinate communication routine,
364 * so that here we can use normal pbc.
366 pbc_null = set_pbc_dd(&pbc, ir->ePBC,
367 DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
375 /* Communicate the coordinates required for the non-local constraints
376 * for LINCS and/or SETTLE.
380 dd_move_x_constraints(cr->dd, box, x, xprime, econq == econqCoord);
384 /* We need to initialize the non-local components of v.
385 * We never actually use these values, but we do increment them,
386 * so we should avoid uninitialized variables and overflows.
388 clear_constraint_quantity_nonlocal(cr->dd, v);
392 if (constr->lincsd != NULL)
394 bOK = constrain_lincs(fplog, bLog, bEner, ir, step, constr->lincsd, md, cr,
396 box, pbc_null, lambda, dvdlambda,
397 invdt, v, vir != NULL, vir_r_m_dr,
399 constr->maxwarn, &constr->warncount_lincs);
400 if (!bOK && constr->maxwarn >= 0)
404 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
405 econstr_names[econtLINCS], gmx_step_str(step, buf));
411 if (constr->nblocks > 0)
416 bOK = bshakef(fplog, constr->shaked,
417 md->invmass, constr->nblocks, constr->sblock,
418 idef, ir, x, xprime, nrnb,
419 constr->lagr, lambda, dvdlambda,
420 invdt, v, vir != NULL, vir_r_m_dr,
421 constr->maxwarn >= 0, econq);
424 bOK = bshakef(fplog, constr->shaked,
425 md->invmass, constr->nblocks, constr->sblock,
426 idef, ir, x, min_proj, nrnb,
427 constr->lagr, lambda, dvdlambda,
428 invdt, NULL, vir != NULL, vir_r_m_dr,
429 constr->maxwarn >= 0, econq);
432 gmx_fatal(FARGS, "Internal error, SHAKE called for constraining something else than coordinates");
436 if (!bOK && constr->maxwarn >= 0)
440 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
441 econstr_names[econtSHAKE], gmx_step_str(step, buf));
449 int calcvir_atom_end;
453 calcvir_atom_end = 0;
457 calcvir_atom_end = md->homenr;
463 #pragma omp parallel for num_threads(nth) schedule(static)
464 for (th = 0; th < nth; th++)
468 int start_th, end_th;
472 clear_mat(constr->vir_r_m_dr_th[th]);
474 constr->settle_error[th] = -1;
477 start_th = (nsettle* th )/nth;
478 end_th = (nsettle*(th+1))/nth;
479 if (start_th >= 0 && end_th - start_th > 0)
481 csettle(constr->settled,
483 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
486 invdt, v ? v[0] : NULL, calcvir_atom_end,
487 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
488 th == 0 ? &settle_error : &constr->settle_error[th]);
491 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
493 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
496 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
500 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
506 case econqForceDispl:
507 #pragma omp parallel for num_threads(nth) schedule(static)
508 for (th = 0; th < nth; th++)
512 int start_th, end_th;
516 clear_mat(constr->vir_r_m_dr_th[th]);
519 start_th = (nsettle* th )/nth;
520 end_th = (nsettle*(th+1))/nth;
522 if (start_th >= 0 && end_th - start_th > 0)
524 settle_proj(constr->settled, econq,
526 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
529 xprime, min_proj, calcvir_atom_end,
530 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th]);
533 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
535 /* This is an overestimate */
536 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
538 case econqDeriv_FlexCon:
539 /* Nothing to do, since the are no flexible constraints in settles */
542 gmx_incons("Unknown constraint quantity for settle");
550 /* Reduce the virial contributions over the threads */
551 for (i = 1; i < nth; i++)
553 m_add(vir_r_m_dr, constr->vir_r_m_dr_th[i], vir_r_m_dr);
557 if (econq == econqCoord)
559 for (i = 1; i < nth; i++)
561 settle_error = std::max(settle_error, constr->settle_error[i]);
564 if (settle_error >= 0)
566 dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
569 "\nstep " "%" GMX_PRId64 ": Water molecule starting at atom %d can not be "
570 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
571 step, ddglatnr(cr->dd, settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
578 /* The normal uses of constrain() pass step_scaling = 1.0.
579 * The call to constrain() for SD1 that passes step_scaling =
580 * 0.5 also passes vir = NULL, so cannot reach this
581 * assertion. This assertion should remain until someone knows
582 * that this path works for their intended purpose, and then
583 * they can use scaled_delta_t instead of ir->delta_t
585 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
589 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
592 vir_fac = 0.5/ir->delta_t;
595 case econqForceDispl:
599 gmx_incons("Unsupported constraint quantity for virial");
604 vir_fac *= 2; /* only constraining over half the distance here */
606 for (i = 0; i < DIM; i++)
608 for (j = 0; j < DIM; j++)
610 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
617 dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
620 if (econq == econqCoord)
622 if (ir->bPull && pull_have_constraint(ir->pull_work))
624 if (EI_DYNAMICS(ir->eI))
626 t = ir->init_t + (step + delta_step)*ir->delta_t;
632 set_pbc(&pbc, ir->ePBC, box);
633 pull_constraint(ir->pull_work, md, &pbc, cr, ir->delta_t, t, x, xprime, v, *vir);
635 if (constr->ed && delta_step > 0)
637 /* apply the essential dynamcs constraints here */
638 do_edsam(ir, step, cr, xprime, v, box, constr->ed);
645 real *constr_rmsd_data(struct gmx_constr *constr)
649 return lincs_rmsd_data(constr->lincsd);
657 real constr_rmsd(struct gmx_constr *constr)
661 return lincs_rmsd(constr->lincsd);
669 static void make_shake_sblock_serial(struct gmx_constr *constr,
670 t_idef *idef, const t_mdatoms *md)
679 /* Since we are processing the local topology,
680 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
682 ncons = idef->il[F_CONSTR].nr/3;
684 init_blocka(&sblocks);
685 gen_sblocks(NULL, 0, md->homenr, idef, &sblocks, FALSE);
688 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
689 nblocks=blocks->multinr[idef->nodeid] - bstart;
692 constr->nblocks = sblocks.nr;
695 fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n",
696 ncons, bstart, constr->nblocks);
699 /* Calculate block number for each atom */
700 inv_sblock = make_invblocka(&sblocks, md->nr);
702 done_blocka(&sblocks);
704 /* Store the block number in temp array and
705 * sort the constraints in order of the sblock number
706 * and the atom numbers, really sorting a segment of the array!
709 pr_idef(fplog, 0, "Before Sort", idef);
711 iatom = idef->il[F_CONSTR].iatoms;
713 for (i = 0; (i < ncons); i++, iatom += 3)
715 for (m = 0; (m < 3); m++)
717 sb[i].iatom[m] = iatom[m];
719 sb[i].blocknr = inv_sblock[iatom[1]];
722 /* Now sort the blocks */
725 pr_sortblock(debug, "Before sorting", ncons, sb);
726 fprintf(debug, "Going to sort constraints\n");
729 qsort(sb, ncons, (size_t)sizeof(*sb), pcomp);
733 pr_sortblock(debug, "After sorting", ncons, sb);
736 iatom = idef->il[F_CONSTR].iatoms;
737 for (i = 0; (i < ncons); i++, iatom += 3)
739 for (m = 0; (m < 3); m++)
741 iatom[m] = sb[i].iatom[m];
745 pr_idef(fplog, 0, "After Sort", idef);
749 snew(constr->sblock, constr->nblocks+1);
751 for (i = 0; (i < ncons); i++)
753 if (sb[i].blocknr != bnr)
756 constr->sblock[j++] = 3*i;
760 constr->sblock[j++] = 3*ncons;
762 if (j != (constr->nblocks+1))
764 fprintf(stderr, "bstart: %d\n", bstart);
765 fprintf(stderr, "j: %d, nblocks: %d, ncons: %d\n",
766 j, constr->nblocks, ncons);
767 for (i = 0; (i < ncons); i++)
769 fprintf(stderr, "i: %5d sb[i].blocknr: %5d\n", i, sb[i].blocknr);
771 for (j = 0; (j <= constr->nblocks); j++)
773 fprintf(stderr, "sblock[%3d]=%5d\n", j, (int)constr->sblock[j]);
775 gmx_fatal(FARGS, "DEATH HORROR: "
776 "sblocks does not match idef->il[F_CONSTR]");
782 static void make_shake_sblock_dd(struct gmx_constr *constr,
783 const t_ilist *ilcon, const t_block *cgs,
784 const gmx_domdec_t *dd)
789 if (dd->ncg_home+1 > constr->sblock_nalloc)
791 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
792 srenew(constr->sblock, constr->sblock_nalloc);
796 iatom = ilcon->iatoms;
799 for (c = 0; c < ncons; c++)
801 if (c == 0 || iatom[1] >= cgs->index[cg+1])
803 constr->sblock[constr->nblocks++] = 3*c;
804 while (iatom[1] >= cgs->index[cg+1])
811 constr->sblock[constr->nblocks] = 3*ncons;
814 t_blocka make_at2con(int start, int natoms,
815 const t_ilist *ilist, const t_iparams *iparams,
816 gmx_bool bDynamics, int *nflexiblecons)
818 int *count, ncon, con, con_tot, nflexcon, ftype, i, a;
825 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
827 ncon = ilist[ftype].nr/3;
828 ia = ilist[ftype].iatoms;
829 for (con = 0; con < ncon; con++)
831 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
832 iparams[ia[0]].constr.dB == 0);
837 if (bDynamics || !bFlexCon)
839 for (i = 1; i < 3; i++)
848 *nflexiblecons = nflexcon;
851 at2con.nalloc_index = at2con.nr+1;
852 snew(at2con.index, at2con.nalloc_index);
854 for (a = 0; a < natoms; a++)
856 at2con.index[a+1] = at2con.index[a] + count[a];
859 at2con.nra = at2con.index[natoms];
860 at2con.nalloc_a = at2con.nra;
861 snew(at2con.a, at2con.nalloc_a);
863 /* The F_CONSTRNC constraints have constraint numbers
864 * that continue after the last F_CONSTR constraint.
867 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
869 ncon = ilist[ftype].nr/3;
870 ia = ilist[ftype].iatoms;
871 for (con = 0; con < ncon; con++)
873 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
874 iparams[ia[0]].constr.dB == 0);
875 if (bDynamics || !bFlexCon)
877 for (i = 1; i < 3; i++)
880 at2con.a[at2con.index[a]+count[a]++] = con_tot;
893 static int *make_at2settle(int natoms, const t_ilist *ilist)
899 /* Set all to no settle */
900 for (a = 0; a < natoms; a++)
905 stride = 1 + NRAL(F_SETTLE);
907 for (s = 0; s < ilist->nr; s += stride)
909 at2s[ilist->iatoms[s+1]] = s/stride;
910 at2s[ilist->iatoms[s+2]] = s/stride;
911 at2s[ilist->iatoms[s+3]] = s/stride;
917 void set_constraints(struct gmx_constr *constr,
918 gmx_localtop_t *top, const t_inputrec *ir,
919 const t_mdatoms *md, t_commrec *cr)
923 const t_ilist *settle;
928 if (constr->ncon_tot > 0)
930 /* We are using the local topology,
931 * so there are only F_CONSTR constraints.
933 ncons = idef->il[F_CONSTR].nr/3;
935 /* With DD we might also need to call LINCS with ncons=0 for
936 * communicating coordinates to other nodes that do have constraints.
938 if (ir->eConstrAlg == econtLINCS)
940 set_lincs(idef, md, EI_DYNAMICS(ir->eI), cr, constr->lincsd);
942 if (ir->eConstrAlg == econtSHAKE)
946 make_shake_sblock_dd(constr, &idef->il[F_CONSTR], &top->cgs, cr->dd);
950 make_shake_sblock_serial(constr, idef, md);
952 if (ncons > constr->lagr_nalloc)
954 constr->lagr_nalloc = over_alloc_dd(ncons);
955 srenew(constr->lagr, constr->lagr_nalloc);
960 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
962 settle = &idef->il[F_SETTLE];
963 iO = settle->iatoms[1];
964 iH = settle->iatoms[2];
966 settle_init(md->massT[iO], md->massT[iH],
967 md->invmass[iO], md->invmass[iH],
968 idef->iparams[settle->iatoms[0]].settle.doh,
969 idef->iparams[settle->iatoms[0]].settle.dhh);
972 /* Make a selection of the local atoms for essential dynamics */
973 if (constr->ed && cr->dd)
975 dd_make_local_ed_indices(cr->dd, constr->ed);
979 static void constr_recur(const t_blocka *at2con,
980 const t_ilist *ilist, const t_iparams *iparams,
982 int at, int depth, int nc, int *path,
983 real r0, real r1, real *r2max,
995 ncon1 = ilist[F_CONSTR].nr/3;
996 ia1 = ilist[F_CONSTR].iatoms;
997 ia2 = ilist[F_CONSTRNC].iatoms;
999 /* Loop over all constraints connected to this atom */
1000 for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
1003 /* Do not walk over already used constraints */
1005 for (a1 = 0; a1 < depth; a1++)
1007 if (con == path[a1])
1014 ia = constr_iatomptr(ncon1, ia1, ia2, con);
1015 /* Flexible constraints currently have length 0, which is incorrect */
1018 len = iparams[ia[0]].constr.dA;
1022 len = iparams[ia[0]].constr.dB;
1024 /* In the worst case the bond directions alternate */
1035 /* Assume angles of 120 degrees between all bonds */
1036 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
1038 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
1041 fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
1042 for (a1 = 0; a1 < depth; a1++)
1044 fprintf(debug, " %d %5.3f",
1046 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
1048 fprintf(debug, " %d %5.3f\n", con, len);
1051 /* Limit the number of recursions to 1000*nc,
1052 * so a call does not take more than a second,
1053 * even for highly connected systems.
1055 if (depth + 1 < nc && *count < 1000*nc)
1067 constr_recur(at2con, ilist, iparams,
1068 bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
1075 static real constr_r_max_moltype(const gmx_moltype_t *molt,
1076 const t_iparams *iparams,
1077 const t_inputrec *ir)
1079 int natoms, nflexcon, *path, at, count;
1082 real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
1084 if (molt->ilist[F_CONSTR].nr == 0 &&
1085 molt->ilist[F_CONSTRNC].nr == 0)
1090 natoms = molt->atoms.nr;
1092 at2con = make_at2con(0, natoms, molt->ilist, iparams,
1093 EI_DYNAMICS(ir->eI), &nflexcon);
1094 snew(path, 1+ir->nProjOrder);
1095 for (at = 0; at < 1+ir->nProjOrder; at++)
1101 for (at = 0; at < natoms; at++)
1107 constr_recur(&at2con, molt->ilist, iparams,
1108 FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
1110 if (ir->efep == efepNO)
1112 rmax = sqrt(r2maxA);
1117 for (at = 0; at < natoms; at++)
1122 constr_recur(&at2con, molt->ilist, iparams,
1123 TRUE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxB, &count);
1125 lam0 = ir->fepvals->init_lambda;
1126 if (EI_DYNAMICS(ir->eI))
1128 lam0 += ir->init_step*ir->fepvals->delta_lambda;
1130 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1131 if (EI_DYNAMICS(ir->eI))
1133 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1134 rmax = std::max(rmax, (1 - lam1)*std::sqrt(r2maxA) + lam1*std::sqrt(r2maxB));
1138 done_blocka(&at2con);
1144 real constr_r_max(FILE *fplog, const gmx_mtop_t *mtop, const t_inputrec *ir)
1150 for (mt = 0; mt < mtop->nmoltype; mt++)
1152 rmax = std::max(rmax,
1153 constr_r_max_moltype(&mtop->moltype[mt],
1154 mtop->ffparams.iparams, ir));
1159 fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
1165 gmx_constr_t init_constraints(FILE *fplog,
1166 const gmx_mtop_t *mtop, const t_inputrec *ir,
1167 gmx_edsam_t ed, t_state *state,
1170 int ncon, nset, nmol, settle_type, i, mt, nflexcon;
1171 struct gmx_constr *constr;
1174 gmx_mtop_ilistloop_t iloop;
1177 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1178 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1179 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
1181 if (ncon+nset == 0 &&
1182 !(ir->bPull && pull_have_constraint(ir->pull_work)) &&
1190 constr->ncon_tot = ncon;
1191 constr->nflexcon = 0;
1194 constr->n_at2con_mt = mtop->nmoltype;
1195 snew(constr->at2con_mt, constr->n_at2con_mt);
1196 for (mt = 0; mt < mtop->nmoltype; mt++)
1198 constr->at2con_mt[mt] = make_at2con(0, mtop->moltype[mt].atoms.nr,
1199 mtop->moltype[mt].ilist,
1200 mtop->ffparams.iparams,
1201 EI_DYNAMICS(ir->eI), &nflexcon);
1202 for (i = 0; i < mtop->nmolblock; i++)
1204 if (mtop->molblock[i].type == mt)
1206 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1211 if (constr->nflexcon > 0)
1215 fprintf(fplog, "There are %d flexible constraints\n",
1217 if (ir->fc_stepsize == 0)
1220 "WARNING: step size for flexible constraining = 0\n"
1221 " All flexible constraints will be rigid.\n"
1222 " Will try to keep all flexible constraints at their original length,\n"
1223 " but the lengths may exhibit some drift.\n\n");
1224 constr->nflexcon = 0;
1227 if (constr->nflexcon > 0)
1229 please_cite(fplog, "Hess2002");
1233 if (ir->eConstrAlg == econtLINCS)
1235 constr->lincsd = init_lincs(fplog, mtop,
1236 constr->nflexcon, constr->at2con_mt,
1237 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1238 ir->nLincsIter, ir->nProjOrder);
1241 if (ir->eConstrAlg == econtSHAKE)
1243 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1245 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1247 if (constr->nflexcon)
1249 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");
1251 please_cite(fplog, "Ryckaert77a");
1254 please_cite(fplog, "Barth95a");
1257 constr->shaked = shake_init();
1263 please_cite(fplog, "Miyamoto92a");
1265 constr->bInterCGsettles = inter_charge_group_settles(mtop);
1267 /* Check that we have only one settle type */
1269 iloop = gmx_mtop_ilistloop_init(mtop);
1270 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
1272 for (i = 0; i < ilist[F_SETTLE].nr; i += 4)
1274 if (settle_type == -1)
1276 settle_type = ilist[F_SETTLE].iatoms[i];
1278 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1281 "The [molecules] section of your topology specifies more than one block of\n"
1282 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1283 "are trying to partition your solvent into different *groups* (e.g. for\n"
1284 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1285 "files specify groups. Otherwise, you may wish to change the least-used\n"
1286 "block of molecules with SETTLE constraints into 3 normal constraints.");
1291 constr->n_at2settle_mt = mtop->nmoltype;
1292 snew(constr->at2settle_mt, constr->n_at2settle_mt);
1293 for (mt = 0; mt < mtop->nmoltype; mt++)
1295 constr->at2settle_mt[mt] =
1296 make_at2settle(mtop->moltype[mt].atoms.nr,
1297 &mtop->moltype[mt].ilist[F_SETTLE]);
1301 if ((ncon + nset) > 0 && ir->epc == epcMTTK)
1303 gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
1306 constr->maxwarn = 999;
1307 env = getenv("GMX_MAXCONSTRWARN");
1310 constr->maxwarn = 0;
1311 sscanf(env, "%8d", &constr->maxwarn);
1315 "Setting the maximum number of constraint warnings to %d\n",
1321 "Setting the maximum number of constraint warnings to %d\n",
1325 if (constr->maxwarn < 0 && fplog)
1327 fprintf(fplog, "maxwarn < 0, will not stop on constraint errors\n");
1329 constr->warncount_lincs = 0;
1330 constr->warncount_settle = 0;
1332 /* Initialize the essential dynamics sampling.
1333 * Put the pointer to the ED struct in constr */
1335 if (ed != NULL || state->edsamstate.nED > 0)
1337 init_edsam(mtop, ir, cr, ed, state->x, state->box, &state->edsamstate);
1340 constr->warn_mtop = mtop;
1345 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1347 return constr->at2con_mt;
1350 const int **atom2settle_moltype(gmx_constr_t constr)
1352 return (const int **)constr->at2settle_mt;
1356 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1358 const gmx_moltype_t *molt;
1362 int *at2cg, cg, a, ftype, i;
1366 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1368 molt = &mtop->moltype[mtop->molblock[mb].type];
1370 if (molt->ilist[F_CONSTR].nr > 0 ||
1371 molt->ilist[F_CONSTRNC].nr > 0 ||
1372 molt->ilist[F_SETTLE].nr > 0)
1375 snew(at2cg, molt->atoms.nr);
1376 for (cg = 0; cg < cgs->nr; cg++)
1378 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1384 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1386 il = &molt->ilist[ftype];
1387 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(ftype))
1389 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1403 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1405 const gmx_moltype_t *molt;
1409 int *at2cg, cg, a, ftype, i;
1413 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1415 molt = &mtop->moltype[mtop->molblock[mb].type];
1417 if (molt->ilist[F_SETTLE].nr > 0)
1420 snew(at2cg, molt->atoms.nr);
1421 for (cg = 0; cg < cgs->nr; cg++)
1423 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1429 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1431 il = &molt->ilist[ftype];
1432 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(F_SETTLE))
1434 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1435 at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1449 /* helper functions for andersen temperature control, because the
1450 * gmx_constr construct is only defined in constr.c. Return the list
1451 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1453 extern int *get_sblock(struct gmx_constr *constr)
1455 return constr->sblock;
1458 extern int get_nblocks(struct gmx_constr *constr)
1460 return constr->nblocks;