Refactor gmx_group_t to SimulationAtomGroups
[gromacs.git] / src / gromacs / mdlib / shellfc.cpp
blob1917fd76755fa15f6a2075f0dfba051398af7b44
1 /*
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.
37 #include "gmxpre.h"
39 #include "shellfc.h"
41 #include <cmath>
42 #include <cstdint>
43 #include <cstdlib>
44 #include <cstring>
46 #include <algorithm>
47 #include <array>
49 #include "gromacs/domdec/dlbtiming.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/gmxlib/chargegroup.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/force.h"
60 #include "gromacs/mdlib/force_flags.h"
61 #include "gromacs/mdlib/gmx_omp_nthreads.h"
62 #include "gromacs/mdlib/mdatoms.h"
63 #include "gromacs/mdlib/vsite.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/enerdata.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/pbcutil/mshift.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/topology/ifunc.h"
73 #include "gromacs/topology/mtop_lookup.h"
74 #include "gromacs/topology/mtop_util.h"
75 #include "gromacs/utility/arrayref.h"
76 #include "gromacs/utility/arraysize.h"
77 #include "gromacs/utility/cstringutil.h"
78 #include "gromacs/utility/fatalerror.h"
79 #include "gromacs/utility/smalloc.h"
81 typedef struct {
82 int nnucl;
83 int shell; /* The shell id */
84 int nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
85 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
86 real k; /* force constant */
87 real k_1; /* 1 over force constant */
88 rvec xold;
89 rvec fold;
90 rvec step;
91 } t_shell;
93 struct gmx_shellfc_t {
94 /* Shell counts, indices, parameters and working data */
95 int nshell_gl; /* The number of shells in the system */
96 t_shell *shell_gl; /* All the shells (for DD only) */
97 int *shell_index_gl; /* Global shell index (for DD only) */
98 gmx_bool bInterCG; /* Are there inter charge-group shells? */
99 int nshell; /* The number of local shells */
100 t_shell *shell; /* The local shells */
101 int shell_nalloc; /* The allocation size of shell */
102 gmx_bool bPredict; /* Predict shell positions */
103 gmx_bool bRequireInit; /* Require initialization of shell positions */
104 int nflexcon; /* The number of flexible constraints */
106 /* Temporary arrays, should be fixed size 2 when fully converted to C++ */
107 PaddedVector<gmx::RVec> *x; /* Array for iterative minimization */
108 PaddedVector<gmx::RVec> *f; /* Array for iterative minimization */
110 /* Flexible constraint working data */
111 rvec *acc_dir; /* Acceleration direction for flexcon */
112 rvec *x_old; /* Old coordinates for flexcon */
113 int flex_nalloc; /* The allocation size of acc_dir and x_old */
114 rvec *adir_xnold; /* Work space for init_adir */
115 rvec *adir_xnew; /* Work space for init_adir */
116 int adir_nalloc; /* Work space for init_adir */
117 std::int64_t numForceEvaluations; /* Total number of force evaluations */
118 int numConvergedIterations; /* Total number of iterations that converged */
122 static void pr_shell(FILE *fplog, int ns, t_shell s[])
124 int i;
126 fprintf(fplog, "SHELL DATA\n");
127 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
128 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
129 for (i = 0; (i < ns); i++)
131 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
132 if (s[i].nnucl == 2)
134 fprintf(fplog, " %5d\n", s[i].nucl2);
136 else if (s[i].nnucl == 3)
138 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
140 else
142 fprintf(fplog, "\n");
147 /* TODO The remain call of this function passes non-NULL mass and NULL
148 * mtop, so this routine can be simplified.
150 * The other code path supported doing prediction before the MD loop
151 * started, but even when called, the prediction was always
152 * over-written by a subsequent call in the MD loop, so has been
153 * removed. */
154 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
155 int ns, t_shell s[],
156 const real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
158 int i, m, s1, n1, n2, n3;
159 real dt_1, fudge, tm, m1, m2, m3;
160 rvec *ptr;
162 /* We introduce a fudge factor for performance reasons: with this choice
163 * the initial force on the shells is about a factor of two lower than
164 * without
166 fudge = 1.0;
168 if (bInit)
170 if (fplog)
172 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
174 ptr = x;
175 dt_1 = 1;
177 else
179 ptr = v;
180 dt_1 = fudge*dt;
183 int molb = 0;
184 for (i = 0; (i < ns); i++)
186 s1 = s[i].shell;
187 if (bInit)
189 clear_rvec(x[s1]);
191 switch (s[i].nnucl)
193 case 1:
194 n1 = s[i].nucl1;
195 for (m = 0; (m < DIM); m++)
197 x[s1][m] += ptr[n1][m]*dt_1;
199 break;
200 case 2:
201 n1 = s[i].nucl1;
202 n2 = s[i].nucl2;
203 if (mass)
205 m1 = mass[n1];
206 m2 = mass[n2];
208 else
210 /* Not the correct masses with FE, but it is just a prediction... */
211 m1 = mtopGetAtomMass(mtop, n1, &molb);
212 m2 = mtopGetAtomMass(mtop, n2, &molb);
214 tm = dt_1/(m1+m2);
215 for (m = 0; (m < DIM); m++)
217 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
219 break;
220 case 3:
221 n1 = s[i].nucl1;
222 n2 = s[i].nucl2;
223 n3 = s[i].nucl3;
224 if (mass)
226 m1 = mass[n1];
227 m2 = mass[n2];
228 m3 = mass[n3];
230 else
232 /* Not the correct masses with FE, but it is just a prediction... */
233 m1 = mtopGetAtomMass(mtop, n1, &molb);
234 m2 = mtopGetAtomMass(mtop, n2, &molb);
235 m3 = mtopGetAtomMass(mtop, n3, &molb);
237 tm = dt_1/(m1+m2+m3);
238 for (m = 0; (m < DIM); m++)
240 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
242 break;
243 default:
244 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
249 /*! \brief Count the different particle types in a system
251 * Routine prints a warning to stderr in case an unknown particle type
252 * is encountered.
253 * \param[in] fplog Print what we have found if not NULL
254 * \param[in] mtop Molecular topology.
255 * \returns Array holding the number of particles of a type
257 static std::array<int, eptNR> countPtypes(FILE *fplog,
258 const gmx_mtop_t *mtop)
260 std::array<int, eptNR> nptype = { { 0 } };
261 /* Count number of shells, and find their indices */
262 for (int i = 0; (i < eptNR); i++)
264 nptype[i] = 0;
267 gmx_mtop_atomloop_block_t aloopb = gmx_mtop_atomloop_block_init(mtop);
268 int nmol;
269 const t_atom *atom;
270 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
272 switch (atom->ptype)
274 case eptAtom:
275 case eptVSite:
276 case eptShell:
277 nptype[atom->ptype] += nmol;
278 break;
279 default:
280 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
281 static_cast<int>(atom->ptype));
284 if (fplog)
286 /* Print the number of each particle type */
287 int n = 0;
288 for (const auto &i : nptype)
290 if (i != 0)
292 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
294 n++;
297 return nptype;
300 gmx_shellfc_t *init_shell_flexcon(FILE *fplog,
301 const gmx_mtop_t *mtop, int nflexcon,
302 int nstcalcenergy,
303 bool usingDomainDecomposition)
305 gmx_shellfc_t *shfc;
306 t_shell *shell;
307 int *shell_index = nullptr, *at2cg;
309 int ns, nshell, nsi;
310 int i, j, type, a_offset, cg, mol, ftype, nra;
311 real qS, alpha;
312 int aS, aN = 0; /* Shell and nucleus */
313 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
314 #define NBT asize(bondtypes)
315 const gmx_ffparams_t *ffparams;
317 std::array<int, eptNR> n = countPtypes(fplog, mtop);
318 nshell = n[eptShell];
320 if (nshell == 0 && nflexcon == 0)
322 /* We're not doing shells or flexible constraints */
323 return nullptr;
326 snew(shfc, 1);
327 shfc->x = new PaddedVector<gmx::RVec>[2] {};
328 shfc->f = new PaddedVector<gmx::RVec>[2] {};
329 shfc->nflexcon = nflexcon;
331 if (nshell == 0)
333 /* Only flexible constraints, no shells.
334 * Note that make_local_shells() does not need to be called.
336 shfc->nshell = 0;
337 shfc->bPredict = FALSE;
339 return shfc;
342 if (nstcalcenergy != 1)
344 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);
346 if (usingDomainDecomposition)
348 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
351 /* We have shells: fill the shell data structure */
353 /* Global system sized array, this should be avoided */
354 snew(shell_index, mtop->natoms);
356 nshell = 0;
357 for (const AtomProxy atomP : AtomRange(*mtop))
359 const t_atom &local = atomP.atom();
360 int i = atomP.globalAtomNumber();
361 if (local.ptype == eptShell)
363 shell_index[i] = nshell++;
367 snew(shell, nshell);
369 /* Initiate the shell structures */
370 for (i = 0; (i < nshell); i++)
372 shell[i].shell = -1;
373 shell[i].nnucl = 0;
374 shell[i].nucl1 = -1;
375 shell[i].nucl2 = -1;
376 shell[i].nucl3 = -1;
377 /* shell[i].bInterCG=FALSE; */
378 shell[i].k_1 = 0;
379 shell[i].k = 0;
382 ffparams = &mtop->ffparams;
384 /* Now fill the structures */
385 shfc->bInterCG = FALSE;
386 ns = 0;
387 a_offset = 0;
388 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
390 const gmx_molblock_t *molb = &mtop->molblock[mb];
391 const gmx_moltype_t *molt = &mtop->moltype[molb->type];
392 const t_block *cgs = &molt->cgs;
394 snew(at2cg, molt->atoms.nr);
395 for (cg = 0; cg < cgs->nr; cg++)
397 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
399 at2cg[i] = cg;
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()); )
411 type = ia[0];
412 ftype = ffparams->functype[type];
413 nra = interaction_function[ftype].nratoms;
415 /* Check whether we have a bond with a shell */
416 aS = -1;
418 switch (bondtypes[j])
420 case F_BONDS:
421 case F_HARMONIC:
422 case F_CUBICBONDS:
423 case F_POLARIZATION:
424 case F_ANHARM_POL:
425 if (atom[ia[1]].ptype == eptShell)
427 aS = ia[1];
428 aN = ia[2];
430 else if (atom[ia[2]].ptype == eptShell)
432 aS = ia[2];
433 aN = ia[1];
435 break;
436 case F_WATER_POL:
437 aN = ia[4]; /* Dummy */
438 aS = ia[5]; /* Shell */
439 break;
440 default:
441 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
444 if (aS != -1)
446 qS = atom[aS].q;
448 /* Check whether one of the particles is a shell... */
449 nsi = shell_index[a_offset+aS];
450 if ((nsi < 0) || (nsi >= nshell))
452 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
453 nsi, nshell, aS);
455 if (shell[nsi].shell == -1)
457 shell[nsi].shell = a_offset + aS;
458 ns++;
460 else if (shell[nsi].shell != a_offset+aS)
462 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
465 if (shell[nsi].nucl1 == -1)
467 shell[nsi].nucl1 = a_offset + aN;
469 else if (shell[nsi].nucl2 == -1)
471 shell[nsi].nucl2 = a_offset + aN;
473 else if (shell[nsi].nucl3 == -1)
475 shell[nsi].nucl3 = a_offset + aN;
477 else
479 if (fplog)
481 pr_shell(fplog, ns, shell);
483 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
485 if (at2cg[aS] != at2cg[aN])
487 /* shell[nsi].bInterCG = TRUE; */
488 shfc->bInterCG = TRUE;
491 switch (bondtypes[j])
493 case F_BONDS:
494 case F_HARMONIC:
495 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
496 break;
497 case F_CUBICBONDS:
498 shell[nsi].k += ffparams->iparams[type].cubic.kb;
499 break;
500 case F_POLARIZATION:
501 case F_ANHARM_POL:
502 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
504 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);
506 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/
507 ffparams->iparams[type].polarize.alpha;
508 break;
509 case F_WATER_POL:
510 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
512 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);
514 alpha = (ffparams->iparams[type].wpol.al_x+
515 ffparams->iparams[type].wpol.al_y+
516 ffparams->iparams[type].wpol.al_z)/3.0;
517 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/alpha;
518 break;
519 default:
520 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
522 shell[nsi].nnucl++;
524 ia += nra+1;
525 i += nra+1;
528 a_offset += molt->atoms.nr;
530 /* Done with this molecule type */
531 sfree(at2cg);
534 /* Verify whether it's all correct */
535 if (ns != nshell)
537 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
540 for (i = 0; (i < ns); i++)
542 shell[i].k_1 = 1.0/shell[i].k;
545 if (debug)
547 pr_shell(debug, ns, shell);
551 shfc->nshell_gl = ns;
552 shfc->shell_gl = shell;
553 shfc->shell_index_gl = shell_index;
555 shfc->bPredict = (getenv("GMX_NOPREDICT") == nullptr);
556 shfc->bRequireInit = FALSE;
557 if (!shfc->bPredict)
559 if (fplog)
561 fprintf(fplog, "\nWill never predict shell positions\n");
564 else
566 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
567 if (shfc->bRequireInit && fplog)
569 fprintf(fplog, "\nWill always initiate shell positions\n");
573 if (shfc->bPredict)
575 if (shfc->bInterCG)
577 if (fplog)
579 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");
581 /* Prediction improves performance, so we should implement either:
582 * 1. communication for the atoms needed for prediction
583 * 2. prediction using the velocities of shells; currently the
584 * shell velocities are zeroed, it's a bit tricky to keep
585 * track of the shell displacements and thus the velocity.
587 shfc->bPredict = FALSE;
591 return shfc;
594 void make_local_shells(const t_commrec *cr,
595 const t_mdatoms *md,
596 gmx_shellfc_t *shfc)
598 t_shell *shell;
599 int a0, a1, *ind, nshell, i;
600 gmx_domdec_t *dd = nullptr;
602 if (DOMAINDECOMP(cr))
604 dd = cr->dd;
605 a0 = 0;
606 a1 = dd_numHomeAtoms(*dd);
608 else
610 /* Single node: we need all shells, just copy the pointer */
611 shfc->nshell = shfc->nshell_gl;
612 shfc->shell = shfc->shell_gl;
614 return;
617 ind = shfc->shell_index_gl;
619 nshell = 0;
620 shell = shfc->shell;
621 for (i = a0; i < a1; i++)
623 if (md->ptype[i] == eptShell)
625 if (nshell+1 > shfc->shell_nalloc)
627 shfc->shell_nalloc = over_alloc_dd(nshell+1);
628 srenew(shell, shfc->shell_nalloc);
630 if (dd)
632 shell[nshell] = shfc->shell_gl[ind[dd->globalAtomIndices[i]]];
634 else
636 shell[nshell] = shfc->shell_gl[ind[i]];
639 /* With inter-cg shells we can no do shell prediction,
640 * so we do not need the nuclei numbers.
642 if (!shfc->bInterCG)
644 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
645 if (shell[nshell].nnucl > 1)
647 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
649 if (shell[nshell].nnucl > 2)
651 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
654 shell[nshell].shell = i;
655 nshell++;
659 shfc->nshell = nshell;
660 shfc->shell = shell;
663 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
665 real xo, yo, zo;
666 real dx, dy, dz;
668 xo = xold[XX];
669 yo = xold[YY];
670 zo = xold[ZZ];
672 dx = f[XX]*step;
673 dy = f[YY]*step;
674 dz = f[ZZ]*step;
676 xnew[XX] = xo+dx;
677 xnew[YY] = yo+dy;
678 xnew[ZZ] = zo+dz;
681 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
683 real xo, yo, zo;
684 real dx, dy, dz;
686 xo = xold[XX];
687 yo = xold[YY];
688 zo = xold[ZZ];
690 dx = f[XX]*step[XX];
691 dy = f[YY]*step[YY];
692 dz = f[ZZ]*step[ZZ];
694 xnew[XX] = xo+dx;
695 xnew[YY] = yo+dy;
696 xnew[ZZ] = zo+dz;
699 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
700 gmx::ArrayRef<gmx::RVec> xnew,
701 const rvec acc_dir[], int homenr, real step)
703 const rvec *xo = as_rvec_array(xold.data());
704 rvec *xn = as_rvec_array(xnew.data());
706 for (int i = 0; i < homenr; i++)
708 do_1pos(xn[i], xo[i], acc_dir[i], step);
712 static void shell_pos_sd(gmx::ArrayRef<const gmx::RVec> xcur,
713 gmx::ArrayRef<gmx::RVec> xnew,
714 gmx::ArrayRef<const gmx::RVec> f,
715 int ns, t_shell s[], int count)
717 const real step_scale_min = 0.8,
718 step_scale_increment = 0.2,
719 step_scale_max = 1.2,
720 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
721 int i, shell, d;
722 real dx, df, k_est;
723 const real zero = 0;
724 #ifdef PRINT_STEP
725 real step_min, step_max;
727 step_min = 1e30;
728 step_max = 0;
729 #endif
730 for (i = 0; (i < ns); i++)
732 shell = s[i].shell;
733 if (count == 1)
735 for (d = 0; d < DIM; d++)
737 s[i].step[d] = s[i].k_1;
738 #ifdef PRINT_STEP
739 step_min = std::min(step_min, s[i].step[d]);
740 step_max = std::max(step_max, s[i].step[d]);
741 #endif
744 else
746 for (d = 0; d < DIM; d++)
748 dx = xcur[shell][d] - s[i].xold[d];
749 df = f[shell][d] - s[i].fold[d];
750 /* -dx/df gets used to generate an interpolated value, but would
751 * cause a NaN if df were binary-equal to zero. Values close to
752 * zero won't cause problems (because of the min() and max()), so
753 * just testing for binary inequality is OK. */
754 if (zero != df)
756 k_est = -dx/df;
757 /* Scale the step size by a factor interpolated from
758 * step_scale_min to step_scale_max, as k_est goes from 0 to
759 * step_scale_multiple * s[i].step[d] */
760 s[i].step[d] =
761 step_scale_min * s[i].step[d] +
762 step_scale_increment * std::min(step_scale_multiple * s[i].step[d], std::max(k_est, zero));
764 else
766 /* Here 0 == df */
767 if (gmx_numzero(dx)) /* 0 == dx */
769 /* Likely this will never happen, but if it does just
770 * don't scale the step. */
772 else /* 0 != dx */
774 s[i].step[d] *= step_scale_max;
777 #ifdef PRINT_STEP
778 step_min = std::min(step_min, s[i].step[d]);
779 step_max = std::max(step_max, s[i].step[d]);
780 #endif
783 copy_rvec(xcur [shell], s[i].xold);
784 copy_rvec(f[shell], s[i].fold);
786 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
788 if (gmx_debug_at)
790 fprintf(debug, "shell[%d] = %d\n", i, shell);
791 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
792 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
793 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
794 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
797 #ifdef PRINT_STEP
798 printf("step %.3e %.3e\n", step_min, step_max);
799 #endif
802 static void decrease_step_size(int nshell, t_shell s[])
804 int i;
806 for (i = 0; i < nshell; i++)
808 svmul(0.8, s[i].step, s[i].step);
812 static void print_epot(FILE *fp, int64_t mdstep, int count, real epot, real df,
813 int ndir, real sf_dir)
815 char buf[22];
817 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
818 gmx_step_str(mdstep, buf), count, epot, df);
819 if (ndir)
821 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir/ndir));
823 else
825 fprintf(fp, "\n");
830 static real rms_force(const t_commrec *cr, gmx::ArrayRef<const gmx::RVec> force, int ns, t_shell s[],
831 int ndir, real *sf_dir, real *Epot)
833 double buf[4];
834 const rvec *f = as_rvec_array(force.data());
836 buf[0] = *sf_dir;
837 for (int i = 0; i < ns; i++)
839 int shell = s[i].shell;
840 buf[0] += norm2(f[shell]);
842 int ntot = ns;
844 if (PAR(cr))
846 buf[1] = ntot;
847 buf[2] = *sf_dir;
848 buf[3] = *Epot;
849 gmx_sumd(4, buf, cr);
850 ntot = gmx::roundToInt(buf[1]);
851 *sf_dir = buf[2];
852 *Epot = buf[3];
854 ntot += ndir;
856 return (ntot ? std::sqrt(buf[0]/ntot) : 0);
859 static void dump_shells(FILE *fp, gmx::ArrayRef<gmx::RVec> f, real ftol, int ns, t_shell s[])
861 int i, shell;
862 real ft2, ff2;
864 ft2 = gmx::square(ftol);
866 for (i = 0; (i < ns); i++)
868 shell = s[i].shell;
869 ff2 = iprod(f[shell], f[shell]);
870 if (ff2 > ft2)
872 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
873 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
878 static void init_adir(gmx_shellfc_t *shfc,
879 gmx::Constraints *constr,
880 const t_inputrec *ir,
881 const t_commrec *cr,
882 int dd_ac1,
883 int64_t step,
884 const t_mdatoms *md,
885 int end,
886 rvec *x_old,
887 rvec *x_init,
888 rvec *x,
889 rvec *f,
890 rvec *acc_dir,
891 matrix box,
892 gmx::ArrayRef<const real> lambda,
893 real *dvdlambda)
895 rvec *xnold, *xnew;
896 double dt, w_dt;
897 int n, d;
898 unsigned short *ptype;
900 if (DOMAINDECOMP(cr))
902 n = dd_ac1;
904 else
906 n = end;
908 if (n > shfc->adir_nalloc)
910 shfc->adir_nalloc = over_alloc_dd(n);
911 srenew(shfc->adir_xnold, shfc->adir_nalloc);
912 srenew(shfc->adir_xnew, shfc->adir_nalloc);
914 xnold = shfc->adir_xnold;
915 xnew = shfc->adir_xnew;
917 ptype = md->ptype;
919 dt = ir->delta_t;
921 /* Does NOT work with freeze or acceleration groups (yet) */
922 for (n = 0; n < end; n++)
924 w_dt = md->invmass[n]*dt;
926 for (d = 0; d < DIM; d++)
928 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
930 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
931 xnew[n][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
933 else
935 xnold[n][d] = x[n][d];
936 xnew[n][d] = x[n][d];
940 constr->apply(FALSE, FALSE, step, 0, 1.0,
941 x, xnold, nullptr, box,
942 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
943 nullptr, nullptr, gmx::ConstraintVariable::Positions);
944 constr->apply(FALSE, FALSE, step, 0, 1.0,
945 x, xnew, nullptr, box,
946 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
947 nullptr, nullptr, gmx::ConstraintVariable::Positions);
949 for (n = 0; n < end; n++)
951 for (d = 0; d < DIM; d++)
953 xnew[n][d] =
954 -(2*x[n][d]-xnold[n][d]-xnew[n][d])/gmx::square(dt)
955 - f[n][d]*md->invmass[n];
957 clear_rvec(acc_dir[n]);
960 /* Project the acceleration on the old bond directions */
961 constr->apply(FALSE, FALSE, step, 0, 1.0,
962 x_old, xnew, acc_dir, box,
963 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
964 nullptr, nullptr, gmx::ConstraintVariable::Deriv_FlexCon);
967 void relax_shell_flexcon(FILE *fplog,
968 const t_commrec *cr,
969 const gmx_multisim_t *ms,
970 gmx_bool bVerbose,
971 gmx_enfrot *enforcedRotation,
972 int64_t mdstep,
973 const t_inputrec *inputrec,
974 gmx_bool bDoNS,
975 int force_flags,
976 gmx_localtop_t *top,
977 gmx::Constraints *constr,
978 gmx_enerdata_t *enerd,
979 t_fcdata *fcd,
980 t_state *state,
981 gmx::ArrayRefWithPadding<gmx::RVec> f,
982 tensor force_vir,
983 const t_mdatoms *md,
984 t_nrnb *nrnb,
985 gmx_wallcycle_t wcycle,
986 t_graph *graph,
987 const SimulationGroups *groups,
988 gmx_shellfc_t *shfc,
989 t_forcerec *fr,
990 gmx::PpForceWorkload *ppForceWorkload,
991 double t,
992 rvec mu_tot,
993 const gmx_vsite_t *vsite,
994 const DDBalanceRegionHandler &ddBalanceRegionHandler)
996 int nshell;
997 t_shell *shell;
998 const t_idef *idef;
999 rvec *acc_dir = nullptr, *x_old = nullptr;
1000 real Epot[2], df[2];
1001 real sf_dir, invdt;
1002 real ftol, dum = 0;
1003 char sbuf[22];
1004 gmx_bool bCont, bInit, bConverged;
1005 int nat, dd_ac0, dd_ac1 = 0, i;
1006 int homenr = md->homenr, end = homenr, cg0, cg1;
1007 int nflexcon, number_steps, d, Min = 0, count = 0;
1008 #define Try (1-Min) /* At start Try = 1 */
1010 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1011 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
1012 ftol = inputrec->em_tol;
1013 number_steps = inputrec->niter;
1014 nshell = shfc->nshell;
1015 shell = shfc->shell;
1016 nflexcon = shfc->nflexcon;
1018 idef = &top->idef;
1020 if (DOMAINDECOMP(cr))
1022 nat = dd_natoms_vsite(cr->dd);
1023 if (nflexcon > 0)
1025 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1026 nat = std::max(nat, dd_ac1);
1029 else
1031 nat = state->natoms;
1034 for (i = 0; (i < 2); i++)
1036 shfc->x[i].resizeWithPadding(nat);
1037 shfc->f[i].resizeWithPadding(nat);
1040 /* Create views that we can swap */
1041 gmx::ArrayRefWithPadding<gmx::RVec> posWithPadding[2];
1042 gmx::ArrayRefWithPadding<gmx::RVec> forceWithPadding[2];
1043 gmx::ArrayRef<gmx::RVec> pos[2];
1044 gmx::ArrayRef<gmx::RVec> force[2];
1045 for (i = 0; (i < 2); i++)
1047 posWithPadding[i] = shfc->x[i].arrayRefWithPadding();
1048 pos[i] = posWithPadding[i].paddedArrayRef();
1049 forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1050 force[i] = forceWithPadding[i].paddedArrayRef();
1053 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1055 /* This is the only time where the coordinates are used
1056 * before do_force is called, which normally puts all
1057 * charge groups in the box.
1059 if (inputrec->cutoff_scheme == ecutsVERLET)
1061 auto xRef = state->x.arrayRefWithPadding().paddedArrayRef();
1062 put_atoms_in_box_omp(fr->ePBC, state->box, xRef.subArray(0, md->homenr), gmx_omp_nthreads_get(emntDefault));
1064 else
1066 cg0 = 0;
1067 cg1 = top->cgs.nr;
1068 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1069 &(top->cgs), state->x.rvec_array(), fr->cg_cm);
1072 if (graph)
1074 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x.rvec_array());
1078 /* After this all coordinate arrays will contain whole charge groups */
1079 if (graph)
1081 shift_self(graph, state->box, state->x.rvec_array());
1084 if (nflexcon)
1086 if (nat > shfc->flex_nalloc)
1088 shfc->flex_nalloc = over_alloc_dd(nat);
1089 srenew(shfc->acc_dir, shfc->flex_nalloc);
1090 srenew(shfc->x_old, shfc->flex_nalloc);
1092 acc_dir = shfc->acc_dir;
1093 x_old = shfc->x_old;
1094 auto x = makeArrayRef(state->x);
1095 auto v = makeArrayRef(state->v);
1096 for (i = 0; i < homenr; i++)
1098 for (d = 0; d < DIM; d++)
1100 shfc->x_old[i][d] =
1101 x[i][d] - v[i][d]*inputrec->delta_t;
1106 /* Do a prediction of the shell positions, when appropriate.
1107 * Without velocities (EM, NM, BD) we only do initial prediction.
1109 if (shfc->bPredict && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1111 predict_shells(fplog, state->x.rvec_array(), state->v.rvec_array(), inputrec->delta_t, nshell, shell,
1112 md->massT, nullptr, bInit);
1115 /* do_force expected the charge groups to be in the box */
1116 if (graph)
1118 unshift_self(graph, state->box, state->x.rvec_array());
1121 /* Calculate the forces first time around */
1122 if (gmx_debug_at)
1124 pr_rvecs(debug, 0, "x b4 do_force", state->x.rvec_array(), homenr);
1126 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1127 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation,
1128 mdstep, nrnb, wcycle, top, groups,
1129 state->box, state->x.arrayRefWithPadding(), &state->hist,
1130 forceWithPadding[Min], force_vir, md, enerd, fcd,
1131 state->lambda, graph,
1132 fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
1133 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags,
1134 ddBalanceRegionHandler);
1136 sf_dir = 0;
1137 if (nflexcon)
1139 init_adir(shfc,
1140 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1141 shfc->x_old, state->x.rvec_array(), state->x.rvec_array(), as_rvec_array(force[Min].data()),
1142 shfc->acc_dir,
1143 state->box, state->lambda, &dum);
1145 for (i = 0; i < end; i++)
1147 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i]);
1150 sum_epot(&(enerd->grpp), enerd->term);
1151 Epot[Min] = enerd->term[F_EPOT];
1153 df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1154 df[Try] = 0;
1155 if (debug)
1157 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1160 if (gmx_debug_at)
1162 pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md->nr);
1165 if (nshell+nflexcon > 0)
1167 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1168 * shell positions are updated, therefore the other particles must
1169 * be set here, in advance.
1171 std::copy(state->x.begin(),
1172 state->x.end(),
1173 posWithPadding[Min].paddedArrayRef().begin());
1174 std::copy(state->x.begin(),
1175 state->x.end(),
1176 posWithPadding[Try].paddedArrayRef().begin());
1179 if (bVerbose && MASTER(cr))
1181 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1184 if (debug)
1186 fprintf(debug, "%17s: %14.10e\n",
1187 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1188 fprintf(debug, "%17s: %14.10e\n",
1189 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1190 fprintf(debug, "%17s: %14.10e\n",
1191 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1192 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1195 /* First check whether we should do shells, or whether the force is
1196 * low enough even without minimization.
1198 bConverged = (df[Min] < ftol);
1200 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1202 if (vsite)
1204 construct_vsites(vsite, as_rvec_array(pos[Min].data()),
1205 inputrec->delta_t, state->v.rvec_array(),
1206 idef->iparams, idef->il,
1207 fr->ePBC, fr->bMolPBC, cr, state->box);
1210 if (nflexcon)
1212 init_adir(shfc,
1213 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1214 x_old, state->x.rvec_array(),
1215 as_rvec_array(pos[Min].data()),
1216 as_rvec_array(force[Min].data()), acc_dir,
1217 state->box, state->lambda, &dum);
1219 directional_sd(pos[Min], pos[Try], acc_dir, end, fr->fc_stepsize);
1222 /* New positions, Steepest descent */
1223 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1225 /* do_force expected the charge groups to be in the box */
1226 if (graph)
1228 unshift_self(graph, state->box, as_rvec_array(pos[Try].data()));
1231 if (gmx_debug_at)
1233 pr_rvecs(debug, 0, "RELAX: pos[Min] ", as_rvec_array(pos[Min].data()), homenr);
1234 pr_rvecs(debug, 0, "RELAX: pos[Try] ", as_rvec_array(pos[Try].data()), homenr);
1236 /* Try the new positions */
1237 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation,
1238 1, nrnb, wcycle,
1239 top, groups, state->box, posWithPadding[Try], &state->hist,
1240 forceWithPadding[Try], force_vir,
1241 md, enerd, fcd, state->lambda, graph,
1242 fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
1243 shellfc_flags,
1244 ddBalanceRegionHandler);
1245 sum_epot(&(enerd->grpp), enerd->term);
1246 if (gmx_debug_at)
1248 pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1249 pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1251 sf_dir = 0;
1252 if (nflexcon)
1254 init_adir(shfc,
1255 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1256 x_old, state->x.rvec_array(),
1257 as_rvec_array(pos[Try].data()),
1258 as_rvec_array(force[Try].data()),
1259 acc_dir, state->box, state->lambda, &dum);
1261 for (i = 0; i < end; i++)
1263 sf_dir += md->massT[i]*norm2(acc_dir[i]);
1267 Epot[Try] = enerd->term[F_EPOT];
1269 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1271 if (debug)
1273 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1276 if (debug)
1278 if (gmx_debug_at)
1280 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1282 if (gmx_debug_at)
1284 fprintf(debug, "SHELL ITER %d\n", count);
1285 dump_shells(debug, force[Try], ftol, nshell, shell);
1289 if (bVerbose && MASTER(cr))
1291 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1294 bConverged = (df[Try] < ftol);
1296 if ((df[Try] < df[Min]))
1298 if (debug)
1300 fprintf(debug, "Swapping Min and Try\n");
1302 if (nflexcon)
1304 /* Correct the velocities for the flexible constraints */
1305 invdt = 1/inputrec->delta_t;
1306 auto v = makeArrayRef(state->v);
1307 for (i = 0; i < end; i++)
1309 for (d = 0; d < DIM; d++)
1311 v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1315 Min = Try;
1317 else
1319 decrease_step_size(nshell, shell);
1322 shfc->numForceEvaluations += count;
1323 if (bConverged)
1325 shfc->numConvergedIterations++;
1327 if (MASTER(cr) && !(bConverged))
1329 /* Note that the energies and virial are incorrect when not converged */
1330 if (fplog)
1332 fprintf(fplog,
1333 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1334 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1336 fprintf(stderr,
1337 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1338 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1341 /* Copy back the coordinates and the forces */
1342 std::copy(pos[Min].begin(), pos[Min].end(), makeArrayRef(state->x).data());
1343 std::copy(force[Min].begin(), force[Min].end(), f.unpaddedArrayRef().begin());
1346 void done_shellfc(FILE *fplog, gmx_shellfc_t *shfc, int64_t numSteps)
1348 if (shfc && fplog && numSteps > 0)
1350 double numStepsAsDouble = static_cast<double>(numSteps);
1351 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1352 (shfc->numConvergedIterations*100.0)/numStepsAsDouble);
1353 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1354 shfc->numForceEvaluations/numStepsAsDouble);
1357 // TODO Deallocate memory in shfc