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,2019, 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/domdec/mdsetup.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/enerdata_utils.h"
60 #include "gromacs/mdlib/force.h"
61 #include "gromacs/mdlib/force_flags.h"
62 #include "gromacs/mdlib/gmx_omp_nthreads.h"
63 #include "gromacs/mdlib/mdatoms.h"
64 #include "gromacs/mdlib/vsite.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/enerdata.h"
67 #include "gromacs/mdtypes/forcerec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/pbcutil/mshift.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/topology/ifunc.h"
74 #include "gromacs/topology/mtop_lookup.h"
75 #include "gromacs/topology/mtop_util.h"
76 #include "gromacs/utility/arrayref.h"
77 #include "gromacs/utility/arraysize.h"
78 #include "gromacs/utility/cstringutil.h"
79 #include "gromacs/utility/fatalerror.h"
80 #include "gromacs/utility/gmxassert.h"
81 #include "gromacs/utility/smalloc.h"
85 int shell; /* The shell id */
86 int nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
87 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
88 real k; /* force constant */
89 real k_1; /* 1 over force constant */
95 struct gmx_shellfc_t {
96 /* Shell counts, indices, parameters and working data */
97 int nshell_gl; /* The number of shells in the system */
98 t_shell *shell_gl; /* All the shells (for DD only) */
99 int *shell_index_gl; /* Global shell index (for DD only) */
100 gmx_bool bInterCG; /* Are there inter charge-group shells? */
101 int nshell; /* The number of local shells */
102 t_shell *shell; /* The local shells */
103 int shell_nalloc; /* The allocation size of shell */
104 gmx_bool bPredict; /* Predict shell positions */
105 gmx_bool bRequireInit; /* Require initialization of shell positions */
106 int nflexcon; /* The number of flexible constraints */
108 /* Temporary arrays, should be fixed size 2 when fully converted to C++ */
109 PaddedHostVector<gmx::RVec> *x; /* Array for iterative minimization */
110 PaddedHostVector<gmx::RVec> *f; /* Array for iterative minimization */
112 /* Flexible constraint working data */
113 rvec *acc_dir; /* Acceleration direction for flexcon */
114 rvec *x_old; /* Old coordinates for flexcon */
115 int flex_nalloc; /* The allocation size of acc_dir and x_old */
116 rvec *adir_xnold; /* Work space for init_adir */
117 rvec *adir_xnew; /* Work space for init_adir */
118 int adir_nalloc; /* Work space for init_adir */
119 std::int64_t numForceEvaluations; /* Total number of force evaluations */
120 int numConvergedIterations; /* Total number of iterations that converged */
124 static void pr_shell(FILE *fplog, int ns, t_shell s[])
128 fprintf(fplog, "SHELL DATA\n");
129 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
130 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
131 for (i = 0; (i < ns); i++)
133 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
136 fprintf(fplog, " %5d\n", s[i].nucl2);
138 else if (s[i].nnucl == 3)
140 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
144 fprintf(fplog, "\n");
149 /* TODO The remain call of this function passes non-NULL mass and NULL
150 * mtop, so this routine can be simplified.
152 * The other code path supported doing prediction before the MD loop
153 * started, but even when called, the prediction was always
154 * over-written by a subsequent call in the MD loop, so has been
156 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
158 const real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
160 int i, m, s1, n1, n2, n3;
161 real dt_1, fudge, tm, m1, m2, m3;
164 GMX_RELEASE_ASSERT(mass || mtop, "Must have masses or a way to look them up");
166 /* We introduce a fudge factor for performance reasons: with this choice
167 * the initial force on the shells is about a factor of two lower than
176 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
188 for (i = 0; (i < ns); i++)
199 for (m = 0; (m < DIM); m++)
201 x[s1][m] += ptr[n1][m]*dt_1;
214 /* Not the correct masses with FE, but it is just a prediction... */
215 m1 = mtopGetAtomMass(mtop, n1, &molb);
216 m2 = mtopGetAtomMass(mtop, n2, &molb);
219 for (m = 0; (m < DIM); m++)
221 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
236 /* Not the correct masses with FE, but it is just a prediction... */
237 m1 = mtopGetAtomMass(mtop, n1, &molb);
238 m2 = mtopGetAtomMass(mtop, n2, &molb);
239 m3 = mtopGetAtomMass(mtop, n3, &molb);
241 tm = dt_1/(m1+m2+m3);
242 for (m = 0; (m < DIM); m++)
244 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
248 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
253 /*! \brief Count the different particle types in a system
255 * Routine prints a warning to stderr in case an unknown particle type
257 * \param[in] fplog Print what we have found if not NULL
258 * \param[in] mtop Molecular topology.
259 * \returns Array holding the number of particles of a type
261 std::array<int, eptNR> countPtypes(FILE *fplog,
262 const gmx_mtop_t *mtop)
264 std::array<int, eptNR> nptype = { { 0 } };
265 /* Count number of shells, and find their indices */
266 for (int i = 0; (i < eptNR); i++)
271 gmx_mtop_atomloop_block_t aloopb = gmx_mtop_atomloop_block_init(mtop);
274 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
281 nptype[atom->ptype] += nmol;
284 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
285 static_cast<int>(atom->ptype));
290 /* Print the number of each particle type */
292 for (const auto &i : nptype)
296 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
304 gmx_shellfc_t *init_shell_flexcon(FILE *fplog,
305 const gmx_mtop_t *mtop, int nflexcon,
307 bool usingDomainDecomposition)
311 int *shell_index = nullptr;
314 int i, j, type, a_offset, mol, ftype, nra;
316 int aS, aN = 0; /* Shell and nucleus */
317 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
318 #define NBT asize(bondtypes)
319 const gmx_ffparams_t *ffparams;
321 std::array<int, eptNR> n = countPtypes(fplog, mtop);
322 nshell = n[eptShell];
324 if (nshell == 0 && nflexcon == 0)
326 /* We're not doing shells or flexible constraints */
331 shfc->x = new PaddedHostVector<gmx::RVec>[2] {};
332 shfc->f = new PaddedHostVector<gmx::RVec>[2] {};
333 shfc->nflexcon = nflexcon;
337 /* Only flexible constraints, no shells.
338 * Note that make_local_shells() does not need to be called.
341 shfc->bPredict = FALSE;
346 if (nstcalcenergy != 1)
348 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);
350 if (usingDomainDecomposition)
352 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
355 /* We have shells: fill the shell data structure */
357 /* Global system sized array, this should be avoided */
358 snew(shell_index, mtop->natoms);
361 for (const AtomProxy atomP : AtomRange(*mtop))
363 const t_atom &local = atomP.atom();
364 int i = atomP.globalAtomNumber();
365 if (local.ptype == eptShell)
367 shell_index[i] = nshell++;
373 /* Initiate the shell structures */
374 for (i = 0; (i < nshell); i++)
381 /* shell[i].bInterCG=FALSE; */
386 ffparams = &mtop->ffparams;
388 /* Now fill the structures */
389 /* TODO: See if we can use update groups that cover shell constructions */
390 shfc->bInterCG = FALSE;
393 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
395 const gmx_molblock_t *molb = &mtop->molblock[mb];
396 const gmx_moltype_t *molt = &mtop->moltype[molb->type];
398 const t_atom *atom = molt->atoms.atom;
399 for (mol = 0; mol < molb->nmol; mol++)
401 for (j = 0; (j < NBT); j++)
403 const int *ia = molt->ilist[bondtypes[j]].iatoms.data();
404 for (i = 0; (i < molt->ilist[bondtypes[j]].size()); )
407 ftype = ffparams->functype[type];
408 nra = interaction_function[ftype].nratoms;
410 /* Check whether we have a bond with a shell */
413 switch (bondtypes[j])
420 if (atom[ia[1]].ptype == eptShell)
425 else if (atom[ia[2]].ptype == eptShell)
432 aN = ia[4]; /* Dummy */
433 aS = ia[5]; /* Shell */
436 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
443 /* Check whether one of the particles is a shell... */
444 nsi = shell_index[a_offset+aS];
445 if ((nsi < 0) || (nsi >= nshell))
447 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
450 if (shell[nsi].shell == -1)
452 shell[nsi].shell = a_offset + aS;
455 else if (shell[nsi].shell != a_offset+aS)
457 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
460 if (shell[nsi].nucl1 == -1)
462 shell[nsi].nucl1 = a_offset + aN;
464 else if (shell[nsi].nucl2 == -1)
466 shell[nsi].nucl2 = a_offset + aN;
468 else if (shell[nsi].nucl3 == -1)
470 shell[nsi].nucl3 = a_offset + aN;
476 pr_shell(fplog, ns, shell);
478 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
482 /* shell[nsi].bInterCG = TRUE; */
483 shfc->bInterCG = TRUE;
486 switch (bondtypes[j])
490 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
493 shell[nsi].k += ffparams->iparams[type].cubic.kb;
497 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
499 gmx_fatal(FARGS, "polarize can not be used with qA(%e) != qB(%e) for atom %d of molecule block %zu", qS, atom[aS].qB, aS+1, mb+1);
501 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/
502 ffparams->iparams[type].polarize.alpha;
505 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
507 gmx_fatal(FARGS, "water_pol can not be used with qA(%e) != qB(%e) for atom %d of molecule block %zu", qS, atom[aS].qB, aS+1, mb+1);
509 alpha = (ffparams->iparams[type].wpol.al_x+
510 ffparams->iparams[type].wpol.al_y+
511 ffparams->iparams[type].wpol.al_z)/3.0;
512 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/alpha;
515 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
523 a_offset += molt->atoms.nr;
525 /* Done with this molecule type */
528 /* Verify whether it's all correct */
531 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
534 for (i = 0; (i < ns); i++)
536 shell[i].k_1 = 1.0/shell[i].k;
541 pr_shell(debug, ns, shell);
545 shfc->nshell_gl = ns;
546 shfc->shell_gl = shell;
547 shfc->shell_index_gl = shell_index;
549 shfc->bPredict = (getenv("GMX_NOPREDICT") == nullptr);
550 shfc->bRequireInit = FALSE;
555 fprintf(fplog, "\nWill never predict shell positions\n");
560 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
561 if (shfc->bRequireInit && fplog)
563 fprintf(fplog, "\nWill always initiate shell positions\n");
573 fprintf(fplog, "\nNOTE: in the current version shell prediction during the crun is disabled\n\n");
575 /* Prediction improves performance, so we should implement either:
576 * 1. communication for the atoms needed for prediction
577 * 2. prediction using the velocities of shells; currently the
578 * shell velocities are zeroed, it's a bit tricky to keep
579 * track of the shell displacements and thus the velocity.
581 shfc->bPredict = FALSE;
588 void make_local_shells(const t_commrec *cr,
593 int a0, a1, *ind, nshell, i;
594 gmx_domdec_t *dd = nullptr;
596 if (DOMAINDECOMP(cr))
600 a1 = dd_numHomeAtoms(*dd);
604 /* Single node: we need all shells, just copy the pointer */
605 shfc->nshell = shfc->nshell_gl;
606 shfc->shell = shfc->shell_gl;
611 ind = shfc->shell_index_gl;
615 for (i = a0; i < a1; i++)
617 if (md->ptype[i] == eptShell)
619 if (nshell+1 > shfc->shell_nalloc)
621 shfc->shell_nalloc = over_alloc_dd(nshell+1);
622 srenew(shell, shfc->shell_nalloc);
626 shell[nshell] = shfc->shell_gl[ind[dd->globalAtomIndices[i]]];
630 shell[nshell] = shfc->shell_gl[ind[i]];
633 /* With inter-cg shells we can no do shell prediction,
634 * so we do not need the nuclei numbers.
638 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
639 if (shell[nshell].nnucl > 1)
641 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
643 if (shell[nshell].nnucl > 2)
645 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
648 shell[nshell].shell = i;
653 shfc->nshell = nshell;
657 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
675 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
693 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
694 gmx::ArrayRef<gmx::RVec> xnew,
695 const rvec acc_dir[], int homenr, real step)
697 const rvec *xo = as_rvec_array(xold.data());
698 rvec *xn = as_rvec_array(xnew.data());
700 for (int i = 0; i < homenr; i++)
702 do_1pos(xn[i], xo[i], acc_dir[i], step);
706 static void shell_pos_sd(gmx::ArrayRef<const gmx::RVec> xcur,
707 gmx::ArrayRef<gmx::RVec> xnew,
708 gmx::ArrayRef<const gmx::RVec> f,
709 int ns, t_shell s[], int count)
711 const real step_scale_min = 0.8,
712 step_scale_increment = 0.2,
713 step_scale_max = 1.2,
714 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
719 real step_min, step_max;
724 for (i = 0; (i < ns); i++)
729 for (d = 0; d < DIM; d++)
731 s[i].step[d] = s[i].k_1;
733 step_min = std::min(step_min, s[i].step[d]);
734 step_max = std::max(step_max, s[i].step[d]);
740 for (d = 0; d < DIM; d++)
742 dx = xcur[shell][d] - s[i].xold[d];
743 df = f[shell][d] - s[i].fold[d];
744 /* -dx/df gets used to generate an interpolated value, but would
745 * cause a NaN if df were binary-equal to zero. Values close to
746 * zero won't cause problems (because of the min() and max()), so
747 * just testing for binary inequality is OK. */
751 /* Scale the step size by a factor interpolated from
752 * step_scale_min to step_scale_max, as k_est goes from 0 to
753 * step_scale_multiple * s[i].step[d] */
755 step_scale_min * s[i].step[d] +
756 step_scale_increment * std::min(step_scale_multiple * s[i].step[d], std::max(k_est, zero));
761 if (gmx_numzero(dx)) /* 0 == dx */
763 /* Likely this will never happen, but if it does just
764 * don't scale the step. */
768 s[i].step[d] *= step_scale_max;
772 step_min = std::min(step_min, s[i].step[d]);
773 step_max = std::max(step_max, s[i].step[d]);
777 copy_rvec(xcur [shell], s[i].xold);
778 copy_rvec(f[shell], s[i].fold);
780 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
784 fprintf(debug, "shell[%d] = %d\n", i, shell);
785 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
786 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
787 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
788 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
792 printf("step %.3e %.3e\n", step_min, step_max);
796 static void decrease_step_size(int nshell, t_shell s[])
800 for (i = 0; i < nshell; i++)
802 svmul(0.8, s[i].step, s[i].step);
806 static void print_epot(FILE *fp, int64_t mdstep, int count, real epot, real df,
807 int ndir, real sf_dir)
811 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
812 gmx_step_str(mdstep, buf), count, epot, df);
815 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir/ndir));
824 static real rms_force(const t_commrec *cr, gmx::ArrayRef<const gmx::RVec> force, int ns, t_shell s[],
825 int ndir, real *sf_dir, real *Epot)
828 const rvec *f = as_rvec_array(force.data());
831 for (int i = 0; i < ns; i++)
833 int shell = s[i].shell;
834 buf[0] += norm2(f[shell]);
843 gmx_sumd(4, buf, cr);
844 ntot = gmx::roundToInt(buf[1]);
850 return (ntot ? std::sqrt(buf[0]/ntot) : 0);
853 static void dump_shells(FILE *fp, gmx::ArrayRef<gmx::RVec> f, real ftol, int ns, t_shell s[])
858 ft2 = gmx::square(ftol);
860 for (i = 0; (i < ns); i++)
863 ff2 = iprod(f[shell], f[shell]);
866 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
867 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
872 static void init_adir(gmx_shellfc_t *shfc,
873 gmx::Constraints *constr,
874 const t_inputrec *ir,
886 gmx::ArrayRef<const real> lambda,
892 unsigned short *ptype;
894 if (DOMAINDECOMP(cr))
902 if (n > shfc->adir_nalloc)
904 shfc->adir_nalloc = over_alloc_dd(n);
905 srenew(shfc->adir_xnold, shfc->adir_nalloc);
906 srenew(shfc->adir_xnew, shfc->adir_nalloc);
908 xnold = shfc->adir_xnold;
909 xnew = shfc->adir_xnew;
915 /* Does NOT work with freeze or acceleration groups (yet) */
916 for (n = 0; n < end; n++)
918 w_dt = md->invmass[n]*dt;
920 for (d = 0; d < DIM; d++)
922 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
924 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
925 xnew[n][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
929 xnold[n][d] = x[n][d];
930 xnew[n][d] = x[n][d];
934 constr->apply(FALSE, FALSE, step, 0, 1.0,
935 x, xnold, nullptr, box,
936 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
937 nullptr, nullptr, gmx::ConstraintVariable::Positions);
938 constr->apply(FALSE, FALSE, step, 0, 1.0,
939 x, xnew, nullptr, box,
940 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
941 nullptr, nullptr, gmx::ConstraintVariable::Positions);
943 for (n = 0; n < end; n++)
945 for (d = 0; d < DIM; d++)
948 -(2*x[n][d]-xnold[n][d]-xnew[n][d])/gmx::square(dt)
949 - f[n][d]*md->invmass[n];
951 clear_rvec(acc_dir[n]);
954 /* Project the acceleration on the old bond directions */
955 constr->apply(FALSE, FALSE, step, 0, 1.0,
956 x_old, xnew, acc_dir, box,
957 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
958 nullptr, nullptr, gmx::ConstraintVariable::Deriv_FlexCon);
961 void relax_shell_flexcon(FILE *fplog,
963 const gmx_multisim_t *ms,
965 gmx_enfrot *enforcedRotation,
967 const t_inputrec *inputrec,
968 gmx::ImdSession *imdSession,
972 const gmx_localtop_t *top,
973 gmx::Constraints *constr,
974 gmx_enerdata_t *enerd,
977 gmx::ArrayRefWithPadding<gmx::RVec> x,
978 gmx::ArrayRefWithPadding<gmx::RVec> v,
980 gmx::ArrayRef<real> lambda,
982 gmx::ArrayRefWithPadding<gmx::RVec> f,
986 gmx_wallcycle_t wcycle,
990 gmx::MdrunScheduleWorkload *runScheduleWork,
993 const gmx_vsite_t *vsite,
994 const DDBalanceRegionHandler &ddBalanceRegionHandler)
996 auto xRvec = as_rvec_array(x.paddedArrayRef().data());
997 auto vRvec = as_rvec_array(v.paddedArrayRef().data());
1002 rvec *acc_dir = nullptr, *x_old = nullptr;
1003 real Epot[2], df[2];
1007 gmx_bool bCont, bInit, bConverged;
1008 int nat, dd_ac0, dd_ac1 = 0, i;
1009 int homenr = md->homenr, end = homenr;
1010 int nflexcon, number_steps, d, Min = 0, count = 0;
1011 #define Try (1-Min) /* At start Try = 1 */
1013 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1014 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
1015 ftol = inputrec->em_tol;
1016 number_steps = inputrec->niter;
1017 nshell = shfc->nshell;
1018 shell = shfc->shell;
1019 nflexcon = shfc->nflexcon;
1023 if (DOMAINDECOMP(cr))
1025 nat = dd_natoms_vsite(cr->dd);
1028 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1029 nat = std::max(nat, dd_ac1);
1037 for (i = 0; (i < 2); i++)
1039 shfc->x[i].resizeWithPadding(nat);
1040 shfc->f[i].resizeWithPadding(nat);
1043 /* Create views that we can swap */
1044 gmx::ArrayRefWithPadding<gmx::RVec> posWithPadding[2];
1045 gmx::ArrayRefWithPadding<gmx::RVec> forceWithPadding[2];
1046 gmx::ArrayRef<gmx::RVec> pos[2];
1047 gmx::ArrayRef<gmx::RVec> force[2];
1048 for (i = 0; (i < 2); i++)
1050 posWithPadding[i] = shfc->x[i].arrayRefWithPadding();
1051 pos[i] = posWithPadding[i].paddedArrayRef();
1052 forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1053 force[i] = forceWithPadding[i].paddedArrayRef();
1056 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1058 /* This is the only time where the coordinates are used
1059 * before do_force is called, which normally puts all
1060 * charge groups in the box.
1062 auto xRef = x.paddedArrayRef();
1063 put_atoms_in_box_omp(fr->ePBC, box, xRef.subArray(0, md->homenr), gmx_omp_nthreads_get(emntDefault));
1067 mk_mshift(fplog, graph, fr->ePBC, box, xRvec);
1071 /* After this all coordinate arrays will contain whole charge groups */
1074 shift_self(graph, box, xRvec);
1079 if (nat > shfc->flex_nalloc)
1081 shfc->flex_nalloc = over_alloc_dd(nat);
1082 srenew(shfc->acc_dir, shfc->flex_nalloc);
1083 srenew(shfc->x_old, shfc->flex_nalloc);
1085 acc_dir = shfc->acc_dir;
1086 x_old = shfc->x_old;
1087 auto xArrayRef = x.paddedArrayRef();
1088 auto vArrayRef = v.paddedArrayRef();
1089 for (i = 0; i < homenr; i++)
1091 for (d = 0; d < DIM; d++)
1094 xArrayRef[i][d] - vArrayRef[i][d]*inputrec->delta_t;
1099 /* Do a prediction of the shell positions, when appropriate.
1100 * Without velocities (EM, NM, BD) we only do initial prediction.
1102 if (shfc->bPredict && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1104 predict_shells(fplog, xRvec, vRvec, inputrec->delta_t, nshell, shell,
1105 md->massT, nullptr, bInit);
1108 /* do_force expected the charge groups to be in the box */
1111 unshift_self(graph, box, xRvec);
1114 /* Calculate the forces first time around */
1117 pr_rvecs(debug, 0, "x b4 do_force", xRvec, homenr);
1119 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1120 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession,
1122 mdstep, nrnb, wcycle, top,
1124 forceWithPadding[Min], force_vir, md, enerd, fcd,
1126 fr, runScheduleWork, vsite, mu_tot, t, nullptr,
1127 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags,
1128 ddBalanceRegionHandler);
1134 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1135 shfc->x_old, xRvec, xRvec, as_rvec_array(force[Min].data()),
1139 for (i = 0; i < end; i++)
1141 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i]);
1144 sum_epot(&(enerd->grpp), enerd->term);
1145 Epot[Min] = enerd->term[F_EPOT];
1147 df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1151 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1156 pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md->nr);
1159 if (nshell+nflexcon > 0)
1161 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1162 * shell positions are updated, therefore the other particles must
1163 * be set here, in advance.
1165 std::copy(x.paddedArrayRef().begin(),
1166 x.paddedArrayRef().end(),
1167 posWithPadding[Min].paddedArrayRef().begin());
1168 std::copy(x.paddedArrayRef().begin(),
1169 x.paddedArrayRef().end(),
1170 posWithPadding[Try].paddedArrayRef().begin());
1173 if (bVerbose && MASTER(cr))
1175 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1180 fprintf(debug, "%17s: %14.10e\n",
1181 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1182 fprintf(debug, "%17s: %14.10e\n",
1183 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1184 fprintf(debug, "%17s: %14.10e\n",
1185 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1186 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1189 /* First check whether we should do shells, or whether the force is
1190 * low enough even without minimization.
1192 bConverged = (df[Min] < ftol);
1194 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1198 construct_vsites(vsite, as_rvec_array(pos[Min].data()),
1199 inputrec->delta_t, vRvec,
1200 idef->iparams, idef->il,
1201 fr->ePBC, fr->bMolPBC, cr, box);
1207 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1209 as_rvec_array(pos[Min].data()),
1210 as_rvec_array(force[Min].data()), acc_dir,
1213 directional_sd(pos[Min], pos[Try], acc_dir, end, fr->fc_stepsize);
1216 /* New positions, Steepest descent */
1217 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1219 /* do_force expected the charge groups to be in the box */
1222 unshift_self(graph, box, as_rvec_array(pos[Try].data()));
1227 pr_rvecs(debug, 0, "RELAX: pos[Min] ", as_rvec_array(pos[Min].data()), homenr);
1228 pr_rvecs(debug, 0, "RELAX: pos[Try] ", as_rvec_array(pos[Try].data()), homenr);
1230 /* Try the new positions */
1231 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession,
1234 top, box, posWithPadding[Try], hist,
1235 forceWithPadding[Try], force_vir,
1236 md, enerd, fcd, lambda, graph,
1237 fr, runScheduleWork, vsite, mu_tot, t, nullptr,
1239 ddBalanceRegionHandler);
1240 sum_epot(&(enerd->grpp), enerd->term);
1243 pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1244 pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1250 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1252 as_rvec_array(pos[Try].data()),
1253 as_rvec_array(force[Try].data()),
1254 acc_dir, box, lambda, &dum);
1256 for (i = 0; i < end; i++)
1258 sf_dir += md->massT[i]*norm2(acc_dir[i]);
1262 Epot[Try] = enerd->term[F_EPOT];
1264 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1268 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1275 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1279 fprintf(debug, "SHELL ITER %d\n", count);
1280 dump_shells(debug, force[Try], ftol, nshell, shell);
1284 if (bVerbose && MASTER(cr))
1286 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1289 bConverged = (df[Try] < ftol);
1291 if ((df[Try] < df[Min]))
1295 fprintf(debug, "Swapping Min and Try\n");
1299 /* Correct the velocities for the flexible constraints */
1300 invdt = 1/inputrec->delta_t;
1301 auto vArrayRef = v.paddedArrayRef();
1302 for (i = 0; i < end; i++)
1304 for (d = 0; d < DIM; d++)
1306 vArrayRef[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1314 decrease_step_size(nshell, shell);
1317 shfc->numForceEvaluations += count;
1320 shfc->numConvergedIterations++;
1322 if (MASTER(cr) && !(bConverged))
1324 /* Note that the energies and virial are incorrect when not converged */
1328 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1329 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1332 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1333 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1336 /* Copy back the coordinates and the forces */
1337 std::copy(pos[Min].begin(), pos[Min].end(), x.paddedArrayRef().data());
1338 std::copy(force[Min].begin(), force[Min].end(), f.unpaddedArrayRef().begin());
1341 void done_shellfc(FILE *fplog, gmx_shellfc_t *shfc, int64_t numSteps)
1343 if (shfc && fplog && numSteps > 0)
1345 double numStepsAsDouble = static_cast<double>(numSteps);
1346 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1347 (shfc->numConvergedIterations*100.0)/numStepsAsDouble);
1348 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1349 shfc->numForceEvaluations/numStepsAsDouble);
1352 // TODO Deallocate memory in shfc