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