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-2008, 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.
45 #include "types/commrec.h"
55 #include "gromacs/math/units.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "chargegroup.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/pbcutil/mshift.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/smalloc.h"
68 atom_id shell; /* The shell id */
69 atom_id nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
70 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
71 real k; /* force constant */
72 real k_1; /* 1 over force constant */
78 typedef struct gmx_shellfc {
79 int nshell_gl; /* The number of shells in the system */
80 t_shell *shell_gl; /* All the shells (for DD only) */
81 int *shell_index_gl; /* Global shell index (for DD only) */
82 gmx_bool bInterCG; /* Are there inter charge-group shells? */
83 int nshell; /* The number of local shells */
84 t_shell *shell; /* The local shells */
85 int shell_nalloc; /* The allocation size of shell */
86 gmx_bool bPredict; /* Predict shell positions */
87 gmx_bool bRequireInit; /* Require initialization of shell positions */
88 int nflexcon; /* The number of flexible constraints */
89 rvec *x[2]; /* Array for iterative minimization */
90 rvec *f[2]; /* Array for iterative minimization */
91 int x_nalloc; /* The allocation size of x and f */
92 rvec *acc_dir; /* Acceleration direction for flexcon */
93 rvec *x_old; /* Old coordinates for flexcon */
94 int flex_nalloc; /* The allocation size of acc_dir and x_old */
95 rvec *adir_xnold; /* Work space for init_adir */
96 rvec *adir_xnew; /* Work space for init_adir */
97 int adir_nalloc; /* Work space for init_adir */
101 static void pr_shell(FILE *fplog, int ns, t_shell s[])
105 fprintf(fplog, "SHELL DATA\n");
106 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
107 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
108 for (i = 0; (i < ns); i++)
110 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
113 fprintf(fplog, " %5d\n", s[i].nucl2);
115 else if (s[i].nnucl == 3)
117 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
121 fprintf(fplog, "\n");
126 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
128 real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
130 int i, m, s1, n1, n2, n3;
131 real dt_1, dt_2, dt_3, fudge, tm, m1, m2, m3;
133 gmx_mtop_atomlookup_t alook = NULL;
138 alook = gmx_mtop_atomlookup_init(mtop);
141 /* We introduce a fudge factor for performance reasons: with this choice
142 * the initial force on the shells is about a factor of two lower than
151 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
162 for (i = 0; (i < ns); i++)
173 for (m = 0; (m < DIM); m++)
175 x[s1][m] += ptr[n1][m]*dt_1;
188 /* Not the correct masses with FE, but it is just a prediction... */
193 for (m = 0; (m < DIM); m++)
195 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
210 /* Not the correct masses with FE, but it is just a prediction... */
211 gmx_mtop_atomnr_to_atom(alook, n1, &atom);
213 gmx_mtop_atomnr_to_atom(alook, n2, &atom);
215 gmx_mtop_atomnr_to_atom(alook, n3, &atom);
218 tm = dt_1/(m1+m2+m3);
219 for (m = 0; (m < DIM); m++)
221 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
225 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
231 gmx_mtop_atomlookup_destroy(alook);
235 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
236 gmx_bool bCutoffSchemeIsVerlet,
237 gmx_mtop_t *mtop, int nflexcon,
240 struct gmx_shellfc *shfc;
242 int *shell_index = NULL, *at2cg;
244 int n[eptNR], ns, nshell, nsi;
245 int i, j, nmol, type, mb, mt, a_offset, cg, mol, ftype, nra;
247 int aS, aN = 0; /* Shell and nucleus */
248 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
249 #define NBT asize(bondtypes)
251 gmx_mtop_atomloop_block_t aloopb;
252 gmx_mtop_atomloop_all_t aloop;
253 gmx_ffparams_t *ffparams;
254 gmx_molblock_t *molb;
258 /* Count number of shells, and find their indices */
259 for (i = 0; (i < eptNR); i++)
264 aloopb = gmx_mtop_atomloop_block_init(mtop);
265 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
267 n[atom->ptype] += nmol;
272 /* Print the number of each particle type */
273 for (i = 0; (i < eptNR); i++)
277 fprintf(fplog, "There are: %d %ss\n", n[i], ptype_str[i]);
282 nshell = n[eptShell];
284 if (nshell == 0 && nflexcon == 0)
286 /* We're not doing shells or flexible constraints */
290 if (bCutoffSchemeIsVerlet)
292 gmx_fatal(FARGS, "The shell code does not work with the Verlet cut-off scheme.\n");
296 shfc->nflexcon = nflexcon;
303 /* We have shells: fill the shell data structure */
305 /* Global system sized array, this should be avoided */
306 snew(shell_index, mtop->natoms);
308 aloop = gmx_mtop_atomloop_all_init(mtop);
310 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
312 if (atom->ptype == eptShell)
314 shell_index[i] = nshell++;
320 /* Initiate the shell structures */
321 for (i = 0; (i < nshell); i++)
323 shell[i].shell = NO_ATID;
325 shell[i].nucl1 = NO_ATID;
326 shell[i].nucl2 = NO_ATID;
327 shell[i].nucl3 = NO_ATID;
328 /* shell[i].bInterCG=FALSE; */
333 ffparams = &mtop->ffparams;
335 /* Now fill the structures */
336 shfc->bInterCG = FALSE;
339 for (mb = 0; mb < mtop->nmolblock; mb++)
341 molb = &mtop->molblock[mb];
342 molt = &mtop->moltype[molb->type];
345 snew(at2cg, molt->atoms.nr);
346 for (cg = 0; cg < cgs->nr; cg++)
348 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
354 atom = molt->atoms.atom;
355 for (mol = 0; mol < molb->nmol; mol++)
357 for (j = 0; (j < NBT); j++)
359 ia = molt->ilist[bondtypes[j]].iatoms;
360 for (i = 0; (i < molt->ilist[bondtypes[j]].nr); )
363 ftype = ffparams->functype[type];
364 nra = interaction_function[ftype].nratoms;
366 /* Check whether we have a bond with a shell */
369 switch (bondtypes[j])
376 if (atom[ia[1]].ptype == eptShell)
381 else if (atom[ia[2]].ptype == eptShell)
388 aN = ia[4]; /* Dummy */
389 aS = ia[5]; /* Shell */
392 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
399 /* Check whether one of the particles is a shell... */
400 nsi = shell_index[a_offset+aS];
401 if ((nsi < 0) || (nsi >= nshell))
403 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
406 if (shell[nsi].shell == NO_ATID)
408 shell[nsi].shell = a_offset + aS;
411 else if (shell[nsi].shell != a_offset+aS)
413 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
416 if (shell[nsi].nucl1 == NO_ATID)
418 shell[nsi].nucl1 = a_offset + aN;
420 else if (shell[nsi].nucl2 == NO_ATID)
422 shell[nsi].nucl2 = a_offset + aN;
424 else if (shell[nsi].nucl3 == NO_ATID)
426 shell[nsi].nucl3 = a_offset + aN;
432 pr_shell(fplog, ns, shell);
434 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
436 if (at2cg[aS] != at2cg[aN])
438 /* shell[nsi].bInterCG = TRUE; */
439 shfc->bInterCG = TRUE;
442 switch (bondtypes[j])
446 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
449 shell[nsi].k += ffparams->iparams[type].cubic.kb;
453 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
455 gmx_fatal(FARGS, "polarize can not be used with qA(%e) != qB(%e) for atom %d of molecule block %d", qS, atom[aS].qB, aS+1, mb+1);
457 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/
458 ffparams->iparams[type].polarize.alpha;
461 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
463 gmx_fatal(FARGS, "water_pol can not be used with qA(%e) != qB(%e) for atom %d of molecule block %d", qS, atom[aS].qB, aS+1, mb+1);
465 alpha = (ffparams->iparams[type].wpol.al_x+
466 ffparams->iparams[type].wpol.al_y+
467 ffparams->iparams[type].wpol.al_z)/3.0;
468 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/alpha;
471 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
479 a_offset += molt->atoms.nr;
481 /* Done with this molecule type */
485 /* Verify whether it's all correct */
488 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
491 for (i = 0; (i < ns); i++)
493 shell[i].k_1 = 1.0/shell[i].k;
498 pr_shell(debug, ns, shell);
502 shfc->nshell_gl = ns;
503 shfc->shell_gl = shell;
504 shfc->shell_index_gl = shell_index;
506 shfc->bPredict = (getenv("GMX_NOPREDICT") == NULL);
507 shfc->bRequireInit = FALSE;
512 fprintf(fplog, "\nWill never predict shell positions\n");
517 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != NULL);
518 if (shfc->bRequireInit && fplog)
520 fprintf(fplog, "\nWill always initiate shell positions\n");
528 predict_shells(fplog, x, NULL, 0, shfc->nshell_gl, shfc->shell_gl,
536 fprintf(fplog, "\nNOTE: there all shells that are connected to particles outside thier own charge group, will not predict shells positions during the run\n\n");
538 shfc->bPredict = FALSE;
545 void make_local_shells(t_commrec *cr, t_mdatoms *md,
546 struct gmx_shellfc *shfc)
549 int a0, a1, *ind, nshell, i;
550 gmx_domdec_t *dd = NULL;
552 if (DOMAINDECOMP(cr))
560 /* Single node: we need all shells, just copy the pointer */
561 shfc->nshell = shfc->nshell_gl;
562 shfc->shell = shfc->shell_gl;
567 ind = shfc->shell_index_gl;
571 for (i = a0; i < a1; i++)
573 if (md->ptype[i] == eptShell)
575 if (nshell+1 > shfc->shell_nalloc)
577 shfc->shell_nalloc = over_alloc_dd(nshell+1);
578 srenew(shell, shfc->shell_nalloc);
582 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
586 shell[nshell] = shfc->shell_gl[ind[i]];
589 /* With inter-cg shells we can no do shell prediction,
590 * so we do not need the nuclei numbers.
594 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
595 if (shell[nshell].nnucl > 1)
597 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
599 if (shell[nshell].nnucl > 2)
601 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
604 shell[nshell].shell = i;
609 shfc->nshell = nshell;
613 static void do_1pos(rvec xnew, rvec xold, rvec f, real step)
631 static void do_1pos3(rvec xnew, rvec xold, rvec f, rvec step)
649 static void directional_sd(rvec xold[], rvec xnew[], rvec acc_dir[],
650 int start, int homenr, real step)
654 for (i = start; i < homenr; i++)
656 do_1pos(xnew[i], xold[i], acc_dir[i], step);
660 static void shell_pos_sd(rvec xcur[], rvec xnew[], rvec f[],
661 int ns, t_shell s[], int count)
663 const real step_scale_min = 0.8,
664 step_scale_increment = 0.2,
665 step_scale_max = 1.2,
666 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
670 real step_min, step_max;
675 for (i = 0; (i < ns); i++)
680 for (d = 0; d < DIM; d++)
682 s[i].step[d] = s[i].k_1;
684 step_min = min(step_min, s[i].step[d]);
685 step_max = max(step_max, s[i].step[d]);
691 for (d = 0; d < DIM; d++)
693 dx = xcur[shell][d] - s[i].xold[d];
694 df = f[shell][d] - s[i].fold[d];
695 /* -dx/df gets used to generate an interpolated value, but would
696 * cause a NaN if df were binary-equal to zero. Values close to
697 * zero won't cause problems (because of the min() and max()), so
698 * just testing for binary inequality is OK. */
702 /* Scale the step size by a factor interpolated from
703 * step_scale_min to step_scale_max, as k_est goes from 0 to
704 * step_scale_multiple * s[i].step[d] */
706 step_scale_min * s[i].step[d] +
707 step_scale_increment * min(step_scale_multiple * s[i].step[d], max(k_est, 0));
712 if (gmx_numzero(dx)) /* 0 == dx */
714 /* Likely this will never happen, but if it does just
715 * don't scale the step. */
719 s[i].step[d] *= step_scale_max;
723 step_min = min(step_min, s[i].step[d]);
724 step_max = max(step_max, s[i].step[d]);
728 copy_rvec(xcur[shell], s[i].xold);
729 copy_rvec(f[shell], s[i].fold);
731 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
735 fprintf(debug, "shell[%d] = %d\n", i, shell);
736 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
737 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
738 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
739 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
743 printf("step %.3e %.3e\n", step_min, step_max);
747 static void decrease_step_size(int nshell, t_shell s[])
751 for (i = 0; i < nshell; i++)
753 svmul(0.8, s[i].step, s[i].step);
757 static void print_epot(FILE *fp, gmx_int64_t mdstep, int count, real epot, real df,
758 int ndir, real sf_dir)
762 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
763 gmx_step_str(mdstep, buf), count, epot, df);
766 fprintf(fp, ", dir. rmsF: %6.2e\n", sqrt(sf_dir/ndir));
775 static real rms_force(t_commrec *cr, rvec f[], int ns, t_shell s[],
776 int ndir, real *sf_dir, real *Epot)
782 for (i = 0; i < ns; i++)
785 buf[0] += norm2(f[shell]);
794 gmx_sumd(4, buf, cr);
795 ntot = (int)(buf[1] + 0.5);
801 return (ntot ? sqrt(buf[0]/ntot) : 0);
804 static void check_pbc(FILE *fp, rvec x[], int shell)
809 for (m = 0; (m < DIM); m++)
811 if (fabs(x[shell][m]-x[now][m]) > 0.3)
813 pr_rvecs(fp, 0, "SHELL-X", x+now, 5);
819 static void dump_shells(FILE *fp, rvec x[], rvec f[], real ftol, int ns, t_shell s[])
826 for (i = 0; (i < ns); i++)
829 ff2 = iprod(f[shell], f[shell]);
832 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
833 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], sqrt(ff2));
835 check_pbc(fp, x, shell);
839 static void init_adir(FILE *log, gmx_shellfc_t shfc,
840 gmx_constr_t constr, t_idef *idef, t_inputrec *ir,
841 t_commrec *cr, int dd_ac1,
842 gmx_int64_t step, t_mdatoms *md, int start, int end,
843 rvec *x_old, rvec *x_init, rvec *x,
844 rvec *f, rvec *acc_dir,
845 gmx_bool bMolPBC, matrix box,
846 real *lambda, real *dvdlambda, t_nrnb *nrnb)
853 unsigned short *ptype;
856 if (DOMAINDECOMP(cr))
864 if (n > shfc->adir_nalloc)
866 shfc->adir_nalloc = over_alloc_dd(n);
867 srenew(shfc->adir_xnold, shfc->adir_nalloc);
868 srenew(shfc->adir_xnew, shfc->adir_nalloc);
870 xnold = shfc->adir_xnold;
871 xnew = shfc->adir_xnew;
877 /* Does NOT work with freeze or acceleration groups (yet) */
878 for (n = start; n < end; n++)
880 w_dt = md->invmass[n]*dt;
882 for (d = 0; d < DIM; d++)
884 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
886 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
887 xnew[n-start][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
891 xnold[n-start][d] = x[n][d];
892 xnew[n-start][d] = x[n][d];
896 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
897 x, xnold-start, NULL, bMolPBC, box,
898 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
899 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
900 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
901 x, xnew-start, NULL, bMolPBC, box,
902 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
903 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
905 for (n = start; n < end; n++)
907 for (d = 0; d < DIM; d++)
910 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
911 - f[n][d]*md->invmass[n];
913 clear_rvec(acc_dir[n]);
916 /* Project the acceleration on the old bond directions */
917 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
918 x_old, xnew-start, acc_dir, bMolPBC, box,
919 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
920 NULL, NULL, nrnb, econqDeriv_FlexCon, FALSE, 0, 0);
923 int relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose,
924 gmx_int64_t mdstep, t_inputrec *inputrec,
925 gmx_bool bDoNS, int force_flags,
928 gmx_enerdata_t *enerd, t_fcdata *fcd,
929 t_state *state, rvec f[],
932 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
934 gmx_groups_t *groups,
935 struct gmx_shellfc *shfc,
938 double t, rvec mu_tot,
939 gmx_bool *bConverged,
946 rvec *pos[2], *force[2], *acc_dir = NULL, *x_old = NULL;
950 real ftol, xiH, xiS, dum = 0;
952 gmx_bool bCont, bInit;
953 int nat, dd_ac0, dd_ac1 = 0, i;
954 int start = 0, homenr = md->homenr, end = start+homenr, cg0, cg1;
955 int nflexcon, g, number_steps, d, Min = 0, count = 0;
956 #define Try (1-Min) /* At start Try = 1 */
958 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
959 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
960 ftol = inputrec->em_tol;
961 number_steps = inputrec->niter;
962 nshell = shfc->nshell;
964 nflexcon = shfc->nflexcon;
968 if (DOMAINDECOMP(cr))
970 nat = dd_natoms_vsite(cr->dd);
973 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
974 nat = max(nat, dd_ac1);
982 if (nat > shfc->x_nalloc)
984 /* Allocate local arrays */
985 shfc->x_nalloc = over_alloc_dd(nat);
986 for (i = 0; (i < 2); i++)
988 srenew(shfc->x[i], shfc->x_nalloc);
989 srenew(shfc->f[i], shfc->x_nalloc);
992 for (i = 0; (i < 2); i++)
995 force[i] = shfc->f[i];
998 /* When we had particle decomposition, this code only worked with
999 * PD when all particles involved with each shell were in the same
1000 * charge group. Not sure if this is still relevant. */
1001 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1003 /* This is the only time where the coordinates are used
1004 * before do_force is called, which normally puts all
1005 * charge groups in the box.
1009 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1010 &(top->cgs), state->x, fr->cg_cm);
1013 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
1017 /* After this all coordinate arrays will contain whole molecules */
1020 shift_self(graph, state->box, state->x);
1025 if (nat > shfc->flex_nalloc)
1027 shfc->flex_nalloc = over_alloc_dd(nat);
1028 srenew(shfc->acc_dir, shfc->flex_nalloc);
1029 srenew(shfc->x_old, shfc->flex_nalloc);
1031 acc_dir = shfc->acc_dir;
1032 x_old = shfc->x_old;
1033 for (i = 0; i < homenr; i++)
1035 for (d = 0; d < DIM; d++)
1038 state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
1043 /* Do a prediction of the shell positions */
1044 if (shfc->bPredict && !bCont)
1046 predict_shells(fplog, state->x, state->v, inputrec->delta_t, nshell, shell,
1047 md->massT, NULL, bInit);
1050 /* do_force expected the charge groups to be in the box */
1053 unshift_self(graph, state->box, state->x);
1056 /* Calculate the forces first time around */
1059 pr_rvecs(debug, 0, "x b4 do_force", state->x + start, homenr);
1061 do_force(fplog, cr, inputrec, mdstep, nrnb, wcycle, top, groups,
1062 state->box, state->x, &state->hist,
1063 force[Min], force_vir, md, enerd, fcd,
1064 state->lambda, graph,
1065 fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1066 (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
1071 init_adir(fplog, shfc,
1072 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1073 shfc->x_old-start, state->x, state->x, force[Min],
1074 shfc->acc_dir-start,
1075 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1077 for (i = start; i < end; i++)
1079 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
1083 Epot[Min] = enerd->term[F_EPOT];
1085 df[Min] = rms_force(cr, shfc->f[Min], nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1089 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1094 pr_rvecs(debug, 0, "force0", force[Min], md->nr);
1097 if (nshell+nflexcon > 0)
1099 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1100 * shell positions are updated, therefore the other particles must
1103 memcpy(pos[Min], state->x, nat*sizeof(state->x[0]));
1104 memcpy(pos[Try], state->x, nat*sizeof(state->x[0]));
1107 if (bVerbose && MASTER(cr))
1109 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1114 fprintf(debug, "%17s: %14.10e\n",
1115 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1116 fprintf(debug, "%17s: %14.10e\n",
1117 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1118 fprintf(debug, "%17s: %14.10e\n",
1119 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1120 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1123 /* First check whether we should do shells, or whether the force is
1124 * low enough even without minimization.
1126 *bConverged = (df[Min] < ftol);
1128 for (count = 1; (!(*bConverged) && (count < number_steps)); count++)
1132 construct_vsites(vsite, pos[Min], inputrec->delta_t, state->v,
1133 idef->iparams, idef->il,
1134 fr->ePBC, fr->bMolPBC, cr, state->box);
1139 init_adir(fplog, shfc,
1140 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1141 x_old-start, state->x, pos[Min], force[Min], acc_dir-start,
1142 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1144 directional_sd(pos[Min], pos[Try], acc_dir-start, start, end,
1148 /* New positions, Steepest descent */
1149 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1151 /* do_force expected the charge groups to be in the box */
1154 unshift_self(graph, state->box, pos[Try]);
1159 pr_rvecs(debug, 0, "RELAX: pos[Min] ", pos[Min] + start, homenr);
1160 pr_rvecs(debug, 0, "RELAX: pos[Try] ", pos[Try] + start, homenr);
1162 /* Try the new positions */
1163 do_force(fplog, cr, inputrec, 1, nrnb, wcycle,
1164 top, groups, state->box, pos[Try], &state->hist,
1165 force[Try], force_vir,
1166 md, enerd, fcd, state->lambda, graph,
1167 fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1172 pr_rvecs(debug, 0, "RELAX: force[Min]", force[Min] + start, homenr);
1173 pr_rvecs(debug, 0, "RELAX: force[Try]", force[Try] + start, homenr);
1178 init_adir(fplog, shfc,
1179 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1180 x_old-start, state->x, pos[Try], force[Try], acc_dir-start,
1181 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1183 for (i = start; i < end; i++)
1185 sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
1189 Epot[Try] = enerd->term[F_EPOT];
1191 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1195 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1202 pr_rvecs(debug, 0, "F na do_force", force[Try] + start, homenr);
1206 fprintf(debug, "SHELL ITER %d\n", count);
1207 dump_shells(debug, pos[Try], force[Try], ftol, nshell, shell);
1211 if (bVerbose && MASTER(cr))
1213 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1216 *bConverged = (df[Try] < ftol);
1218 if ((df[Try] < df[Min]))
1222 fprintf(debug, "Swapping Min and Try\n");
1226 /* Correct the velocities for the flexible constraints */
1227 invdt = 1/inputrec->delta_t;
1228 for (i = start; i < end; i++)
1230 for (d = 0; d < DIM; d++)
1232 state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1240 decrease_step_size(nshell, shell);
1243 if (MASTER(cr) && !(*bConverged))
1245 /* Note that the energies and virial are incorrect when not converged */
1249 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1250 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1253 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1254 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1257 /* Copy back the coordinates and the forces */
1258 memcpy(state->x, pos[Min], nat*sizeof(state->x[0]));
1259 memcpy(f, force[Min], nat*sizeof(f[0]));