Turn SystemAtomIterator into proper iterator
[gromacs.git] / src / gromacs / mdlib / shellfc.cpp
blob02e4c87c95f48aab58474cd22b40b8c473a7fcbb
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/mdrun.h"
62 #include "gromacs/mdlib/sim_util.h"
63 #include "gromacs/mdlib/vsite.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/forcerec.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/state.h"
69 #include "gromacs/pbcutil/mshift.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/topology/ifunc.h"
72 #include "gromacs/topology/mtop_lookup.h"
73 #include "gromacs/topology/mtop_util.h"
74 #include "gromacs/utility/arrayref.h"
75 #include "gromacs/utility/arraysize.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/smalloc.h"
80 typedef struct {
81 int nnucl;
82 int shell; /* The shell id */
83 int nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
84 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
85 real k; /* force constant */
86 real k_1; /* 1 over force constant */
87 rvec xold;
88 rvec fold;
89 rvec step;
90 } t_shell;
92 struct gmx_shellfc_t {
93 /* Shell counts, indices, parameters and working data */
94 int nshell_gl; /* The number of shells in the system */
95 t_shell *shell_gl; /* All the shells (for DD only) */
96 int *shell_index_gl; /* Global shell index (for DD only) */
97 gmx_bool bInterCG; /* Are there inter charge-group shells? */
98 int nshell; /* The number of local shells */
99 t_shell *shell; /* The local shells */
100 int shell_nalloc; /* The allocation size of shell */
101 gmx_bool bPredict; /* Predict shell positions */
102 gmx_bool bRequireInit; /* Require initialization of shell positions */
103 int nflexcon; /* The number of flexible constraints */
105 /* Temporary arrays, should be fixed size 2 when fully converted to C++ */
106 PaddedVector<gmx::RVec> *x; /* Array for iterative minimization */
107 PaddedVector<gmx::RVec> *f; /* Array for iterative minimization */
109 /* Flexible constraint working data */
110 rvec *acc_dir; /* Acceleration direction for flexcon */
111 rvec *x_old; /* Old coordinates for flexcon */
112 int flex_nalloc; /* The allocation size of acc_dir and x_old */
113 rvec *adir_xnold; /* Work space for init_adir */
114 rvec *adir_xnew; /* Work space for init_adir */
115 int adir_nalloc; /* Work space for init_adir */
116 std::int64_t numForceEvaluations; /* Total number of force evaluations */
117 int numConvergedIterations; /* Total number of iterations that converged */
121 static void pr_shell(FILE *fplog, int ns, t_shell s[])
123 int i;
125 fprintf(fplog, "SHELL DATA\n");
126 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
127 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
128 for (i = 0; (i < ns); i++)
130 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
131 if (s[i].nnucl == 2)
133 fprintf(fplog, " %5d\n", s[i].nucl2);
135 else if (s[i].nnucl == 3)
137 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
139 else
141 fprintf(fplog, "\n");
146 /* TODO The remain call of this function passes non-NULL mass and NULL
147 * mtop, so this routine can be simplified.
149 * The other code path supported doing prediction before the MD loop
150 * started, but even when called, the prediction was always
151 * over-written by a subsequent call in the MD loop, so has been
152 * removed. */
153 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
154 int ns, t_shell s[],
155 const real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
157 int i, m, s1, n1, n2, n3;
158 real dt_1, fudge, tm, m1, m2, m3;
159 rvec *ptr;
161 /* We introduce a fudge factor for performance reasons: with this choice
162 * the initial force on the shells is about a factor of two lower than
163 * without
165 fudge = 1.0;
167 if (bInit)
169 if (fplog)
171 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
173 ptr = x;
174 dt_1 = 1;
176 else
178 ptr = v;
179 dt_1 = fudge*dt;
182 int molb = 0;
183 for (i = 0; (i < ns); i++)
185 s1 = s[i].shell;
186 if (bInit)
188 clear_rvec(x[s1]);
190 switch (s[i].nnucl)
192 case 1:
193 n1 = s[i].nucl1;
194 for (m = 0; (m < DIM); m++)
196 x[s1][m] += ptr[n1][m]*dt_1;
198 break;
199 case 2:
200 n1 = s[i].nucl1;
201 n2 = s[i].nucl2;
202 if (mass)
204 m1 = mass[n1];
205 m2 = mass[n2];
207 else
209 /* Not the correct masses with FE, but it is just a prediction... */
210 m1 = mtopGetAtomMass(mtop, n1, &molb);
211 m2 = mtopGetAtomMass(mtop, n2, &molb);
213 tm = dt_1/(m1+m2);
214 for (m = 0; (m < DIM); m++)
216 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
218 break;
219 case 3:
220 n1 = s[i].nucl1;
221 n2 = s[i].nucl2;
222 n3 = s[i].nucl3;
223 if (mass)
225 m1 = mass[n1];
226 m2 = mass[n2];
227 m3 = mass[n3];
229 else
231 /* Not the correct masses with FE, but it is just a prediction... */
232 m1 = mtopGetAtomMass(mtop, n1, &molb);
233 m2 = mtopGetAtomMass(mtop, n2, &molb);
234 m3 = mtopGetAtomMass(mtop, n3, &molb);
236 tm = dt_1/(m1+m2+m3);
237 for (m = 0; (m < DIM); m++)
239 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
241 break;
242 default:
243 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
248 /*! \brief Count the different particle types in a system
250 * Routine prints a warning to stderr in case an unknown particle type
251 * is encountered.
252 * \param[in] fplog Print what we have found if not NULL
253 * \param[in] mtop Molecular topology.
254 * \returns Array holding the number of particles of a type
256 static std::array<int, eptNR> countPtypes(FILE *fplog,
257 const gmx_mtop_t *mtop)
259 std::array<int, eptNR> nptype = { { 0 } };
260 /* Count number of shells, and find their indices */
261 for (int i = 0; (i < eptNR); i++)
263 nptype[i] = 0;
266 gmx_mtop_atomloop_block_t aloopb = gmx_mtop_atomloop_block_init(mtop);
267 int nmol;
268 const t_atom *atom;
269 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
271 switch (atom->ptype)
273 case eptAtom:
274 case eptVSite:
275 case eptShell:
276 nptype[atom->ptype] += nmol;
277 break;
278 default:
279 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
280 static_cast<int>(atom->ptype));
283 if (fplog)
285 /* Print the number of each particle type */
286 int n = 0;
287 for (const auto &i : nptype)
289 if (i != 0)
291 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
293 n++;
296 return nptype;
299 gmx_shellfc_t *init_shell_flexcon(FILE *fplog,
300 const gmx_mtop_t *mtop, int nflexcon,
301 int nstcalcenergy,
302 bool usingDomainDecomposition)
304 gmx_shellfc_t *shfc;
305 t_shell *shell;
306 int *shell_index = nullptr, *at2cg;
308 int ns, nshell, nsi;
309 int i, j, type, a_offset, cg, mol, ftype, nra;
310 real qS, alpha;
311 int aS, aN = 0; /* Shell and nucleus */
312 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
313 #define NBT asize(bondtypes)
314 const gmx_ffparams_t *ffparams;
316 std::array<int, eptNR> n = countPtypes(fplog, mtop);
317 nshell = n[eptShell];
319 if (nshell == 0 && nflexcon == 0)
321 /* We're not doing shells or flexible constraints */
322 return nullptr;
325 snew(shfc, 1);
326 shfc->x = new PaddedVector<gmx::RVec>[2] {};
327 shfc->f = new PaddedVector<gmx::RVec>[2] {};
328 shfc->nflexcon = nflexcon;
330 if (nshell == 0)
332 /* Only flexible constraints, no shells.
333 * Note that make_local_shells() does not need to be called.
335 shfc->nshell = 0;
336 shfc->bPredict = FALSE;
338 return shfc;
341 if (nstcalcenergy != 1)
343 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);
345 if (usingDomainDecomposition)
347 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
350 /* We have shells: fill the shell data structure */
352 /* Global system sized array, this should be avoided */
353 snew(shell_index, mtop->natoms);
355 nshell = 0;
356 for (const AtomProxy &atomP : AtomRange(*mtop))
358 const t_atom &local = atomP.atom();
359 int i = atomP.globalAtomNumber();
360 if (local.ptype == eptShell)
362 shell_index[i] = nshell++;
366 snew(shell, nshell);
368 /* Initiate the shell structures */
369 for (i = 0; (i < nshell); i++)
371 shell[i].shell = -1;
372 shell[i].nnucl = 0;
373 shell[i].nucl1 = -1;
374 shell[i].nucl2 = -1;
375 shell[i].nucl3 = -1;
376 /* shell[i].bInterCG=FALSE; */
377 shell[i].k_1 = 0;
378 shell[i].k = 0;
381 ffparams = &mtop->ffparams;
383 /* Now fill the structures */
384 shfc->bInterCG = FALSE;
385 ns = 0;
386 a_offset = 0;
387 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
389 const gmx_molblock_t *molb = &mtop->molblock[mb];
390 const gmx_moltype_t *molt = &mtop->moltype[molb->type];
391 const t_block *cgs = &molt->cgs;
393 snew(at2cg, molt->atoms.nr);
394 for (cg = 0; cg < cgs->nr; cg++)
396 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
398 at2cg[i] = cg;
402 const t_atom *atom = molt->atoms.atom;
403 for (mol = 0; mol < molb->nmol; mol++)
405 for (j = 0; (j < NBT); j++)
407 const int *ia = molt->ilist[bondtypes[j]].iatoms.data();
408 for (i = 0; (i < molt->ilist[bondtypes[j]].size()); )
410 type = ia[0];
411 ftype = ffparams->functype[type];
412 nra = interaction_function[ftype].nratoms;
414 /* Check whether we have a bond with a shell */
415 aS = -1;
417 switch (bondtypes[j])
419 case F_BONDS:
420 case F_HARMONIC:
421 case F_CUBICBONDS:
422 case F_POLARIZATION:
423 case F_ANHARM_POL:
424 if (atom[ia[1]].ptype == eptShell)
426 aS = ia[1];
427 aN = ia[2];
429 else if (atom[ia[2]].ptype == eptShell)
431 aS = ia[2];
432 aN = ia[1];
434 break;
435 case F_WATER_POL:
436 aN = ia[4]; /* Dummy */
437 aS = ia[5]; /* Shell */
438 break;
439 default:
440 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
443 if (aS != -1)
445 qS = atom[aS].q;
447 /* Check whether one of the particles is a shell... */
448 nsi = shell_index[a_offset+aS];
449 if ((nsi < 0) || (nsi >= nshell))
451 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
452 nsi, nshell, aS);
454 if (shell[nsi].shell == -1)
456 shell[nsi].shell = a_offset + aS;
457 ns++;
459 else if (shell[nsi].shell != a_offset+aS)
461 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
464 if (shell[nsi].nucl1 == -1)
466 shell[nsi].nucl1 = a_offset + aN;
468 else if (shell[nsi].nucl2 == -1)
470 shell[nsi].nucl2 = a_offset + aN;
472 else if (shell[nsi].nucl3 == -1)
474 shell[nsi].nucl3 = a_offset + aN;
476 else
478 if (fplog)
480 pr_shell(fplog, ns, shell);
482 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
484 if (at2cg[aS] != at2cg[aN])
486 /* shell[nsi].bInterCG = TRUE; */
487 shfc->bInterCG = TRUE;
490 switch (bondtypes[j])
492 case F_BONDS:
493 case F_HARMONIC:
494 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
495 break;
496 case F_CUBICBONDS:
497 shell[nsi].k += ffparams->iparams[type].cubic.kb;
498 break;
499 case F_POLARIZATION:
500 case F_ANHARM_POL:
501 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
503 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);
505 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/
506 ffparams->iparams[type].polarize.alpha;
507 break;
508 case F_WATER_POL:
509 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
511 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);
513 alpha = (ffparams->iparams[type].wpol.al_x+
514 ffparams->iparams[type].wpol.al_y+
515 ffparams->iparams[type].wpol.al_z)/3.0;
516 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/alpha;
517 break;
518 default:
519 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
521 shell[nsi].nnucl++;
523 ia += nra+1;
524 i += nra+1;
527 a_offset += molt->atoms.nr;
529 /* Done with this molecule type */
530 sfree(at2cg);
533 /* Verify whether it's all correct */
534 if (ns != nshell)
536 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
539 for (i = 0; (i < ns); i++)
541 shell[i].k_1 = 1.0/shell[i].k;
544 if (debug)
546 pr_shell(debug, ns, shell);
550 shfc->nshell_gl = ns;
551 shfc->shell_gl = shell;
552 shfc->shell_index_gl = shell_index;
554 shfc->bPredict = (getenv("GMX_NOPREDICT") == nullptr);
555 shfc->bRequireInit = FALSE;
556 if (!shfc->bPredict)
558 if (fplog)
560 fprintf(fplog, "\nWill never predict shell positions\n");
563 else
565 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
566 if (shfc->bRequireInit && fplog)
568 fprintf(fplog, "\nWill always initiate shell positions\n");
572 if (shfc->bPredict)
574 if (shfc->bInterCG)
576 if (fplog)
578 fprintf(fplog, "\nNOTE: there all shells that are connected to particles outside thier own charge group, will not predict shells positions during the run\n\n");
580 /* Prediction improves performance, so we should implement either:
581 * 1. communication for the atoms needed for prediction
582 * 2. prediction using the velocities of shells; currently the
583 * shell velocities are zeroed, it's a bit tricky to keep
584 * track of the shell displacements and thus the velocity.
586 shfc->bPredict = FALSE;
590 return shfc;
593 void make_local_shells(const t_commrec *cr,
594 const t_mdatoms *md,
595 gmx_shellfc_t *shfc)
597 t_shell *shell;
598 int a0, a1, *ind, nshell, i;
599 gmx_domdec_t *dd = nullptr;
601 if (DOMAINDECOMP(cr))
603 dd = cr->dd;
604 a0 = 0;
605 a1 = dd_numHomeAtoms(*dd);
607 else
609 /* Single node: we need all shells, just copy the pointer */
610 shfc->nshell = shfc->nshell_gl;
611 shfc->shell = shfc->shell_gl;
613 return;
616 ind = shfc->shell_index_gl;
618 nshell = 0;
619 shell = shfc->shell;
620 for (i = a0; i < a1; i++)
622 if (md->ptype[i] == eptShell)
624 if (nshell+1 > shfc->shell_nalloc)
626 shfc->shell_nalloc = over_alloc_dd(nshell+1);
627 srenew(shell, shfc->shell_nalloc);
629 if (dd)
631 shell[nshell] = shfc->shell_gl[ind[dd->globalAtomIndices[i]]];
633 else
635 shell[nshell] = shfc->shell_gl[ind[i]];
638 /* With inter-cg shells we can no do shell prediction,
639 * so we do not need the nuclei numbers.
641 if (!shfc->bInterCG)
643 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
644 if (shell[nshell].nnucl > 1)
646 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
648 if (shell[nshell].nnucl > 2)
650 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
653 shell[nshell].shell = i;
654 nshell++;
658 shfc->nshell = nshell;
659 shfc->shell = shell;
662 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
664 real xo, yo, zo;
665 real dx, dy, dz;
667 xo = xold[XX];
668 yo = xold[YY];
669 zo = xold[ZZ];
671 dx = f[XX]*step;
672 dy = f[YY]*step;
673 dz = f[ZZ]*step;
675 xnew[XX] = xo+dx;
676 xnew[YY] = yo+dy;
677 xnew[ZZ] = zo+dz;
680 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
682 real xo, yo, zo;
683 real dx, dy, dz;
685 xo = xold[XX];
686 yo = xold[YY];
687 zo = xold[ZZ];
689 dx = f[XX]*step[XX];
690 dy = f[YY]*step[YY];
691 dz = f[ZZ]*step[ZZ];
693 xnew[XX] = xo+dx;
694 xnew[YY] = yo+dy;
695 xnew[ZZ] = zo+dz;
698 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
699 gmx::ArrayRef<gmx::RVec> xnew,
700 const rvec acc_dir[], int homenr, real step)
702 const rvec *xo = as_rvec_array(xold.data());
703 rvec *xn = as_rvec_array(xnew.data());
705 for (int i = 0; i < homenr; i++)
707 do_1pos(xn[i], xo[i], acc_dir[i], step);
711 static void shell_pos_sd(gmx::ArrayRef<const gmx::RVec> xcur,
712 gmx::ArrayRef<gmx::RVec> xnew,
713 gmx::ArrayRef<const gmx::RVec> f,
714 int ns, t_shell s[], int count)
716 const real step_scale_min = 0.8,
717 step_scale_increment = 0.2,
718 step_scale_max = 1.2,
719 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
720 int i, shell, d;
721 real dx, df, k_est;
722 const real zero = 0;
723 #ifdef PRINT_STEP
724 real step_min, step_max;
726 step_min = 1e30;
727 step_max = 0;
728 #endif
729 for (i = 0; (i < ns); i++)
731 shell = s[i].shell;
732 if (count == 1)
734 for (d = 0; d < DIM; d++)
736 s[i].step[d] = s[i].k_1;
737 #ifdef PRINT_STEP
738 step_min = std::min(step_min, s[i].step[d]);
739 step_max = std::max(step_max, s[i].step[d]);
740 #endif
743 else
745 for (d = 0; d < DIM; d++)
747 dx = xcur[shell][d] - s[i].xold[d];
748 df = f[shell][d] - s[i].fold[d];
749 /* -dx/df gets used to generate an interpolated value, but would
750 * cause a NaN if df were binary-equal to zero. Values close to
751 * zero won't cause problems (because of the min() and max()), so
752 * just testing for binary inequality is OK. */
753 if (zero != df)
755 k_est = -dx/df;
756 /* Scale the step size by a factor interpolated from
757 * step_scale_min to step_scale_max, as k_est goes from 0 to
758 * step_scale_multiple * s[i].step[d] */
759 s[i].step[d] =
760 step_scale_min * s[i].step[d] +
761 step_scale_increment * std::min(step_scale_multiple * s[i].step[d], std::max(k_est, zero));
763 else
765 /* Here 0 == df */
766 if (gmx_numzero(dx)) /* 0 == dx */
768 /* Likely this will never happen, but if it does just
769 * don't scale the step. */
771 else /* 0 != dx */
773 s[i].step[d] *= step_scale_max;
776 #ifdef PRINT_STEP
777 step_min = std::min(step_min, s[i].step[d]);
778 step_max = std::max(step_max, s[i].step[d]);
779 #endif
782 copy_rvec(xcur [shell], s[i].xold);
783 copy_rvec(f[shell], s[i].fold);
785 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
787 if (gmx_debug_at)
789 fprintf(debug, "shell[%d] = %d\n", i, shell);
790 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
791 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
792 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
793 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
796 #ifdef PRINT_STEP
797 printf("step %.3e %.3e\n", step_min, step_max);
798 #endif
801 static void decrease_step_size(int nshell, t_shell s[])
803 int i;
805 for (i = 0; i < nshell; i++)
807 svmul(0.8, s[i].step, s[i].step);
811 static void print_epot(FILE *fp, int64_t mdstep, int count, real epot, real df,
812 int ndir, real sf_dir)
814 char buf[22];
816 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
817 gmx_step_str(mdstep, buf), count, epot, df);
818 if (ndir)
820 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir/ndir));
822 else
824 fprintf(fp, "\n");
829 static real rms_force(const t_commrec *cr, gmx::ArrayRef<const gmx::RVec> force, int ns, t_shell s[],
830 int ndir, real *sf_dir, real *Epot)
832 double buf[4];
833 const rvec *f = as_rvec_array(force.data());
835 buf[0] = *sf_dir;
836 for (int i = 0; i < ns; i++)
838 int shell = s[i].shell;
839 buf[0] += norm2(f[shell]);
841 int ntot = ns;
843 if (PAR(cr))
845 buf[1] = ntot;
846 buf[2] = *sf_dir;
847 buf[3] = *Epot;
848 gmx_sumd(4, buf, cr);
849 ntot = gmx::roundToInt(buf[1]);
850 *sf_dir = buf[2];
851 *Epot = buf[3];
853 ntot += ndir;
855 return (ntot ? std::sqrt(buf[0]/ntot) : 0);
858 static void check_pbc(FILE *fp, gmx::ArrayRef<gmx::RVec> x, int shell)
860 int m, now;
862 now = shell-4;
863 for (m = 0; (m < DIM); m++)
865 if (std::fabs(x[shell][m]-x[now][m]) > 0.3)
867 pr_rvecs(fp, 0, "SHELL-X", as_rvec_array(x.data())+now, 5);
868 break;
873 static void dump_shells(FILE *fp, gmx::ArrayRef<gmx::RVec> x, gmx::ArrayRef<gmx::RVec> f, real ftol, int ns, t_shell s[])
875 int i, shell;
876 real ft2, ff2;
878 ft2 = gmx::square(ftol);
880 for (i = 0; (i < ns); i++)
882 shell = s[i].shell;
883 ff2 = iprod(f[shell], f[shell]);
884 if (ff2 > ft2)
886 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
887 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
889 check_pbc(fp, x, shell);
893 static void init_adir(gmx_shellfc_t *shfc,
894 gmx::Constraints *constr,
895 const t_inputrec *ir,
896 const t_commrec *cr,
897 int dd_ac1,
898 int64_t step,
899 const t_mdatoms *md,
900 int end,
901 rvec *x_old,
902 rvec *x_init,
903 rvec *x,
904 rvec *f,
905 rvec *acc_dir,
906 matrix box,
907 gmx::ArrayRef<const real> lambda,
908 real *dvdlambda)
910 rvec *xnold, *xnew;
911 double dt, w_dt;
912 int n, d;
913 unsigned short *ptype;
915 if (DOMAINDECOMP(cr))
917 n = dd_ac1;
919 else
921 n = end;
923 if (n > shfc->adir_nalloc)
925 shfc->adir_nalloc = over_alloc_dd(n);
926 srenew(shfc->adir_xnold, shfc->adir_nalloc);
927 srenew(shfc->adir_xnew, shfc->adir_nalloc);
929 xnold = shfc->adir_xnold;
930 xnew = shfc->adir_xnew;
932 ptype = md->ptype;
934 dt = ir->delta_t;
936 /* Does NOT work with freeze or acceleration groups (yet) */
937 for (n = 0; n < end; n++)
939 w_dt = md->invmass[n]*dt;
941 for (d = 0; d < DIM; d++)
943 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
945 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
946 xnew[n][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
948 else
950 xnold[n][d] = x[n][d];
951 xnew[n][d] = x[n][d];
955 constr->apply(FALSE, FALSE, step, 0, 1.0,
956 x, xnold, nullptr, box,
957 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
958 nullptr, nullptr, gmx::ConstraintVariable::Positions);
959 constr->apply(FALSE, FALSE, step, 0, 1.0,
960 x, xnew, nullptr, box,
961 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
962 nullptr, nullptr, gmx::ConstraintVariable::Positions);
964 for (n = 0; n < end; n++)
966 for (d = 0; d < DIM; d++)
968 xnew[n][d] =
969 -(2*x[n][d]-xnold[n][d]-xnew[n][d])/gmx::square(dt)
970 - f[n][d]*md->invmass[n];
972 clear_rvec(acc_dir[n]);
975 /* Project the acceleration on the old bond directions */
976 constr->apply(FALSE, FALSE, step, 0, 1.0,
977 x_old, xnew, acc_dir, box,
978 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
979 nullptr, nullptr, gmx::ConstraintVariable::Deriv_FlexCon);
982 void relax_shell_flexcon(FILE *fplog,
983 const t_commrec *cr,
984 const gmx_multisim_t *ms,
985 gmx_bool bVerbose,
986 gmx_enfrot *enforcedRotation,
987 int64_t mdstep,
988 const t_inputrec *inputrec,
989 gmx_bool bDoNS,
990 int force_flags,
991 gmx_localtop_t *top,
992 gmx::Constraints *constr,
993 gmx_enerdata_t *enerd,
994 t_fcdata *fcd,
995 t_state *state,
996 gmx::ArrayRefWithPadding<gmx::RVec> f,
997 tensor force_vir,
998 const t_mdatoms *md,
999 t_nrnb *nrnb,
1000 gmx_wallcycle_t wcycle,
1001 t_graph *graph,
1002 const gmx_groups_t *groups,
1003 gmx_shellfc_t *shfc,
1004 t_forcerec *fr,
1005 gmx::PpForceWorkload *ppForceWorkload,
1006 double t,
1007 rvec mu_tot,
1008 const gmx_vsite_t *vsite,
1009 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion,
1010 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion)
1012 int nshell;
1013 t_shell *shell;
1014 const t_idef *idef;
1015 rvec *acc_dir = nullptr, *x_old = nullptr;
1016 real Epot[2], df[2];
1017 real sf_dir, invdt;
1018 real ftol, dum = 0;
1019 char sbuf[22];
1020 gmx_bool bCont, bInit, bConverged;
1021 int nat, dd_ac0, dd_ac1 = 0, i;
1022 int homenr = md->homenr, end = homenr, cg0, cg1;
1023 int nflexcon, number_steps, d, Min = 0, count = 0;
1024 #define Try (1-Min) /* At start Try = 1 */
1026 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1027 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
1028 ftol = inputrec->em_tol;
1029 number_steps = inputrec->niter;
1030 nshell = shfc->nshell;
1031 shell = shfc->shell;
1032 nflexcon = shfc->nflexcon;
1034 idef = &top->idef;
1036 if (DOMAINDECOMP(cr))
1038 nat = dd_natoms_vsite(cr->dd);
1039 if (nflexcon > 0)
1041 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1042 nat = std::max(nat, dd_ac1);
1045 else
1047 nat = state->natoms;
1050 for (i = 0; (i < 2); i++)
1052 shfc->x[i].resizeWithPadding(nat);
1053 shfc->f[i].resizeWithPadding(nat);
1056 /* Create views that we can swap */
1057 gmx::ArrayRefWithPadding<gmx::RVec> posWithPadding[2];
1058 gmx::ArrayRefWithPadding<gmx::RVec> forceWithPadding[2];
1059 gmx::ArrayRef<gmx::RVec> pos[2];
1060 gmx::ArrayRef<gmx::RVec> force[2];
1061 for (i = 0; (i < 2); i++)
1063 posWithPadding[i] = shfc->x[i].arrayRefWithPadding();
1064 pos[i] = posWithPadding[i].paddedArrayRef();
1065 forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1066 force[i] = forceWithPadding[i].paddedArrayRef();
1069 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1071 /* This is the only time where the coordinates are used
1072 * before do_force is called, which normally puts all
1073 * charge groups in the box.
1075 if (inputrec->cutoff_scheme == ecutsVERLET)
1077 auto xRef = state->x.arrayRefWithPadding().paddedArrayRef();
1078 put_atoms_in_box_omp(fr->ePBC, state->box, xRef.subArray(0, md->homenr));
1080 else
1082 cg0 = 0;
1083 cg1 = top->cgs.nr;
1084 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1085 &(top->cgs), state->x.rvec_array(), fr->cg_cm);
1088 if (graph)
1090 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x.rvec_array());
1094 /* After this all coordinate arrays will contain whole charge groups */
1095 if (graph)
1097 shift_self(graph, state->box, state->x.rvec_array());
1100 if (nflexcon)
1102 if (nat > shfc->flex_nalloc)
1104 shfc->flex_nalloc = over_alloc_dd(nat);
1105 srenew(shfc->acc_dir, shfc->flex_nalloc);
1106 srenew(shfc->x_old, shfc->flex_nalloc);
1108 acc_dir = shfc->acc_dir;
1109 x_old = shfc->x_old;
1110 auto x = makeArrayRef(state->x);
1111 auto v = makeArrayRef(state->v);
1112 for (i = 0; i < homenr; i++)
1114 for (d = 0; d < DIM; d++)
1116 shfc->x_old[i][d] =
1117 x[i][d] - v[i][d]*inputrec->delta_t;
1122 /* Do a prediction of the shell positions, when appropriate.
1123 * Without velocities (EM, NM, BD) we only do initial prediction.
1125 if (shfc->bPredict && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1127 predict_shells(fplog, state->x.rvec_array(), state->v.rvec_array(), inputrec->delta_t, nshell, shell,
1128 md->massT, nullptr, bInit);
1131 /* do_force expected the charge groups to be in the box */
1132 if (graph)
1134 unshift_self(graph, state->box, state->x.rvec_array());
1137 /* Calculate the forces first time around */
1138 if (gmx_debug_at)
1140 pr_rvecs(debug, 0, "x b4 do_force", state->x.rvec_array(), homenr);
1142 int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1143 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation,
1144 mdstep, nrnb, wcycle, top, groups,
1145 state->box, state->x.arrayRefWithPadding(), &state->hist,
1146 forceWithPadding[Min], force_vir, md, enerd, fcd,
1147 state->lambda, graph,
1148 fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
1149 (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags,
1150 ddOpenBalanceRegion, ddCloseBalanceRegion);
1152 sf_dir = 0;
1153 if (nflexcon)
1155 init_adir(shfc,
1156 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1157 shfc->x_old, state->x.rvec_array(), state->x.rvec_array(), as_rvec_array(force[Min].data()),
1158 shfc->acc_dir,
1159 state->box, state->lambda, &dum);
1161 for (i = 0; i < end; i++)
1163 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i]);
1166 sum_epot(&(enerd->grpp), enerd->term);
1167 Epot[Min] = enerd->term[F_EPOT];
1169 df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1170 df[Try] = 0;
1171 if (debug)
1173 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1176 if (gmx_debug_at)
1178 pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md->nr);
1181 if (nshell+nflexcon > 0)
1183 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1184 * shell positions are updated, therefore the other particles must
1185 * be set here, in advance.
1187 std::copy(state->x.begin(),
1188 state->x.end(),
1189 posWithPadding[Min].paddedArrayRef().begin());
1190 std::copy(state->x.begin(),
1191 state->x.end(),
1192 posWithPadding[Try].paddedArrayRef().begin());
1195 if (bVerbose && MASTER(cr))
1197 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1200 if (debug)
1202 fprintf(debug, "%17s: %14.10e\n",
1203 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1204 fprintf(debug, "%17s: %14.10e\n",
1205 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1206 fprintf(debug, "%17s: %14.10e\n",
1207 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1208 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1211 /* First check whether we should do shells, or whether the force is
1212 * low enough even without minimization.
1214 bConverged = (df[Min] < ftol);
1216 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1218 if (vsite)
1220 construct_vsites(vsite, as_rvec_array(pos[Min].data()),
1221 inputrec->delta_t, state->v.rvec_array(),
1222 idef->iparams, idef->il,
1223 fr->ePBC, fr->bMolPBC, cr, state->box);
1226 if (nflexcon)
1228 init_adir(shfc,
1229 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1230 x_old, state->x.rvec_array(),
1231 as_rvec_array(pos[Min].data()),
1232 as_rvec_array(force[Min].data()), acc_dir,
1233 state->box, state->lambda, &dum);
1235 directional_sd(pos[Min], pos[Try], acc_dir, end, fr->fc_stepsize);
1238 /* New positions, Steepest descent */
1239 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1241 /* do_force expected the charge groups to be in the box */
1242 if (graph)
1244 unshift_self(graph, state->box, as_rvec_array(pos[Try].data()));
1247 if (gmx_debug_at)
1249 pr_rvecs(debug, 0, "RELAX: pos[Min] ", as_rvec_array(pos[Min].data()), homenr);
1250 pr_rvecs(debug, 0, "RELAX: pos[Try] ", as_rvec_array(pos[Try].data()), homenr);
1252 /* Try the new positions */
1253 do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation,
1254 1, nrnb, wcycle,
1255 top, groups, state->box, posWithPadding[Try], &state->hist,
1256 forceWithPadding[Try], force_vir,
1257 md, enerd, fcd, state->lambda, graph,
1258 fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
1259 shellfc_flags,
1260 ddOpenBalanceRegion, ddCloseBalanceRegion);
1261 sum_epot(&(enerd->grpp), enerd->term);
1262 if (gmx_debug_at)
1264 pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1265 pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1267 sf_dir = 0;
1268 if (nflexcon)
1270 init_adir(shfc,
1271 constr, inputrec, cr, dd_ac1, mdstep, md, end,
1272 x_old, state->x.rvec_array(),
1273 as_rvec_array(pos[Try].data()),
1274 as_rvec_array(force[Try].data()),
1275 acc_dir, state->box, state->lambda, &dum);
1277 for (i = 0; i < end; i++)
1279 sf_dir += md->massT[i]*norm2(acc_dir[i]);
1283 Epot[Try] = enerd->term[F_EPOT];
1285 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1287 if (debug)
1289 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1292 if (debug)
1294 if (gmx_debug_at)
1296 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1298 if (gmx_debug_at)
1300 fprintf(debug, "SHELL ITER %d\n", count);
1301 dump_shells(debug, pos[Try], force[Try], ftol, nshell, shell);
1305 if (bVerbose && MASTER(cr))
1307 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1310 bConverged = (df[Try] < ftol);
1312 if ((df[Try] < df[Min]))
1314 if (debug)
1316 fprintf(debug, "Swapping Min and Try\n");
1318 if (nflexcon)
1320 /* Correct the velocities for the flexible constraints */
1321 invdt = 1/inputrec->delta_t;
1322 auto v = makeArrayRef(state->v);
1323 for (i = 0; i < end; i++)
1325 for (d = 0; d < DIM; d++)
1327 v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1331 Min = Try;
1333 else
1335 decrease_step_size(nshell, shell);
1338 shfc->numForceEvaluations += count;
1339 if (bConverged)
1341 shfc->numConvergedIterations++;
1343 if (MASTER(cr) && !(bConverged))
1345 /* Note that the energies and virial are incorrect when not converged */
1346 if (fplog)
1348 fprintf(fplog,
1349 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1350 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1352 fprintf(stderr,
1353 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1354 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1357 /* Copy back the coordinates and the forces */
1358 std::copy(pos[Min].begin(), pos[Min].end(), makeArrayRef(state->x).data());
1359 std::copy(force[Min].begin(), force[Min].end(), f.unpaddedArrayRef().begin());
1362 void done_shellfc(FILE *fplog, gmx_shellfc_t *shfc, int64_t numSteps)
1364 if (shfc && fplog && numSteps > 0)
1366 double numStepsAsDouble = static_cast<double>(numSteps);
1367 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1368 (shfc->numConvergedIterations*100.0)/numStepsAsDouble);
1369 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1370 shfc->numForceEvaluations/numStepsAsDouble);
1373 // TODO Deallocate memory in shfc