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/chargegroup.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"
82 #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 */
96 struct gmx_shellfc_t {
97 /* Shell counts, indices, parameters and working data */
98 int nshell_gl; /* The number of shells in the system */
99 t_shell *shell_gl; /* All the shells (for DD only) */
100 int *shell_index_gl; /* Global shell index (for DD only) */
101 gmx_bool bInterCG; /* Are there inter charge-group shells? */
102 int nshell; /* The number of local shells */
103 t_shell *shell; /* The local shells */
104 int shell_nalloc; /* The allocation size of shell */
105 gmx_bool bPredict; /* Predict shell positions */
106 gmx_bool bRequireInit; /* Require initialization of shell positions */
107 int nflexcon; /* The number of flexible constraints */
109 /* Temporary arrays, should be fixed size 2 when fully converted to C++ */
110 PaddedVector<gmx::RVec> *x; /* Array for iterative minimization */
111 PaddedVector<gmx::RVec> *f; /* Array for iterative minimization */
113 /* Flexible constraint working data */
114 rvec *acc_dir; /* Acceleration direction for flexcon */
115 rvec *x_old; /* Old coordinates for flexcon */
116 int flex_nalloc; /* The allocation size of acc_dir and x_old */
117 rvec *adir_xnold; /* Work space for init_adir */
118 rvec *adir_xnew; /* Work space for init_adir */
119 int adir_nalloc; /* 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, int ns, t_shell s[])
129 fprintf(fplog, "SHELL DATA\n");
130 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
131 "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, rvec x[], rvec v[], real dt,
159 const real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
161 int i, m, s1, n1, n2, n3;
162 real dt_1, fudge, tm, m1, m2, m3;
165 GMX_RELEASE_ASSERT(mass || mtop, "Must have masses or a way to look them up");
167 /* We introduce a fudge factor for performance reasons: with this choice
168 * 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");
189 for (i = 0; (i < ns); i++)
200 for (m = 0; (m < DIM); m++)
202 x[s1][m] += ptr[n1][m]*dt_1;
215 /* Not the correct masses with FE, but it is just a prediction... */
216 m1 = mtopGetAtomMass(mtop, n1, &molb);
217 m2 = mtopGetAtomMass(mtop, n2, &molb);
220 for (m = 0; (m < DIM); m++)
222 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
237 /* Not the correct masses with FE, but it is just a prediction... */
238 m1 = mtopGetAtomMass(mtop, n1, &molb);
239 m2 = mtopGetAtomMass(mtop, n2, &molb);
240 m3 = mtopGetAtomMass(mtop, n3, &molb);
242 tm = dt_1/(m1+m2+m3);
243 for (m = 0; (m < DIM); m++)
245 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
249 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
254 /*! \brief Count the different particle types in a system
256 * Routine prints a warning to stderr in case an unknown particle type
258 * \param[in] fplog Print what we have found if not NULL
259 * \param[in] mtop Molecular topology.
260 * \returns Array holding the number of particles of a type
262 std::array<int, eptNR> countPtypes(FILE *fplog,
263 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))
282 nptype[atom->ptype] += nmol;
285 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
286 static_cast<int>(atom->ptype));
291 /* Print the number of each particle type */
293 for (const auto &i : nptype)
297 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
305 gmx_shellfc_t *init_shell_flexcon(FILE *fplog,
306 const gmx_mtop_t *mtop, int nflexcon,
308 bool usingDomainDecomposition)
312 int *shell_index = nullptr, *at2cg;
315 int i, j, type, a_offset, cg, mol, ftype, nra;
317 int aS, aN = 0; /* Shell and nucleus */
318 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
319 #define NBT asize(bondtypes)
320 const gmx_ffparams_t *ffparams;
322 std::array<int, eptNR> n = countPtypes(fplog, mtop);
323 nshell = n[eptShell];
325 if (nshell == 0 && nflexcon == 0)
327 /* We're not doing shells or flexible constraints */
332 shfc->x = new PaddedVector<gmx::RVec>[2] {};
333 shfc->f = new PaddedVector<gmx::RVec>[2] {};
334 shfc->nflexcon = nflexcon;
338 /* Only flexible constraints, no shells.
339 * Note that make_local_shells() does not need to be called.
342 shfc->bPredict = FALSE;
347 if (nstcalcenergy != 1)
349 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);
351 if (usingDomainDecomposition)
353 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
356 /* We have shells: fill the shell data structure */
358 /* Global system sized array, this should be avoided */
359 snew(shell_index, mtop->natoms);
362 for (const AtomProxy atomP : AtomRange(*mtop))
364 const t_atom &local = atomP.atom();
365 int i = atomP.globalAtomNumber();
366 if (local.ptype == eptShell)
368 shell_index[i] = nshell++;
374 /* Initiate the shell structures */
375 for (i = 0; (i < nshell); i++)
382 /* shell[i].bInterCG=FALSE; */
387 ffparams = &mtop->ffparams;
389 /* Now fill the structures */
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];
397 const t_block *cgs = &molt->cgs;
399 snew(at2cg, molt->atoms.nr);
400 for (cg = 0; cg < cgs->nr; cg++)
402 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
408 const t_atom *atom = molt->atoms.atom;
409 for (mol = 0; mol < molb->nmol; mol++)
411 for (j = 0; (j < NBT); j++)
413 const int *ia = molt->ilist[bondtypes[j]].iatoms.data();
414 for (i = 0; (i < molt->ilist[bondtypes[j]].size()); )
417 ftype = ffparams->functype[type];
418 nra = interaction_function[ftype].nratoms;
420 /* Check whether we have a bond with a shell */
423 switch (bondtypes[j])
430 if (atom[ia[1]].ptype == eptShell)
435 else if (atom[ia[2]].ptype == eptShell)
442 aN = ia[4]; /* Dummy */
443 aS = ia[5]; /* Shell */
446 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
453 /* Check whether one of the particles is a shell... */
454 nsi = shell_index[a_offset+aS];
455 if ((nsi < 0) || (nsi >= nshell))
457 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
460 if (shell[nsi].shell == -1)
462 shell[nsi].shell = a_offset + aS;
465 else if (shell[nsi].shell != a_offset+aS)
467 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
470 if (shell[nsi].nucl1 == -1)
472 shell[nsi].nucl1 = a_offset + aN;
474 else if (shell[nsi].nucl2 == -1)
476 shell[nsi].nucl2 = a_offset + aN;
478 else if (shell[nsi].nucl3 == -1)
480 shell[nsi].nucl3 = a_offset + aN;
486 pr_shell(fplog, ns, shell);
488 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
490 if (at2cg[aS] != at2cg[aN])
492 /* shell[nsi].bInterCG = TRUE; */
493 shfc->bInterCG = TRUE;
496 switch (bondtypes[j])
500 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
503 shell[nsi].k += ffparams->iparams[type].cubic.kb;
507 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
509 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);
511 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/
512 ffparams->iparams[type].polarize.alpha;
515 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
517 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);
519 alpha = (ffparams->iparams[type].wpol.al_x+
520 ffparams->iparams[type].wpol.al_y+
521 ffparams->iparams[type].wpol.al_z)/3.0;
522 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/alpha;
525 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
533 a_offset += molt->atoms.nr;
535 /* Done with this molecule type */
539 /* Verify whether it's all correct */
542 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
545 for (i = 0; (i < ns); i++)
547 shell[i].k_1 = 1.0/shell[i].k;
552 pr_shell(debug, ns, shell);
556 shfc->nshell_gl = ns;
557 shfc->shell_gl = shell;
558 shfc->shell_index_gl = shell_index;
560 shfc->bPredict = (getenv("GMX_NOPREDICT") == nullptr);
561 shfc->bRequireInit = FALSE;
566 fprintf(fplog, "\nWill never predict shell positions\n");
571 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
572 if (shfc->bRequireInit && fplog)
574 fprintf(fplog, "\nWill always initiate shell positions\n");
584 fprintf(fplog, "\nNOTE: there are shells that are connected to particles outside their own charge group, will not predict shells positions during the run\n\n");
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,
604 int a0, a1, *ind, nshell, i;
605 gmx_domdec_t *dd = nullptr;
607 if (DOMAINDECOMP(cr))
611 a1 = dd_numHomeAtoms(*dd);
615 /* Single node: we need all shells, just copy the pointer */
616 shfc->nshell = shfc->nshell_gl;
617 shfc->shell = shfc->shell_gl;
622 ind = shfc->shell_index_gl;
626 for (i = a0; i < a1; i++)
628 if (md->ptype[i] == eptShell)
630 if (nshell+1 > shfc->shell_nalloc)
632 shfc->shell_nalloc = over_alloc_dd(nshell+1);
633 srenew(shell, shfc->shell_nalloc);
637 shell[nshell] = shfc->shell_gl[ind[dd->globalAtomIndices[i]]];
641 shell[nshell] = shfc->shell_gl[ind[i]];
644 /* With inter-cg shells we can no do shell prediction,
645 * so we do not need the nuclei numbers.
649 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
650 if (shell[nshell].nnucl > 1)
652 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
654 if (shell[nshell].nnucl > 2)
656 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
659 shell[nshell].shell = i;
664 shfc->nshell = nshell;
668 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
686 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
704 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
705 gmx::ArrayRef<gmx::RVec> xnew,
706 const rvec acc_dir[], int homenr, real step)
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,
720 int ns, t_shell s[], int count)
722 const real step_scale_min = 0.8,
723 step_scale_increment = 0.2,
724 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] */
766 step_scale_min * s[i].step[d] +
767 step_scale_increment * std::min(step_scale_multiple * s[i].step[d], std::max(k_est, zero));
772 if (gmx_numzero(dx)) /* 0 == dx */
774 /* Likely this will never happen, but if it does just
775 * don't scale the step. */
779 s[i].step[d] *= step_scale_max;
783 step_min = std::min(step_min, s[i].step[d]);
784 step_max = std::max(step_max, s[i].step[d]);
788 copy_rvec(xcur [shell], s[i].xold);
789 copy_rvec(f[shell], s[i].fold);
791 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
795 fprintf(debug, "shell[%d] = %d\n", i, shell);
796 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
797 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
798 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
799 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
803 printf("step %.3e %.3e\n", step_min, step_max);
807 static void decrease_step_size(int nshell, t_shell s[])
811 for (i = 0; i < nshell; i++)
813 svmul(0.8, s[i].step, s[i].step);
817 static void print_epot(FILE *fp, int64_t mdstep, int count, real epot, real df,
818 int ndir, real sf_dir)
822 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
823 gmx_step_str(mdstep, buf), count, epot, df);
826 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir/ndir));
835 static real rms_force(const t_commrec *cr, gmx::ArrayRef<const gmx::RVec> force, int ns, t_shell s[],
836 int ndir, real *sf_dir, real *Epot)
839 const rvec *f = as_rvec_array(force.data());
842 for (int i = 0; i < ns; i++)
844 int shell = s[i].shell;
845 buf[0] += norm2(f[shell]);
854 gmx_sumd(4, buf, cr);
855 ntot = gmx::roundToInt(buf[1]);
861 return (ntot ? std::sqrt(buf[0]/ntot) : 0);
864 static void dump_shells(FILE *fp, gmx::ArrayRef<gmx::RVec> f, real ftol, int ns, t_shell s[])
869 ft2 = gmx::square(ftol);
871 for (i = 0; (i < ns); i++)
874 ff2 = iprod(f[shell], f[shell]);
877 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
878 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
883 static void init_adir(gmx_shellfc_t *shfc,
884 gmx::Constraints *constr,
885 const t_inputrec *ir,
897 gmx::ArrayRef<const real> lambda,
903 unsigned short *ptype;
905 if (DOMAINDECOMP(cr))
913 if (n > shfc->adir_nalloc)
915 shfc->adir_nalloc = over_alloc_dd(n);
916 srenew(shfc->adir_xnold, shfc->adir_nalloc);
917 srenew(shfc->adir_xnew, shfc->adir_nalloc);
919 xnold = shfc->adir_xnold;
920 xnew = shfc->adir_xnew;
926 /* Does NOT work with freeze or acceleration groups (yet) */
927 for (n = 0; n < end; n++)
929 w_dt = md->invmass[n]*dt;
931 for (d = 0; d < DIM; d++)
933 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
935 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
936 xnew[n][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
940 xnold[n][d] = x[n][d];
941 xnew[n][d] = x[n][d];
945 constr->apply(FALSE, FALSE, step, 0, 1.0,
946 x, xnold, nullptr, box,
947 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
948 nullptr, nullptr, gmx::ConstraintVariable::Positions);
949 constr->apply(FALSE, FALSE, step, 0, 1.0,
950 x, xnew, nullptr, box,
951 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
952 nullptr, nullptr, gmx::ConstraintVariable::Positions);
954 for (n = 0; n < end; n++)
956 for (d = 0; d < DIM; d++)
959 -(2*x[n][d]-xnold[n][d]-xnew[n][d])/gmx::square(dt)
960 - f[n][d]*md->invmass[n];
962 clear_rvec(acc_dir[n]);
965 /* Project the acceleration on the old bond directions */
966 constr->apply(FALSE, FALSE, step, 0, 1.0,
967 x_old, xnew, acc_dir, box,
968 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
969 nullptr, nullptr, gmx::ConstraintVariable::Deriv_FlexCon);
972 void relax_shell_flexcon(FILE *fplog,
974 const gmx_multisim_t *ms,
976 gmx_enfrot *enforcedRotation,
978 const t_inputrec *inputrec,
979 gmx::ImdSession *imdSession,
983 const gmx_localtop_t *top,
984 gmx::Constraints *constr,
985 gmx_enerdata_t *enerd,
988 gmx::ArrayRefWithPadding<gmx::RVec> x,
989 gmx::ArrayRefWithPadding<gmx::RVec> v,
991 gmx::ArrayRef<real> lambda,
993 gmx::ArrayRefWithPadding<gmx::RVec> f,
997 gmx_wallcycle_t wcycle,
1001 gmx::MdScheduleWorkload *mdScheduleWork,
1004 const gmx_vsite_t *vsite,
1005 const DDBalanceRegionHandler &ddBalanceRegionHandler)
1007 auto xRvec = as_rvec_array(x.paddedArrayRef().data());
1008 auto vRvec = as_rvec_array(v.paddedArrayRef().data());
1013 rvec *acc_dir = nullptr, *x_old = nullptr;
1014 real Epot[2], df[2];
1018 gmx_bool bCont, bInit, bConverged;
1019 int nat, dd_ac0, dd_ac1 = 0, i;
1020 int homenr = md->homenr, end = homenr;
1021 int nflexcon, number_steps, d, Min = 0, count = 0;
1022 #define Try (1-Min) /* At start Try = 1 */
1024 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1025 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
1026 ftol = inputrec->em_tol;
1027 number_steps = inputrec->niter;
1028 nshell = shfc->nshell;
1029 shell = shfc->shell;
1030 nflexcon = shfc->nflexcon;
1034 if (DOMAINDECOMP(cr))
1036 nat = dd_natoms_vsite(cr->dd);
1039 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1040 nat = std::max(nat, dd_ac1);
1048 for (i = 0; (i < 2); i++)
1050 shfc->x[i].resizeWithPadding(nat);
1051 shfc->f[i].resizeWithPadding(nat);
1054 /* Create views that we can swap */
1055 gmx::ArrayRefWithPadding<gmx::RVec> posWithPadding[2];
1056 gmx::ArrayRefWithPadding<gmx::RVec> forceWithPadding[2];
1057 gmx::ArrayRef<gmx::RVec> pos[2];
1058 gmx::ArrayRef<gmx::RVec> force[2];
1059 for (i = 0; (i < 2); i++)
1061 posWithPadding[i] = shfc->x[i].arrayRefWithPadding();
1062 pos[i] = posWithPadding[i].paddedArrayRef();
1063 forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1064 force[i] = forceWithPadding[i].paddedArrayRef();
1067 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1069 /* This is the only time where the coordinates are used
1070 * before do_force is called, which normally puts all
1071 * charge groups in the box.
1073 auto xRef = x.paddedArrayRef();
1074 put_atoms_in_box_omp(fr->ePBC, box, xRef.subArray(0, md->homenr), gmx_omp_nthreads_get(emntDefault));
1078 mk_mshift(fplog, graph, fr->ePBC, box, xRvec);
1082 /* After this all coordinate arrays will contain whole charge groups */
1085 shift_self(graph, box, xRvec);
1090 if (nat > shfc->flex_nalloc)
1092 shfc->flex_nalloc = over_alloc_dd(nat);
1093 srenew(shfc->acc_dir, shfc->flex_nalloc);
1094 srenew(shfc->x_old, shfc->flex_nalloc);
1096 acc_dir = shfc->acc_dir;
1097 x_old = shfc->x_old;
1098 auto xArrayRef = x.paddedArrayRef();
1099 auto vArrayRef = v.paddedArrayRef();
1100 for (i = 0; i < homenr; i++)
1102 for (d = 0; d < DIM; d++)
1105 xArrayRef[i][d] - vArrayRef[i][d]*inputrec->delta_t;
1110 /* Do a prediction of the shell positions, when appropriate.
1111 * Without velocities (EM, NM, BD) we only do initial prediction.
1113 if (shfc->bPredict && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1115 predict_shells(fplog, xRvec, vRvec, inputrec->delta_t, nshell, shell,
1116 md->massT, nullptr, bInit);
1119 /* do_force expected the charge groups to be in the box */
1122 unshift_self(graph, box, xRvec);
1125 /* Calculate the forces first time around */
1128 pr_rvecs(debug, 0, "x b4 do_force", xRvec, homenr);
1130 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1131 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession,
1133 mdstep, nrnb, wcycle, top,
1135 forceWithPadding[Min], force_vir, md, enerd, fcd,
1137 fr, mdScheduleWork, vsite, mu_tot, t, nullptr,
1138 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags,
1139 ddBalanceRegionHandler);
1145 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1146 shfc->x_old, xRvec, xRvec, as_rvec_array(force[Min].data()),
1150 for (i = 0; i < end; i++)
1152 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i]);
1155 sum_epot(&(enerd->grpp), enerd->term);
1156 Epot[Min] = enerd->term[F_EPOT];
1158 df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1162 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1167 pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md->nr);
1170 if (nshell+nflexcon > 0)
1172 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1173 * shell positions are updated, therefore the other particles must
1174 * be set here, in advance.
1176 std::copy(x.paddedArrayRef().begin(),
1177 x.paddedArrayRef().end(),
1178 posWithPadding[Min].paddedArrayRef().begin());
1179 std::copy(x.paddedArrayRef().begin(),
1180 x.paddedArrayRef().end(),
1181 posWithPadding[Try].paddedArrayRef().begin());
1184 if (bVerbose && MASTER(cr))
1186 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1191 fprintf(debug, "%17s: %14.10e\n",
1192 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1193 fprintf(debug, "%17s: %14.10e\n",
1194 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1195 fprintf(debug, "%17s: %14.10e\n",
1196 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1197 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1200 /* First check whether we should do shells, or whether the force is
1201 * low enough even without minimization.
1203 bConverged = (df[Min] < ftol);
1205 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1209 construct_vsites(vsite, as_rvec_array(pos[Min].data()),
1210 inputrec->delta_t, vRvec,
1211 idef->iparams, idef->il,
1212 fr->ePBC, fr->bMolPBC, cr, box);
1218 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1220 as_rvec_array(pos[Min].data()),
1221 as_rvec_array(force[Min].data()), acc_dir,
1224 directional_sd(pos[Min], pos[Try], acc_dir, end, fr->fc_stepsize);
1227 /* New positions, Steepest descent */
1228 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1230 /* do_force expected the charge groups to be in the box */
1233 unshift_self(graph, box, as_rvec_array(pos[Try].data()));
1238 pr_rvecs(debug, 0, "RELAX: pos[Min] ", as_rvec_array(pos[Min].data()), homenr);
1239 pr_rvecs(debug, 0, "RELAX: pos[Try] ", as_rvec_array(pos[Try].data()), homenr);
1241 /* Try the new positions */
1242 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession,
1245 top, box, posWithPadding[Try], hist,
1246 forceWithPadding[Try], force_vir,
1247 md, enerd, fcd, lambda, graph,
1248 fr, mdScheduleWork, vsite, mu_tot, t, nullptr,
1250 ddBalanceRegionHandler);
1251 sum_epot(&(enerd->grpp), enerd->term);
1254 pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1255 pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1261 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1263 as_rvec_array(pos[Try].data()),
1264 as_rvec_array(force[Try].data()),
1265 acc_dir, box, lambda, &dum);
1267 for (i = 0; i < end; i++)
1269 sf_dir += md->massT[i]*norm2(acc_dir[i]);
1273 Epot[Try] = enerd->term[F_EPOT];
1275 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1279 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1286 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1290 fprintf(debug, "SHELL ITER %d\n", count);
1291 dump_shells(debug, force[Try], ftol, nshell, shell);
1295 if (bVerbose && MASTER(cr))
1297 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1300 bConverged = (df[Try] < ftol);
1302 if ((df[Try] < df[Min]))
1306 fprintf(debug, "Swapping Min and Try\n");
1310 /* Correct the velocities for the flexible constraints */
1311 invdt = 1/inputrec->delta_t;
1312 auto vArrayRef = v.paddedArrayRef();
1313 for (i = 0; i < end; i++)
1315 for (d = 0; d < DIM; d++)
1317 vArrayRef[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1325 decrease_step_size(nshell, shell);
1328 shfc->numForceEvaluations += count;
1331 shfc->numConvergedIterations++;
1333 if (MASTER(cr) && !(bConverged))
1335 /* Note that the energies and virial are incorrect when not converged */
1339 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1340 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1343 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1344 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1347 /* Copy back the coordinates and the forces */
1348 std::copy(pos[Min].begin(), pos[Min].end(), x.paddedArrayRef().data());
1349 std::copy(force[Min].begin(), force[Min].end(), f.unpaddedArrayRef().begin());
1352 void done_shellfc(FILE *fplog, gmx_shellfc_t *shfc, int64_t numSteps)
1354 if (shfc && fplog && numSteps > 0)
1356 double numStepsAsDouble = static_cast<double>(numSteps);
1357 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1358 (shfc->numConvergedIterations*100.0)/numStepsAsDouble);
1359 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1360 shfc->numForceEvaluations/numStepsAsDouble);
1363 // TODO Deallocate memory in shfc