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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
50 #include "gromacs/domdec/dlbtiming.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/mdsetup.h"
54 #include "gromacs/gmxlib/network.h"
55 #include "gromacs/math/functions.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/math/vecdump.h"
59 #include "gromacs/mdlib/constr.h"
60 #include "gromacs/mdlib/enerdata_utils.h"
61 #include "gromacs/mdlib/force.h"
62 #include "gromacs/mdlib/force_flags.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdlib/mdatoms.h"
65 #include "gromacs/mdlib/vsite.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/enerdata.h"
68 #include "gromacs/mdtypes/forcerec.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/mdtypes/state.h"
72 #include "gromacs/pbcutil/mshift.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/topology/ifunc.h"
75 #include "gromacs/topology/mtop_lookup.h"
76 #include "gromacs/topology/mtop_util.h"
77 #include "gromacs/utility/arrayref.h"
78 #include "gromacs/utility/arraysize.h"
79 #include "gromacs/utility/cstringutil.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/gmxassert.h"
84 using gmx::ArrayRefWithPadding;
89 int nnucl = 0; /* The number of nuclei */
90 int shellIndex = -1; /* The shell index */
91 int nucl1 = -1; /* The first nuclei connected to the shell */
92 int nucl2 = -1; /* The second nuclei connected to the shell */
93 int nucl3 = -1; /* The third nuclei connected to the shell */
94 real k = 0; /* force constant */
95 real k_1 = 0; /* 1 over force constant */
96 rvec xold; /* The old shell coordinates */
97 rvec fold; /* The old force on the shell */
98 rvec step; /* Step size for steepest descents */
103 /* Shell counts, indices, parameters and working data */
104 std::vector<t_shell> shell_gl; /* All the shells (for DD only) */
105 std::vector<int> shell_index_gl; /* Global shell index (for DD only) */
106 gmx_bool bInterCG; /* Are there inter charge-group shells? */
107 std::vector<t_shell> shells; /* The local shells */
108 bool predictShells = false; /* Predict shell positions */
109 bool requireInit = false; /* Require initialization of shell positions */
110 int nflexcon = 0; /* The number of flexible constraints */
112 std::array<PaddedHostVector<gmx::RVec>, 2> x; /* Coordinate buffers for iterative minimization */
113 std::array<PaddedHostVector<gmx::RVec>, 2> f; /* Force buffers for iterative minimization */
115 /* Flexible constraint working data */
116 std::vector<RVec> acc_dir; /* Acceleration direction for flexcon */
117 gmx::PaddedVector<RVec> x_old; /* Old coordinates for flexcon */
118 gmx::PaddedVector<RVec> adir_xnold; /* Work space for init_adir */
119 gmx::PaddedVector<RVec> adir_xnew; /* Work space for init_adir */
120 std::int64_t numForceEvaluations; /* Total number of force evaluations */
121 int numConvergedIterations; /* Total number of iterations that converged */
125 static void pr_shell(FILE* fplog, ArrayRef<const t_shell> shells)
127 fprintf(fplog, "SHELL DATA\n");
128 fprintf(fplog, "%5s %8s %5s %5s %5s\n", "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
129 for (const t_shell& shell : shells)
131 fprintf(fplog, "%5d %8.3f %5d", shell.shellIndex, 1.0 / shell.k_1, shell.nucl1);
132 if (shell.nnucl == 2)
134 fprintf(fplog, " %5d\n", shell.nucl2);
136 else if (shell.nnucl == 3)
138 fprintf(fplog, " %5d %5d\n", shell.nucl2, shell.nucl3);
142 fprintf(fplog, "\n");
147 /* TODO The remain call of this function passes non-NULL mass and NULL
148 * mtop, so this routine can be simplified.
150 * The other code path supported doing prediction before the MD loop
151 * started, but even when called, the prediction was always
152 * over-written by a subsequent call in the MD loop, so has been
154 static void predict_shells(FILE* fplog,
158 ArrayRef<const t_shell> shells,
164 real dt_1, fudge, tm, m1, m2, m3;
166 GMX_RELEASE_ASSERT(mass || mtop, "Must have masses or a way to look them up");
168 /* We introduce a fudge factor for performance reasons: with this choice
169 * the initial force on the shells is about a factor of two lower than
179 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
191 for (const t_shell& shell : shells)
193 const int s1 = shell.shellIndex;
202 for (m = 0; (m < DIM); m++)
204 x[s1][m] += xOrV[n1][m] * dt_1;
217 /* Not the correct masses with FE, but it is just a prediction... */
218 m1 = mtopGetAtomMass(mtop, n1, &molb);
219 m2 = mtopGetAtomMass(mtop, n2, &molb);
221 tm = dt_1 / (m1 + m2);
222 for (m = 0; (m < DIM); m++)
224 x[s1][m] += (m1 * xOrV[n1][m] + m2 * xOrV[n2][m]) * tm;
239 /* Not the correct masses with FE, but it is just a prediction... */
240 m1 = mtopGetAtomMass(mtop, n1, &molb);
241 m2 = mtopGetAtomMass(mtop, n2, &molb);
242 m3 = mtopGetAtomMass(mtop, n3, &molb);
244 tm = dt_1 / (m1 + m2 + m3);
245 for (m = 0; (m < DIM); m++)
247 x[s1][m] += (m1 * xOrV[n1][m] + m2 * xOrV[n2][m] + m3 * xOrV[n3][m]) * tm;
250 default: gmx_fatal(FARGS, "Shell %d has %d nuclei!", s1, shell.nnucl);
255 /*! \brief Count the different particle types in a system
257 * Routine prints a warning to stderr in case an unknown particle type
259 * \param[in] fplog Print what we have found if not NULL
260 * \param[in] mtop Molecular topology.
261 * \returns Array holding the number of particles of a type
263 std::array<int, eptNR> countPtypes(FILE* fplog, const gmx_mtop_t* mtop)
265 std::array<int, eptNR> nptype = { { 0 } };
266 /* Count number of shells, and find their indices */
267 for (int i = 0; (i < eptNR); i++)
272 gmx_mtop_atomloop_block_t aloopb = gmx_mtop_atomloop_block_init(mtop);
275 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
281 case eptShell: nptype[atom->ptype] += nmol; break;
283 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
284 static_cast<int>(atom->ptype));
289 /* Print the number of each particle type */
291 for (const auto& i : nptype)
295 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
303 gmx_shellfc_t* init_shell_flexcon(FILE* fplog, const gmx_mtop_t* mtop, int nflexcon, int nstcalcenergy, bool usingDomainDecomposition)
308 int i, j, type, a_offset, mol, ftype, nra;
310 int aS, aN = 0; /* Shell and nucleus */
311 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
312 #define NBT asize(bondtypes)
313 const gmx_ffparams_t* ffparams;
315 std::array<int, eptNR> n = countPtypes(fplog, mtop);
316 nshell = n[eptShell];
318 if (nshell == 0 && nflexcon == 0)
320 /* We're not doing shells or flexible constraints */
324 shfc = new gmx_shellfc_t;
325 shfc->nflexcon = nflexcon;
329 /* Only flexible constraints, no shells.
330 * Note that make_local_shells() does not need to be called.
335 if (nstcalcenergy != 1)
338 "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is "
339 "not supported in combination with shell particles.\nPlease make a new tpr file.",
342 if (usingDomainDecomposition)
346 "Shell particles are not implemented with domain decomposition, use a single rank");
349 /* We have shells: fill the shell data structure */
351 /* Global system sized array, this should be avoided */
352 std::vector<int> shell_index(mtop->natoms);
355 for (const AtomProxy atomP : AtomRange(*mtop))
357 const t_atom& local = atomP.atom();
358 int i = atomP.globalAtomNumber();
359 if (local.ptype == eptShell)
361 shell_index[i] = nshell++;
365 std::vector<t_shell> shell(nshell);
367 ffparams = &mtop->ffparams;
369 /* Now fill the structures */
370 /* TODO: See if we can use update groups that cover shell constructions */
371 shfc->bInterCG = FALSE;
374 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
376 const gmx_molblock_t* molb = &mtop->molblock[mb];
377 const gmx_moltype_t* molt = &mtop->moltype[molb->type];
379 const t_atom* atom = molt->atoms.atom;
380 for (mol = 0; mol < molb->nmol; mol++)
382 for (j = 0; (j < NBT); j++)
384 const int* ia = molt->ilist[bondtypes[j]].iatoms.data();
385 for (i = 0; (i < molt->ilist[bondtypes[j]].size());)
388 ftype = ffparams->functype[type];
389 nra = interaction_function[ftype].nratoms;
391 /* Check whether we have a bond with a shell */
394 switch (bondtypes[j])
401 if (atom[ia[1]].ptype == eptShell)
406 else if (atom[ia[2]].ptype == eptShell)
413 aN = ia[4]; /* Dummy */
414 aS = ia[5]; /* Shell */
416 default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
423 /* Check whether one of the particles is a shell... */
424 nsi = shell_index[a_offset + aS];
425 if ((nsi < 0) || (nsi >= nshell))
427 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d", nsi, nshell, aS);
429 if (shell[nsi].shellIndex == -1)
431 shell[nsi].shellIndex = a_offset + aS;
434 else if (shell[nsi].shellIndex != a_offset + aS)
436 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
439 if (shell[nsi].nucl1 == -1)
441 shell[nsi].nucl1 = a_offset + aN;
443 else if (shell[nsi].nucl2 == -1)
445 shell[nsi].nucl2 = a_offset + aN;
447 else if (shell[nsi].nucl3 == -1)
449 shell[nsi].nucl3 = a_offset + aN;
455 pr_shell(fplog, shell);
457 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
461 /* shell[nsi].bInterCG = TRUE; */
462 shfc->bInterCG = TRUE;
465 switch (bondtypes[j])
469 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
472 shell[nsi].k += ffparams->iparams[type].cubic.kb;
476 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS * 10))
479 "polarize can not be used with qA(%e) != qB(%e) for "
480 "atom %d of molecule block %zu",
481 qS, atom[aS].qB, aS + 1, mb + 1);
483 shell[nsi].k += gmx::square(qS) * ONE_4PI_EPS0
484 / ffparams->iparams[type].polarize.alpha;
487 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS * 10))
490 "water_pol can not be used with qA(%e) != qB(%e) for "
491 "atom %d of molecule block %zu",
492 qS, atom[aS].qB, aS + 1, mb + 1);
494 alpha = (ffparams->iparams[type].wpol.al_x
495 + ffparams->iparams[type].wpol.al_y
496 + ffparams->iparams[type].wpol.al_z)
498 shell[nsi].k += gmx::square(qS) * ONE_4PI_EPS0 / alpha;
500 default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
508 a_offset += molt->atoms.nr;
510 /* Done with this molecule type */
513 /* Verify whether it's all correct */
516 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
519 for (i = 0; (i < ns); i++)
521 shell[i].k_1 = 1.0 / shell[i].k;
526 pr_shell(debug, shell);
530 shfc->shell_gl = shell;
531 shfc->shell_index_gl = shell_index;
533 shfc->predictShells = (getenv("GMX_NOPREDICT") == nullptr);
534 shfc->requireInit = false;
535 if (!shfc->predictShells)
539 fprintf(fplog, "\nWill never predict shell positions\n");
544 shfc->requireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
545 if (shfc->requireInit && fplog)
547 fprintf(fplog, "\nWill always initiate shell positions\n");
551 if (shfc->predictShells)
558 "\nNOTE: in the current version shell prediction during the crun is "
561 /* Prediction improves performance, so we should implement either:
562 * 1. communication for the atoms needed for prediction
563 * 2. prediction using the velocities of shells; currently the
564 * shell velocities are zeroed, it's a bit tricky to keep
565 * track of the shell displacements and thus the velocity.
567 shfc->predictShells = false;
574 void make_local_shells(const t_commrec* cr, const t_mdatoms* md, gmx_shellfc_t* shfc)
577 gmx_domdec_t* dd = nullptr;
579 if (DOMAINDECOMP(cr))
583 a1 = dd_numHomeAtoms(*dd);
587 /* Single node: we need all shells, copy them */
588 shfc->shells = shfc->shell_gl;
593 ArrayRef<const int> ind = shfc->shell_index_gl;
595 std::vector<t_shell>& shells = shfc->shells;
597 for (int i = a0; i < a1; i++)
599 if (md->ptype[i] == eptShell)
603 shells.push_back(shfc->shell_gl[ind[dd->globalAtomIndices[i]]]);
607 shells.push_back(shfc->shell_gl[ind[i]]);
609 t_shell& shell = shells.back();
611 /* With inter-cg shells we can no do shell prediction,
612 * so we do not need the nuclei numbers.
616 shell.nucl1 = i + shell.nucl1 - shell.shellIndex;
619 shell.nucl2 = i + shell.nucl2 - shell.shellIndex;
623 shell.nucl3 = i + shell.nucl3 - shell.shellIndex;
626 shell.shellIndex = i;
631 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
649 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
658 dx = f[XX] * step[XX];
659 dy = f[YY] * step[YY];
660 dz = f[ZZ] * step[ZZ];
667 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
668 gmx::ArrayRef<gmx::RVec> xnew,
669 ArrayRef<const gmx::RVec> acc_dir,
673 const rvec* xo = as_rvec_array(xold.data());
674 rvec* xn = as_rvec_array(xnew.data());
676 for (int i = 0; i < homenr; i++)
678 do_1pos(xn[i], xo[i], acc_dir[i], step);
682 static void shell_pos_sd(gmx::ArrayRef<const gmx::RVec> xcur,
683 gmx::ArrayRef<gmx::RVec> xnew,
684 gmx::ArrayRef<const gmx::RVec> f,
685 ArrayRef<t_shell> shells,
688 const real step_scale_min = 0.8, step_scale_increment = 0.2, step_scale_max = 1.2,
689 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
694 real step_min, step_max;
699 for (t_shell& shell : shells)
701 const int ind = shell.shellIndex;
704 for (d = 0; d < DIM; d++)
706 shell.step[d] = shell.k_1;
708 step_min = std::min(step_min, shell.step[d]);
709 step_max = std::max(step_max, shell.step[d]);
715 for (d = 0; d < DIM; d++)
717 dx = xcur[ind][d] - shell.xold[d];
718 df = f[ind][d] - shell.fold[d];
719 /* -dx/df gets used to generate an interpolated value, but would
720 * cause a NaN if df were binary-equal to zero. Values close to
721 * zero won't cause problems (because of the min() and max()), so
722 * just testing for binary inequality is OK. */
726 /* Scale the step size by a factor interpolated from
727 * step_scale_min to step_scale_max, as k_est goes from 0 to
728 * step_scale_multiple * shell.step[d] */
729 shell.step[d] = step_scale_min * shell.step[d]
730 + step_scale_increment
731 * std::min(step_scale_multiple * shell.step[d],
732 std::max(k_est, zero));
737 if (gmx_numzero(dx)) /* 0 == dx */
739 /* Likely this will never happen, but if it does just
740 * don't scale the step. */
744 shell.step[d] *= step_scale_max;
748 step_min = std::min(step_min, shell.step[d]);
749 step_max = std::max(step_max, shell.step[d]);
753 copy_rvec(xcur[ind], shell.xold);
754 copy_rvec(f[ind], shell.fold);
756 do_1pos3(xnew[ind], xcur[ind], f[ind], shell.step);
760 fprintf(debug, "shell = %d\n", ind);
761 pr_rvec(debug, 0, "fshell", f[ind], DIM, TRUE);
762 pr_rvec(debug, 0, "xold", xcur[ind], DIM, TRUE);
763 pr_rvec(debug, 0, "step", shell.step, DIM, TRUE);
764 pr_rvec(debug, 0, "xnew", xnew[ind], DIM, TRUE);
768 printf("step %.3e %.3e\n", step_min, step_max);
772 static void decrease_step_size(ArrayRef<t_shell> shells)
774 for (t_shell& shell : shells)
776 svmul(0.8, shell.step, shell.step);
780 static void print_epot(FILE* fp, int64_t mdstep, int count, real epot, real df, int ndir, real sf_dir)
784 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e", gmx_step_str(mdstep, buf), count, epot, df);
787 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir / ndir));
796 static real rms_force(const t_commrec* cr,
797 gmx::ArrayRef<const gmx::RVec> force,
798 ArrayRef<const t_shell> shells,
804 const rvec* f = as_rvec_array(force.data());
807 for (const t_shell& shell : shells)
809 buf[0] += norm2(f[shell.shellIndex]);
811 int ntot = shells.ssize();
818 gmx_sumd(4, buf, cr);
819 ntot = gmx::roundToInt(buf[1]);
825 return (ntot ? std::sqrt(buf[0] / ntot) : 0);
828 static void dump_shells(FILE* fp, gmx::ArrayRef<gmx::RVec> f, real ftol, ArrayRef<const t_shell> shells)
832 ft2 = gmx::square(ftol);
834 for (const t_shell& shell : shells)
836 const int ind = shell.shellIndex;
837 ff2 = iprod(f[ind], f[ind]);
840 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n", ind, f[ind][XX],
841 f[ind][YY], f[ind][ZZ], std::sqrt(ff2));
846 static void init_adir(gmx_shellfc_t* shfc,
847 gmx::Constraints* constr,
848 const t_inputrec* ir,
854 ArrayRefWithPadding<RVec> xOld,
855 ArrayRef<RVec> x_init,
856 ArrayRefWithPadding<RVec> xCurrent,
858 ArrayRef<RVec> acc_dir,
860 gmx::ArrayRef<const real> lambda,
865 unsigned short* ptype;
867 if (DOMAINDECOMP(cr))
875 shfc->adir_xnold.resizeWithPadding(n);
876 shfc->adir_xnew.resizeWithPadding(n);
877 rvec* xnold = as_rvec_array(shfc->adir_xnold.data());
878 rvec* xnew = as_rvec_array(shfc->adir_xnew.data());
879 rvec* x_old = as_rvec_array(xOld.paddedArrayRef().data());
880 rvec* x = as_rvec_array(xCurrent.paddedArrayRef().data());
886 /* Does NOT work with freeze or acceleration groups (yet) */
887 for (n = 0; n < end; n++)
889 w_dt = md->invmass[n] * dt;
891 for (d = 0; d < DIM; d++)
893 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
895 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
896 xnew[n][d] = 2 * x[n][d] - x_old[n][d] + f[n][d] * w_dt * dt;
900 xnold[n][d] = x[n][d];
901 xnew[n][d] = x[n][d];
905 constr->apply(FALSE, FALSE, step, 0, 1.0, xCurrent, shfc->adir_xnold.arrayRefWithPadding(),
906 ArrayRef<RVec>(), box, lambda[efptBONDED], &(dvdlambda[efptBONDED]),
907 ArrayRefWithPadding<RVec>(), nullptr, gmx::ConstraintVariable::Positions);
908 constr->apply(FALSE, FALSE, step, 0, 1.0, xCurrent, shfc->adir_xnew.arrayRefWithPadding(),
909 ArrayRef<RVec>(), box, lambda[efptBONDED], &(dvdlambda[efptBONDED]),
910 ArrayRefWithPadding<RVec>(), nullptr, gmx::ConstraintVariable::Positions);
912 for (n = 0; n < end; n++)
914 for (d = 0; d < DIM; d++)
916 xnew[n][d] = -(2 * x[n][d] - xnold[n][d] - xnew[n][d]) / gmx::square(dt)
917 - f[n][d] * md->invmass[n];
919 clear_rvec(acc_dir[n]);
922 /* Project the acceleration on the old bond directions */
923 constr->apply(FALSE, FALSE, step, 0, 1.0, xOld, shfc->adir_xnew.arrayRefWithPadding(), acc_dir,
924 box, lambda[efptBONDED], &(dvdlambda[efptBONDED]), ArrayRefWithPadding<RVec>(),
925 nullptr, gmx::ConstraintVariable::Deriv_FlexCon);
928 void relax_shell_flexcon(FILE* fplog,
930 const gmx_multisim_t* ms,
932 gmx_enfrot* enforcedRotation,
934 const t_inputrec* inputrec,
935 gmx::ImdSession* imdSession,
939 const gmx_localtop_t* top,
940 gmx::Constraints* constr,
941 gmx_enerdata_t* enerd,
944 gmx::ArrayRefWithPadding<gmx::RVec> xPadded,
945 gmx::ArrayRefWithPadding<gmx::RVec> vPadded,
947 gmx::ArrayRef<real> lambda,
949 gmx::ArrayRefWithPadding<gmx::RVec> f,
953 gmx_wallcycle_t wcycle,
957 gmx::MdrunScheduleWorkload* runScheduleWork,
960 const gmx_vsite_t* vsite,
961 const DDBalanceRegionHandler& ddBalanceRegionHandler)
967 int nat, dd_ac0, dd_ac1 = 0, i;
968 int homenr = md->homenr, end = homenr;
969 int d, Min = 0, count = 0;
970 #define Try (1 - Min) /* At start Try = 1 */
972 const bool bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
973 const bool bInit = (mdstep == inputrec->init_step) || shfc->requireInit;
974 const real ftol = inputrec->em_tol;
975 const int number_steps = inputrec->niter;
976 ArrayRef<t_shell> shells = shfc->shells;
977 const int nflexcon = shfc->nflexcon;
979 const t_idef* idef = &top->idef;
981 if (DOMAINDECOMP(cr))
983 nat = dd_natoms_vsite(cr->dd);
986 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
987 nat = std::max(nat, dd_ac1);
995 for (i = 0; (i < 2); i++)
997 shfc->x[i].resizeWithPadding(nat);
998 shfc->f[i].resizeWithPadding(nat);
1001 /* Create views that we can swap for trail and minimum for positions and forces */
1002 gmx::ArrayRefWithPadding<gmx::RVec> posWithPadding[2];
1003 gmx::ArrayRefWithPadding<gmx::RVec> forceWithPadding[2];
1004 gmx::ArrayRef<gmx::RVec> pos[2];
1005 gmx::ArrayRef<gmx::RVec> force[2];
1006 for (i = 0; (i < 2); i++)
1008 posWithPadding[i] = shfc->x[i].arrayRefWithPadding();
1009 pos[i] = posWithPadding[i].paddedArrayRef();
1010 forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1011 force[i] = forceWithPadding[i].paddedArrayRef();
1014 ArrayRef<RVec> x = xPadded.unpaddedArrayRef();
1015 ArrayRef<RVec> v = vPadded.unpaddedArrayRef();
1017 if (bDoNS && inputrec->pbcType != PbcType::No && !DOMAINDECOMP(cr))
1019 /* This is the only time where the coordinates are used
1020 * before do_force is called, which normally puts all
1021 * charge groups in the box.
1023 put_atoms_in_box_omp(fr->pbcType, box, x.subArray(0, md->homenr),
1024 gmx_omp_nthreads_get(emntDefault));
1028 mk_mshift(fplog, graph, fr->pbcType, box, as_rvec_array(x.data()));
1032 /* After this all coordinate arrays will contain whole charge groups */
1035 shift_self(graph, box, as_rvec_array(x.data()));
1040 shfc->acc_dir.resize(nat);
1041 shfc->x_old.resizeWithPadding(nat);
1042 ArrayRef<RVec> x_old = shfc->x_old.arrayRefWithPadding().unpaddedArrayRef();
1043 for (i = 0; i < homenr; i++)
1045 for (d = 0; d < DIM; d++)
1047 x_old[i][d] = x[i][d] - v[i][d] * inputrec->delta_t;
1052 /* Do a prediction of the shell positions, when appropriate.
1053 * Without velocities (EM, NM, BD) we only do initial prediction.
1055 if (shfc->predictShells && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1057 predict_shells(fplog, x, v, inputrec->delta_t, shells, md->massT, nullptr, bInit);
1060 /* do_force expected the charge groups to be in the box */
1063 unshift_self(graph, box, as_rvec_array(x.data()));
1066 /* Calculate the forces first time around */
1069 pr_rvecs(debug, 0, "x b4 do_force", as_rvec_array(x.data()), homenr);
1071 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1072 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession, pull_work, mdstep,
1073 nrnb, wcycle, top, box, xPadded, hist, forceWithPadding[Min], force_vir, md, enerd,
1074 fcd, lambda, graph, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
1075 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags, ddBalanceRegionHandler);
1080 init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end, shfc->x_old.arrayRefWithPadding(),
1081 x, xPadded, force[Min], shfc->acc_dir, box, lambda, &dum);
1083 for (i = 0; i < end; i++)
1085 sf_dir += md->massT[i] * norm2(shfc->acc_dir[i]);
1088 sum_epot(&(enerd->grpp), enerd->term);
1089 Epot[Min] = enerd->term[F_EPOT];
1091 df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), shells, nflexcon, &sf_dir, &Epot[Min]);
1095 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1100 pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md->nr);
1103 if (!shells.empty() || nflexcon > 0)
1105 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1106 * shell positions are updated, therefore the other particles must
1107 * be set here, in advance.
1109 std::copy(xPadded.paddedArrayRef().begin(), xPadded.paddedArrayRef().end(),
1110 posWithPadding[Min].paddedArrayRef().begin());
1111 std::copy(xPadded.paddedArrayRef().begin(), xPadded.paddedArrayRef().end(),
1112 posWithPadding[Try].paddedArrayRef().begin());
1115 if (bVerbose && MASTER(cr))
1117 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1122 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1123 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1124 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1125 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1128 /* First check whether we should do shells, or whether the force is
1129 * low enough even without minimization.
1131 bool bConverged = (df[Min] < ftol);
1133 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1137 construct_vsites(vsite, as_rvec_array(pos[Min].data()), inputrec->delta_t,
1138 as_rvec_array(v.data()), idef->iparams, idef->il, fr->pbcType,
1139 fr->bMolPBC, cr, box);
1144 init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end,
1145 shfc->x_old.arrayRefWithPadding(), x, posWithPadding[Min], force[Min],
1146 shfc->acc_dir, box, lambda, &dum);
1148 directional_sd(pos[Min], pos[Try], shfc->acc_dir, end, fr->fc_stepsize);
1151 /* New positions, Steepest descent */
1152 shell_pos_sd(pos[Min], pos[Try], force[Min], shells, count);
1154 /* do_force expected the charge groups to be in the box */
1157 unshift_self(graph, box, as_rvec_array(pos[Try].data()));
1162 pr_rvecs(debug, 0, "RELAX: pos[Min] ", as_rvec_array(pos[Min].data()), homenr);
1163 pr_rvecs(debug, 0, "RELAX: pos[Try] ", as_rvec_array(pos[Try].data()), homenr);
1165 /* Try the new positions */
1166 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession, pull_work, 1, nrnb,
1167 wcycle, top, box, posWithPadding[Try], hist, forceWithPadding[Try], force_vir, md,
1168 enerd, fcd, lambda, graph, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
1169 shellfc_flags, ddBalanceRegionHandler);
1170 sum_epot(&(enerd->grpp), enerd->term);
1173 pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1174 pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1179 init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end,
1180 shfc->x_old.arrayRefWithPadding(), x, posWithPadding[Try], force[Try],
1181 shfc->acc_dir, box, lambda, &dum);
1183 ArrayRef<const gmx::RVec> acc_dir = shfc->acc_dir;
1184 for (i = 0; i < end; i++)
1186 sf_dir += md->massT[i] * norm2(acc_dir[i]);
1190 Epot[Try] = enerd->term[F_EPOT];
1192 df[Try] = rms_force(cr, force[Try], shells, nflexcon, &sf_dir, &Epot[Try]);
1196 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1203 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1207 fprintf(debug, "SHELL ITER %d\n", count);
1208 dump_shells(debug, force[Try], ftol, shells);
1212 if (bVerbose && MASTER(cr))
1214 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1217 bConverged = (df[Try] < ftol);
1219 if ((df[Try] < df[Min]))
1223 fprintf(debug, "Swapping Min and Try\n");
1227 /* Correct the velocities for the flexible constraints */
1228 invdt = 1 / inputrec->delta_t;
1229 for (i = 0; i < end; i++)
1231 for (d = 0; d < DIM; d++)
1233 v[i][d] += (pos[Try][i][d] - pos[Min][i][d]) * invdt;
1241 decrease_step_size(shells);
1244 shfc->numForceEvaluations += count;
1247 shfc->numConvergedIterations++;
1249 if (MASTER(cr) && !(bConverged))
1251 /* Note that the energies and virial are incorrect when not converged */
1254 fprintf(fplog, "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1255 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1257 fprintf(stderr, "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1258 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1261 /* Copy back the coordinates and the forces */
1262 std::copy(pos[Min].begin(), pos[Min].end(), x.data());
1263 std::copy(force[Min].begin(), force[Min].end(), f.unpaddedArrayRef().begin());
1266 void done_shellfc(FILE* fplog, gmx_shellfc_t* shfc, int64_t numSteps)
1268 if (shfc && fplog && numSteps > 0)
1270 double numStepsAsDouble = static_cast<double>(numSteps);
1271 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1272 (shfc->numConvergedIterations * 100.0) / numStepsAsDouble);
1273 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1274 shfc->numForceEvaluations / numStepsAsDouble);