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.
44 #include "gromacs/legacyheaders/typedefs.h"
45 #include "gromacs/legacyheaders/types/commrec.h"
46 #include "gromacs/legacyheaders/txtdump.h"
47 #include "gromacs/legacyheaders/force.h"
48 #include "gromacs/legacyheaders/mdrun.h"
49 #include "gromacs/legacyheaders/mdatoms.h"
50 #include "gromacs/legacyheaders/vsite.h"
51 #include "gromacs/legacyheaders/network.h"
52 #include "gromacs/legacyheaders/names.h"
53 #include "gromacs/legacyheaders/constr.h"
54 #include "gromacs/legacyheaders/domdec.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/legacyheaders/shellfc.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/legacyheaders/chargegroup.h"
59 #include "gromacs/legacyheaders/macros.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/pbcutil/mshift.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/smalloc.h"
69 atom_id shell; /* The shell id */
70 atom_id nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
71 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
72 real k; /* force constant */
73 real k_1; /* 1 over force constant */
79 typedef struct gmx_shellfc {
80 int nshell_gl; /* The number of shells in the system */
81 t_shell *shell_gl; /* All the shells (for DD only) */
82 int *shell_index_gl; /* Global shell index (for DD only) */
83 gmx_bool bInterCG; /* Are there inter charge-group shells? */
84 int nshell; /* The number of local shells */
85 t_shell *shell; /* The local shells */
86 int shell_nalloc; /* The allocation size of shell */
87 gmx_bool bPredict; /* Predict shell positions */
88 gmx_bool bRequireInit; /* Require initialization of shell positions */
89 int nflexcon; /* The number of flexible constraints */
90 rvec *x[2]; /* Array for iterative minimization */
91 rvec *f[2]; /* Array for iterative minimization */
92 int x_nalloc; /* The allocation size of x and f */
93 rvec *acc_dir; /* Acceleration direction for flexcon */
94 rvec *x_old; /* Old coordinates for flexcon */
95 int flex_nalloc; /* The allocation size of acc_dir and x_old */
96 rvec *adir_xnold; /* Work space for init_adir */
97 rvec *adir_xnew; /* Work space for init_adir */
98 int adir_nalloc; /* Work space for init_adir */
102 static void pr_shell(FILE *fplog, int ns, t_shell s[])
106 fprintf(fplog, "SHELL DATA\n");
107 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
108 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
109 for (i = 0; (i < ns); i++)
111 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
114 fprintf(fplog, " %5d\n", s[i].nucl2);
116 else if (s[i].nnucl == 3)
118 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
122 fprintf(fplog, "\n");
127 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
129 real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
131 int i, m, s1, n1, n2, n3;
132 real dt_1, dt_2, dt_3, fudge, tm, m1, m2, m3;
134 gmx_mtop_atomlookup_t alook = NULL;
139 alook = gmx_mtop_atomlookup_init(mtop);
142 /* We introduce a fudge factor for performance reasons: with this choice
143 * the initial force on the shells is about a factor of two lower than
152 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
163 for (i = 0; (i < ns); i++)
174 for (m = 0; (m < DIM); m++)
176 x[s1][m] += ptr[n1][m]*dt_1;
189 /* Not the correct masses with FE, but it is just a prediction... */
194 for (m = 0; (m < DIM); m++)
196 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
211 /* Not the correct masses with FE, but it is just a prediction... */
212 gmx_mtop_atomnr_to_atom(alook, n1, &atom);
214 gmx_mtop_atomnr_to_atom(alook, n2, &atom);
216 gmx_mtop_atomnr_to_atom(alook, n3, &atom);
219 tm = dt_1/(m1+m2+m3);
220 for (m = 0; (m < DIM); m++)
222 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
226 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
232 gmx_mtop_atomlookup_destroy(alook);
236 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
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 */
291 shfc->nflexcon = nflexcon;
298 /* We have shells: fill the shell data structure */
300 /* Global system sized array, this should be avoided */
301 snew(shell_index, mtop->natoms);
303 aloop = gmx_mtop_atomloop_all_init(mtop);
305 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
307 if (atom->ptype == eptShell)
309 shell_index[i] = nshell++;
315 /* Initiate the shell structures */
316 for (i = 0; (i < nshell); i++)
318 shell[i].shell = NO_ATID;
320 shell[i].nucl1 = NO_ATID;
321 shell[i].nucl2 = NO_ATID;
322 shell[i].nucl3 = NO_ATID;
323 /* shell[i].bInterCG=FALSE; */
328 ffparams = &mtop->ffparams;
330 /* Now fill the structures */
331 shfc->bInterCG = FALSE;
334 for (mb = 0; mb < mtop->nmolblock; mb++)
336 molb = &mtop->molblock[mb];
337 molt = &mtop->moltype[molb->type];
340 snew(at2cg, molt->atoms.nr);
341 for (cg = 0; cg < cgs->nr; cg++)
343 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
349 atom = molt->atoms.atom;
350 for (mol = 0; mol < molb->nmol; mol++)
352 for (j = 0; (j < NBT); j++)
354 ia = molt->ilist[bondtypes[j]].iatoms;
355 for (i = 0; (i < molt->ilist[bondtypes[j]].nr); )
358 ftype = ffparams->functype[type];
359 nra = interaction_function[ftype].nratoms;
361 /* Check whether we have a bond with a shell */
364 switch (bondtypes[j])
371 if (atom[ia[1]].ptype == eptShell)
376 else if (atom[ia[2]].ptype == eptShell)
383 aN = ia[4]; /* Dummy */
384 aS = ia[5]; /* Shell */
387 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
394 /* Check whether one of the particles is a shell... */
395 nsi = shell_index[a_offset+aS];
396 if ((nsi < 0) || (nsi >= nshell))
398 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
401 if (shell[nsi].shell == NO_ATID)
403 shell[nsi].shell = a_offset + aS;
406 else if (shell[nsi].shell != a_offset+aS)
408 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
411 if (shell[nsi].nucl1 == NO_ATID)
413 shell[nsi].nucl1 = a_offset + aN;
415 else if (shell[nsi].nucl2 == NO_ATID)
417 shell[nsi].nucl2 = a_offset + aN;
419 else if (shell[nsi].nucl3 == NO_ATID)
421 shell[nsi].nucl3 = a_offset + aN;
427 pr_shell(fplog, ns, shell);
429 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
431 if (at2cg[aS] != at2cg[aN])
433 /* shell[nsi].bInterCG = TRUE; */
434 shfc->bInterCG = TRUE;
437 switch (bondtypes[j])
441 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
444 shell[nsi].k += ffparams->iparams[type].cubic.kb;
448 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
450 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);
452 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/
453 ffparams->iparams[type].polarize.alpha;
456 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
458 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);
460 alpha = (ffparams->iparams[type].wpol.al_x+
461 ffparams->iparams[type].wpol.al_y+
462 ffparams->iparams[type].wpol.al_z)/3.0;
463 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/alpha;
466 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
474 a_offset += molt->atoms.nr;
476 /* Done with this molecule type */
480 /* Verify whether it's all correct */
483 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
486 for (i = 0; (i < ns); i++)
488 shell[i].k_1 = 1.0/shell[i].k;
493 pr_shell(debug, ns, shell);
497 shfc->nshell_gl = ns;
498 shfc->shell_gl = shell;
499 shfc->shell_index_gl = shell_index;
501 shfc->bPredict = (getenv("GMX_NOPREDICT") == NULL);
502 shfc->bRequireInit = FALSE;
507 fprintf(fplog, "\nWill never predict shell positions\n");
512 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != NULL);
513 if (shfc->bRequireInit && fplog)
515 fprintf(fplog, "\nWill always initiate shell positions\n");
523 predict_shells(fplog, x, NULL, 0, shfc->nshell_gl, shfc->shell_gl,
531 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");
533 /* Prediction improves performance, so we should implement either:
534 * 1. communication for the atoms needed for prediction
535 * 2. prediction using the velocities of shells; currently the
536 * shell velocities are zeroed, it's a bit tricky to keep
537 * track of the shell displacements and thus the velocity.
539 shfc->bPredict = FALSE;
546 void make_local_shells(t_commrec *cr, t_mdatoms *md,
547 struct gmx_shellfc *shfc)
550 int a0, a1, *ind, nshell, i;
551 gmx_domdec_t *dd = NULL;
553 if (DOMAINDECOMP(cr))
561 /* Single node: we need all shells, just copy the pointer */
562 shfc->nshell = shfc->nshell_gl;
563 shfc->shell = shfc->shell_gl;
568 ind = shfc->shell_index_gl;
572 for (i = a0; i < a1; i++)
574 if (md->ptype[i] == eptShell)
576 if (nshell+1 > shfc->shell_nalloc)
578 shfc->shell_nalloc = over_alloc_dd(nshell+1);
579 srenew(shell, shfc->shell_nalloc);
583 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
587 shell[nshell] = shfc->shell_gl[ind[i]];
590 /* With inter-cg shells we can no do shell prediction,
591 * so we do not need the nuclei numbers.
595 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
596 if (shell[nshell].nnucl > 1)
598 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
600 if (shell[nshell].nnucl > 2)
602 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
605 shell[nshell].shell = i;
610 shfc->nshell = nshell;
614 static void do_1pos(rvec xnew, rvec xold, rvec f, real step)
632 static void do_1pos3(rvec xnew, rvec xold, rvec f, rvec step)
650 static void directional_sd(rvec xold[], rvec xnew[], rvec acc_dir[],
651 int start, int homenr, real step)
655 for (i = start; i < homenr; i++)
657 do_1pos(xnew[i], xold[i], acc_dir[i], step);
661 static void shell_pos_sd(rvec xcur[], rvec xnew[], rvec f[],
662 int ns, t_shell s[], int count)
664 const real step_scale_min = 0.8,
665 step_scale_increment = 0.2,
666 step_scale_max = 1.2,
667 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
671 real step_min, step_max;
676 for (i = 0; (i < ns); i++)
681 for (d = 0; d < DIM; d++)
683 s[i].step[d] = s[i].k_1;
685 step_min = min(step_min, s[i].step[d]);
686 step_max = max(step_max, s[i].step[d]);
692 for (d = 0; d < DIM; d++)
694 dx = xcur[shell][d] - s[i].xold[d];
695 df = f[shell][d] - s[i].fold[d];
696 /* -dx/df gets used to generate an interpolated value, but would
697 * cause a NaN if df were binary-equal to zero. Values close to
698 * zero won't cause problems (because of the min() and max()), so
699 * just testing for binary inequality is OK. */
703 /* Scale the step size by a factor interpolated from
704 * step_scale_min to step_scale_max, as k_est goes from 0 to
705 * step_scale_multiple * s[i].step[d] */
707 step_scale_min * s[i].step[d] +
708 step_scale_increment * min(step_scale_multiple * s[i].step[d], max(k_est, 0));
713 if (gmx_numzero(dx)) /* 0 == dx */
715 /* Likely this will never happen, but if it does just
716 * don't scale the step. */
720 s[i].step[d] *= step_scale_max;
724 step_min = min(step_min, s[i].step[d]);
725 step_max = max(step_max, s[i].step[d]);
729 copy_rvec(xcur[shell], s[i].xold);
730 copy_rvec(f[shell], s[i].fold);
732 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
736 fprintf(debug, "shell[%d] = %d\n", i, shell);
737 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
738 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
739 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
740 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
744 printf("step %.3e %.3e\n", step_min, step_max);
748 static void decrease_step_size(int nshell, t_shell s[])
752 for (i = 0; i < nshell; i++)
754 svmul(0.8, s[i].step, s[i].step);
758 static void print_epot(FILE *fp, gmx_int64_t mdstep, int count, real epot, real df,
759 int ndir, real sf_dir)
763 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
764 gmx_step_str(mdstep, buf), count, epot, df);
767 fprintf(fp, ", dir. rmsF: %6.2e\n", sqrt(sf_dir/ndir));
776 static real rms_force(t_commrec *cr, rvec f[], int ns, t_shell s[],
777 int ndir, real *sf_dir, real *Epot)
783 for (i = 0; i < ns; i++)
786 buf[0] += norm2(f[shell]);
795 gmx_sumd(4, buf, cr);
796 ntot = (int)(buf[1] + 0.5);
802 return (ntot ? sqrt(buf[0]/ntot) : 0);
805 static void check_pbc(FILE *fp, rvec x[], int shell)
810 for (m = 0; (m < DIM); m++)
812 if (fabs(x[shell][m]-x[now][m]) > 0.3)
814 pr_rvecs(fp, 0, "SHELL-X", x+now, 5);
820 static void dump_shells(FILE *fp, rvec x[], rvec f[], real ftol, int ns, t_shell s[])
827 for (i = 0; (i < ns); i++)
830 ff2 = iprod(f[shell], f[shell]);
833 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
834 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], sqrt(ff2));
836 check_pbc(fp, x, shell);
840 static void init_adir(FILE *log, gmx_shellfc_t shfc,
841 gmx_constr_t constr, t_idef *idef, t_inputrec *ir,
842 t_commrec *cr, int dd_ac1,
843 gmx_int64_t step, t_mdatoms *md, int start, int end,
844 rvec *x_old, rvec *x_init, rvec *x,
845 rvec *f, rvec *acc_dir,
846 gmx_bool bMolPBC, matrix box,
847 real *lambda, real *dvdlambda, t_nrnb *nrnb)
854 unsigned short *ptype;
857 if (DOMAINDECOMP(cr))
865 if (n > shfc->adir_nalloc)
867 shfc->adir_nalloc = over_alloc_dd(n);
868 srenew(shfc->adir_xnold, shfc->adir_nalloc);
869 srenew(shfc->adir_xnew, shfc->adir_nalloc);
871 xnold = shfc->adir_xnold;
872 xnew = shfc->adir_xnew;
878 /* Does NOT work with freeze or acceleration groups (yet) */
879 for (n = start; n < end; n++)
881 w_dt = md->invmass[n]*dt;
883 for (d = 0; d < DIM; d++)
885 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
887 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
888 xnew[n-start][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
892 xnold[n-start][d] = x[n][d];
893 xnew[n-start][d] = x[n][d];
897 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
898 x, xnold-start, NULL, bMolPBC, box,
899 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
900 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
901 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
902 x, xnew-start, NULL, bMolPBC, box,
903 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
904 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
906 for (n = start; n < end; n++)
908 for (d = 0; d < DIM; d++)
911 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
912 - f[n][d]*md->invmass[n];
914 clear_rvec(acc_dir[n]);
917 /* Project the acceleration on the old bond directions */
918 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
919 x_old, xnew-start, acc_dir, bMolPBC, box,
920 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
921 NULL, NULL, nrnb, econqDeriv_FlexCon, FALSE, 0, 0);
924 int relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose,
925 gmx_int64_t mdstep, t_inputrec *inputrec,
926 gmx_bool bDoNS, int force_flags,
929 gmx_enerdata_t *enerd, t_fcdata *fcd,
930 t_state *state, rvec f[],
933 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
935 gmx_groups_t *groups,
936 struct gmx_shellfc *shfc,
939 double t, rvec mu_tot,
940 gmx_bool *bConverged,
947 rvec *pos[2], *force[2], *acc_dir = NULL, *x_old = NULL;
951 real ftol, xiH, xiS, dum = 0;
953 gmx_bool bCont, bInit;
954 int nat, dd_ac0, dd_ac1 = 0, i;
955 int start = 0, homenr = md->homenr, end = start+homenr, cg0, cg1;
956 int nflexcon, g, number_steps, d, Min = 0, count = 0;
957 #define Try (1-Min) /* At start Try = 1 */
959 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
960 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
961 ftol = inputrec->em_tol;
962 number_steps = inputrec->niter;
963 nshell = shfc->nshell;
965 nflexcon = shfc->nflexcon;
969 if (DOMAINDECOMP(cr))
971 nat = dd_natoms_vsite(cr->dd);
974 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
975 nat = max(nat, dd_ac1);
983 if (nat > shfc->x_nalloc)
985 /* Allocate local arrays */
986 shfc->x_nalloc = over_alloc_dd(nat);
987 for (i = 0; (i < 2); i++)
989 srenew(shfc->x[i], shfc->x_nalloc);
990 srenew(shfc->f[i], shfc->x_nalloc);
993 for (i = 0; (i < 2); i++)
996 force[i] = shfc->f[i];
999 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1001 /* This is the only time where the coordinates are used
1002 * before do_force is called, which normally puts all
1003 * charge groups in the box.
1005 if (inputrec->cutoff_scheme == ecutsVERLET)
1007 put_atoms_in_box_omp(fr->ePBC, state->box, md->homenr, state->x);
1013 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1014 &(top->cgs), state->x, fr->cg_cm);
1019 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
1023 /* After this all coordinate arrays will contain whole charge groups */
1026 shift_self(graph, state->box, state->x);
1031 if (nat > shfc->flex_nalloc)
1033 shfc->flex_nalloc = over_alloc_dd(nat);
1034 srenew(shfc->acc_dir, shfc->flex_nalloc);
1035 srenew(shfc->x_old, shfc->flex_nalloc);
1037 acc_dir = shfc->acc_dir;
1038 x_old = shfc->x_old;
1039 for (i = 0; i < homenr; i++)
1041 for (d = 0; d < DIM; d++)
1044 state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
1049 /* Do a prediction of the shell positions */
1050 if (shfc->bPredict && !bCont)
1052 predict_shells(fplog, state->x, state->v, inputrec->delta_t, nshell, shell,
1053 md->massT, NULL, bInit);
1056 /* do_force expected the charge groups to be in the box */
1059 unshift_self(graph, state->box, state->x);
1062 /* Calculate the forces first time around */
1065 pr_rvecs(debug, 0, "x b4 do_force", state->x + start, homenr);
1067 do_force(fplog, cr, inputrec, mdstep, nrnb, wcycle, top, groups,
1068 state->box, state->x, &state->hist,
1069 force[Min], force_vir, md, enerd, fcd,
1070 state->lambda, graph,
1071 fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1072 (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
1077 init_adir(fplog, shfc,
1078 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1079 shfc->x_old-start, state->x, state->x, force[Min],
1080 shfc->acc_dir-start,
1081 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1083 for (i = start; i < end; i++)
1085 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
1089 Epot[Min] = enerd->term[F_EPOT];
1091 df[Min] = rms_force(cr, shfc->f[Min], nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1095 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1100 pr_rvecs(debug, 0, "force0", force[Min], md->nr);
1103 if (nshell+nflexcon > 0)
1105 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1106 * shell positions are updated, therefore the other particles must
1109 memcpy(pos[Min], state->x, nat*sizeof(state->x[0]));
1110 memcpy(pos[Try], state->x, nat*sizeof(state->x[0]));
1113 if (bVerbose && MASTER(cr))
1115 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1120 fprintf(debug, "%17s: %14.10e\n",
1121 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1122 fprintf(debug, "%17s: %14.10e\n",
1123 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1124 fprintf(debug, "%17s: %14.10e\n",
1125 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1126 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1129 /* First check whether we should do shells, or whether the force is
1130 * low enough even without minimization.
1132 *bConverged = (df[Min] < ftol);
1134 for (count = 1; (!(*bConverged) && (count < number_steps)); count++)
1138 construct_vsites(vsite, pos[Min], inputrec->delta_t, state->v,
1139 idef->iparams, idef->il,
1140 fr->ePBC, fr->bMolPBC, cr, state->box);
1145 init_adir(fplog, shfc,
1146 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1147 x_old-start, state->x, pos[Min], force[Min], acc_dir-start,
1148 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1150 directional_sd(pos[Min], pos[Try], acc_dir-start, start, end,
1154 /* New positions, Steepest descent */
1155 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1157 /* do_force expected the charge groups to be in the box */
1160 unshift_self(graph, state->box, pos[Try]);
1165 pr_rvecs(debug, 0, "RELAX: pos[Min] ", pos[Min] + start, homenr);
1166 pr_rvecs(debug, 0, "RELAX: pos[Try] ", pos[Try] + start, homenr);
1168 /* Try the new positions */
1169 do_force(fplog, cr, inputrec, 1, nrnb, wcycle,
1170 top, groups, state->box, pos[Try], &state->hist,
1171 force[Try], force_vir,
1172 md, enerd, fcd, state->lambda, graph,
1173 fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1178 pr_rvecs(debug, 0, "RELAX: force[Min]", force[Min] + start, homenr);
1179 pr_rvecs(debug, 0, "RELAX: force[Try]", force[Try] + start, homenr);
1184 init_adir(fplog, shfc,
1185 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1186 x_old-start, state->x, pos[Try], force[Try], acc_dir-start,
1187 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1189 for (i = start; i < end; i++)
1191 sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
1195 Epot[Try] = enerd->term[F_EPOT];
1197 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1201 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1208 pr_rvecs(debug, 0, "F na do_force", force[Try] + start, homenr);
1212 fprintf(debug, "SHELL ITER %d\n", count);
1213 dump_shells(debug, pos[Try], force[Try], ftol, nshell, shell);
1217 if (bVerbose && MASTER(cr))
1219 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1222 *bConverged = (df[Try] < ftol);
1224 if ((df[Try] < df[Min]))
1228 fprintf(debug, "Swapping Min and Try\n");
1232 /* Correct the velocities for the flexible constraints */
1233 invdt = 1/inputrec->delta_t;
1234 for (i = start; i < end; i++)
1236 for (d = 0; d < DIM; d++)
1238 state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1246 decrease_step_size(nshell, shell);
1249 if (MASTER(cr) && !(*bConverged))
1251 /* Note that the energies and virial are incorrect when not converged */
1255 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1256 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1259 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1260 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1263 /* Copy back the coordinates and the forces */
1264 memcpy(state->x, pos[Min], nat*sizeof(state->x[0]));
1265 memcpy(f, force[Min], nat*sizeof(f[0]));