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.
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/legacyheaders/types/commrec.h"
46 #include "gromacs/legacyheaders/constr.h"
47 #include "gromacs/legacyheaders/copyrite.h"
48 #include "gromacs/legacyheaders/mdrun.h"
49 #include "gromacs/legacyheaders/nrnb.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/txtdump.h"
53 #include "gromacs/legacyheaders/domdec.h"
54 #include "gromacs/fileio/pdbio.h"
55 #include "gromacs/legacyheaders/splitter.h"
56 #include "gromacs/topology/mtop_util.h"
57 #include "gromacs/fileio/gmxfio.h"
58 #include "gromacs/legacyheaders/macros.h"
59 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
60 #include "gromacs/essentialdynamics/edsam.h"
61 #include "gromacs/pulling/pull.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/topology/block.h"
65 #include "gromacs/topology/invblock.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/smalloc.h"
69 typedef struct gmx_constr {
70 int ncon_tot; /* The total number of constraints */
71 int nflexcon; /* The number of flexible constraints */
72 int n_at2con_mt; /* The size of at2con = #moltypes */
73 t_blocka *at2con_mt; /* A list of atoms to constraints */
74 int n_at2settle_mt; /* The size of at2settle = #moltypes */
75 int **at2settle_mt; /* A list of atoms to settles */
76 gmx_bool bInterCGsettles;
77 gmx_lincsdata_t lincsd; /* LINCS data */
78 gmx_shakedata_t shaked; /* SHAKE data */
79 gmx_settledata_t settled; /* SETTLE data */
80 int nblocks; /* The number of SHAKE blocks */
81 int *sblock; /* The SHAKE blocks */
82 int sblock_nalloc; /* The allocation size of sblock */
83 real *lagr; /* Lagrange multipliers for SHAKE */
84 int lagr_nalloc; /* The allocation size of lagr */
85 int maxwarn; /* The maximum number of warnings */
88 gmx_edsam_t ed; /* The essential dynamics data */
90 tensor *vir_r_m_dr_th; /* Thread local working data */
91 int *settle_error; /* Thread local working data */
93 gmx_mtop_t *warn_mtop; /* Only used for printing warnings */
101 /* delta_t should be used instead of ir->delta_t, to permit the time
102 step to be scaled by the calling code */
103 static void *init_vetavars(t_vetavars *vars,
104 gmx_bool constr_deriv,
105 real veta, real vetanew,
106 t_inputrec *ir, real delta_t,
107 gmx_ekindata_t *ekind, gmx_bool bPscal)
112 /* first, set the alpha integrator variable */
113 if ((ir->opts.nrdf[0] > 0) && bPscal)
115 vars->alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
121 g = 0.5*veta*delta_t;
122 vars->rscale = exp(g)*series_sinhx(g);
123 g = -0.25*vars->alpha*veta*delta_t;
124 vars->vscale = exp(g)*series_sinhx(g);
125 vars->rvscale = vars->vscale*vars->rscale;
126 vars->veta = vetanew;
130 snew(vars->vscale_nhc, ir->opts.ngtc);
131 if ((ekind == NULL) || (!bPscal))
133 for (i = 0; i < ir->opts.ngtc; i++)
135 vars->vscale_nhc[i] = 1;
140 for (i = 0; i < ir->opts.ngtc; i++)
142 vars->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
148 vars->vscale_nhc = NULL;
154 static void free_vetavars(t_vetavars *vars)
156 if (vars->vscale_nhc != NULL)
158 sfree(vars->vscale_nhc);
162 static int pcomp(const void *p1, const void *p2)
165 atom_id min1, min2, max1, max2;
166 t_sortblock *a1 = (t_sortblock *)p1;
167 t_sortblock *a2 = (t_sortblock *)p2;
169 db = a1->blocknr-a2->blocknr;
176 min1 = min(a1->iatom[1], a1->iatom[2]);
177 max1 = max(a1->iatom[1], a1->iatom[2]);
178 min2 = min(a2->iatom[1], a2->iatom[2]);
179 max2 = max(a2->iatom[1], a2->iatom[2]);
191 int n_flexible_constraints(struct gmx_constr *constr)
197 nflexcon = constr->nflexcon;
207 void too_many_constraint_warnings(int eConstrAlg, int warncount)
209 const char *abort = "- aborting to avoid logfile runaway.\n"
210 "This normally happens when your system is not sufficiently equilibrated,"
211 "or if you are changing lambda too fast in free energy simulations.\n";
214 "Too many %s warnings (%d)\n"
215 "If you know what you are doing you can %s"
216 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
217 "but normally it is better to fix the problem",
218 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
219 (eConstrAlg == econtLINCS) ?
220 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
223 static void write_constr_pdb(const char *fn, const char *title,
225 int start, int homenr, t_commrec *cr,
226 rvec x[], matrix box)
230 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
235 if (DOMAINDECOMP(cr))
238 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
245 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
249 sprintf(fname, "%s.pdb", fn);
252 out = gmx_fio_fopen(fname, "w");
254 fprintf(out, "TITLE %s\n", title);
255 gmx_write_pdb_box(out, -1, box);
256 for (i = start; i < start+homenr; i++)
260 if (i >= dd->nat_home && i < dd_ac0)
264 ii = dd->gatindex[i];
270 gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
271 gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, anm, ' ', resnm, ' ', resnr, ' ',
272 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 0.0, "");
274 fprintf(out, "TER\n");
279 static void dump_confs(FILE *fplog, gmx_int64_t step, gmx_mtop_t *mtop,
280 int start, int homenr, t_commrec *cr,
281 rvec x[], rvec xprime[], matrix box)
283 char buf[256], buf2[22];
285 char *env = getenv("GMX_SUPPRESS_DUMP");
291 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
292 write_constr_pdb(buf, "initial coordinates",
293 mtop, start, homenr, cr, x, box);
294 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
295 write_constr_pdb(buf, "coordinates after constraining",
296 mtop, start, homenr, cr, xprime, box);
299 fprintf(fplog, "Wrote pdb files with previous and current coordinates\n");
301 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
304 static void pr_sortblock(FILE *fp, const char *title, int nsb, t_sortblock sb[])
308 fprintf(fp, "%s\n", title);
309 for (i = 0; (i < nsb); i++)
311 fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
312 i, sb[i].iatom[0], sb[i].iatom[1], sb[i].iatom[2],
317 gmx_bool constrain(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
318 struct gmx_constr *constr,
319 t_idef *idef, t_inputrec *ir, gmx_ekindata_t *ekind,
321 gmx_int64_t step, int delta_step,
324 rvec *x, rvec *xprime, rvec *min_proj,
325 gmx_bool bMolPBC, matrix box,
326 real lambda, real *dvdlambda,
327 rvec *v, tensor *vir,
328 t_nrnb *nrnb, int econq, gmx_bool bPscal,
329 real veta, real vetanew)
332 int start, homenr, nrend;
334 int ncons, settle_error;
338 real invdt, vir_fac, t;
341 t_pbc pbc, *pbc_null;
346 if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
348 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");
356 nrend = start+homenr;
358 scaled_delta_t = step_scaling * ir->delta_t;
360 /* set constants for pressure control integration */
361 init_vetavars(&vetavar, econq != econqCoord,
362 veta, vetanew, ir, scaled_delta_t, ekind, bPscal);
364 /* Prepare time step for use in constraint implementations, and
365 avoid generating inf when ir->delta_t = 0. */
366 if (ir->delta_t == 0)
372 invdt = 1.0/scaled_delta_t;
375 if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
377 /* Set the constraint lengths for the step at which this configuration
378 * is meant to be. The invmasses should not be changed.
380 lambda += delta_step*ir->fepvals->delta_lambda;
385 clear_mat(vir_r_m_dr);
390 settle = &idef->il[F_SETTLE];
391 nsettle = settle->nr/(1+NRAL(F_SETTLE));
395 nth = gmx_omp_nthreads_get(emntSETTLE);
402 if (nth > 1 && constr->vir_r_m_dr_th == NULL)
404 snew(constr->vir_r_m_dr_th, nth);
405 snew(constr->settle_error, nth);
410 /* We do not need full pbc when constraints do not cross charge groups,
411 * i.e. when dd->constraint_comm==NULL.
412 * Note that PBC for constraints is different from PBC for bondeds.
413 * For constraints there is both forward and backward communication.
415 if (ir->ePBC != epbcNONE &&
416 (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm == NULL))
418 /* With pbc=screw the screw has been changed to a shift
419 * by the constraint coordinate communication routine,
420 * so that here we can use normal pbc.
422 pbc_null = set_pbc_dd(&pbc, ir->ePBC, cr->dd, FALSE, box);
429 /* Communicate the coordinates required for the non-local constraints
430 * for LINCS and/or SETTLE.
434 dd_move_x_constraints(cr->dd, box, x, xprime, econq == econqCoord);
437 if (constr->lincsd != NULL)
439 bOK = constrain_lincs(fplog, bLog, bEner, ir, step, constr->lincsd, md, cr,
441 box, pbc_null, lambda, dvdlambda,
442 invdt, v, vir != NULL, vir_r_m_dr,
444 constr->maxwarn, &constr->warncount_lincs);
445 if (!bOK && constr->maxwarn >= 0)
449 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
450 econstr_names[econtLINCS], gmx_step_str(step, buf));
456 if (constr->nblocks > 0)
461 bOK = bshakef(fplog, constr->shaked,
462 md->invmass, constr->nblocks, constr->sblock,
463 idef, ir, x, xprime, nrnb,
464 constr->lagr, lambda, dvdlambda,
465 invdt, v, vir != NULL, vir_r_m_dr,
466 constr->maxwarn >= 0, econq, &vetavar);
469 bOK = bshakef(fplog, constr->shaked,
470 md->invmass, constr->nblocks, constr->sblock,
471 idef, ir, x, min_proj, nrnb,
472 constr->lagr, lambda, dvdlambda,
473 invdt, NULL, vir != NULL, vir_r_m_dr,
474 constr->maxwarn >= 0, econq, &vetavar);
477 gmx_fatal(FARGS, "Internal error, SHAKE called for constraining something else than coordinates");
481 if (!bOK && constr->maxwarn >= 0)
485 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
486 econstr_names[econtSHAKE], gmx_step_str(step, buf));
494 int calcvir_atom_end;
498 calcvir_atom_end = 0;
502 calcvir_atom_end = md->homenr;
508 #pragma omp parallel for num_threads(nth) schedule(static)
509 for (th = 0; th < nth; th++)
511 int start_th, end_th;
515 clear_mat(constr->vir_r_m_dr_th[th]);
518 start_th = (nsettle* th )/nth;
519 end_th = (nsettle*(th+1))/nth;
520 if (start_th >= 0 && end_th - start_th > 0)
522 csettle(constr->settled,
524 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
527 invdt, v ? v[0] : NULL, calcvir_atom_end,
528 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
529 th == 0 ? &settle_error : &constr->settle_error[th],
533 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
536 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
540 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
546 case econqForceDispl:
547 #pragma omp parallel for num_threads(nth) schedule(static)
548 for (th = 0; th < nth; th++)
550 int start_th, end_th;
554 clear_mat(constr->vir_r_m_dr_th[th]);
557 start_th = (nsettle* th )/nth;
558 end_th = (nsettle*(th+1))/nth;
560 if (start_th >= 0 && end_th - start_th > 0)
562 settle_proj(constr->settled, econq,
564 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
567 xprime, min_proj, calcvir_atom_end,
568 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
572 /* This is an overestimate */
573 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
575 case econqDeriv_FlexCon:
576 /* Nothing to do, since the are no flexible constraints in settles */
579 gmx_incons("Unknown constraint quantity for settle");
585 /* Combine virial and error info of the other threads */
586 for (i = 1; i < nth; i++)
588 m_add(vir_r_m_dr, constr->vir_r_m_dr_th[i], vir_r_m_dr);
589 settle_error = constr->settle_error[i];
592 if (econq == econqCoord && settle_error >= 0)
595 if (constr->maxwarn >= 0)
599 "\nstep " "%"GMX_PRId64 ": Water molecule starting at atom %d can not be "
600 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
601 step, ddglatnr(cr->dd, settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
604 fprintf(fplog, "%s", buf);
606 fprintf(stderr, "%s", buf);
607 constr->warncount_settle++;
608 if (constr->warncount_settle > constr->maxwarn)
610 too_many_constraint_warnings(-1, constr->warncount_settle);
617 free_vetavars(&vetavar);
621 /* The normal uses of constrain() pass step_scaling = 1.0.
622 * The call to constrain() for SD1 that passes step_scaling =
623 * 0.5 also passes vir = NULL, so cannot reach this
624 * assertion. This assertion should remain until someone knows
625 * that this path works for their intended purpose, and then
626 * they can use scaled_delta_t instead of ir->delta_t
628 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
632 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
635 vir_fac = 0.5/ir->delta_t;
638 case econqForceDispl:
643 gmx_incons("Unsupported constraint quantity for virial");
648 vir_fac *= 2; /* only constraining over half the distance here */
650 for (i = 0; i < DIM; i++)
652 for (j = 0; j < DIM; j++)
654 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
661 dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
664 if (econq == econqCoord)
666 if (ir->ePull == epullCONSTRAINT)
668 if (EI_DYNAMICS(ir->eI))
670 t = ir->init_t + (step + delta_step)*ir->delta_t;
676 set_pbc(&pbc, ir->ePBC, box);
677 pull_constraint(ir->pull, md, &pbc, cr, ir->delta_t, t, x, xprime, v, *vir);
679 if (constr->ed && delta_step > 0)
681 /* apply the essential dynamcs constraints here */
682 do_edsam(ir, step, cr, xprime, v, box, constr->ed);
689 real *constr_rmsd_data(struct gmx_constr *constr)
693 return lincs_rmsd_data(constr->lincsd);
701 real constr_rmsd(struct gmx_constr *constr, gmx_bool bSD2)
705 return lincs_rmsd(constr->lincsd, bSD2);
713 static void make_shake_sblock_serial(struct gmx_constr *constr,
714 t_idef *idef, t_mdatoms *md)
723 /* Since we are processing the local topology,
724 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
726 ncons = idef->il[F_CONSTR].nr/3;
728 init_blocka(&sblocks);
729 gen_sblocks(NULL, 0, md->homenr, idef, &sblocks, FALSE);
732 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
733 nblocks=blocks->multinr[idef->nodeid] - bstart;
736 constr->nblocks = sblocks.nr;
739 fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n",
740 ncons, bstart, constr->nblocks);
743 /* Calculate block number for each atom */
744 inv_sblock = make_invblocka(&sblocks, md->nr);
746 done_blocka(&sblocks);
748 /* Store the block number in temp array and
749 * sort the constraints in order of the sblock number
750 * and the atom numbers, really sorting a segment of the array!
753 pr_idef(fplog, 0, "Before Sort", idef);
755 iatom = idef->il[F_CONSTR].iatoms;
757 for (i = 0; (i < ncons); i++, iatom += 3)
759 for (m = 0; (m < 3); m++)
761 sb[i].iatom[m] = iatom[m];
763 sb[i].blocknr = inv_sblock[iatom[1]];
766 /* Now sort the blocks */
769 pr_sortblock(debug, "Before sorting", ncons, sb);
770 fprintf(debug, "Going to sort constraints\n");
773 qsort(sb, ncons, (size_t)sizeof(*sb), pcomp);
777 pr_sortblock(debug, "After sorting", ncons, sb);
780 iatom = idef->il[F_CONSTR].iatoms;
781 for (i = 0; (i < ncons); i++, iatom += 3)
783 for (m = 0; (m < 3); m++)
785 iatom[m] = sb[i].iatom[m];
789 pr_idef(fplog, 0, "After Sort", idef);
793 snew(constr->sblock, constr->nblocks+1);
795 for (i = 0; (i < ncons); i++)
797 if (sb[i].blocknr != bnr)
800 constr->sblock[j++] = 3*i;
804 constr->sblock[j++] = 3*ncons;
806 if (j != (constr->nblocks+1))
808 fprintf(stderr, "bstart: %d\n", bstart);
809 fprintf(stderr, "j: %d, nblocks: %d, ncons: %d\n",
810 j, constr->nblocks, ncons);
811 for (i = 0; (i < ncons); i++)
813 fprintf(stderr, "i: %5d sb[i].blocknr: %5d\n", i, sb[i].blocknr);
815 for (j = 0; (j <= constr->nblocks); j++)
817 fprintf(stderr, "sblock[%3d]=%5d\n", j, (int)constr->sblock[j]);
819 gmx_fatal(FARGS, "DEATH HORROR: "
820 "sblocks does not match idef->il[F_CONSTR]");
826 static void make_shake_sblock_dd(struct gmx_constr *constr,
827 t_ilist *ilcon, t_block *cgs,
833 if (dd->ncg_home+1 > constr->sblock_nalloc)
835 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
836 srenew(constr->sblock, constr->sblock_nalloc);
840 iatom = ilcon->iatoms;
843 for (c = 0; c < ncons; c++)
845 if (c == 0 || iatom[1] >= cgs->index[cg+1])
847 constr->sblock[constr->nblocks++] = 3*c;
848 while (iatom[1] >= cgs->index[cg+1])
855 constr->sblock[constr->nblocks] = 3*ncons;
858 t_blocka make_at2con(int start, int natoms,
859 t_ilist *ilist, t_iparams *iparams,
860 gmx_bool bDynamics, int *nflexiblecons)
862 int *count, ncon, con, con_tot, nflexcon, ftype, i, a;
869 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
871 ncon = ilist[ftype].nr/3;
872 ia = ilist[ftype].iatoms;
873 for (con = 0; con < ncon; con++)
875 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
876 iparams[ia[0]].constr.dB == 0);
881 if (bDynamics || !bFlexCon)
883 for (i = 1; i < 3; i++)
892 *nflexiblecons = nflexcon;
895 at2con.nalloc_index = at2con.nr+1;
896 snew(at2con.index, at2con.nalloc_index);
898 for (a = 0; a < natoms; a++)
900 at2con.index[a+1] = at2con.index[a] + count[a];
903 at2con.nra = at2con.index[natoms];
904 at2con.nalloc_a = at2con.nra;
905 snew(at2con.a, at2con.nalloc_a);
907 /* The F_CONSTRNC constraints have constraint numbers
908 * that continue after the last F_CONSTR constraint.
911 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
913 ncon = ilist[ftype].nr/3;
914 ia = ilist[ftype].iatoms;
915 for (con = 0; con < ncon; con++)
917 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
918 iparams[ia[0]].constr.dB == 0);
919 if (bDynamics || !bFlexCon)
921 for (i = 1; i < 3; i++)
924 at2con.a[at2con.index[a]+count[a]++] = con_tot;
937 static int *make_at2settle(int natoms, const t_ilist *ilist)
943 /* Set all to no settle */
944 for (a = 0; a < natoms; a++)
949 stride = 1 + NRAL(F_SETTLE);
951 for (s = 0; s < ilist->nr; s += stride)
953 at2s[ilist->iatoms[s+1]] = s/stride;
954 at2s[ilist->iatoms[s+2]] = s/stride;
955 at2s[ilist->iatoms[s+3]] = s/stride;
961 void set_constraints(struct gmx_constr *constr,
962 gmx_localtop_t *top, t_inputrec *ir,
963 t_mdatoms *md, t_commrec *cr)
972 if (constr->ncon_tot > 0)
974 /* We are using the local topology,
975 * so there are only F_CONSTR constraints.
977 ncons = idef->il[F_CONSTR].nr/3;
979 /* With DD we might also need to call LINCS with ncons=0 for
980 * communicating coordinates to other nodes that do have constraints.
982 if (ir->eConstrAlg == econtLINCS)
984 set_lincs(idef, md, EI_DYNAMICS(ir->eI), cr, constr->lincsd);
986 if (ir->eConstrAlg == econtSHAKE)
990 make_shake_sblock_dd(constr, &idef->il[F_CONSTR], &top->cgs, cr->dd);
994 make_shake_sblock_serial(constr, idef, md);
996 if (ncons > constr->lagr_nalloc)
998 constr->lagr_nalloc = over_alloc_dd(ncons);
999 srenew(constr->lagr, constr->lagr_nalloc);
1004 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
1006 settle = &idef->il[F_SETTLE];
1007 iO = settle->iatoms[1];
1008 iH = settle->iatoms[2];
1010 settle_init(md->massT[iO], md->massT[iH],
1011 md->invmass[iO], md->invmass[iH],
1012 idef->iparams[settle->iatoms[0]].settle.doh,
1013 idef->iparams[settle->iatoms[0]].settle.dhh);
1016 /* Make a selection of the local atoms for essential dynamics */
1017 if (constr->ed && cr->dd)
1019 dd_make_local_ed_indices(cr->dd, constr->ed);
1023 static void constr_recur(t_blocka *at2con,
1024 t_ilist *ilist, t_iparams *iparams, gmx_bool bTopB,
1025 int at, int depth, int nc, int *path,
1026 real r0, real r1, real *r2max,
1038 ncon1 = ilist[F_CONSTR].nr/3;
1039 ia1 = ilist[F_CONSTR].iatoms;
1040 ia2 = ilist[F_CONSTRNC].iatoms;
1042 /* Loop over all constraints connected to this atom */
1043 for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
1046 /* Do not walk over already used constraints */
1048 for (a1 = 0; a1 < depth; a1++)
1050 if (con == path[a1])
1057 ia = constr_iatomptr(ncon1, ia1, ia2, con);
1058 /* Flexible constraints currently have length 0, which is incorrect */
1061 len = iparams[ia[0]].constr.dA;
1065 len = iparams[ia[0]].constr.dB;
1067 /* In the worst case the bond directions alternate */
1078 /* Assume angles of 120 degrees between all bonds */
1079 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
1081 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
1084 fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
1085 for (a1 = 0; a1 < depth; a1++)
1087 fprintf(debug, " %d %5.3f",
1089 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
1091 fprintf(debug, " %d %5.3f\n", con, len);
1094 /* Limit the number of recursions to 1000*nc,
1095 * so a call does not take more than a second,
1096 * even for highly connected systems.
1098 if (depth + 1 < nc && *count < 1000*nc)
1110 constr_recur(at2con, ilist, iparams,
1111 bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
1118 static real constr_r_max_moltype(gmx_moltype_t *molt, t_iparams *iparams,
1121 int natoms, nflexcon, *path, at, count;
1124 real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
1126 if (molt->ilist[F_CONSTR].nr == 0 &&
1127 molt->ilist[F_CONSTRNC].nr == 0)
1132 natoms = molt->atoms.nr;
1134 at2con = make_at2con(0, natoms, molt->ilist, iparams,
1135 EI_DYNAMICS(ir->eI), &nflexcon);
1136 snew(path, 1+ir->nProjOrder);
1137 for (at = 0; at < 1+ir->nProjOrder; at++)
1143 for (at = 0; at < natoms; at++)
1149 constr_recur(&at2con, molt->ilist, iparams,
1150 FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
1152 if (ir->efep == efepNO)
1154 rmax = sqrt(r2maxA);
1159 for (at = 0; at < natoms; at++)
1164 constr_recur(&at2con, molt->ilist, iparams,
1165 TRUE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxB, &count);
1167 lam0 = ir->fepvals->init_lambda;
1168 if (EI_DYNAMICS(ir->eI))
1170 lam0 += ir->init_step*ir->fepvals->delta_lambda;
1172 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1173 if (EI_DYNAMICS(ir->eI))
1175 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1176 rmax = max(rmax, (1 - lam1)*sqrt(r2maxA) + lam1*sqrt(r2maxB));
1180 done_blocka(&at2con);
1186 real constr_r_max(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir)
1192 for (mt = 0; mt < mtop->nmoltype; mt++)
1195 constr_r_max_moltype(&mtop->moltype[mt],
1196 mtop->ffparams.iparams, ir));
1201 fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
1207 gmx_constr_t init_constraints(FILE *fplog,
1208 gmx_mtop_t *mtop, t_inputrec *ir,
1209 gmx_edsam_t ed, t_state *state,
1212 int ncon, nset, nmol, settle_type, i, natoms, mt, nflexcon;
1213 struct gmx_constr *constr;
1216 gmx_mtop_ilistloop_t iloop;
1219 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1220 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1221 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
1223 if (ncon+nset == 0 && ir->ePull != epullCONSTRAINT && ed == NULL)
1230 constr->ncon_tot = ncon;
1231 constr->nflexcon = 0;
1234 constr->n_at2con_mt = mtop->nmoltype;
1235 snew(constr->at2con_mt, constr->n_at2con_mt);
1236 for (mt = 0; mt < mtop->nmoltype; mt++)
1238 constr->at2con_mt[mt] = make_at2con(0, mtop->moltype[mt].atoms.nr,
1239 mtop->moltype[mt].ilist,
1240 mtop->ffparams.iparams,
1241 EI_DYNAMICS(ir->eI), &nflexcon);
1242 for (i = 0; i < mtop->nmolblock; i++)
1244 if (mtop->molblock[i].type == mt)
1246 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1251 if (constr->nflexcon > 0)
1255 fprintf(fplog, "There are %d flexible constraints\n",
1257 if (ir->fc_stepsize == 0)
1260 "WARNING: step size for flexible constraining = 0\n"
1261 " All flexible constraints will be rigid.\n"
1262 " Will try to keep all flexible constraints at their original length,\n"
1263 " but the lengths may exhibit some drift.\n\n");
1264 constr->nflexcon = 0;
1267 if (constr->nflexcon > 0)
1269 please_cite(fplog, "Hess2002");
1273 if (ir->eConstrAlg == econtLINCS)
1275 constr->lincsd = init_lincs(fplog, mtop,
1276 constr->nflexcon, constr->at2con_mt,
1277 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1278 ir->nLincsIter, ir->nProjOrder);
1281 if (ir->eConstrAlg == econtSHAKE)
1283 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1285 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1287 if (constr->nflexcon)
1289 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");
1291 please_cite(fplog, "Ryckaert77a");
1294 please_cite(fplog, "Barth95a");
1297 constr->shaked = shake_init();
1303 please_cite(fplog, "Miyamoto92a");
1305 constr->bInterCGsettles = inter_charge_group_settles(mtop);
1307 /* Check that we have only one settle type */
1309 iloop = gmx_mtop_ilistloop_init(mtop);
1310 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
1312 for (i = 0; i < ilist[F_SETTLE].nr; i += 4)
1314 if (settle_type == -1)
1316 settle_type = ilist[F_SETTLE].iatoms[i];
1318 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1321 "The [molecules] section of your topology specifies more than one block of\n"
1322 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1323 "are trying to partition your solvent into different *groups* (e.g. for\n"
1324 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1325 "files specify groups. Otherwise, you may wish to change the least-used\n"
1326 "block of molecules with SETTLE constraints into 3 normal constraints.");
1331 constr->n_at2settle_mt = mtop->nmoltype;
1332 snew(constr->at2settle_mt, constr->n_at2settle_mt);
1333 for (mt = 0; mt < mtop->nmoltype; mt++)
1335 constr->at2settle_mt[mt] =
1336 make_at2settle(mtop->moltype[mt].atoms.nr,
1337 &mtop->moltype[mt].ilist[F_SETTLE]);
1341 constr->maxwarn = 999;
1342 env = getenv("GMX_MAXCONSTRWARN");
1345 constr->maxwarn = 0;
1346 sscanf(env, "%d", &constr->maxwarn);
1350 "Setting the maximum number of constraint warnings to %d\n",
1356 "Setting the maximum number of constraint warnings to %d\n",
1360 if (constr->maxwarn < 0 && fplog)
1362 fprintf(fplog, "maxwarn < 0, will not stop on constraint errors\n");
1364 constr->warncount_lincs = 0;
1365 constr->warncount_settle = 0;
1367 /* Initialize the essential dynamics sampling.
1368 * Put the pointer to the ED struct in constr */
1370 if (ed != NULL || state->edsamstate.nED > 0)
1372 init_edsam(mtop, ir, cr, ed, state->x, state->box, &state->edsamstate);
1375 constr->warn_mtop = mtop;
1380 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1382 return constr->at2con_mt;
1385 const int **atom2settle_moltype(gmx_constr_t constr)
1387 return (const int **)constr->at2settle_mt;
1391 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1393 const gmx_moltype_t *molt;
1397 int nat, *at2cg, cg, a, ftype, i;
1401 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1403 molt = &mtop->moltype[mtop->molblock[mb].type];
1405 if (molt->ilist[F_CONSTR].nr > 0 ||
1406 molt->ilist[F_CONSTRNC].nr > 0 ||
1407 molt->ilist[F_SETTLE].nr > 0)
1410 snew(at2cg, molt->atoms.nr);
1411 for (cg = 0; cg < cgs->nr; cg++)
1413 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1419 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1421 il = &molt->ilist[ftype];
1422 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(ftype))
1424 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1438 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1440 const gmx_moltype_t *molt;
1444 int nat, *at2cg, cg, a, ftype, i;
1448 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1450 molt = &mtop->moltype[mtop->molblock[mb].type];
1452 if (molt->ilist[F_SETTLE].nr > 0)
1455 snew(at2cg, molt->atoms.nr);
1456 for (cg = 0; cg < cgs->nr; cg++)
1458 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1464 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1466 il = &molt->ilist[ftype];
1467 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(F_SETTLE))
1469 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1470 at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1484 /* helper functions for andersen temperature control, because the
1485 * gmx_constr construct is only defined in constr.c. Return the list
1486 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1488 extern int *get_sblock(struct gmx_constr *constr)
1490 return constr->sblock;
1493 extern int get_nblocks(struct gmx_constr *constr)
1495 return constr->nblocks;