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/smalloc.h"
85 int shell; /* The shell id */
86 int nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
87 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
88 real k; /* force constant */
89 real k_1; /* 1 over force constant */
95 struct gmx_shellfc_t {
96 /* Shell counts, indices, parameters and working data */
97 int nshell_gl; /* The number of shells in the system */
98 t_shell *shell_gl; /* All the shells (for DD only) */
99 int *shell_index_gl; /* Global shell index (for DD only) */
100 gmx_bool bInterCG; /* Are there inter charge-group shells? */
101 int nshell; /* The number of local shells */
102 t_shell *shell; /* The local shells */
103 int shell_nalloc; /* The allocation size of shell */
104 gmx_bool bPredict; /* Predict shell positions */
105 gmx_bool bRequireInit; /* Require initialization of shell positions */
106 int nflexcon; /* The number of flexible constraints */
108 /* Temporary arrays, should be fixed size 2 when fully converted to C++ */
109 PaddedVector<gmx::RVec> *x; /* Array for iterative minimization */
110 PaddedVector<gmx::RVec> *f; /* Array for iterative minimization */
112 /* Flexible constraint working data */
113 rvec *acc_dir; /* Acceleration direction for flexcon */
114 rvec *x_old; /* Old coordinates for flexcon */
115 int flex_nalloc; /* The allocation size of acc_dir and x_old */
116 rvec *adir_xnold; /* Work space for init_adir */
117 rvec *adir_xnew; /* Work space for init_adir */
118 int adir_nalloc; /* Work space for init_adir */
119 std::int64_t numForceEvaluations; /* Total number of force evaluations */
120 int numConvergedIterations; /* Total number of iterations that converged */
124 static void pr_shell(FILE *fplog, int ns, t_shell s[])
128 fprintf(fplog, "SHELL DATA\n");
129 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
130 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
131 for (i = 0; (i < ns); i++)
133 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
136 fprintf(fplog, " %5d\n", s[i].nucl2);
138 else if (s[i].nnucl == 3)
140 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
144 fprintf(fplog, "\n");
149 /* TODO The remain call of this function passes non-NULL mass and NULL
150 * mtop, so this routine can be simplified.
152 * The other code path supported doing prediction before the MD loop
153 * started, but even when called, the prediction was always
154 * over-written by a subsequent call in the MD loop, so has been
156 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
158 const real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
160 int i, m, s1, n1, n2, n3;
161 real dt_1, fudge, tm, m1, m2, m3;
164 /* We introduce a fudge factor for performance reasons: with this choice
165 * the initial force on the shells is about a factor of two lower than
174 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
186 for (i = 0; (i < ns); i++)
197 for (m = 0; (m < DIM); m++)
199 x[s1][m] += ptr[n1][m]*dt_1;
212 /* Not the correct masses with FE, but it is just a prediction... */
213 m1 = mtopGetAtomMass(mtop, n1, &molb);
214 m2 = mtopGetAtomMass(mtop, n2, &molb);
217 for (m = 0; (m < DIM); m++)
219 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
234 /* Not the correct masses with FE, but it is just a prediction... */
235 m1 = mtopGetAtomMass(mtop, n1, &molb);
236 m2 = mtopGetAtomMass(mtop, n2, &molb);
237 m3 = mtopGetAtomMass(mtop, n3, &molb);
239 tm = dt_1/(m1+m2+m3);
240 for (m = 0; (m < DIM); m++)
242 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
246 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
251 /*! \brief Count the different particle types in a system
253 * Routine prints a warning to stderr in case an unknown particle type
255 * \param[in] fplog Print what we have found if not NULL
256 * \param[in] mtop Molecular topology.
257 * \returns Array holding the number of particles of a type
259 static std::array<int, eptNR> countPtypes(FILE *fplog,
260 const gmx_mtop_t *mtop)
262 std::array<int, eptNR> nptype = { { 0 } };
263 /* Count number of shells, and find their indices */
264 for (int i = 0; (i < eptNR); i++)
269 gmx_mtop_atomloop_block_t aloopb = gmx_mtop_atomloop_block_init(mtop);
272 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
279 nptype[atom->ptype] += nmol;
282 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
283 static_cast<int>(atom->ptype));
288 /* Print the number of each particle type */
290 for (const auto &i : nptype)
294 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
302 gmx_shellfc_t *init_shell_flexcon(FILE *fplog,
303 const gmx_mtop_t *mtop, int nflexcon,
305 bool usingDomainDecomposition)
309 int *shell_index = nullptr, *at2cg;
312 int i, j, type, a_offset, cg, mol, ftype, nra;
314 int aS, aN = 0; /* Shell and nucleus */
315 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
316 #define NBT asize(bondtypes)
317 const gmx_ffparams_t *ffparams;
319 std::array<int, eptNR> n = countPtypes(fplog, mtop);
320 nshell = n[eptShell];
322 if (nshell == 0 && nflexcon == 0)
324 /* We're not doing shells or flexible constraints */
329 shfc->x = new PaddedVector<gmx::RVec>[2] {};
330 shfc->f = new PaddedVector<gmx::RVec>[2] {};
331 shfc->nflexcon = nflexcon;
335 /* Only flexible constraints, no shells.
336 * Note that make_local_shells() does not need to be called.
339 shfc->bPredict = FALSE;
344 if (nstcalcenergy != 1)
346 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);
348 if (usingDomainDecomposition)
350 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
353 /* We have shells: fill the shell data structure */
355 /* Global system sized array, this should be avoided */
356 snew(shell_index, mtop->natoms);
359 for (const AtomProxy atomP : AtomRange(*mtop))
361 const t_atom &local = atomP.atom();
362 int i = atomP.globalAtomNumber();
363 if (local.ptype == eptShell)
365 shell_index[i] = nshell++;
371 /* Initiate the shell structures */
372 for (i = 0; (i < nshell); i++)
379 /* shell[i].bInterCG=FALSE; */
384 ffparams = &mtop->ffparams;
386 /* Now fill the structures */
387 shfc->bInterCG = FALSE;
390 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
392 const gmx_molblock_t *molb = &mtop->molblock[mb];
393 const gmx_moltype_t *molt = &mtop->moltype[molb->type];
394 const t_block *cgs = &molt->cgs;
396 snew(at2cg, molt->atoms.nr);
397 for (cg = 0; cg < cgs->nr; cg++)
399 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
405 const t_atom *atom = molt->atoms.atom;
406 for (mol = 0; mol < molb->nmol; mol++)
408 for (j = 0; (j < NBT); j++)
410 const int *ia = molt->ilist[bondtypes[j]].iatoms.data();
411 for (i = 0; (i < molt->ilist[bondtypes[j]].size()); )
414 ftype = ffparams->functype[type];
415 nra = interaction_function[ftype].nratoms;
417 /* Check whether we have a bond with a shell */
420 switch (bondtypes[j])
427 if (atom[ia[1]].ptype == eptShell)
432 else if (atom[ia[2]].ptype == eptShell)
439 aN = ia[4]; /* Dummy */
440 aS = ia[5]; /* Shell */
443 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
450 /* Check whether one of the particles is a shell... */
451 nsi = shell_index[a_offset+aS];
452 if ((nsi < 0) || (nsi >= nshell))
454 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
457 if (shell[nsi].shell == -1)
459 shell[nsi].shell = a_offset + aS;
462 else if (shell[nsi].shell != a_offset+aS)
464 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
467 if (shell[nsi].nucl1 == -1)
469 shell[nsi].nucl1 = a_offset + aN;
471 else if (shell[nsi].nucl2 == -1)
473 shell[nsi].nucl2 = a_offset + aN;
475 else if (shell[nsi].nucl3 == -1)
477 shell[nsi].nucl3 = a_offset + aN;
483 pr_shell(fplog, ns, shell);
485 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
487 if (at2cg[aS] != at2cg[aN])
489 /* shell[nsi].bInterCG = TRUE; */
490 shfc->bInterCG = TRUE;
493 switch (bondtypes[j])
497 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
500 shell[nsi].k += ffparams->iparams[type].cubic.kb;
504 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
506 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);
508 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/
509 ffparams->iparams[type].polarize.alpha;
512 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
514 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);
516 alpha = (ffparams->iparams[type].wpol.al_x+
517 ffparams->iparams[type].wpol.al_y+
518 ffparams->iparams[type].wpol.al_z)/3.0;
519 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/alpha;
522 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
530 a_offset += molt->atoms.nr;
532 /* Done with this molecule type */
536 /* Verify whether it's all correct */
539 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
542 for (i = 0; (i < ns); i++)
544 shell[i].k_1 = 1.0/shell[i].k;
549 pr_shell(debug, ns, shell);
553 shfc->nshell_gl = ns;
554 shfc->shell_gl = shell;
555 shfc->shell_index_gl = shell_index;
557 shfc->bPredict = (getenv("GMX_NOPREDICT") == nullptr);
558 shfc->bRequireInit = FALSE;
563 fprintf(fplog, "\nWill never predict shell positions\n");
568 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
569 if (shfc->bRequireInit && fplog)
571 fprintf(fplog, "\nWill always initiate shell positions\n");
581 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");
583 /* Prediction improves performance, so we should implement either:
584 * 1. communication for the atoms needed for prediction
585 * 2. prediction using the velocities of shells; currently the
586 * shell velocities are zeroed, it's a bit tricky to keep
587 * track of the shell displacements and thus the velocity.
589 shfc->bPredict = FALSE;
596 void make_local_shells(const t_commrec *cr,
601 int a0, a1, *ind, nshell, i;
602 gmx_domdec_t *dd = nullptr;
604 if (DOMAINDECOMP(cr))
608 a1 = dd_numHomeAtoms(*dd);
612 /* Single node: we need all shells, just copy the pointer */
613 shfc->nshell = shfc->nshell_gl;
614 shfc->shell = shfc->shell_gl;
619 ind = shfc->shell_index_gl;
623 for (i = a0; i < a1; i++)
625 if (md->ptype[i] == eptShell)
627 if (nshell+1 > shfc->shell_nalloc)
629 shfc->shell_nalloc = over_alloc_dd(nshell+1);
630 srenew(shell, shfc->shell_nalloc);
634 shell[nshell] = shfc->shell_gl[ind[dd->globalAtomIndices[i]]];
638 shell[nshell] = shfc->shell_gl[ind[i]];
641 /* With inter-cg shells we can no do shell prediction,
642 * so we do not need the nuclei numbers.
646 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
647 if (shell[nshell].nnucl > 1)
649 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
651 if (shell[nshell].nnucl > 2)
653 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
656 shell[nshell].shell = i;
661 shfc->nshell = nshell;
665 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
683 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
701 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
702 gmx::ArrayRef<gmx::RVec> xnew,
703 const rvec acc_dir[], int homenr, real step)
705 const rvec *xo = as_rvec_array(xold.data());
706 rvec *xn = as_rvec_array(xnew.data());
708 for (int i = 0; i < homenr; i++)
710 do_1pos(xn[i], xo[i], acc_dir[i], step);
714 static void shell_pos_sd(gmx::ArrayRef<const gmx::RVec> xcur,
715 gmx::ArrayRef<gmx::RVec> xnew,
716 gmx::ArrayRef<const gmx::RVec> f,
717 int ns, t_shell s[], int count)
719 const real step_scale_min = 0.8,
720 step_scale_increment = 0.2,
721 step_scale_max = 1.2,
722 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
727 real step_min, step_max;
732 for (i = 0; (i < ns); i++)
737 for (d = 0; d < DIM; d++)
739 s[i].step[d] = s[i].k_1;
741 step_min = std::min(step_min, s[i].step[d]);
742 step_max = std::max(step_max, s[i].step[d]);
748 for (d = 0; d < DIM; d++)
750 dx = xcur[shell][d] - s[i].xold[d];
751 df = f[shell][d] - s[i].fold[d];
752 /* -dx/df gets used to generate an interpolated value, but would
753 * cause a NaN if df were binary-equal to zero. Values close to
754 * zero won't cause problems (because of the min() and max()), so
755 * just testing for binary inequality is OK. */
759 /* Scale the step size by a factor interpolated from
760 * step_scale_min to step_scale_max, as k_est goes from 0 to
761 * step_scale_multiple * s[i].step[d] */
763 step_scale_min * s[i].step[d] +
764 step_scale_increment * std::min(step_scale_multiple * s[i].step[d], std::max(k_est, zero));
769 if (gmx_numzero(dx)) /* 0 == dx */
771 /* Likely this will never happen, but if it does just
772 * don't scale the step. */
776 s[i].step[d] *= step_scale_max;
780 step_min = std::min(step_min, s[i].step[d]);
781 step_max = std::max(step_max, s[i].step[d]);
785 copy_rvec(xcur [shell], s[i].xold);
786 copy_rvec(f[shell], s[i].fold);
788 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
792 fprintf(debug, "shell[%d] = %d\n", i, shell);
793 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
794 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
795 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
796 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
800 printf("step %.3e %.3e\n", step_min, step_max);
804 static void decrease_step_size(int nshell, t_shell s[])
808 for (i = 0; i < nshell; i++)
810 svmul(0.8, s[i].step, s[i].step);
814 static void print_epot(FILE *fp, int64_t mdstep, int count, real epot, real df,
815 int ndir, real sf_dir)
819 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
820 gmx_step_str(mdstep, buf), count, epot, df);
823 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir/ndir));
832 static real rms_force(const t_commrec *cr, gmx::ArrayRef<const gmx::RVec> force, int ns, t_shell s[],
833 int ndir, real *sf_dir, real *Epot)
836 const rvec *f = as_rvec_array(force.data());
839 for (int i = 0; i < ns; i++)
841 int shell = s[i].shell;
842 buf[0] += norm2(f[shell]);
851 gmx_sumd(4, buf, cr);
852 ntot = gmx::roundToInt(buf[1]);
858 return (ntot ? std::sqrt(buf[0]/ntot) : 0);
861 static void dump_shells(FILE *fp, gmx::ArrayRef<gmx::RVec> f, real ftol, int ns, t_shell s[])
866 ft2 = gmx::square(ftol);
868 for (i = 0; (i < ns); i++)
871 ff2 = iprod(f[shell], f[shell]);
874 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
875 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
880 static void init_adir(gmx_shellfc_t *shfc,
881 gmx::Constraints *constr,
882 const t_inputrec *ir,
894 gmx::ArrayRef<const real> lambda,
900 unsigned short *ptype;
902 if (DOMAINDECOMP(cr))
910 if (n > shfc->adir_nalloc)
912 shfc->adir_nalloc = over_alloc_dd(n);
913 srenew(shfc->adir_xnold, shfc->adir_nalloc);
914 srenew(shfc->adir_xnew, shfc->adir_nalloc);
916 xnold = shfc->adir_xnold;
917 xnew = shfc->adir_xnew;
923 /* Does NOT work with freeze or acceleration groups (yet) */
924 for (n = 0; n < end; n++)
926 w_dt = md->invmass[n]*dt;
928 for (d = 0; d < DIM; d++)
930 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
932 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
933 xnew[n][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
937 xnold[n][d] = x[n][d];
938 xnew[n][d] = x[n][d];
942 constr->apply(FALSE, FALSE, step, 0, 1.0,
943 x, xnold, nullptr, box,
944 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
945 nullptr, nullptr, gmx::ConstraintVariable::Positions);
946 constr->apply(FALSE, FALSE, step, 0, 1.0,
947 x, xnew, nullptr, box,
948 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
949 nullptr, nullptr, gmx::ConstraintVariable::Positions);
951 for (n = 0; n < end; n++)
953 for (d = 0; d < DIM; d++)
956 -(2*x[n][d]-xnold[n][d]-xnew[n][d])/gmx::square(dt)
957 - f[n][d]*md->invmass[n];
959 clear_rvec(acc_dir[n]);
962 /* Project the acceleration on the old bond directions */
963 constr->apply(FALSE, FALSE, step, 0, 1.0,
964 x_old, xnew, acc_dir, box,
965 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
966 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,
981 gmx::Constraints *constr,
982 gmx_enerdata_t *enerd,
985 gmx::ArrayRefWithPadding<gmx::RVec> f,
989 gmx_wallcycle_t wcycle,
993 gmx::PpForceWorkload *ppForceWorkload,
996 const gmx_vsite_t *vsite,
997 const DDBalanceRegionHandler &ddBalanceRegionHandler)
1002 rvec *acc_dir = nullptr, *x_old = nullptr;
1003 real Epot[2], df[2];
1007 gmx_bool bCont, bInit, bConverged;
1008 int nat, dd_ac0, dd_ac1 = 0, i;
1009 int homenr = md->homenr, end = homenr, cg0, cg1;
1010 int nflexcon, number_steps, d, Min = 0, count = 0;
1011 #define Try (1-Min) /* At start Try = 1 */
1013 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1014 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
1015 ftol = inputrec->em_tol;
1016 number_steps = inputrec->niter;
1017 nshell = shfc->nshell;
1018 shell = shfc->shell;
1019 nflexcon = shfc->nflexcon;
1023 if (DOMAINDECOMP(cr))
1025 nat = dd_natoms_vsite(cr->dd);
1028 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1029 nat = std::max(nat, dd_ac1);
1034 nat = state->natoms;
1037 for (i = 0; (i < 2); i++)
1039 shfc->x[i].resizeWithPadding(nat);
1040 shfc->f[i].resizeWithPadding(nat);
1043 /* Create views that we can swap */
1044 gmx::ArrayRefWithPadding<gmx::RVec> posWithPadding[2];
1045 gmx::ArrayRefWithPadding<gmx::RVec> forceWithPadding[2];
1046 gmx::ArrayRef<gmx::RVec> pos[2];
1047 gmx::ArrayRef<gmx::RVec> force[2];
1048 for (i = 0; (i < 2); i++)
1050 posWithPadding[i] = shfc->x[i].arrayRefWithPadding();
1051 pos[i] = posWithPadding[i].paddedArrayRef();
1052 forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1053 force[i] = forceWithPadding[i].paddedArrayRef();
1056 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1058 /* This is the only time where the coordinates are used
1059 * before do_force is called, which normally puts all
1060 * charge groups in the box.
1062 if (inputrec->cutoff_scheme == ecutsVERLET)
1064 auto xRef = state->x.arrayRefWithPadding().paddedArrayRef();
1065 put_atoms_in_box_omp(fr->ePBC, state->box, xRef.subArray(0, md->homenr), gmx_omp_nthreads_get(emntDefault));
1071 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1072 &(top->cgs), state->x.rvec_array(), fr->cg_cm);
1077 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x.rvec_array());
1081 /* After this all coordinate arrays will contain whole charge groups */
1084 shift_self(graph, state->box, state->x.rvec_array());
1089 if (nat > shfc->flex_nalloc)
1091 shfc->flex_nalloc = over_alloc_dd(nat);
1092 srenew(shfc->acc_dir, shfc->flex_nalloc);
1093 srenew(shfc->x_old, shfc->flex_nalloc);
1095 acc_dir = shfc->acc_dir;
1096 x_old = shfc->x_old;
1097 auto x = makeArrayRef(state->x);
1098 auto v = makeArrayRef(state->v);
1099 for (i = 0; i < homenr; i++)
1101 for (d = 0; d < DIM; d++)
1104 x[i][d] - v[i][d]*inputrec->delta_t;
1109 /* Do a prediction of the shell positions, when appropriate.
1110 * Without velocities (EM, NM, BD) we only do initial prediction.
1112 if (shfc->bPredict && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1114 predict_shells(fplog, state->x.rvec_array(), state->v.rvec_array(), inputrec->delta_t, nshell, shell,
1115 md->massT, nullptr, bInit);
1118 /* do_force expected the charge groups to be in the box */
1121 unshift_self(graph, state->box, state->x.rvec_array());
1124 /* Calculate the forces first time around */
1127 pr_rvecs(debug, 0, "x b4 do_force", state->x.rvec_array(), homenr);
1129 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1130 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession,
1132 mdstep, nrnb, wcycle, top,
1133 state->box, state->x.arrayRefWithPadding(), &state->hist,
1134 forceWithPadding[Min], force_vir, md, enerd, fcd,
1135 state->lambda, graph,
1136 fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
1137 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags,
1138 ddBalanceRegionHandler);
1144 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1145 shfc->x_old, state->x.rvec_array(), state->x.rvec_array(), as_rvec_array(force[Min].data()),
1147 state->box, state->lambda, &dum);
1149 for (i = 0; i < end; i++)
1151 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i]);
1154 sum_epot(&(enerd->grpp), enerd->term);
1155 Epot[Min] = enerd->term[F_EPOT];
1157 df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1161 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1166 pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md->nr);
1169 if (nshell+nflexcon > 0)
1171 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1172 * shell positions are updated, therefore the other particles must
1173 * be set here, in advance.
1175 std::copy(state->x.begin(),
1177 posWithPadding[Min].paddedArrayRef().begin());
1178 std::copy(state->x.begin(),
1180 posWithPadding[Try].paddedArrayRef().begin());
1183 if (bVerbose && MASTER(cr))
1185 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1190 fprintf(debug, "%17s: %14.10e\n",
1191 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1192 fprintf(debug, "%17s: %14.10e\n",
1193 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1194 fprintf(debug, "%17s: %14.10e\n",
1195 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1196 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1199 /* First check whether we should do shells, or whether the force is
1200 * low enough even without minimization.
1202 bConverged = (df[Min] < ftol);
1204 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1208 construct_vsites(vsite, as_rvec_array(pos[Min].data()),
1209 inputrec->delta_t, state->v.rvec_array(),
1210 idef->iparams, idef->il,
1211 fr->ePBC, fr->bMolPBC, cr, state->box);
1217 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1218 x_old, state->x.rvec_array(),
1219 as_rvec_array(pos[Min].data()),
1220 as_rvec_array(force[Min].data()), acc_dir,
1221 state->box, state->lambda, &dum);
1223 directional_sd(pos[Min], pos[Try], acc_dir, end, fr->fc_stepsize);
1226 /* New positions, Steepest descent */
1227 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1229 /* do_force expected the charge groups to be in the box */
1232 unshift_self(graph, state->box, as_rvec_array(pos[Try].data()));
1237 pr_rvecs(debug, 0, "RELAX: pos[Min] ", as_rvec_array(pos[Min].data()), homenr);
1238 pr_rvecs(debug, 0, "RELAX: pos[Try] ", as_rvec_array(pos[Try].data()), homenr);
1240 /* Try the new positions */
1241 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession,
1244 top, state->box, posWithPadding[Try], &state->hist,
1245 forceWithPadding[Try], force_vir,
1246 md, enerd, fcd, state->lambda, graph,
1247 fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
1249 ddBalanceRegionHandler);
1250 sum_epot(&(enerd->grpp), enerd->term);
1253 pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1254 pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1260 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1261 x_old, state->x.rvec_array(),
1262 as_rvec_array(pos[Try].data()),
1263 as_rvec_array(force[Try].data()),
1264 acc_dir, state->box, state->lambda, &dum);
1266 for (i = 0; i < end; i++)
1268 sf_dir += md->massT[i]*norm2(acc_dir[i]);
1272 Epot[Try] = enerd->term[F_EPOT];
1274 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1278 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1285 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1289 fprintf(debug, "SHELL ITER %d\n", count);
1290 dump_shells(debug, force[Try], ftol, nshell, shell);
1294 if (bVerbose && MASTER(cr))
1296 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1299 bConverged = (df[Try] < ftol);
1301 if ((df[Try] < df[Min]))
1305 fprintf(debug, "Swapping Min and Try\n");
1309 /* Correct the velocities for the flexible constraints */
1310 invdt = 1/inputrec->delta_t;
1311 auto v = makeArrayRef(state->v);
1312 for (i = 0; i < end; i++)
1314 for (d = 0; d < DIM; d++)
1316 v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1324 decrease_step_size(nshell, shell);
1327 shfc->numForceEvaluations += count;
1330 shfc->numConvergedIterations++;
1332 if (MASTER(cr) && !(bConverged))
1334 /* Note that the energies and virial are incorrect when not converged */
1338 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1339 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1342 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1343 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1346 /* Copy back the coordinates and the forces */
1347 std::copy(pos[Min].begin(), pos[Min].end(), makeArrayRef(state->x).data());
1348 std::copy(force[Min].begin(), force[Min].end(), f.unpaddedArrayRef().begin());
1351 void done_shellfc(FILE *fplog, gmx_shellfc_t *shfc, int64_t numSteps)
1353 if (shfc && fplog && numSteps > 0)
1355 double numStepsAsDouble = static_cast<double>(numSteps);
1356 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1357 (shfc->numConvergedIterations*100.0)/numStepsAsDouble);
1358 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1359 shfc->numForceEvaluations/numStepsAsDouble);
1362 // TODO Deallocate memory in shfc