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