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,2021, 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/utilities.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/math/vecdump.h"
60 #include "gromacs/mdlib/constr.h"
61 #include "gromacs/mdlib/enerdata_utils.h"
62 #include "gromacs/mdlib/force.h"
63 #include "gromacs/mdlib/force_flags.h"
64 #include "gromacs/mdlib/gmx_omp_nthreads.h"
65 #include "gromacs/mdlib/vsite.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/enerdata.h"
68 #include "gromacs/mdtypes/forcebuffers.h"
69 #include "gromacs/mdtypes/forcerec.h"
70 #include "gromacs/mdtypes/inputrec.h"
71 #include "gromacs/mdtypes/md_enums.h"
72 #include "gromacs/mdtypes/mdatom.h"
73 #include "gromacs/mdtypes/state.h"
74 #include "gromacs/pbcutil/pbc.h"
75 #include "gromacs/topology/ifunc.h"
76 #include "gromacs/topology/mtop_lookup.h"
77 #include "gromacs/topology/mtop_util.h"
78 #include "gromacs/utility/arrayref.h"
79 #include "gromacs/utility/arraysize.h"
80 #include "gromacs/utility/cstringutil.h"
81 #include "gromacs/utility/fatalerror.h"
82 #include "gromacs/utility/gmxassert.h"
85 using gmx::ArrayRefWithPadding;
90 int nnucl = 0; /* The number of nuclei */
91 int shellIndex = -1; /* The shell index */
92 int nucl1 = -1; /* The first nuclei connected to the shell */
93 int nucl2 = -1; /* The second nuclei connected to the shell */
94 int nucl3 = -1; /* The third nuclei connected to the shell */
95 real k = 0; /* force constant */
96 real k_1 = 0; /* 1 over force constant */
97 rvec xold; /* The old shell coordinates */
98 rvec fold; /* The old force on the shell */
99 rvec step; /* Step size for steepest descents */
104 /* Shell counts, indices, parameters and working data */
105 std::vector<t_shell> shell_gl; /* All the shells (for DD only) */
106 std::vector<int> shell_index_gl; /* Global shell index (for DD only) */
107 gmx_bool bInterCG; /* Are there inter charge-group shells? */
108 std::vector<t_shell> shells; /* The local shells */
109 bool predictShells = false; /* Predict shell positions */
110 bool requireInit = false; /* Require initialization of shell positions */
111 int nflexcon = 0; /* The number of flexible constraints */
113 std::array<PaddedHostVector<RVec>, 2> x; /* Coordinate buffers for iterative minimization */
114 std::array<PaddedHostVector<RVec>, 2> f; /* Force buffers for iterative minimization */
116 /* Flexible constraint working data */
117 std::vector<RVec> acc_dir; /* Acceleration direction for flexcon */
118 gmx::PaddedVector<RVec> x_old; /* Old coordinates for flexcon */
119 gmx::PaddedVector<RVec> adir_xnold; /* Work space for init_adir */
120 gmx::PaddedVector<RVec> adir_xnew; /* 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, ArrayRef<const t_shell> shells)
128 fprintf(fplog, "SHELL DATA\n");
129 fprintf(fplog, "%5s %8s %5s %5s %5s\n", "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
130 for (const t_shell& shell : shells)
132 fprintf(fplog, "%5d %8.3f %5d", shell.shellIndex, 1.0 / shell.k_1, shell.nucl1);
133 if (shell.nnucl == 2)
135 fprintf(fplog, " %5d\n", shell.nucl2);
137 else if (shell.nnucl == 3)
139 fprintf(fplog, " %5d %5d\n", shell.nucl2, shell.nucl3);
143 fprintf(fplog, "\n");
148 /* TODO The remain call of this function passes non-NULL mass and NULL
149 * mtop, so this routine can be simplified.
151 * The other code path supported doing prediction before the MD loop
152 * started, but even when called, the prediction was always
153 * over-written by a subsequent call in the MD loop, so has been
155 static void predict_shells(FILE* fplog,
159 ArrayRef<const t_shell> shells,
160 gmx::ArrayRef<const real> mass,
164 real dt_1, fudge, tm, m1, m2, m3;
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
177 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
188 for (const t_shell& shell : shells)
190 const int s1 = shell.shellIndex;
199 for (m = 0; (m < DIM); m++)
201 x[s1][m] += xOrV[n1][m] * dt_1;
209 tm = dt_1 / (m1 + m2);
210 for (m = 0; (m < DIM); m++)
212 x[s1][m] += (m1 * xOrV[n1][m] + m2 * xOrV[n2][m]) * tm;
222 tm = dt_1 / (m1 + m2 + m3);
223 for (m = 0; (m < DIM); m++)
225 x[s1][m] += (m1 * xOrV[n1][m] + m2 * xOrV[n2][m] + m3 * xOrV[n3][m]) * tm;
228 default: gmx_fatal(FARGS, "Shell %d has %d nuclei!", s1, shell.nnucl);
233 gmx_shellfc_t* init_shell_flexcon(FILE* fplog,
234 const gmx_mtop_t& mtop,
237 bool usingDomainDecomposition,
243 int i, j, type, a_offset, mol, ftype, nra;
245 int aS, aN = 0; /* Shell and nucleus */
246 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
247 #define NBT asize(bondtypes)
248 const gmx_ffparams_t* ffparams;
250 const gmx::EnumerationArray<ParticleType, int> numParticles = gmx_mtop_particletype_count(mtop);
253 /* Print the number of each particle type */
254 for (const auto entry : gmx::keysOf(numParticles))
256 const int number = numParticles[entry];
259 fprintf(fplog, "There are: %d %ss\n", number, enumValueToString(entry));
264 nshell = numParticles[ParticleType::Shell];
266 if (nshell == 0 && nflexcon == 0)
268 /* We're not doing shells or flexible constraints */
272 shfc = new gmx_shellfc_t;
273 shfc->nflexcon = nflexcon;
277 /* Only flexible constraints, no shells.
278 * Note that make_local_shells() does not need to be called.
283 if (nstcalcenergy != 1)
286 "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is "
287 "not supported in combination with shell particles.\nPlease make a new tpr file.",
290 if (nshell > 0 && usingDomainDecomposition)
294 "Shell particles are not implemented with domain decomposition, use a single rank");
297 /* We have shells: fill the shell data structure */
299 /* Global system sized array, this should be avoided */
300 std::vector<int> shell_index(mtop.natoms);
303 for (const AtomProxy atomP : AtomRange(mtop))
305 const t_atom& local = atomP.atom();
306 int i = atomP.globalAtomNumber();
307 if (local.ptype == ParticleType::Shell)
309 shell_index[i] = nshell++;
313 std::vector<t_shell> shell(nshell);
315 ffparams = &mtop.ffparams;
317 /* Now fill the structures */
318 /* TODO: See if we can use update groups that cover shell constructions */
319 shfc->bInterCG = FALSE;
322 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
324 const gmx_molblock_t& molb = mtop.molblock[mb];
325 const gmx_moltype_t& molt = mtop.moltype[molb.type];
327 const t_atom* atom = molt.atoms.atom;
328 for (mol = 0; mol < molb.nmol; mol++)
330 for (j = 0; (j < NBT); j++)
332 const int* ia = molt.ilist[bondtypes[j]].iatoms.data();
333 for (i = 0; (i < molt.ilist[bondtypes[j]].size());)
336 ftype = ffparams->functype[type];
337 nra = interaction_function[ftype].nratoms;
339 /* Check whether we have a bond with a shell */
342 switch (bondtypes[j])
349 if (atom[ia[1]].ptype == ParticleType::Shell)
354 else if (atom[ia[2]].ptype == ParticleType::Shell)
361 aN = ia[4]; /* Dummy */
362 aS = ia[5]; /* Shell */
364 default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
371 /* Check whether one of the particles is a shell... */
372 nsi = shell_index[a_offset + aS];
373 if ((nsi < 0) || (nsi >= nshell))
375 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d", nsi, nshell, aS);
377 if (shell[nsi].shellIndex == -1)
379 shell[nsi].shellIndex = a_offset + aS;
382 else if (shell[nsi].shellIndex != a_offset + aS)
384 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
387 if (shell[nsi].nucl1 == -1)
389 shell[nsi].nucl1 = a_offset + aN;
391 else if (shell[nsi].nucl2 == -1)
393 shell[nsi].nucl2 = a_offset + aN;
395 else if (shell[nsi].nucl3 == -1)
397 shell[nsi].nucl3 = a_offset + aN;
403 pr_shell(fplog, shell);
405 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
409 /* shell[nsi].bInterCG = TRUE; */
410 shfc->bInterCG = TRUE;
413 switch (bondtypes[j])
417 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
420 shell[nsi].k += ffparams->iparams[type].cubic.kb;
424 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS * 10))
427 "polarize can not be used with qA(%e) != qB(%e) for "
428 "atom %d of molecule block %zu",
434 shell[nsi].k += gmx::square(qS) * gmx::c_one4PiEps0
435 / ffparams->iparams[type].polarize.alpha;
438 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS * 10))
441 "water_pol can not be used with qA(%e) != qB(%e) for "
442 "atom %d of molecule block %zu",
448 alpha = (ffparams->iparams[type].wpol.al_x
449 + ffparams->iparams[type].wpol.al_y
450 + ffparams->iparams[type].wpol.al_z)
452 shell[nsi].k += gmx::square(qS) * gmx::c_one4PiEps0 / alpha;
454 default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
462 a_offset += molt.atoms.nr;
464 /* Done with this molecule type */
467 /* Verify whether it's all correct */
470 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
473 for (i = 0; (i < ns); i++)
475 shell[i].k_1 = 1.0 / shell[i].k;
480 pr_shell(debug, shell);
484 shfc->shell_gl = shell;
485 shfc->shell_index_gl = shell_index;
487 shfc->predictShells = (getenv("GMX_NOPREDICT") == nullptr);
488 shfc->requireInit = false;
489 if (!shfc->predictShells)
493 fprintf(fplog, "\nWill never predict shell positions\n");
498 shfc->requireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
499 if (shfc->requireInit && fplog)
501 fprintf(fplog, "\nWill always initiate shell positions\n");
505 if (shfc->predictShells)
512 "\nNOTE: in the current version shell prediction during the crun is "
515 /* Prediction improves performance, so we should implement either:
516 * 1. communication for the atoms needed for prediction
517 * 2. prediction using the velocities of shells; currently the
518 * shell velocities are zeroed, it's a bit tricky to keep
519 * track of the shell displacements and thus the velocity.
521 shfc->predictShells = false;
525 /* shfc->x is used as a coordinate buffer for the sim_util's `do_force` function, and
526 * when using PME it must be pinned. */
529 for (i = 0; i < 2; i++)
531 changePinningPolicy(&shfc->x[i], gmx::PinningPolicy::PinnedIfSupported);
538 void gmx::make_local_shells(const t_commrec* cr, const t_mdatoms& md, gmx_shellfc_t* shfc)
541 gmx_domdec_t* dd = nullptr;
543 if (haveDDAtomOrdering(*cr))
547 a1 = dd_numHomeAtoms(*dd);
551 /* Single node: we need all shells, copy them */
552 shfc->shells = shfc->shell_gl;
557 ArrayRef<const int> ind = shfc->shell_index_gl;
559 std::vector<t_shell>& shells = shfc->shells;
561 auto* ptype = md.ptype;
562 for (int i = a0; i < a1; i++)
564 if (ptype[i] == ParticleType::Shell)
568 shells.push_back(shfc->shell_gl[ind[dd->globalAtomIndices[i]]]);
572 shells.push_back(shfc->shell_gl[ind[i]]);
574 t_shell& shell = shells.back();
576 /* With inter-cg shells we can no do shell prediction,
577 * so we do not need the nuclei numbers.
581 shell.nucl1 = i + shell.nucl1 - shell.shellIndex;
584 shell.nucl2 = i + shell.nucl2 - shell.shellIndex;
588 shell.nucl3 = i + shell.nucl3 - shell.shellIndex;
591 shell.shellIndex = i;
596 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
614 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
623 dx = f[XX] * step[XX];
624 dy = f[YY] * step[YY];
625 dz = f[ZZ] * step[ZZ];
632 static void directional_sd(ArrayRef<const RVec> xold,
634 ArrayRef<const RVec> acc_dir,
638 const rvec* xo = as_rvec_array(xold.data());
639 rvec* xn = as_rvec_array(xnew.data());
641 for (int i = 0; i < homenr; i++)
643 do_1pos(xn[i], xo[i], acc_dir[i], step);
647 static void shell_pos_sd(ArrayRef<const RVec> xcur,
649 ArrayRef<const RVec> f,
650 ArrayRef<t_shell> shells,
653 const real step_scale_min = 0.8, step_scale_increment = 0.2, step_scale_max = 1.2,
654 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
659 real step_min, step_max;
664 for (t_shell& shell : shells)
666 const int ind = shell.shellIndex;
669 for (d = 0; d < DIM; d++)
671 shell.step[d] = shell.k_1;
673 step_min = std::min(step_min, shell.step[d]);
674 step_max = std::max(step_max, shell.step[d]);
680 for (d = 0; d < DIM; d++)
682 dx = xcur[ind][d] - shell.xold[d];
683 df = f[ind][d] - shell.fold[d];
684 /* -dx/df gets used to generate an interpolated value, but would
685 * cause a NaN if df were binary-equal to zero. Values close to
686 * zero won't cause problems (because of the min() and max()), so
687 * just testing for binary inequality is OK. */
691 /* Scale the step size by a factor interpolated from
692 * step_scale_min to step_scale_max, as k_est goes from 0 to
693 * step_scale_multiple * shell.step[d] */
694 shell.step[d] = step_scale_min * shell.step[d]
695 + step_scale_increment
696 * std::min(step_scale_multiple * shell.step[d],
697 std::max(k_est, zero));
702 if (gmx_numzero(dx)) /* 0 == dx */
704 /* Likely this will never happen, but if it does just
705 * don't scale the step. */
709 shell.step[d] *= step_scale_max;
713 step_min = std::min(step_min, shell.step[d]);
714 step_max = std::max(step_max, shell.step[d]);
718 copy_rvec(xcur[ind], shell.xold);
719 copy_rvec(f[ind], shell.fold);
721 do_1pos3(xnew[ind], xcur[ind], f[ind], shell.step);
725 fprintf(debug, "shell = %d\n", ind);
726 pr_rvec(debug, 0, "fshell", f[ind], DIM, TRUE);
727 pr_rvec(debug, 0, "xold", xcur[ind], DIM, TRUE);
728 pr_rvec(debug, 0, "step", shell.step, DIM, TRUE);
729 pr_rvec(debug, 0, "xnew", xnew[ind], DIM, TRUE);
733 printf("step %.3e %.3e\n", step_min, step_max);
737 static void decrease_step_size(ArrayRef<t_shell> shells)
739 for (t_shell& shell : shells)
741 svmul(0.8, shell.step, shell.step);
745 static void print_epot(FILE* fp, int64_t mdstep, int count, real epot, real df, int ndir, real sf_dir)
749 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e", gmx_step_str(mdstep, buf), count, epot, df);
752 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir / ndir));
761 static real rms_force(const t_commrec* cr,
762 ArrayRef<const RVec> force,
763 ArrayRef<const t_shell> shells,
769 const rvec* f = as_rvec_array(force.data());
772 for (const t_shell& shell : shells)
774 buf[0] += norm2(f[shell.shellIndex]);
776 int ntot = shells.ssize();
783 gmx_sumd(4, buf, cr);
784 ntot = gmx::roundToInt(buf[1]);
790 return (ntot ? std::sqrt(buf[0] / ntot) : 0);
793 static void dump_shells(FILE* fp, ArrayRef<RVec> f, real ftol, ArrayRef<const t_shell> shells)
797 ft2 = gmx::square(ftol);
799 for (const t_shell& shell : shells)
801 const int ind = shell.shellIndex;
802 ff2 = iprod(f[ind], f[ind]);
806 "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
816 static void init_adir(gmx_shellfc_t* shfc,
817 gmx::Constraints* constr,
818 const t_inputrec* ir,
824 ArrayRefWithPadding<RVec> xOld,
825 ArrayRef<RVec> x_init,
826 ArrayRefWithPadding<RVec> xCurrent,
828 ArrayRef<RVec> acc_dir,
830 ArrayRef<const real> lambda,
836 if (haveDDAtomOrdering(*cr))
844 shfc->adir_xnold.resizeWithPadding(n);
845 shfc->adir_xnew.resizeWithPadding(n);
846 rvec* xnold = as_rvec_array(shfc->adir_xnold.data());
847 rvec* xnew = as_rvec_array(shfc->adir_xnew.data());
848 rvec* x_old = as_rvec_array(xOld.paddedArrayRef().data());
849 rvec* x = as_rvec_array(xCurrent.paddedArrayRef().data());
851 auto* ptype = md.ptype;
852 auto invmass = gmx::arrayRefFromArray(md.invmass, md.nr);
855 /* Does NOT work with freeze groups (yet) */
856 for (n = 0; n < end; n++)
858 w_dt = invmass[n] * dt;
860 for (d = 0; d < DIM; d++)
862 if ((ptype[n] != ParticleType::VSite) && (ptype[n] != ParticleType::Shell))
864 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
865 xnew[n][d] = 2 * x[n][d] - x_old[n][d] + f[n][d] * w_dt * dt;
869 xnold[n][d] = x[n][d];
870 xnew[n][d] = x[n][d];
874 bool needsLogging = false;
875 bool computeEnergy = false;
876 bool computeVirial = false;
877 constr->apply(needsLogging,
883 shfc->adir_xnold.arrayRefWithPadding(),
886 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Bonded)],
887 &(dvdlambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Bonded)]),
891 gmx::ConstraintVariable::Positions);
892 constr->apply(needsLogging,
898 shfc->adir_xnew.arrayRefWithPadding(),
901 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Bonded)],
902 &(dvdlambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Bonded)]),
906 gmx::ConstraintVariable::Positions);
908 for (n = 0; n < end; n++)
910 for (d = 0; d < DIM; d++)
912 xnew[n][d] = -(2 * x[n][d] - xnold[n][d] - xnew[n][d]) / gmx::square(dt)
913 - f[n][d] * invmass[n];
915 clear_rvec(acc_dir[n]);
918 /* Project the acceleration on the old bond directions */
919 constr->apply(needsLogging,
925 shfc->adir_xnew.arrayRefWithPadding(),
928 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Bonded)],
929 &(dvdlambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Bonded)]),
933 gmx::ConstraintVariable::Deriv_FlexCon);
936 void relax_shell_flexcon(FILE* fplog,
938 const gmx_multisim_t* ms,
940 gmx_enfrot* enforcedRotation,
942 const t_inputrec* inputrec,
943 gmx::ImdSession* imdSession,
947 const gmx_localtop_t* top,
948 gmx::Constraints* constr,
949 gmx_enerdata_t* enerd,
951 ArrayRefWithPadding<RVec> xPadded,
952 ArrayRefWithPadding<RVec> vPadded,
954 ArrayRef<real> lambda,
955 const history_t* hist,
956 gmx::ForceBuffersView* f,
959 CpuPpLongRangeNonbondeds* longRangeNonbondeds,
961 gmx_wallcycle* wcycle,
964 gmx::MdrunScheduleWorkload* runScheduleWork,
967 gmx::VirtualSitesHandler* vsite,
968 const DDBalanceRegionHandler& ddBalanceRegionHandler)
974 int nat, dd_ac0, dd_ac1 = 0, i;
975 int homenr = md.homenr, end = homenr;
976 int d, Min = 0, count = 0;
977 #define Try (1 - Min) /* At start Try = 1 */
979 const bool bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
980 const bool bInit = (mdstep == inputrec->init_step) || shfc->requireInit;
981 const real ftol = inputrec->em_tol;
982 const int number_steps = inputrec->niter;
983 ArrayRef<t_shell> shells = shfc->shells;
984 const int nflexcon = shfc->nflexcon;
986 if (haveDDAtomOrdering(*cr))
988 nat = dd_natoms_vsite(*cr->dd);
991 dd_get_constraint_range(*cr->dd, &dd_ac0, &dd_ac1);
992 nat = std::max(nat, dd_ac1);
1000 for (i = 0; (i < 2); i++)
1002 shfc->x[i].resizeWithPadding(nat);
1003 shfc->f[i].resizeWithPadding(nat);
1006 /* Create views that we can swap for trail and minimum for positions and forces */
1007 ArrayRefWithPadding<RVec> posWithPadding[2];
1008 ArrayRefWithPadding<RVec> forceWithPadding[2];
1009 ArrayRef<RVec> pos[2];
1010 ArrayRef<RVec> force[2];
1011 for (i = 0; (i < 2); i++)
1013 posWithPadding[i] = shfc->x[i].arrayRefWithPadding();
1014 pos[i] = posWithPadding[i].paddedArrayRef();
1015 forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1016 force[i] = forceWithPadding[i].paddedArrayRef();
1019 ArrayRef<RVec> x = xPadded.unpaddedArrayRef();
1020 ArrayRef<RVec> v = vPadded.unpaddedArrayRef();
1022 if (bDoNS && inputrec->pbcType != PbcType::No && !haveDDAtomOrdering(*cr))
1024 /* This is the only time where the coordinates are used
1025 * before do_force is called, which normally puts all
1026 * charge groups in the box.
1028 put_atoms_in_box_omp(
1029 fr->pbcType, box, x.subArray(0, md.homenr), gmx_omp_nthreads_get(ModuleMultiThread::Default));
1034 shfc->acc_dir.resize(nat);
1035 shfc->x_old.resizeWithPadding(nat);
1036 ArrayRef<RVec> x_old = shfc->x_old.arrayRefWithPadding().unpaddedArrayRef();
1037 for (i = 0; i < homenr; i++)
1039 for (d = 0; d < DIM; d++)
1041 x_old[i][d] = x[i][d] - v[i][d] * inputrec->delta_t;
1046 auto massT = gmx::arrayRefFromArray(md.massT, md.nr);
1047 /* Do a prediction of the shell positions, when appropriate.
1048 * Without velocities (EM, NM, BD) we only do initial prediction.
1050 if (shfc->predictShells && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1052 predict_shells(fplog, x, v, inputrec->delta_t, shells, massT, bInit);
1055 /* Calculate the forces first time around */
1058 pr_rvecs(debug, 0, "x b4 do_force", as_rvec_array(x.data()), homenr);
1060 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1061 gmx::ForceBuffersView forceViewInit = gmx::ForceBuffersView(forceWithPadding[Min], {}, false);
1088 longRangeNonbondeds,
1089 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags,
1090 ddBalanceRegionHandler);
1103 shfc->x_old.arrayRefWithPadding(),
1112 for (i = 0; i < end; i++)
1114 sf_dir += massT[i] * norm2(shfc->acc_dir[i]);
1117 accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals.get());
1118 Epot[Min] = enerd->term[F_EPOT];
1120 df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), shells, nflexcon, &sf_dir, &Epot[Min]);
1124 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1129 pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md.nr);
1132 if (!shells.empty() || nflexcon > 0)
1134 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1135 * shell positions are updated, therefore the other particles must
1136 * be set here, in advance.
1138 std::copy(xPadded.paddedArrayRef().begin(),
1139 xPadded.paddedArrayRef().end(),
1140 posWithPadding[Min].paddedArrayRef().begin());
1141 std::copy(xPadded.paddedArrayRef().begin(),
1142 xPadded.paddedArrayRef().end(),
1143 posWithPadding[Try].paddedArrayRef().begin());
1146 if (bVerbose && MASTER(cr))
1148 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1153 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1154 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1155 fprintf(debug, "%17s: %14.10e\n", interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1156 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1159 /* First check whether we should do shells, or whether the force is
1160 * low enough even without minimization.
1162 bool bConverged = (df[Min] < ftol);
1164 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1168 vsite->construct(pos[Min], v, box, gmx::VSiteOperation::PositionsAndVelocities);
1181 shfc->x_old.arrayRefWithPadding(),
1183 posWithPadding[Min],
1190 directional_sd(pos[Min], pos[Try], shfc->acc_dir, end, fr->fc_stepsize);
1193 /* New positions, Steepest descent */
1194 shell_pos_sd(pos[Min], pos[Try], force[Min], shells, count);
1198 pr_rvecs(debug, 0, "RELAX: pos[Min] ", as_rvec_array(pos[Min].data()), homenr);
1199 pr_rvecs(debug, 0, "RELAX: pos[Try] ", as_rvec_array(pos[Try].data()), homenr);
1201 /* Try the new positions */
1202 gmx::ForceBuffersView forceViewTry = gmx::ForceBuffersView(forceWithPadding[Try], {}, false);
1216 posWithPadding[Try],
1229 longRangeNonbondeds,
1231 ddBalanceRegionHandler);
1232 accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals.get());
1235 pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1236 pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1249 shfc->x_old.arrayRefWithPadding(),
1251 posWithPadding[Try],
1258 ArrayRef<const RVec> acc_dir = shfc->acc_dir;
1259 for (i = 0; i < end; i++)
1261 sf_dir += massT[i] * norm2(acc_dir[i]);
1265 Epot[Try] = enerd->term[F_EPOT];
1267 df[Try] = rms_force(cr, force[Try], shells, nflexcon, &sf_dir, &Epot[Try]);
1271 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1278 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1282 fprintf(debug, "SHELL ITER %d\n", count);
1283 dump_shells(debug, force[Try], ftol, shells);
1287 if (bVerbose && MASTER(cr))
1289 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1292 bConverged = (df[Try] < ftol);
1294 if ((df[Try] < df[Min]))
1298 fprintf(debug, "Swapping Min and Try\n");
1302 /* Correct the velocities for the flexible constraints */
1303 invdt = 1 / inputrec->delta_t;
1304 for (i = 0; i < end; i++)
1306 for (d = 0; d < DIM; d++)
1308 v[i][d] += (pos[Try][i][d] - pos[Min][i][d]) * invdt;
1316 decrease_step_size(shells);
1319 shfc->numForceEvaluations += count;
1322 shfc->numConvergedIterations++;
1324 if (MASTER(cr) && !(bConverged))
1326 /* Note that the energies and virial are incorrect when not converged */
1330 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1331 gmx_step_str(mdstep, sbuf),
1336 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1337 gmx_step_str(mdstep, sbuf),
1342 /* Copy back the coordinates and the forces */
1343 std::copy(pos[Min].begin(), pos[Min].end(), x.data());
1344 std::copy(force[Min].begin(), force[Min].end(), f->force().begin());
1347 void done_shellfc(FILE* fplog, gmx_shellfc_t* shfc, int64_t numSteps)
1349 if (shfc && fplog && numSteps > 0)
1351 double numStepsAsDouble = static_cast<double>(numSteps);
1353 "Fraction of iterations that converged: %.2f %%\n",
1354 (shfc->numConvergedIterations * 100.0) / numStepsAsDouble);
1356 "Average number of force evaluations per MD step: %.2f\n\n",
1357 shfc->numForceEvaluations / numStepsAsDouble);