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,2015,2016,2017,2018, 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.
49 #include "gromacs/domdec/dlbtiming.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/gmxlib/chargegroup.h"
53 #include "gromacs/gmxlib/network.h"
54 #include "gromacs/math/functions.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/math/vecdump.h"
58 #include "gromacs/mdlib/constr.h"
59 #include "gromacs/mdlib/force.h"
60 #include "gromacs/mdlib/mdrun.h"
61 #include "gromacs/mdlib/sim_util.h"
62 #include "gromacs/mdlib/vsite.h"
63 #include "gromacs/mdtypes/commrec.h"
64 #include "gromacs/mdtypes/inputrec.h"
65 #include "gromacs/mdtypes/md_enums.h"
66 #include "gromacs/mdtypes/state.h"
67 #include "gromacs/pbcutil/mshift.h"
68 #include "gromacs/pbcutil/pbc.h"
69 #include "gromacs/topology/mtop_lookup.h"
70 #include "gromacs/topology/mtop_util.h"
71 #include "gromacs/utility/arrayref.h"
72 #include "gromacs/utility/arraysize.h"
73 #include "gromacs/utility/cstringutil.h"
74 #include "gromacs/utility/fatalerror.h"
75 #include "gromacs/utility/smalloc.h"
79 int shell; /* The shell id */
80 int nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
81 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
82 real k; /* force constant */
83 real k_1; /* 1 over force constant */
89 struct gmx_shellfc_t {
90 /* Shell counts, indices, parameters and working data */
91 int nshell_gl; /* The number of shells in the system */
92 t_shell *shell_gl; /* All the shells (for DD only) */
93 int *shell_index_gl; /* Global shell index (for DD only) */
94 gmx_bool bInterCG; /* Are there inter charge-group shells? */
95 int nshell; /* The number of local shells */
96 t_shell *shell; /* The local shells */
97 int shell_nalloc; /* The allocation size of shell */
98 gmx_bool bPredict; /* Predict shell positions */
99 gmx_bool bRequireInit; /* Require initialization of shell positions */
100 int nflexcon; /* The number of flexible constraints */
102 /* Temporary arrays, should be fixed size 2 when fully converted to C++ */
103 PaddedRVecVector *x; /* Array for iterative minimization */
104 PaddedRVecVector *f; /* Array for iterative minimization */
106 /* Flexible constraint working data */
107 rvec *acc_dir; /* Acceleration direction for flexcon */
108 rvec *x_old; /* Old coordinates for flexcon */
109 int flex_nalloc; /* The allocation size of acc_dir and x_old */
110 rvec *adir_xnold; /* Work space for init_adir */
111 rvec *adir_xnew; /* Work space for init_adir */
112 int adir_nalloc; /* Work space for init_adir */
113 std::int64_t numForceEvaluations; /* Total number of force evaluations */
114 int numConvergedIterations; /* Total number of iterations that converged */
118 static void pr_shell(FILE *fplog, int ns, t_shell s[])
122 fprintf(fplog, "SHELL DATA\n");
123 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
124 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
125 for (i = 0; (i < ns); i++)
127 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
130 fprintf(fplog, " %5d\n", s[i].nucl2);
132 else if (s[i].nnucl == 3)
134 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
138 fprintf(fplog, "\n");
143 /* TODO The remain call of this function passes non-NULL mass and NULL
144 * mtop, so this routine can be simplified.
146 * The other code path supported doing prediction before the MD loop
147 * started, but even when called, the prediction was always
148 * over-written by a subsequent call in the MD loop, so has been
150 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
152 real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
154 int i, m, s1, n1, n2, n3;
155 real dt_1, fudge, tm, m1, m2, m3;
158 /* We introduce a fudge factor for performance reasons: with this choice
159 * the initial force on the shells is about a factor of two lower than
168 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
180 for (i = 0; (i < ns); i++)
191 for (m = 0; (m < DIM); m++)
193 x[s1][m] += ptr[n1][m]*dt_1;
206 /* Not the correct masses with FE, but it is just a prediction... */
207 m1 = mtopGetAtomMass(mtop, n1, &molb);
208 m2 = mtopGetAtomMass(mtop, n2, &molb);
211 for (m = 0; (m < DIM); m++)
213 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
228 /* Not the correct masses with FE, but it is just a prediction... */
229 m1 = mtopGetAtomMass(mtop, n1, &molb);
230 m2 = mtopGetAtomMass(mtop, n2, &molb);
231 m3 = mtopGetAtomMass(mtop, n3, &molb);
233 tm = dt_1/(m1+m2+m3);
234 for (m = 0; (m < DIM); m++)
236 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
240 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
245 /*! \brief Count the different particle types in a system
247 * Routine prints a warning to stderr in case an unknown particle type
249 * \param[in] fplog Print what we have found if not NULL
250 * \param[in] mtop Molecular topology.
251 * \returns Array holding the number of particles of a type
253 static std::array<int, eptNR> countPtypes(FILE *fplog,
256 std::array<int, eptNR> nptype = { { 0 } };
257 /* Count number of shells, and find their indices */
258 for (int i = 0; (i < eptNR); i++)
263 gmx_mtop_atomloop_block_t aloopb = gmx_mtop_atomloop_block_init(mtop);
266 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
273 nptype[atom->ptype] += nmol;
276 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
277 static_cast<int>(atom->ptype));
282 /* Print the number of each particle type */
284 for (const auto &i : nptype)
288 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
296 gmx_shellfc_t *init_shell_flexcon(FILE *fplog,
297 gmx_mtop_t *mtop, int nflexcon,
299 bool usingDomainDecomposition)
303 int *shell_index = nullptr, *at2cg;
307 int i, j, type, mb, a_offset, cg, mol, ftype, nra;
309 int aS, aN = 0; /* Shell and nucleus */
310 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
311 #define NBT asize(bondtypes)
313 gmx_mtop_atomloop_all_t aloop;
314 gmx_ffparams_t *ffparams;
315 gmx_molblock_t *molb;
319 std::array<int, eptNR> n = countPtypes(fplog, mtop);
320 nshell = n[eptShell];
322 if (nshell == 0 && nflexcon == 0)
324 /* We're not doing shells or flexible constraints */
329 shfc->x = new PaddedRVecVector[2] {};
330 shfc->f = new PaddedRVecVector[2] {};
331 shfc->nflexcon = nflexcon;
335 /* Only flexible constraints, no shells.
336 * Note that make_local_shells() does not need to be called.
339 shfc->bPredict = FALSE;
344 if (nstcalcenergy != 1)
346 gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combination with shell particles.\nPlease make a new tpr file.", nstcalcenergy);
348 if (usingDomainDecomposition)
350 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
353 /* We have shells: fill the shell data structure */
355 /* Global system sized array, this should be avoided */
356 snew(shell_index, mtop->natoms);
358 aloop = gmx_mtop_atomloop_all_init(mtop);
360 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
362 if (atom->ptype == eptShell)
364 shell_index[i] = nshell++;
370 /* Initiate the shell structures */
371 for (i = 0; (i < nshell); i++)
378 /* shell[i].bInterCG=FALSE; */
383 ffparams = &mtop->ffparams;
385 /* Now fill the structures */
386 shfc->bInterCG = FALSE;
389 for (mb = 0; mb < mtop->nmolblock; mb++)
391 molb = &mtop->molblock[mb];
392 molt = &mtop->moltype[molb->type];
395 snew(at2cg, molt->atoms.nr);
396 for (cg = 0; cg < cgs->nr; cg++)
398 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
404 atom = molt->atoms.atom;
405 for (mol = 0; mol < molb->nmol; mol++)
407 for (j = 0; (j < NBT); j++)
409 ia = molt->ilist[bondtypes[j]].iatoms;
410 for (i = 0; (i < molt->ilist[bondtypes[j]].nr); )
413 ftype = ffparams->functype[type];
414 nra = interaction_function[ftype].nratoms;
416 /* Check whether we have a bond with a shell */
419 switch (bondtypes[j])
426 if (atom[ia[1]].ptype == eptShell)
431 else if (atom[ia[2]].ptype == eptShell)
438 aN = ia[4]; /* Dummy */
439 aS = ia[5]; /* Shell */
442 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
449 /* Check whether one of the particles is a shell... */
450 nsi = shell_index[a_offset+aS];
451 if ((nsi < 0) || (nsi >= nshell))
453 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
456 if (shell[nsi].shell == -1)
458 shell[nsi].shell = a_offset + aS;
461 else if (shell[nsi].shell != a_offset+aS)
463 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
466 if (shell[nsi].nucl1 == -1)
468 shell[nsi].nucl1 = a_offset + aN;
470 else if (shell[nsi].nucl2 == -1)
472 shell[nsi].nucl2 = a_offset + aN;
474 else if (shell[nsi].nucl3 == -1)
476 shell[nsi].nucl3 = a_offset + aN;
482 pr_shell(fplog, ns, shell);
484 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
486 if (at2cg[aS] != at2cg[aN])
488 /* shell[nsi].bInterCG = TRUE; */
489 shfc->bInterCG = TRUE;
492 switch (bondtypes[j])
496 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
499 shell[nsi].k += ffparams->iparams[type].cubic.kb;
503 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
505 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);
507 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/
508 ffparams->iparams[type].polarize.alpha;
511 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
513 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);
515 alpha = (ffparams->iparams[type].wpol.al_x+
516 ffparams->iparams[type].wpol.al_y+
517 ffparams->iparams[type].wpol.al_z)/3.0;
518 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/alpha;
521 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
529 a_offset += molt->atoms.nr;
531 /* Done with this molecule type */
535 /* Verify whether it's all correct */
538 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
541 for (i = 0; (i < ns); i++)
543 shell[i].k_1 = 1.0/shell[i].k;
548 pr_shell(debug, ns, shell);
552 shfc->nshell_gl = ns;
553 shfc->shell_gl = shell;
554 shfc->shell_index_gl = shell_index;
556 shfc->bPredict = (getenv("GMX_NOPREDICT") == nullptr);
557 shfc->bRequireInit = FALSE;
562 fprintf(fplog, "\nWill never predict shell positions\n");
567 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
568 if (shfc->bRequireInit && fplog)
570 fprintf(fplog, "\nWill always initiate shell positions\n");
580 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");
582 /* Prediction improves performance, so we should implement either:
583 * 1. communication for the atoms needed for prediction
584 * 2. prediction using the velocities of shells; currently the
585 * shell velocities are zeroed, it's a bit tricky to keep
586 * track of the shell displacements and thus the velocity.
588 shfc->bPredict = FALSE;
595 void make_local_shells(t_commrec *cr, t_mdatoms *md,
599 int a0, a1, *ind, nshell, i;
600 gmx_domdec_t *dd = nullptr;
602 if (DOMAINDECOMP(cr))
610 /* Single node: we need all shells, just copy the pointer */
611 shfc->nshell = shfc->nshell_gl;
612 shfc->shell = shfc->shell_gl;
617 ind = shfc->shell_index_gl;
621 for (i = a0; i < a1; i++)
623 if (md->ptype[i] == eptShell)
625 if (nshell+1 > shfc->shell_nalloc)
627 shfc->shell_nalloc = over_alloc_dd(nshell+1);
628 srenew(shell, shfc->shell_nalloc);
632 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
636 shell[nshell] = shfc->shell_gl[ind[i]];
639 /* With inter-cg shells we can no do shell prediction,
640 * so we do not need the nuclei numbers.
644 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
645 if (shell[nshell].nnucl > 1)
647 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
649 if (shell[nshell].nnucl > 2)
651 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
654 shell[nshell].shell = i;
659 shfc->nshell = nshell;
663 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
681 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
699 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
700 gmx::ArrayRef<gmx::RVec> xnew,
701 const rvec acc_dir[], int homenr, real step)
703 const rvec *xo = as_rvec_array(xold.data());
704 rvec *xn = as_rvec_array(xnew.data());
706 for (int i = 0; i < homenr; i++)
708 do_1pos(xn[i], xo[i], acc_dir[i], step);
712 static void shell_pos_sd(gmx::ArrayRef<const gmx::RVec> xcur,
713 gmx::ArrayRef<gmx::RVec> xnew,
714 gmx::ArrayRef<const gmx::RVec> f,
715 int ns, t_shell s[], int count)
717 const real step_scale_min = 0.8,
718 step_scale_increment = 0.2,
719 step_scale_max = 1.2,
720 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
725 real step_min, step_max;
730 for (i = 0; (i < ns); i++)
735 for (d = 0; d < DIM; d++)
737 s[i].step[d] = s[i].k_1;
739 step_min = std::min(step_min, s[i].step[d]);
740 step_max = std::max(step_max, s[i].step[d]);
746 for (d = 0; d < DIM; d++)
748 dx = xcur[shell][d] - s[i].xold[d];
749 df = f[shell][d] - s[i].fold[d];
750 /* -dx/df gets used to generate an interpolated value, but would
751 * cause a NaN if df were binary-equal to zero. Values close to
752 * zero won't cause problems (because of the min() and max()), so
753 * just testing for binary inequality is OK. */
757 /* Scale the step size by a factor interpolated from
758 * step_scale_min to step_scale_max, as k_est goes from 0 to
759 * step_scale_multiple * s[i].step[d] */
761 step_scale_min * s[i].step[d] +
762 step_scale_increment * std::min(step_scale_multiple * s[i].step[d], std::max(k_est, zero));
767 if (gmx_numzero(dx)) /* 0 == dx */
769 /* Likely this will never happen, but if it does just
770 * don't scale the step. */
774 s[i].step[d] *= step_scale_max;
778 step_min = std::min(step_min, s[i].step[d]);
779 step_max = std::max(step_max, s[i].step[d]);
783 copy_rvec(xcur [shell], s[i].xold);
784 copy_rvec(f[shell], s[i].fold);
786 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
790 fprintf(debug, "shell[%d] = %d\n", i, shell);
791 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
792 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
793 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
794 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
798 printf("step %.3e %.3e\n", step_min, step_max);
802 static void decrease_step_size(int nshell, t_shell s[])
806 for (i = 0; i < nshell; i++)
808 svmul(0.8, s[i].step, s[i].step);
812 static void print_epot(FILE *fp, gmx_int64_t mdstep, int count, real epot, real df,
813 int ndir, real sf_dir)
817 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
818 gmx_step_str(mdstep, buf), count, epot, df);
821 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir/ndir));
830 static real rms_force(t_commrec *cr, gmx::ArrayRef<const gmx::RVec> force, int ns, t_shell s[],
831 int ndir, real *sf_dir, real *Epot)
834 const rvec *f = as_rvec_array(force.data());
837 for (int i = 0; i < ns; i++)
839 int shell = s[i].shell;
840 buf[0] += norm2(f[shell]);
849 gmx_sumd(4, buf, cr);
850 ntot = (int)(buf[1] + 0.5);
856 return (ntot ? std::sqrt(buf[0]/ntot) : 0);
859 static void check_pbc(FILE *fp, gmx::ArrayRef<gmx::RVec> x, int shell)
864 for (m = 0; (m < DIM); m++)
866 if (fabs(x[shell][m]-x[now][m]) > 0.3)
868 pr_rvecs(fp, 0, "SHELL-X", as_rvec_array(x.data())+now, 5);
874 static void dump_shells(FILE *fp, gmx::ArrayRef<gmx::RVec> x, gmx::ArrayRef<gmx::RVec> f, real ftol, int ns, t_shell s[])
879 ft2 = gmx::square(ftol);
881 for (i = 0; (i < ns); i++)
884 ff2 = iprod(f[shell], f[shell]);
887 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
888 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
890 check_pbc(fp, x, shell);
894 static void init_adir(FILE *log, gmx_shellfc_t *shfc,
895 gmx_constr_t constr, t_idef *idef, t_inputrec *ir,
896 t_commrec *cr, int dd_ac1,
897 gmx_int64_t step, t_mdatoms *md, int end,
898 rvec *x_old, rvec *x_init, rvec *x,
899 rvec *f, rvec *acc_dir,
900 gmx_bool bMolPBC, matrix box,
901 gmx::ArrayRef<const real> lambda, real *dvdlambda,
907 unsigned short *ptype;
909 if (DOMAINDECOMP(cr))
917 if (n > shfc->adir_nalloc)
919 shfc->adir_nalloc = over_alloc_dd(n);
920 srenew(shfc->adir_xnold, shfc->adir_nalloc);
921 srenew(shfc->adir_xnew, shfc->adir_nalloc);
923 xnold = shfc->adir_xnold;
924 xnew = shfc->adir_xnew;
930 /* Does NOT work with freeze or acceleration groups (yet) */
931 for (n = 0; n < end; n++)
933 w_dt = md->invmass[n]*dt;
935 for (d = 0; d < DIM; d++)
937 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
939 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
940 xnew[n][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
944 xnold[n][d] = x[n][d];
945 xnew[n][d] = x[n][d];
949 constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
950 x, xnold, nullptr, bMolPBC, box,
951 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
952 nullptr, nullptr, nrnb, econqCoord);
953 constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
954 x, xnew, nullptr, bMolPBC, box,
955 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
956 nullptr, nullptr, nrnb, econqCoord);
958 for (n = 0; n < end; n++)
960 for (d = 0; d < DIM; d++)
963 -(2*x[n][d]-xnold[n][d]-xnew[n][d])/gmx::square(dt)
964 - f[n][d]*md->invmass[n];
966 clear_rvec(acc_dir[n]);
969 /* Project the acceleration on the old bond directions */
970 constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
971 x_old, xnew, acc_dir, bMolPBC, box,
972 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
973 nullptr, nullptr, nrnb, econqDeriv_FlexCon);
976 void relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose,
977 gmx_int64_t mdstep, t_inputrec *inputrec,
978 gmx_bool bDoNS, int force_flags,
981 gmx_enerdata_t *enerd, t_fcdata *fcd,
982 t_state *state, PaddedRVecVector *f,
985 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
987 gmx_groups_t *groups,
991 double t, rvec mu_tot,
993 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion,
994 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion)
999 rvec *acc_dir = nullptr, *x_old = nullptr;
1000 real Epot[2], df[2];
1004 gmx_bool bCont, bInit, bConverged;
1005 int nat, dd_ac0, dd_ac1 = 0, i;
1006 int homenr = md->homenr, end = homenr, cg0, cg1;
1007 int nflexcon, number_steps, d, Min = 0, count = 0;
1008 #define Try (1-Min) /* At start Try = 1 */
1010 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1011 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
1012 ftol = inputrec->em_tol;
1013 number_steps = inputrec->niter;
1014 nshell = shfc->nshell;
1015 shell = shfc->shell;
1016 nflexcon = shfc->nflexcon;
1020 if (DOMAINDECOMP(cr))
1022 nat = dd_natoms_vsite(cr->dd);
1025 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1026 nat = std::max(nat, dd_ac1);
1031 nat = state->natoms;
1034 for (i = 0; (i < 2); i++)
1036 shfc->f[i].resize(gmx::paddedRVecVectorSize(nat));
1039 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1041 /* This is the only time where the coordinates are used
1042 * before do_force is called, which normally puts all
1043 * charge groups in the box.
1045 if (inputrec->cutoff_scheme == ecutsVERLET)
1047 auto xRef = makeArrayRef(state->x);
1048 put_atoms_in_box_omp(fr->ePBC, state->box, xRef.subArray(0, md->homenr));
1054 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1055 &(top->cgs), as_rvec_array(state->x.data()), fr->cg_cm);
1060 mk_mshift(fplog, graph, fr->ePBC, state->box, as_rvec_array(state->x.data()));
1064 /* After this all coordinate arrays will contain whole charge groups */
1067 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1072 if (nat > shfc->flex_nalloc)
1074 shfc->flex_nalloc = over_alloc_dd(nat);
1075 srenew(shfc->acc_dir, shfc->flex_nalloc);
1076 srenew(shfc->x_old, shfc->flex_nalloc);
1078 acc_dir = shfc->acc_dir;
1079 x_old = shfc->x_old;
1080 for (i = 0; i < homenr; i++)
1082 for (d = 0; d < DIM; d++)
1085 state->x[i][d] - state->v[i][d]*inputrec->delta_t;
1090 /* Do a prediction of the shell positions, when appropriate.
1091 * Without velocities (EM, NM, BD) we only do initial prediction.
1093 if (shfc->bPredict && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1095 predict_shells(fplog, as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), inputrec->delta_t, nshell, shell,
1096 md->massT, nullptr, bInit);
1099 /* do_force expected the charge groups to be in the box */
1102 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1105 /* Calculate the forces first time around */
1108 pr_rvecs(debug, 0, "x b4 do_force", as_rvec_array(state->x.data()), homenr);
1110 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1111 do_force(fplog, cr, inputrec, mdstep, nrnb, wcycle, top, groups,
1112 state->box, state->x, &state->hist,
1113 shfc->f[Min], force_vir, md, enerd, fcd,
1114 state->lambda, graph,
1115 fr, vsite, mu_tot, t, nullptr, bBornRadii,
1116 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags,
1117 ddOpenBalanceRegion, ddCloseBalanceRegion);
1122 init_adir(fplog, shfc,
1123 constr, idef, inputrec, cr, dd_ac1, mdstep, md, end,
1124 shfc->x_old, as_rvec_array(state->x.data()), as_rvec_array(state->x.data()), as_rvec_array(shfc->f[Min].data()),
1126 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1128 for (i = 0; i < end; i++)
1130 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i]);
1133 sum_epot(&(enerd->grpp), enerd->term);
1134 Epot[Min] = enerd->term[F_EPOT];
1135 df[Min] = rms_force(cr, shfc->f[Min], nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1139 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1144 pr_rvecs(debug, 0, "force0", as_rvec_array(shfc->f[Min].data()), md->nr);
1147 if (nshell+nflexcon > 0)
1149 /* Copy x to shfc->x[Min] & shfc->x[Try]: during minimization only the
1150 * shell positions are updated, therefore the other particles must
1151 * be set here, in advance.
1153 shfc->x[Min] = PaddedRVecVector(std::begin(state->x), std::end(state->x));
1154 shfc->x[Try] = PaddedRVecVector(std::begin(state->x), std::end(state->x));
1157 if (bVerbose && MASTER(cr))
1159 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1164 fprintf(debug, "%17s: %14.10e\n",
1165 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1166 fprintf(debug, "%17s: %14.10e\n",
1167 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1168 fprintf(debug, "%17s: %14.10e\n",
1169 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1170 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1173 /* First check whether we should do shells, or whether the force is
1174 * low enough even without minimization.
1176 bConverged = (df[Min] < ftol);
1178 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1182 construct_vsites(vsite, as_rvec_array(shfc->x[Min].data()),
1183 inputrec->delta_t, as_rvec_array(state->v.data()),
1184 idef->iparams, idef->il,
1185 fr->ePBC, fr->bMolPBC, cr, state->box);
1190 init_adir(fplog, shfc,
1191 constr, idef, inputrec, cr, dd_ac1, mdstep, md, end,
1192 x_old, as_rvec_array(state->x.data()), as_rvec_array(shfc->x[Min].data()), as_rvec_array(shfc->f[Min].data()), acc_dir,
1193 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1195 directional_sd(shfc->x[Min], shfc->x[Try], acc_dir, end, fr->fc_stepsize);
1198 /* New positions, Steepest descent */
1199 shell_pos_sd(shfc->x[Min], shfc->x[Try], shfc->f[Min], nshell, shell, count);
1201 /* do_force expected the charge groups to be in the box */
1204 unshift_self(graph, state->box, as_rvec_array(shfc->x[Try].data()));
1209 pr_rvecs(debug, 0, "RELAX: shfc->x[Min] ", as_rvec_array(shfc->x[Min].data()), homenr);
1210 pr_rvecs(debug, 0, "RELAX: shfc->x[Try] ", as_rvec_array(shfc->x[Try].data()), homenr);
1212 /* Try the new positions */
1213 do_force(fplog, cr, inputrec, 1, nrnb, wcycle,
1214 top, groups, state->box, shfc->x[Try], &state->hist,
1215 shfc->f[Try], force_vir,
1216 md, enerd, fcd, state->lambda, graph,
1217 fr, vsite, mu_tot, t, nullptr, bBornRadii,
1219 ddOpenBalanceRegion, ddCloseBalanceRegion);
1220 sum_epot(&(enerd->grpp), enerd->term);
1223 pr_rvecs(debug, 0, "RELAX: shfc->f[Min]", as_rvec_array(shfc->f[Min].data()), homenr);
1224 pr_rvecs(debug, 0, "RELAX: shfc->f[Try]", as_rvec_array(shfc->f[Try].data()), homenr);
1229 init_adir(fplog, shfc,
1230 constr, idef, inputrec, cr, dd_ac1, mdstep, md, end,
1231 x_old, as_rvec_array(state->x.data()), as_rvec_array(shfc->x[Try].data()), as_rvec_array(shfc->f[Try].data()), acc_dir,
1232 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1234 for (i = 0; i < end; i++)
1236 sf_dir += md->massT[i]*norm2(acc_dir[i]);
1240 Epot[Try] = enerd->term[F_EPOT];
1242 df[Try] = rms_force(cr, shfc->f[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1246 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1253 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(shfc->f[Try].data()), homenr);
1257 fprintf(debug, "SHELL ITER %d\n", count);
1258 dump_shells(debug, shfc->x[Try], shfc->f[Try], ftol, nshell, shell);
1262 if (bVerbose && MASTER(cr))
1264 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1267 bConverged = (df[Try] < ftol);
1269 if ((df[Try] < df[Min]))
1273 fprintf(debug, "Swapping Min and Try\n");
1277 /* Correct the velocities for the flexible constraints */
1278 invdt = 1/inputrec->delta_t;
1279 for (i = 0; i < end; i++)
1281 for (d = 0; d < DIM; d++)
1283 state->v[i][d] += (shfc->x[Try][i][d] - shfc->x[Min][i][d])*invdt;
1291 decrease_step_size(nshell, shell);
1294 shfc->numForceEvaluations += count;
1297 shfc->numConvergedIterations++;
1299 if (MASTER(cr) && !(bConverged))
1301 /* Note that the energies and virial are incorrect when not converged */
1305 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1306 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1309 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1310 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1313 /* Copy back the coordinates and the forces */
1314 std::copy(shfc->x[Min].begin(), shfc->x[Min].end(), state->x.begin());
1315 std::copy(shfc->f[Min].begin(), shfc->f[Min].end(), f->begin());
1318 void done_shellfc(FILE *fplog, gmx_shellfc_t *shfc, gmx_int64_t numSteps)
1320 if (shfc && fplog && numSteps > 0)
1322 double numStepsAsDouble = static_cast<double>(numSteps);
1323 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1324 (shfc->numConvergedIterations*100.0)/numStepsAsDouble);
1325 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1326 shfc->numForceEvaluations/numStepsAsDouble);
1329 // TODO Deallocate memory in shfc