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"
86 int shell; /* The shell id */
87 int nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
88 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
89 real k; /* force constant */
90 real k_1; /* 1 over force constant */
98 /* Shell counts, indices, parameters and working data */
99 int nshell_gl; /* The number of shells in the system */
100 t_shell* shell_gl; /* All the shells (for DD only) */
101 int* shell_index_gl; /* Global shell index (for DD only) */
102 gmx_bool bInterCG; /* Are there inter charge-group shells? */
103 int nshell; /* The number of local shells */
104 t_shell* shell; /* The local shells */
105 int shell_nalloc; /* The allocation size of shell */
106 gmx_bool bPredict; /* Predict shell positions */
107 gmx_bool bRequireInit; /* Require initialization of shell positions */
108 int nflexcon; /* The number of flexible constraints */
110 /* Temporary arrays, should be fixed size 2 when fully converted to C++ */
111 PaddedHostVector<gmx::RVec>* x; /* Array for iterative minimization */
112 PaddedHostVector<gmx::RVec>* f; /* Array for iterative minimization */
114 /* Flexible constraint working data */
115 rvec* acc_dir; /* Acceleration direction for flexcon */
116 rvec* x_old; /* Old coordinates for flexcon */
117 int flex_nalloc; /* The allocation size of acc_dir and x_old */
118 rvec* adir_xnold; /* Work space for init_adir */
119 rvec* adir_xnew; /* Work space for init_adir */
120 int adir_nalloc; /* Work space for init_adir */
121 std::int64_t numForceEvaluations; /* Total number of force evaluations */
122 int numConvergedIterations; /* Total number of iterations that converged */
126 static void pr_shell(FILE* fplog, int ns, t_shell s[])
130 fprintf(fplog, "SHELL DATA\n");
131 fprintf(fplog, "%5s %8s %5s %5s %5s\n", "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
132 for (i = 0; (i < ns); i++)
134 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0 / s[i].k_1, s[i].nucl1);
137 fprintf(fplog, " %5d\n", s[i].nucl2);
139 else if (s[i].nnucl == 3)
141 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
145 fprintf(fplog, "\n");
150 /* TODO The remain call of this function passes non-NULL mass and NULL
151 * mtop, so this routine can be simplified.
153 * The other code path supported doing prediction before the MD loop
154 * started, but even when called, the prediction was always
155 * over-written by a subsequent call in the MD loop, so has been
157 static void predict_shells(FILE* fplog,
167 int i, m, s1, n1, n2, n3;
168 real dt_1, fudge, tm, m1, m2, m3;
171 GMX_RELEASE_ASSERT(mass || mtop, "Must have masses or a way to look them up");
173 /* We introduce a fudge factor for performance reasons: with this choice
174 * the initial force on the shells is about a factor of two lower than
183 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
195 for (i = 0; (i < ns); i++)
206 for (m = 0; (m < DIM); m++)
208 x[s1][m] += ptr[n1][m] * dt_1;
221 /* Not the correct masses with FE, but it is just a prediction... */
222 m1 = mtopGetAtomMass(mtop, n1, &molb);
223 m2 = mtopGetAtomMass(mtop, n2, &molb);
225 tm = dt_1 / (m1 + m2);
226 for (m = 0; (m < DIM); m++)
228 x[s1][m] += (m1 * ptr[n1][m] + m2 * ptr[n2][m]) * tm;
243 /* Not the correct masses with FE, but it is just a prediction... */
244 m1 = mtopGetAtomMass(mtop, n1, &molb);
245 m2 = mtopGetAtomMass(mtop, n2, &molb);
246 m3 = mtopGetAtomMass(mtop, n3, &molb);
248 tm = dt_1 / (m1 + m2 + m3);
249 for (m = 0; (m < DIM); m++)
251 x[s1][m] += (m1 * ptr[n1][m] + m2 * ptr[n2][m] + m3 * ptr[n3][m]) * tm;
254 default: gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
259 /*! \brief Count the different particle types in a system
261 * Routine prints a warning to stderr in case an unknown particle type
263 * \param[in] fplog Print what we have found if not NULL
264 * \param[in] mtop Molecular topology.
265 * \returns Array holding the number of particles of a type
267 std::array<int, eptNR> countPtypes(FILE* fplog, const gmx_mtop_t* mtop)
269 std::array<int, eptNR> nptype = { { 0 } };
270 /* Count number of shells, and find their indices */
271 for (int i = 0; (i < eptNR); i++)
276 gmx_mtop_atomloop_block_t aloopb = gmx_mtop_atomloop_block_init(mtop);
279 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
285 case eptShell: nptype[atom->ptype] += nmol; break;
287 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
288 static_cast<int>(atom->ptype));
293 /* Print the number of each particle type */
295 for (const auto& i : nptype)
299 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
307 gmx_shellfc_t* init_shell_flexcon(FILE* fplog, const gmx_mtop_t* mtop, int nflexcon, int nstcalcenergy, 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)
349 "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is "
350 "not supported in combination with shell particles.\nPlease make a new tpr file.",
353 if (usingDomainDecomposition)
357 "Shell particles are not implemented with domain decomposition, use a single rank");
360 /* We have shells: fill the shell data structure */
362 /* Global system sized array, this should be avoided */
363 snew(shell_index, mtop->natoms);
366 for (const AtomProxy atomP : AtomRange(*mtop))
368 const t_atom& local = atomP.atom();
369 int i = atomP.globalAtomNumber();
370 if (local.ptype == eptShell)
372 shell_index[i] = nshell++;
378 /* Initiate the shell structures */
379 for (i = 0; (i < nshell); i++)
386 /* shell[i].bInterCG=FALSE; */
391 ffparams = &mtop->ffparams;
393 /* Now fill the structures */
394 /* TODO: See if we can use update groups that cover shell constructions */
395 shfc->bInterCG = FALSE;
398 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
400 const gmx_molblock_t* molb = &mtop->molblock[mb];
401 const gmx_moltype_t* molt = &mtop->moltype[molb->type];
403 const t_atom* atom = molt->atoms.atom;
404 for (mol = 0; mol < molb->nmol; mol++)
406 for (j = 0; (j < NBT); j++)
408 const int* ia = molt->ilist[bondtypes[j]].iatoms.data();
409 for (i = 0; (i < molt->ilist[bondtypes[j]].size());)
412 ftype = ffparams->functype[type];
413 nra = interaction_function[ftype].nratoms;
415 /* Check whether we have a bond with a shell */
418 switch (bondtypes[j])
425 if (atom[ia[1]].ptype == eptShell)
430 else if (atom[ia[2]].ptype == eptShell)
437 aN = ia[4]; /* Dummy */
438 aS = ia[5]; /* Shell */
440 default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
447 /* Check whether one of the particles is a shell... */
448 nsi = shell_index[a_offset + aS];
449 if ((nsi < 0) || (nsi >= nshell))
451 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d", nsi, nshell, aS);
453 if (shell[nsi].shell == -1)
455 shell[nsi].shell = a_offset + aS;
458 else if (shell[nsi].shell != a_offset + aS)
460 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
463 if (shell[nsi].nucl1 == -1)
465 shell[nsi].nucl1 = a_offset + aN;
467 else if (shell[nsi].nucl2 == -1)
469 shell[nsi].nucl2 = a_offset + aN;
471 else if (shell[nsi].nucl3 == -1)
473 shell[nsi].nucl3 = a_offset + aN;
479 pr_shell(fplog, ns, shell);
481 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
485 /* shell[nsi].bInterCG = TRUE; */
486 shfc->bInterCG = TRUE;
489 switch (bondtypes[j])
493 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
496 shell[nsi].k += ffparams->iparams[type].cubic.kb;
500 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS * 10))
503 "polarize can not be used with qA(%e) != qB(%e) for "
504 "atom %d of molecule block %zu",
505 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))
514 "water_pol can not be used with qA(%e) != qB(%e) for "
515 "atom %d of molecule block %zu",
516 qS, atom[aS].qB, aS + 1, mb + 1);
518 alpha = (ffparams->iparams[type].wpol.al_x
519 + ffparams->iparams[type].wpol.al_y
520 + ffparams->iparams[type].wpol.al_z)
522 shell[nsi].k += gmx::square(qS) * ONE_4PI_EPS0 / alpha;
524 default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
532 a_offset += molt->atoms.nr;
534 /* Done with this molecule type */
537 /* Verify whether it's all correct */
540 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
543 for (i = 0; (i < ns); i++)
545 shell[i].k_1 = 1.0 / shell[i].k;
550 pr_shell(debug, ns, shell);
554 shfc->nshell_gl = ns;
555 shfc->shell_gl = shell;
556 shfc->shell_index_gl = shell_index;
558 shfc->bPredict = (getenv("GMX_NOPREDICT") == nullptr);
559 shfc->bRequireInit = FALSE;
564 fprintf(fplog, "\nWill never predict shell positions\n");
569 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
570 if (shfc->bRequireInit && fplog)
572 fprintf(fplog, "\nWill always initiate shell positions\n");
583 "\nNOTE: in the current version shell prediction during the crun is "
586 /* Prediction improves performance, so we should implement either:
587 * 1. communication for the atoms needed for prediction
588 * 2. prediction using the velocities of shells; currently the
589 * shell velocities are zeroed, it's a bit tricky to keep
590 * track of the shell displacements and thus the velocity.
592 shfc->bPredict = FALSE;
599 void make_local_shells(const t_commrec* cr, const t_mdatoms* md, gmx_shellfc_t* shfc)
602 int a0, a1, *ind, nshell, i;
603 gmx_domdec_t* dd = nullptr;
605 if (DOMAINDECOMP(cr))
609 a1 = dd_numHomeAtoms(*dd);
613 /* Single node: we need all shells, just copy the pointer */
614 shfc->nshell = shfc->nshell_gl;
615 shfc->shell = shfc->shell_gl;
620 ind = shfc->shell_index_gl;
624 for (i = a0; i < a1; i++)
626 if (md->ptype[i] == eptShell)
628 if (nshell + 1 > shfc->shell_nalloc)
630 shfc->shell_nalloc = over_alloc_dd(nshell + 1);
631 srenew(shell, shfc->shell_nalloc);
635 shell[nshell] = shfc->shell_gl[ind[dd->globalAtomIndices[i]]];
639 shell[nshell] = shfc->shell_gl[ind[i]];
642 /* With inter-cg shells we can no do shell prediction,
643 * so we do not need the nuclei numbers.
647 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
648 if (shell[nshell].nnucl > 1)
650 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
652 if (shell[nshell].nnucl > 2)
654 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
657 shell[nshell].shell = i;
662 shfc->nshell = nshell;
666 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
684 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
693 dx = f[XX] * step[XX];
694 dy = f[YY] * step[YY];
695 dz = f[ZZ] * step[ZZ];
702 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
703 gmx::ArrayRef<gmx::RVec> xnew,
704 const rvec acc_dir[],
708 const rvec* xo = as_rvec_array(xold.data());
709 rvec* xn = as_rvec_array(xnew.data());
711 for (int i = 0; i < homenr; i++)
713 do_1pos(xn[i], xo[i], acc_dir[i], step);
717 static void shell_pos_sd(gmx::ArrayRef<const gmx::RVec> xcur,
718 gmx::ArrayRef<gmx::RVec> xnew,
719 gmx::ArrayRef<const gmx::RVec> f,
724 const real step_scale_min = 0.8, step_scale_increment = 0.2, step_scale_max = 1.2,
725 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
730 real step_min, step_max;
735 for (i = 0; (i < ns); i++)
740 for (d = 0; d < DIM; d++)
742 s[i].step[d] = s[i].k_1;
744 step_min = std::min(step_min, s[i].step[d]);
745 step_max = std::max(step_max, s[i].step[d]);
751 for (d = 0; d < DIM; d++)
753 dx = xcur[shell][d] - s[i].xold[d];
754 df = f[shell][d] - s[i].fold[d];
755 /* -dx/df gets used to generate an interpolated value, but would
756 * cause a NaN if df were binary-equal to zero. Values close to
757 * zero won't cause problems (because of the min() and max()), so
758 * just testing for binary inequality is OK. */
762 /* Scale the step size by a factor interpolated from
763 * step_scale_min to step_scale_max, as k_est goes from 0 to
764 * step_scale_multiple * s[i].step[d] */
765 s[i].step[d] = step_scale_min * s[i].step[d]
766 + step_scale_increment
767 * std::min(step_scale_multiple * s[i].step[d],
768 std::max(k_est, zero));
773 if (gmx_numzero(dx)) /* 0 == dx */
775 /* Likely this will never happen, but if it does just
776 * don't scale the step. */
780 s[i].step[d] *= step_scale_max;
784 step_min = std::min(step_min, s[i].step[d]);
785 step_max = std::max(step_max, s[i].step[d]);
789 copy_rvec(xcur[shell], s[i].xold);
790 copy_rvec(f[shell], s[i].fold);
792 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
796 fprintf(debug, "shell[%d] = %d\n", i, shell);
797 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
798 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
799 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
800 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
804 printf("step %.3e %.3e\n", step_min, step_max);
808 static void decrease_step_size(int nshell, t_shell s[])
812 for (i = 0; i < nshell; i++)
814 svmul(0.8, s[i].step, s[i].step);
818 static void print_epot(FILE* fp, int64_t mdstep, int count, real epot, real df, int ndir, real sf_dir)
822 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e", gmx_step_str(mdstep, buf), count, epot, df);
825 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir / ndir));
834 static real rms_force(const t_commrec* cr,
835 gmx::ArrayRef<const gmx::RVec> force,
843 const rvec* f = as_rvec_array(force.data());
846 for (int i = 0; i < ns; i++)
848 int shell = s[i].shell;
849 buf[0] += norm2(f[shell]);
858 gmx_sumd(4, buf, cr);
859 ntot = gmx::roundToInt(buf[1]);
865 return (ntot ? std::sqrt(buf[0] / ntot) : 0);
868 static void dump_shells(FILE* fp, gmx::ArrayRef<gmx::RVec> f, real ftol, int ns, t_shell s[])
873 ft2 = gmx::square(ftol);
875 for (i = 0; (i < ns); i++)
878 ff2 = iprod(f[shell], f[shell]);
881 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n", shell,
882 f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
887 static void init_adir(gmx_shellfc_t* shfc,
888 gmx::Constraints* constr,
889 const t_inputrec* ir,
901 gmx::ArrayRef<const real> lambda,
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 constr->apply(FALSE, FALSE, step, 0, 1.0, x, xnold, nullptr, box, lambda[efptBONDED],
950 &(dvdlambda[efptBONDED]), nullptr, nullptr, gmx::ConstraintVariable::Positions);
951 constr->apply(FALSE, FALSE, step, 0, 1.0, x, xnew, nullptr, box, lambda[efptBONDED],
952 &(dvdlambda[efptBONDED]), nullptr, nullptr, gmx::ConstraintVariable::Positions);
954 for (n = 0; n < end; n++)
956 for (d = 0; d < DIM; d++)
958 xnew[n][d] = -(2 * x[n][d] - xnold[n][d] - xnew[n][d]) / gmx::square(dt)
959 - f[n][d] * md->invmass[n];
961 clear_rvec(acc_dir[n]);
964 /* Project the acceleration on the old bond directions */
965 constr->apply(FALSE, FALSE, step, 0, 1.0, x_old, xnew, acc_dir, box, lambda[efptBONDED],
966 &(dvdlambda[efptBONDED]), nullptr, nullptr, gmx::ConstraintVariable::Deriv_FlexCon);
969 void relax_shell_flexcon(FILE* fplog,
971 const gmx_multisim_t* ms,
973 gmx_enfrot* enforcedRotation,
975 const t_inputrec* inputrec,
976 gmx::ImdSession* imdSession,
980 const gmx_localtop_t* top,
981 gmx::Constraints* constr,
982 gmx_enerdata_t* enerd,
985 gmx::ArrayRefWithPadding<gmx::RVec> x,
986 gmx::ArrayRefWithPadding<gmx::RVec> v,
988 gmx::ArrayRef<real> lambda,
990 gmx::ArrayRefWithPadding<gmx::RVec> f,
994 gmx_wallcycle_t wcycle,
998 gmx::MdrunScheduleWorkload* runScheduleWork,
1001 const gmx_vsite_t* vsite,
1002 const DDBalanceRegionHandler& ddBalanceRegionHandler)
1004 auto xRvec = as_rvec_array(x.paddedArrayRef().data());
1005 auto vRvec = as_rvec_array(v.paddedArrayRef().data());
1010 rvec * acc_dir = nullptr, *x_old = nullptr;
1011 real Epot[2], df[2];
1015 gmx_bool bCont, bInit, bConverged;
1016 int nat, dd_ac0, dd_ac1 = 0, i;
1017 int homenr = md->homenr, end = homenr;
1018 int nflexcon, number_steps, d, Min = 0, count = 0;
1019 #define Try (1 - Min) /* At start Try = 1 */
1021 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1022 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
1023 ftol = inputrec->em_tol;
1024 number_steps = inputrec->niter;
1025 nshell = shfc->nshell;
1026 shell = shfc->shell;
1027 nflexcon = shfc->nflexcon;
1031 if (DOMAINDECOMP(cr))
1033 nat = dd_natoms_vsite(cr->dd);
1036 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1037 nat = std::max(nat, dd_ac1);
1045 for (i = 0; (i < 2); i++)
1047 shfc->x[i].resizeWithPadding(nat);
1048 shfc->f[i].resizeWithPadding(nat);
1051 /* Create views that we can swap */
1052 gmx::ArrayRefWithPadding<gmx::RVec> posWithPadding[2];
1053 gmx::ArrayRefWithPadding<gmx::RVec> forceWithPadding[2];
1054 gmx::ArrayRef<gmx::RVec> pos[2];
1055 gmx::ArrayRef<gmx::RVec> force[2];
1056 for (i = 0; (i < 2); i++)
1058 posWithPadding[i] = shfc->x[i].arrayRefWithPadding();
1059 pos[i] = posWithPadding[i].paddedArrayRef();
1060 forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1061 force[i] = forceWithPadding[i].paddedArrayRef();
1064 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1066 /* This is the only time where the coordinates are used
1067 * before do_force is called, which normally puts all
1068 * charge groups in the box.
1070 auto xRef = x.paddedArrayRef();
1071 put_atoms_in_box_omp(fr->ePBC, box, xRef.subArray(0, md->homenr),
1072 gmx_omp_nthreads_get(emntDefault));
1076 mk_mshift(fplog, graph, fr->ePBC, box, xRvec);
1080 /* After this all coordinate arrays will contain whole charge groups */
1083 shift_self(graph, box, xRvec);
1088 if (nat > shfc->flex_nalloc)
1090 shfc->flex_nalloc = over_alloc_dd(nat);
1091 srenew(shfc->acc_dir, shfc->flex_nalloc);
1092 srenew(shfc->x_old, shfc->flex_nalloc);
1094 acc_dir = shfc->acc_dir;
1095 x_old = shfc->x_old;
1096 auto xArrayRef = x.paddedArrayRef();
1097 auto vArrayRef = v.paddedArrayRef();
1098 for (i = 0; i < homenr; i++)
1100 for (d = 0; d < DIM; d++)
1102 shfc->x_old[i][d] = xArrayRef[i][d] - vArrayRef[i][d] * inputrec->delta_t;
1107 /* Do a prediction of the shell positions, when appropriate.
1108 * Without velocities (EM, NM, BD) we only do initial prediction.
1110 if (shfc->bPredict && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1112 predict_shells(fplog, xRvec, vRvec, inputrec->delta_t, nshell, shell, md->massT, nullptr, bInit);
1115 /* do_force expected the charge groups to be in the box */
1118 unshift_self(graph, box, xRvec);
1121 /* Calculate the forces first time around */
1124 pr_rvecs(debug, 0, "x b4 do_force", xRvec, homenr);
1126 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1127 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession, pull_work, mdstep,
1128 nrnb, wcycle, top, box, x, hist, forceWithPadding[Min], force_vir, md, enerd, fcd,
1129 lambda, graph, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
1130 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags, ddBalanceRegionHandler);
1135 init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end, shfc->x_old, xRvec, xRvec,
1136 as_rvec_array(force[Min].data()), shfc->acc_dir, box, lambda, &dum);
1138 for (i = 0; i < end; i++)
1140 sf_dir += md->massT[i] * norm2(shfc->acc_dir[i]);
1143 sum_epot(&(enerd->grpp), enerd->term);
1144 Epot[Min] = enerd->term[F_EPOT];
1146 df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), nshell, shell, nflexcon,
1147 &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(), x.paddedArrayRef().end(),
1166 posWithPadding[Min].paddedArrayRef().begin());
1167 std::copy(x.paddedArrayRef().begin(), x.paddedArrayRef().end(),
1168 posWithPadding[Try].paddedArrayRef().begin());
1171 if (bVerbose && MASTER(cr))
1173 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1178 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1179 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1180 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1181 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1184 /* First check whether we should do shells, or whether the force is
1185 * low enough even without minimization.
1187 bConverged = (df[Min] < ftol);
1189 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1193 construct_vsites(vsite, as_rvec_array(pos[Min].data()), inputrec->delta_t, vRvec,
1194 idef->iparams, idef->il, fr->ePBC, fr->bMolPBC, cr, box);
1199 init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end, x_old, xRvec,
1200 as_rvec_array(pos[Min].data()), as_rvec_array(force[Min].data()), acc_dir,
1203 directional_sd(pos[Min], pos[Try], acc_dir, end, fr->fc_stepsize);
1206 /* New positions, Steepest descent */
1207 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1209 /* do_force expected the charge groups to be in the box */
1212 unshift_self(graph, box, as_rvec_array(pos[Try].data()));
1217 pr_rvecs(debug, 0, "RELAX: pos[Min] ", as_rvec_array(pos[Min].data()), homenr);
1218 pr_rvecs(debug, 0, "RELAX: pos[Try] ", as_rvec_array(pos[Try].data()), homenr);
1220 /* Try the new positions */
1221 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession, pull_work, 1, nrnb,
1222 wcycle, top, box, posWithPadding[Try], hist, forceWithPadding[Try], force_vir, md,
1223 enerd, fcd, lambda, graph, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
1224 shellfc_flags, ddBalanceRegionHandler);
1225 sum_epot(&(enerd->grpp), enerd->term);
1228 pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1229 pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1234 init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end, x_old, xRvec,
1235 as_rvec_array(pos[Try].data()), as_rvec_array(force[Try].data()), acc_dir,
1238 for (i = 0; i < end; i++)
1240 sf_dir += md->massT[i] * norm2(acc_dir[i]);
1244 Epot[Try] = enerd->term[F_EPOT];
1246 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1250 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1257 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1261 fprintf(debug, "SHELL ITER %d\n", count);
1262 dump_shells(debug, force[Try], ftol, nshell, shell);
1266 if (bVerbose && MASTER(cr))
1268 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1271 bConverged = (df[Try] < ftol);
1273 if ((df[Try] < df[Min]))
1277 fprintf(debug, "Swapping Min and Try\n");
1281 /* Correct the velocities for the flexible constraints */
1282 invdt = 1 / inputrec->delta_t;
1283 auto vArrayRef = v.paddedArrayRef();
1284 for (i = 0; i < end; i++)
1286 for (d = 0; d < DIM; d++)
1288 vArrayRef[i][d] += (pos[Try][i][d] - pos[Min][i][d]) * invdt;
1296 decrease_step_size(nshell, shell);
1299 shfc->numForceEvaluations += count;
1302 shfc->numConvergedIterations++;
1304 if (MASTER(cr) && !(bConverged))
1306 /* Note that the energies and virial are incorrect when not converged */
1309 fprintf(fplog, "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1310 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1312 fprintf(stderr, "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1313 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1316 /* Copy back the coordinates and the forces */
1317 std::copy(pos[Min].begin(), pos[Min].end(), x.paddedArrayRef().data());
1318 std::copy(force[Min].begin(), force[Min].end(), f.unpaddedArrayRef().begin());
1321 void done_shellfc(FILE* fplog, gmx_shellfc_t* shfc, int64_t numSteps)
1323 if (shfc && fplog && numSteps > 0)
1325 double numStepsAsDouble = static_cast<double>(numSteps);
1326 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1327 (shfc->numConvergedIterations * 100.0) / numStepsAsDouble);
1328 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1329 shfc->numForceEvaluations / numStepsAsDouble);
1332 // TODO Deallocate memory in shfc