Add check for shells and GPU update
[alexxy/gromacs.git] / src / gromacs / mdrun / shellfc.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2008, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "shellfc.h"
41
42 #include <cmath>
43 #include <cstdint>
44 #include <cstdlib>
45 #include <cstring>
46
47 #include <algorithm>
48 #include <array>
49
50 #include "gromacs/domdec/dlbtiming.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/mdsetup.h"
54 #include "gromacs/gmxlib/network.h"
55 #include "gromacs/math/functions.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/math/vecdump.h"
59 #include "gromacs/mdlib/constr.h"
60 #include "gromacs/mdlib/enerdata_utils.h"
61 #include "gromacs/mdlib/force.h"
62 #include "gromacs/mdlib/force_flags.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdlib/mdatoms.h"
65 #include "gromacs/mdlib/vsite.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/enerdata.h"
68 #include "gromacs/mdtypes/forcerec.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/mdtypes/state.h"
72 #include "gromacs/pbcutil/mshift.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/topology/ifunc.h"
75 #include "gromacs/topology/mtop_lookup.h"
76 #include "gromacs/topology/mtop_util.h"
77 #include "gromacs/utility/arrayref.h"
78 #include "gromacs/utility/arraysize.h"
79 #include "gromacs/utility/cstringutil.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/gmxassert.h"
82 #include "gromacs/utility/smalloc.h"
83
84 typedef struct
85 {
86     int nnucl;
87     int shell;                  /* The shell id                         */
88     int nucl1, nucl2, nucl3;    /* The nuclei connected to the shell    */
89     /* gmx_bool    bInterCG; */ /* Coupled to nuclei outside cg?        */
90     real k;                     /* force constant                       */
91     real k_1;                   /* 1 over force constant                */
92     rvec xold;
93     rvec fold;
94     rvec step;
95 } t_shell;
96
97 struct gmx_shellfc_t
98 {
99     /* Shell counts, indices, parameters and working data */
100     int      nshell_gl;      /* The number of shells in the system        */
101     t_shell* shell_gl;       /* All the shells (for DD only)              */
102     int*     shell_index_gl; /* Global shell index (for DD only)          */
103     gmx_bool bInterCG;       /* Are there inter charge-group shells?      */
104     int      nshell;         /* The number of local shells                */
105     t_shell* shell;          /* The local shells                          */
106     int      shell_nalloc;   /* The allocation size of shell              */
107     gmx_bool bPredict;       /* Predict shell positions                   */
108     gmx_bool bRequireInit;   /* Require initialization of shell positions */
109     int      nflexcon;       /* The number of flexible constraints        */
110
111     /* Temporary arrays, should be fixed size 2 when fully converted to C++ */
112     PaddedHostVector<gmx::RVec>* x; /* Array for iterative minimization          */
113     PaddedHostVector<gmx::RVec>* f; /* Array for iterative minimization          */
114
115     /* Flexible constraint working data */
116     rvec*        acc_dir;                /* Acceleration direction for flexcon        */
117     rvec*        x_old;                  /* Old coordinates for flexcon               */
118     int          flex_nalloc;            /* The allocation size of acc_dir and x_old  */
119     rvec*        adir_xnold;             /* Work space for init_adir                  */
120     rvec*        adir_xnew;              /* Work space for init_adir                  */
121     int          adir_nalloc;            /* Work space for init_adir                  */
122     std::int64_t numForceEvaluations;    /* Total number of force evaluations         */
123     int          numConvergedIterations; /* Total number of iterations that converged */
124 };
125
126
127 static void pr_shell(FILE* fplog, int ns, t_shell s[])
128 {
129     int i;
130
131     fprintf(fplog, "SHELL DATA\n");
132     fprintf(fplog, "%5s  %8s  %5s  %5s  %5s\n", "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
133     for (i = 0; (i < ns); i++)
134     {
135         fprintf(fplog, "%5d  %8.3f  %5d", s[i].shell, 1.0 / s[i].k_1, s[i].nucl1);
136         if (s[i].nnucl == 2)
137         {
138             fprintf(fplog, "  %5d\n", s[i].nucl2);
139         }
140         else if (s[i].nnucl == 3)
141         {
142             fprintf(fplog, "  %5d  %5d\n", s[i].nucl2, s[i].nucl3);
143         }
144         else
145         {
146             fprintf(fplog, "\n");
147         }
148     }
149 }
150
151 /* TODO The remain call of this function passes non-NULL mass and NULL
152  * mtop, so this routine can be simplified.
153  *
154  * The other code path supported doing prediction before the MD loop
155  * started, but even when called, the prediction was always
156  * over-written by a subsequent call in the MD loop, so has been
157  * removed. */
158 static void predict_shells(FILE*       fplog,
159                            rvec        x[],
160                            rvec        v[],
161                            real        dt,
162                            int         ns,
163                            t_shell     s[],
164                            const real  mass[],
165                            gmx_mtop_t* mtop,
166                            gmx_bool    bInit)
167 {
168     int   i, m, s1, n1, n2, n3;
169     real  dt_1, fudge, tm, m1, m2, m3;
170     rvec* ptr;
171
172     GMX_RELEASE_ASSERT(mass || mtop, "Must have masses or a way to look them up");
173
174     /* We introduce a fudge factor for performance reasons: with this choice
175      * the initial force on the shells is about a factor of two lower than
176      * without
177      */
178     fudge = 1.0;
179
180     if (bInit)
181     {
182         if (fplog)
183         {
184             fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
185         }
186         ptr  = x;
187         dt_1 = 1;
188     }
189     else
190     {
191         ptr  = v;
192         dt_1 = fudge * dt;
193     }
194
195     int molb = 0;
196     for (i = 0; (i < ns); i++)
197     {
198         s1 = s[i].shell;
199         if (bInit)
200         {
201             clear_rvec(x[s1]);
202         }
203         switch (s[i].nnucl)
204         {
205             case 1:
206                 n1 = s[i].nucl1;
207                 for (m = 0; (m < DIM); m++)
208                 {
209                     x[s1][m] += ptr[n1][m] * dt_1;
210                 }
211                 break;
212             case 2:
213                 n1 = s[i].nucl1;
214                 n2 = s[i].nucl2;
215                 if (mass)
216                 {
217                     m1 = mass[n1];
218                     m2 = mass[n2];
219                 }
220                 else
221                 {
222                     /* Not the correct masses with FE, but it is just a prediction... */
223                     m1 = mtopGetAtomMass(mtop, n1, &molb);
224                     m2 = mtopGetAtomMass(mtop, n2, &molb);
225                 }
226                 tm = dt_1 / (m1 + m2);
227                 for (m = 0; (m < DIM); m++)
228                 {
229                     x[s1][m] += (m1 * ptr[n1][m] + m2 * ptr[n2][m]) * tm;
230                 }
231                 break;
232             case 3:
233                 n1 = s[i].nucl1;
234                 n2 = s[i].nucl2;
235                 n3 = s[i].nucl3;
236                 if (mass)
237                 {
238                     m1 = mass[n1];
239                     m2 = mass[n2];
240                     m3 = mass[n3];
241                 }
242                 else
243                 {
244                     /* Not the correct masses with FE, but it is just a prediction... */
245                     m1 = mtopGetAtomMass(mtop, n1, &molb);
246                     m2 = mtopGetAtomMass(mtop, n2, &molb);
247                     m3 = mtopGetAtomMass(mtop, n3, &molb);
248                 }
249                 tm = dt_1 / (m1 + m2 + m3);
250                 for (m = 0; (m < DIM); m++)
251                 {
252                     x[s1][m] += (m1 * ptr[n1][m] + m2 * ptr[n2][m] + m3 * ptr[n3][m]) * tm;
253                 }
254                 break;
255             default: gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
256         }
257     }
258 }
259
260 gmx_shellfc_t* init_shell_flexcon(FILE* fplog, const gmx_mtop_t* mtop, int nflexcon, int nstcalcenergy, bool usingDomainDecomposition)
261 {
262     gmx_shellfc_t* shfc;
263     t_shell*       shell;
264     int*           shell_index = nullptr;
265
266     int  ns, nshell, nsi;
267     int  i, j, type, a_offset, mol, ftype, nra;
268     real qS, alpha;
269     int  aS, aN = 0; /* Shell and nucleus */
270     int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
271 #define NBT asize(bondtypes)
272     const gmx_ffparams_t* ffparams;
273
274     const std::array<int, eptNR> numParticles = gmx_mtop_particletype_count(*mtop);
275     if (fplog)
276     {
277         /* Print the number of each particle type */
278         int pType = 0;
279         for (const auto& n : numParticles)
280         {
281             if (n != 0)
282             {
283                 fprintf(fplog, "There are: %d %ss\n", n, ptype_str[pType]);
284             }
285             pType++;
286         }
287     }
288
289     nshell = numParticles[eptShell];
290
291     if (nshell == 0 && nflexcon == 0)
292     {
293         /* We're not doing shells or flexible constraints */
294         return nullptr;
295     }
296
297     snew(shfc, 1);
298     shfc->x        = new PaddedHostVector<gmx::RVec>[2] {};
299     shfc->f        = new PaddedHostVector<gmx::RVec>[2] {};
300     shfc->nflexcon = nflexcon;
301
302     if (nshell == 0)
303     {
304         /* Only flexible constraints, no shells.
305          * Note that make_local_shells() does not need to be called.
306          */
307         shfc->nshell   = 0;
308         shfc->bPredict = FALSE;
309
310         return shfc;
311     }
312
313     if (nstcalcenergy != 1)
314     {
315         gmx_fatal(FARGS,
316                   "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is "
317                   "not supported in combination with shell particles.\nPlease make a new tpr file.",
318                   nstcalcenergy);
319     }
320     if (usingDomainDecomposition)
321     {
322         gmx_fatal(
323                 FARGS,
324                 "Shell particles are not implemented with domain decomposition, use a single rank");
325     }
326
327     /* We have shells: fill the shell data structure */
328
329     /* Global system sized array, this should be avoided */
330     snew(shell_index, mtop->natoms);
331
332     nshell = 0;
333     for (const AtomProxy atomP : AtomRange(*mtop))
334     {
335         const t_atom& local = atomP.atom();
336         int           i     = atomP.globalAtomNumber();
337         if (local.ptype == eptShell)
338         {
339             shell_index[i] = nshell++;
340         }
341     }
342
343     snew(shell, nshell);
344
345     /* Initiate the shell structures */
346     for (i = 0; (i < nshell); i++)
347     {
348         shell[i].shell = -1;
349         shell[i].nnucl = 0;
350         shell[i].nucl1 = -1;
351         shell[i].nucl2 = -1;
352         shell[i].nucl3 = -1;
353         /* shell[i].bInterCG=FALSE; */
354         shell[i].k_1 = 0;
355         shell[i].k   = 0;
356     }
357
358     ffparams = &mtop->ffparams;
359
360     /* Now fill the structures */
361     /* TODO: See if we can use update groups that cover shell constructions */
362     shfc->bInterCG = FALSE;
363     ns             = 0;
364     a_offset       = 0;
365     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
366     {
367         const gmx_molblock_t* molb = &mtop->molblock[mb];
368         const gmx_moltype_t*  molt = &mtop->moltype[molb->type];
369
370         const t_atom* atom = molt->atoms.atom;
371         for (mol = 0; mol < molb->nmol; mol++)
372         {
373             for (j = 0; (j < NBT); j++)
374             {
375                 const int* ia = molt->ilist[bondtypes[j]].iatoms.data();
376                 for (i = 0; (i < molt->ilist[bondtypes[j]].size());)
377                 {
378                     type  = ia[0];
379                     ftype = ffparams->functype[type];
380                     nra   = interaction_function[ftype].nratoms;
381
382                     /* Check whether we have a bond with a shell */
383                     aS = -1;
384
385                     switch (bondtypes[j])
386                     {
387                         case F_BONDS:
388                         case F_HARMONIC:
389                         case F_CUBICBONDS:
390                         case F_POLARIZATION:
391                         case F_ANHARM_POL:
392                             if (atom[ia[1]].ptype == eptShell)
393                             {
394                                 aS = ia[1];
395                                 aN = ia[2];
396                             }
397                             else if (atom[ia[2]].ptype == eptShell)
398                             {
399                                 aS = ia[2];
400                                 aN = ia[1];
401                             }
402                             break;
403                         case F_WATER_POL:
404                             aN = ia[4]; /* Dummy */
405                             aS = ia[5]; /* Shell */
406                             break;
407                         default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
408                     }
409
410                     if (aS != -1)
411                     {
412                         qS = atom[aS].q;
413
414                         /* Check whether one of the particles is a shell... */
415                         nsi = shell_index[a_offset + aS];
416                         if ((nsi < 0) || (nsi >= nshell))
417                         {
418                             gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d", nsi, nshell, aS);
419                         }
420                         if (shell[nsi].shell == -1)
421                         {
422                             shell[nsi].shell = a_offset + aS;
423                             ns++;
424                         }
425                         else if (shell[nsi].shell != a_offset + aS)
426                         {
427                             gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
428                         }
429
430                         if (shell[nsi].nucl1 == -1)
431                         {
432                             shell[nsi].nucl1 = a_offset + aN;
433                         }
434                         else if (shell[nsi].nucl2 == -1)
435                         {
436                             shell[nsi].nucl2 = a_offset + aN;
437                         }
438                         else if (shell[nsi].nucl3 == -1)
439                         {
440                             shell[nsi].nucl3 = a_offset + aN;
441                         }
442                         else
443                         {
444                             if (fplog)
445                             {
446                                 pr_shell(fplog, ns, shell);
447                             }
448                             gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
449                         }
450                         if (aS != aN)
451                         {
452                             /* shell[nsi].bInterCG = TRUE; */
453                             shfc->bInterCG = TRUE;
454                         }
455
456                         switch (bondtypes[j])
457                         {
458                             case F_BONDS:
459                             case F_HARMONIC:
460                                 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
461                                 break;
462                             case F_CUBICBONDS:
463                                 shell[nsi].k += ffparams->iparams[type].cubic.kb;
464                                 break;
465                             case F_POLARIZATION:
466                             case F_ANHARM_POL:
467                                 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS * 10))
468                                 {
469                                     gmx_fatal(FARGS,
470                                               "polarize can not be used with qA(%e) != qB(%e) for "
471                                               "atom %d of molecule block %zu",
472                                               qS, atom[aS].qB, aS + 1, mb + 1);
473                                 }
474                                 shell[nsi].k += gmx::square(qS) * ONE_4PI_EPS0
475                                                 / ffparams->iparams[type].polarize.alpha;
476                                 break;
477                             case F_WATER_POL:
478                                 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS * 10))
479                                 {
480                                     gmx_fatal(FARGS,
481                                               "water_pol can not be used with qA(%e) != qB(%e) for "
482                                               "atom %d of molecule block %zu",
483                                               qS, atom[aS].qB, aS + 1, mb + 1);
484                                 }
485                                 alpha = (ffparams->iparams[type].wpol.al_x
486                                          + ffparams->iparams[type].wpol.al_y
487                                          + ffparams->iparams[type].wpol.al_z)
488                                         / 3.0;
489                                 shell[nsi].k += gmx::square(qS) * ONE_4PI_EPS0 / alpha;
490                                 break;
491                             default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
492                         }
493                         shell[nsi].nnucl++;
494                     }
495                     ia += nra + 1;
496                     i += nra + 1;
497                 }
498             }
499             a_offset += molt->atoms.nr;
500         }
501         /* Done with this molecule type */
502     }
503
504     /* Verify whether it's all correct */
505     if (ns != nshell)
506     {
507         gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
508     }
509
510     for (i = 0; (i < ns); i++)
511     {
512         shell[i].k_1 = 1.0 / shell[i].k;
513     }
514
515     if (debug)
516     {
517         pr_shell(debug, ns, shell);
518     }
519
520
521     shfc->nshell_gl      = ns;
522     shfc->shell_gl       = shell;
523     shfc->shell_index_gl = shell_index;
524
525     shfc->bPredict     = (getenv("GMX_NOPREDICT") == nullptr);
526     shfc->bRequireInit = FALSE;
527     if (!shfc->bPredict)
528     {
529         if (fplog)
530         {
531             fprintf(fplog, "\nWill never predict shell positions\n");
532         }
533     }
534     else
535     {
536         shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
537         if (shfc->bRequireInit && fplog)
538         {
539             fprintf(fplog, "\nWill always initiate shell positions\n");
540         }
541     }
542
543     if (shfc->bPredict)
544     {
545         if (shfc->bInterCG)
546         {
547             if (fplog)
548             {
549                 fprintf(fplog,
550                         "\nNOTE: in the current version shell prediction during the crun is "
551                         "disabled\n\n");
552             }
553             /* Prediction improves performance, so we should implement either:
554              * 1. communication for the atoms needed for prediction
555              * 2. prediction using the velocities of shells; currently the
556              *    shell velocities are zeroed, it's a bit tricky to keep
557              *    track of the shell displacements and thus the velocity.
558              */
559             shfc->bPredict = FALSE;
560         }
561     }
562
563     return shfc;
564 }
565
566 void make_local_shells(const t_commrec* cr, const t_mdatoms* md, gmx_shellfc_t* shfc)
567 {
568     t_shell*      shell;
569     int           a0, a1, *ind, nshell, i;
570     gmx_domdec_t* dd = nullptr;
571
572     if (DOMAINDECOMP(cr))
573     {
574         dd = cr->dd;
575         a0 = 0;
576         a1 = dd_numHomeAtoms(*dd);
577     }
578     else
579     {
580         /* Single node: we need all shells, just copy the pointer */
581         shfc->nshell = shfc->nshell_gl;
582         shfc->shell  = shfc->shell_gl;
583
584         return;
585     }
586
587     ind = shfc->shell_index_gl;
588
589     nshell = 0;
590     shell  = shfc->shell;
591     for (i = a0; i < a1; i++)
592     {
593         if (md->ptype[i] == eptShell)
594         {
595             if (nshell + 1 > shfc->shell_nalloc)
596             {
597                 shfc->shell_nalloc = over_alloc_dd(nshell + 1);
598                 srenew(shell, shfc->shell_nalloc);
599             }
600             if (dd)
601             {
602                 shell[nshell] = shfc->shell_gl[ind[dd->globalAtomIndices[i]]];
603             }
604             else
605             {
606                 shell[nshell] = shfc->shell_gl[ind[i]];
607             }
608
609             /* With inter-cg shells we can no do shell prediction,
610              * so we do not need the nuclei numbers.
611              */
612             if (!shfc->bInterCG)
613             {
614                 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
615                 if (shell[nshell].nnucl > 1)
616                 {
617                     shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
618                 }
619                 if (shell[nshell].nnucl > 2)
620                 {
621                     shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
622                 }
623             }
624             shell[nshell].shell = i;
625             nshell++;
626         }
627     }
628
629     shfc->nshell = nshell;
630     shfc->shell  = shell;
631 }
632
633 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
634 {
635     real xo, yo, zo;
636     real dx, dy, dz;
637
638     xo = xold[XX];
639     yo = xold[YY];
640     zo = xold[ZZ];
641
642     dx = f[XX] * step;
643     dy = f[YY] * step;
644     dz = f[ZZ] * step;
645
646     xnew[XX] = xo + dx;
647     xnew[YY] = yo + dy;
648     xnew[ZZ] = zo + dz;
649 }
650
651 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
652 {
653     real xo, yo, zo;
654     real dx, dy, dz;
655
656     xo = xold[XX];
657     yo = xold[YY];
658     zo = xold[ZZ];
659
660     dx = f[XX] * step[XX];
661     dy = f[YY] * step[YY];
662     dz = f[ZZ] * step[ZZ];
663
664     xnew[XX] = xo + dx;
665     xnew[YY] = yo + dy;
666     xnew[ZZ] = zo + dz;
667 }
668
669 static void directional_sd(gmx::ArrayRef<const gmx::RVec> xold,
670                            gmx::ArrayRef<gmx::RVec>       xnew,
671                            const rvec                     acc_dir[],
672                            int                            homenr,
673                            real                           step)
674 {
675     const rvec* xo = as_rvec_array(xold.data());
676     rvec*       xn = as_rvec_array(xnew.data());
677
678     for (int i = 0; i < homenr; i++)
679     {
680         do_1pos(xn[i], xo[i], acc_dir[i], step);
681     }
682 }
683
684 static void shell_pos_sd(gmx::ArrayRef<const gmx::RVec> xcur,
685                          gmx::ArrayRef<gmx::RVec>       xnew,
686                          gmx::ArrayRef<const gmx::RVec> f,
687                          int                            ns,
688                          t_shell                        s[],
689                          int                            count)
690 {
691     const real step_scale_min = 0.8, step_scale_increment = 0.2, step_scale_max = 1.2,
692                step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
693     int        i, shell, d;
694     real       dx, df, k_est;
695     const real zero = 0;
696 #ifdef PRINT_STEP
697     real step_min, step_max;
698
699     step_min = 1e30;
700     step_max = 0;
701 #endif
702     for (i = 0; (i < ns); i++)
703     {
704         shell = s[i].shell;
705         if (count == 1)
706         {
707             for (d = 0; d < DIM; d++)
708             {
709                 s[i].step[d] = s[i].k_1;
710 #ifdef PRINT_STEP
711                 step_min = std::min(step_min, s[i].step[d]);
712                 step_max = std::max(step_max, s[i].step[d]);
713 #endif
714             }
715         }
716         else
717         {
718             for (d = 0; d < DIM; d++)
719             {
720                 dx = xcur[shell][d] - s[i].xold[d];
721                 df = f[shell][d] - s[i].fold[d];
722                 /* -dx/df gets used to generate an interpolated value, but would
723                  * cause a NaN if df were binary-equal to zero. Values close to
724                  * zero won't cause problems (because of the min() and max()), so
725                  * just testing for binary inequality is OK. */
726                 if (zero != df)
727                 {
728                     k_est = -dx / df;
729                     /* Scale the step size by a factor interpolated from
730                      * step_scale_min to step_scale_max, as k_est goes from 0 to
731                      * step_scale_multiple * s[i].step[d] */
732                     s[i].step[d] = step_scale_min * s[i].step[d]
733                                    + step_scale_increment
734                                              * std::min(step_scale_multiple * s[i].step[d],
735                                                         std::max(k_est, zero));
736                 }
737                 else
738                 {
739                     /* Here 0 == df */
740                     if (gmx_numzero(dx)) /* 0 == dx */
741                     {
742                         /* Likely this will never happen, but if it does just
743                          * don't scale the step. */
744                     }
745                     else /* 0 != dx */
746                     {
747                         s[i].step[d] *= step_scale_max;
748                     }
749                 }
750 #ifdef PRINT_STEP
751                 step_min = std::min(step_min, s[i].step[d]);
752                 step_max = std::max(step_max, s[i].step[d]);
753 #endif
754             }
755         }
756         copy_rvec(xcur[shell], s[i].xold);
757         copy_rvec(f[shell], s[i].fold);
758
759         do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
760
761         if (gmx_debug_at)
762         {
763             fprintf(debug, "shell[%d] = %d\n", i, shell);
764             pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
765             pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
766             pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
767             pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
768         }
769     }
770 #ifdef PRINT_STEP
771     printf("step %.3e %.3e\n", step_min, step_max);
772 #endif
773 }
774
775 static void decrease_step_size(int nshell, t_shell s[])
776 {
777     int i;
778
779     for (i = 0; i < nshell; i++)
780     {
781         svmul(0.8, s[i].step, s[i].step);
782     }
783 }
784
785 static void print_epot(FILE* fp, int64_t mdstep, int count, real epot, real df, int ndir, real sf_dir)
786 {
787     char buf[22];
788
789     fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e", gmx_step_str(mdstep, buf), count, epot, df);
790     if (ndir)
791     {
792         fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir / ndir));
793     }
794     else
795     {
796         fprintf(fp, "\n");
797     }
798 }
799
800
801 static real rms_force(const t_commrec*               cr,
802                       gmx::ArrayRef<const gmx::RVec> force,
803                       int                            ns,
804                       t_shell                        s[],
805                       int                            ndir,
806                       real*                          sf_dir,
807                       real*                          Epot)
808 {
809     double      buf[4];
810     const rvec* f = as_rvec_array(force.data());
811
812     buf[0] = *sf_dir;
813     for (int i = 0; i < ns; i++)
814     {
815         int shell = s[i].shell;
816         buf[0] += norm2(f[shell]);
817     }
818     int ntot = ns;
819
820     if (PAR(cr))
821     {
822         buf[1] = ntot;
823         buf[2] = *sf_dir;
824         buf[3] = *Epot;
825         gmx_sumd(4, buf, cr);
826         ntot    = gmx::roundToInt(buf[1]);
827         *sf_dir = buf[2];
828         *Epot   = buf[3];
829     }
830     ntot += ndir;
831
832     return (ntot ? std::sqrt(buf[0] / ntot) : 0);
833 }
834
835 static void dump_shells(FILE* fp, gmx::ArrayRef<gmx::RVec> f, real ftol, int ns, t_shell s[])
836 {
837     int  i, shell;
838     real ft2, ff2;
839
840     ft2 = gmx::square(ftol);
841
842     for (i = 0; (i < ns); i++)
843     {
844         shell = s[i].shell;
845         ff2   = iprod(f[shell], f[shell]);
846         if (ff2 > ft2)
847         {
848             fprintf(fp, "SHELL %5d, force %10.5f  %10.5f  %10.5f, |f| %10.5f\n", shell,
849                     f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
850         }
851     }
852 }
853
854 static void init_adir(gmx_shellfc_t*            shfc,
855                       gmx::Constraints*         constr,
856                       const t_inputrec*         ir,
857                       const t_commrec*          cr,
858                       int                       dd_ac1,
859                       int64_t                   step,
860                       const t_mdatoms*          md,
861                       int                       end,
862                       rvec*                     x_old,
863                       rvec*                     x_init,
864                       rvec*                     x,
865                       rvec*                     f,
866                       rvec*                     acc_dir,
867                       const matrix              box,
868                       gmx::ArrayRef<const real> lambda,
869                       real*                     dvdlambda)
870 {
871     rvec *          xnold, *xnew;
872     double          dt, w_dt;
873     int             n, d;
874     unsigned short* ptype;
875
876     if (DOMAINDECOMP(cr))
877     {
878         n = dd_ac1;
879     }
880     else
881     {
882         n = end;
883     }
884     if (n > shfc->adir_nalloc)
885     {
886         shfc->adir_nalloc = over_alloc_dd(n);
887         srenew(shfc->adir_xnold, shfc->adir_nalloc);
888         srenew(shfc->adir_xnew, shfc->adir_nalloc);
889     }
890     xnold = shfc->adir_xnold;
891     xnew  = shfc->adir_xnew;
892
893     ptype = md->ptype;
894
895     dt = ir->delta_t;
896
897     /* Does NOT work with freeze or acceleration groups (yet) */
898     for (n = 0; n < end; n++)
899     {
900         w_dt = md->invmass[n] * dt;
901
902         for (d = 0; d < DIM; d++)
903         {
904             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
905             {
906                 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
907                 xnew[n][d]  = 2 * x[n][d] - x_old[n][d] + f[n][d] * w_dt * dt;
908             }
909             else
910             {
911                 xnold[n][d] = x[n][d];
912                 xnew[n][d]  = x[n][d];
913             }
914         }
915     }
916     constr->apply(FALSE, FALSE, step, 0, 1.0, x, xnold, nullptr, box, lambda[efptBONDED],
917                   &(dvdlambda[efptBONDED]), nullptr, nullptr, gmx::ConstraintVariable::Positions);
918     constr->apply(FALSE, FALSE, step, 0, 1.0, x, xnew, nullptr, box, lambda[efptBONDED],
919                   &(dvdlambda[efptBONDED]), nullptr, nullptr, gmx::ConstraintVariable::Positions);
920
921     for (n = 0; n < end; n++)
922     {
923         for (d = 0; d < DIM; d++)
924         {
925             xnew[n][d] = -(2 * x[n][d] - xnold[n][d] - xnew[n][d]) / gmx::square(dt)
926                          - f[n][d] * md->invmass[n];
927         }
928         clear_rvec(acc_dir[n]);
929     }
930
931     /* Project the acceleration on the old bond directions */
932     constr->apply(FALSE, FALSE, step, 0, 1.0, x_old, xnew, acc_dir, box, lambda[efptBONDED],
933                   &(dvdlambda[efptBONDED]), nullptr, nullptr, gmx::ConstraintVariable::Deriv_FlexCon);
934 }
935
936 void relax_shell_flexcon(FILE*                               fplog,
937                          const t_commrec*                    cr,
938                          const gmx_multisim_t*               ms,
939                          gmx_bool                            bVerbose,
940                          gmx_enfrot*                         enforcedRotation,
941                          int64_t                             mdstep,
942                          const t_inputrec*                   inputrec,
943                          gmx::ImdSession*                    imdSession,
944                          pull_t*                             pull_work,
945                          gmx_bool                            bDoNS,
946                          int                                 force_flags,
947                          const gmx_localtop_t*               top,
948                          gmx::Constraints*                   constr,
949                          gmx_enerdata_t*                     enerd,
950                          t_fcdata*                           fcd,
951                          int                                 natoms,
952                          gmx::ArrayRefWithPadding<gmx::RVec> x,
953                          gmx::ArrayRefWithPadding<gmx::RVec> v,
954                          const matrix                        box,
955                          gmx::ArrayRef<real>                 lambda,
956                          history_t*                          hist,
957                          gmx::ArrayRefWithPadding<gmx::RVec> f,
958                          tensor                              force_vir,
959                          const t_mdatoms*                    md,
960                          t_nrnb*                             nrnb,
961                          gmx_wallcycle_t                     wcycle,
962                          t_graph*                            graph,
963                          gmx_shellfc_t*                      shfc,
964                          t_forcerec*                         fr,
965                          gmx::MdrunScheduleWorkload*         runScheduleWork,
966                          double                              t,
967                          rvec                                mu_tot,
968                          const gmx_vsite_t*                  vsite,
969                          const DDBalanceRegionHandler&       ddBalanceRegionHandler)
970 {
971     auto xRvec = as_rvec_array(x.paddedArrayRef().data());
972     auto vRvec = as_rvec_array(v.paddedArrayRef().data());
973
974     int           nshell;
975     t_shell*      shell;
976     const t_idef* idef;
977     rvec *        acc_dir = nullptr, *x_old = nullptr;
978     real          Epot[2], df[2];
979     real          sf_dir, invdt;
980     real          ftol, dum = 0;
981     char          sbuf[22];
982     gmx_bool      bCont, bInit, bConverged;
983     int           nat, dd_ac0, dd_ac1 = 0, i;
984     int           homenr = md->homenr, end = homenr;
985     int           nflexcon, number_steps, d, Min = 0, count = 0;
986 #define Try (1 - Min) /* At start Try = 1 */
987
988     bCont        = (mdstep == inputrec->init_step) && inputrec->bContinuation;
989     bInit        = (mdstep == inputrec->init_step) || shfc->bRequireInit;
990     ftol         = inputrec->em_tol;
991     number_steps = inputrec->niter;
992     nshell       = shfc->nshell;
993     shell        = shfc->shell;
994     nflexcon     = shfc->nflexcon;
995
996     idef = &top->idef;
997
998     if (DOMAINDECOMP(cr))
999     {
1000         nat = dd_natoms_vsite(cr->dd);
1001         if (nflexcon > 0)
1002         {
1003             dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1004             nat = std::max(nat, dd_ac1);
1005         }
1006     }
1007     else
1008     {
1009         nat = natoms;
1010     }
1011
1012     for (i = 0; (i < 2); i++)
1013     {
1014         shfc->x[i].resizeWithPadding(nat);
1015         shfc->f[i].resizeWithPadding(nat);
1016     }
1017
1018     /* Create views that we can swap */
1019     gmx::ArrayRefWithPadding<gmx::RVec> posWithPadding[2];
1020     gmx::ArrayRefWithPadding<gmx::RVec> forceWithPadding[2];
1021     gmx::ArrayRef<gmx::RVec>            pos[2];
1022     gmx::ArrayRef<gmx::RVec>            force[2];
1023     for (i = 0; (i < 2); i++)
1024     {
1025         posWithPadding[i]   = shfc->x[i].arrayRefWithPadding();
1026         pos[i]              = posWithPadding[i].paddedArrayRef();
1027         forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1028         force[i]            = forceWithPadding[i].paddedArrayRef();
1029     }
1030
1031     if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1032     {
1033         /* This is the only time where the coordinates are used
1034          * before do_force is called, which normally puts all
1035          * charge groups in the box.
1036          */
1037         auto xRef = x.paddedArrayRef();
1038         put_atoms_in_box_omp(fr->ePBC, box, xRef.subArray(0, md->homenr),
1039                              gmx_omp_nthreads_get(emntDefault));
1040
1041         if (graph)
1042         {
1043             mk_mshift(fplog, graph, fr->ePBC, box, xRvec);
1044         }
1045     }
1046
1047     /* After this all coordinate arrays will contain whole charge groups */
1048     if (graph)
1049     {
1050         shift_self(graph, box, xRvec);
1051     }
1052
1053     if (nflexcon)
1054     {
1055         if (nat > shfc->flex_nalloc)
1056         {
1057             shfc->flex_nalloc = over_alloc_dd(nat);
1058             srenew(shfc->acc_dir, shfc->flex_nalloc);
1059             srenew(shfc->x_old, shfc->flex_nalloc);
1060         }
1061         acc_dir        = shfc->acc_dir;
1062         x_old          = shfc->x_old;
1063         auto xArrayRef = x.paddedArrayRef();
1064         auto vArrayRef = v.paddedArrayRef();
1065         for (i = 0; i < homenr; i++)
1066         {
1067             for (d = 0; d < DIM; d++)
1068             {
1069                 shfc->x_old[i][d] = xArrayRef[i][d] - vArrayRef[i][d] * inputrec->delta_t;
1070             }
1071         }
1072     }
1073
1074     /* Do a prediction of the shell positions, when appropriate.
1075      * Without velocities (EM, NM, BD) we only do initial prediction.
1076      */
1077     if (shfc->bPredict && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1078     {
1079         predict_shells(fplog, xRvec, vRvec, inputrec->delta_t, nshell, shell, md->massT, nullptr, bInit);
1080     }
1081
1082     /* do_force expected the charge groups to be in the box */
1083     if (graph)
1084     {
1085         unshift_self(graph, box, xRvec);
1086     }
1087
1088     /* Calculate the forces first time around */
1089     if (gmx_debug_at)
1090     {
1091         pr_rvecs(debug, 0, "x b4 do_force", xRvec, homenr);
1092     }
1093     int shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1094     do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession, pull_work, mdstep,
1095              nrnb, wcycle, top, box, x, hist, forceWithPadding[Min], force_vir, md, enerd, fcd,
1096              lambda, graph, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
1097              (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags, ddBalanceRegionHandler);
1098
1099     sf_dir = 0;
1100     if (nflexcon)
1101     {
1102         init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end, shfc->x_old, xRvec, xRvec,
1103                   as_rvec_array(force[Min].data()), shfc->acc_dir, box, lambda, &dum);
1104
1105         for (i = 0; i < end; i++)
1106         {
1107             sf_dir += md->massT[i] * norm2(shfc->acc_dir[i]);
1108         }
1109     }
1110     sum_epot(&(enerd->grpp), enerd->term);
1111     Epot[Min] = enerd->term[F_EPOT];
1112
1113     df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), nshell, shell, nflexcon,
1114                         &sf_dir, &Epot[Min]);
1115     df[Try] = 0;
1116     if (debug)
1117     {
1118         fprintf(debug, "df = %g  %g\n", df[Min], df[Try]);
1119     }
1120
1121     if (gmx_debug_at)
1122     {
1123         pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md->nr);
1124     }
1125
1126     if (nshell + nflexcon > 0)
1127     {
1128         /* Copy x to pos[Min] & pos[Try]: during minimization only the
1129          * shell positions are updated, therefore the other particles must
1130          * be set here, in advance.
1131          */
1132         std::copy(x.paddedArrayRef().begin(), x.paddedArrayRef().end(),
1133                   posWithPadding[Min].paddedArrayRef().begin());
1134         std::copy(x.paddedArrayRef().begin(), x.paddedArrayRef().end(),
1135                   posWithPadding[Try].paddedArrayRef().begin());
1136     }
1137
1138     if (bVerbose && MASTER(cr))
1139     {
1140         print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1141     }
1142
1143     if (debug)
1144     {
1145         fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1146         fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1147         fprintf(debug, "%17s: %14.10e\n", interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1148         fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1149     }
1150
1151     /* First check whether we should do shells, or whether the force is
1152      * low enough even without minimization.
1153      */
1154     bConverged = (df[Min] < ftol);
1155
1156     for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1157     {
1158         if (vsite)
1159         {
1160             construct_vsites(vsite, as_rvec_array(pos[Min].data()), inputrec->delta_t, vRvec,
1161                              idef->iparams, idef->il, fr->ePBC, fr->bMolPBC, cr, box);
1162         }
1163
1164         if (nflexcon)
1165         {
1166             init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end, x_old, xRvec,
1167                       as_rvec_array(pos[Min].data()), as_rvec_array(force[Min].data()), acc_dir,
1168                       box, lambda, &dum);
1169
1170             directional_sd(pos[Min], pos[Try], acc_dir, end, fr->fc_stepsize);
1171         }
1172
1173         /* New positions, Steepest descent */
1174         shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1175
1176         /* do_force expected the charge groups to be in the box */
1177         if (graph)
1178         {
1179             unshift_self(graph, box, as_rvec_array(pos[Try].data()));
1180         }
1181
1182         if (gmx_debug_at)
1183         {
1184             pr_rvecs(debug, 0, "RELAX: pos[Min]  ", as_rvec_array(pos[Min].data()), homenr);
1185             pr_rvecs(debug, 0, "RELAX: pos[Try]  ", as_rvec_array(pos[Try].data()), homenr);
1186         }
1187         /* Try the new positions */
1188         do_force(fplog, cr, ms, inputrec, nullptr, enforcedRotation, imdSession, pull_work, 1, nrnb,
1189                  wcycle, top, box, posWithPadding[Try], hist, forceWithPadding[Try], force_vir, md,
1190                  enerd, fcd, lambda, graph, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
1191                  shellfc_flags, ddBalanceRegionHandler);
1192         sum_epot(&(enerd->grpp), enerd->term);
1193         if (gmx_debug_at)
1194         {
1195             pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1196             pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1197         }
1198         sf_dir = 0;
1199         if (nflexcon)
1200         {
1201             init_adir(shfc, constr, inputrec, cr, dd_ac1, mdstep, md, end, x_old, xRvec,
1202                       as_rvec_array(pos[Try].data()), as_rvec_array(force[Try].data()), acc_dir,
1203                       box, lambda, &dum);
1204
1205             for (i = 0; i < end; i++)
1206             {
1207                 sf_dir += md->massT[i] * norm2(acc_dir[i]);
1208             }
1209         }
1210
1211         Epot[Try] = enerd->term[F_EPOT];
1212
1213         df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1214
1215         if (debug)
1216         {
1217             fprintf(debug, "df = %g  %g\n", df[Min], df[Try]);
1218         }
1219
1220         if (debug)
1221         {
1222             if (gmx_debug_at)
1223             {
1224                 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1225             }
1226             if (gmx_debug_at)
1227             {
1228                 fprintf(debug, "SHELL ITER %d\n", count);
1229                 dump_shells(debug, force[Try], ftol, nshell, shell);
1230             }
1231         }
1232
1233         if (bVerbose && MASTER(cr))
1234         {
1235             print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1236         }
1237
1238         bConverged = (df[Try] < ftol);
1239
1240         if ((df[Try] < df[Min]))
1241         {
1242             if (debug)
1243             {
1244                 fprintf(debug, "Swapping Min and Try\n");
1245             }
1246             if (nflexcon)
1247             {
1248                 /* Correct the velocities for the flexible constraints */
1249                 invdt          = 1 / inputrec->delta_t;
1250                 auto vArrayRef = v.paddedArrayRef();
1251                 for (i = 0; i < end; i++)
1252                 {
1253                     for (d = 0; d < DIM; d++)
1254                     {
1255                         vArrayRef[i][d] += (pos[Try][i][d] - pos[Min][i][d]) * invdt;
1256                     }
1257                 }
1258             }
1259             Min = Try;
1260         }
1261         else
1262         {
1263             decrease_step_size(nshell, shell);
1264         }
1265     }
1266     shfc->numForceEvaluations += count;
1267     if (bConverged)
1268     {
1269         shfc->numConvergedIterations++;
1270     }
1271     if (MASTER(cr) && !(bConverged))
1272     {
1273         /* Note that the energies and virial are incorrect when not converged */
1274         if (fplog)
1275         {
1276             fprintf(fplog, "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1277                     gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1278         }
1279         fprintf(stderr, "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1280                 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1281     }
1282
1283     /* Copy back the coordinates and the forces */
1284     std::copy(pos[Min].begin(), pos[Min].end(), x.paddedArrayRef().data());
1285     std::copy(force[Min].begin(), force[Min].end(), f.unpaddedArrayRef().begin());
1286 }
1287
1288 void done_shellfc(FILE* fplog, gmx_shellfc_t* shfc, int64_t numSteps)
1289 {
1290     if (shfc && fplog && numSteps > 0)
1291     {
1292         double numStepsAsDouble = static_cast<double>(numSteps);
1293         fprintf(fplog, "Fraction of iterations that converged:           %.2f %%\n",
1294                 (shfc->numConvergedIterations * 100.0) / numStepsAsDouble);
1295         fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1296                 shfc->numForceEvaluations / numStepsAsDouble);
1297     }
1298
1299     // TODO Deallocate memory in shfc
1300 }