38b68d157ace0627e76648f2f1b7a8f7cba06352
[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,2021, 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,
454                                               atom[aS].qB,
455                                               aS + 1,
456                                               mb + 1);
457                                 }
458                                 shell[nsi].k += gmx::square(qS) * ONE_4PI_EPS0
459                                                 / ffparams->iparams[type].polarize.alpha;
460                                 break;
461                             case F_WATER_POL:
462                                 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS * 10))
463                                 {
464                                     gmx_fatal(FARGS,
465                                               "water_pol can not be used with qA(%e) != qB(%e) for "
466                                               "atom %d of molecule block %zu",
467                                               qS,
468                                               atom[aS].qB,
469                                               aS + 1,
470                                               mb + 1);
471                                 }
472                                 alpha = (ffparams->iparams[type].wpol.al_x
473                                          + ffparams->iparams[type].wpol.al_y
474                                          + ffparams->iparams[type].wpol.al_z)
475                                         / 3.0;
476                                 shell[nsi].k += gmx::square(qS) * ONE_4PI_EPS0 / alpha;
477                                 break;
478                             default: gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
479                         }
480                         shell[nsi].nnucl++;
481                     }
482                     ia += nra + 1;
483                     i += nra + 1;
484                 }
485             }
486             a_offset += molt->atoms.nr;
487         }
488         /* Done with this molecule type */
489     }
490
491     /* Verify whether it's all correct */
492     if (ns != nshell)
493     {
494         gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
495     }
496
497     for (i = 0; (i < ns); i++)
498     {
499         shell[i].k_1 = 1.0 / shell[i].k;
500     }
501
502     if (debug)
503     {
504         pr_shell(debug, shell);
505     }
506
507
508     shfc->shell_gl       = shell;
509     shfc->shell_index_gl = shell_index;
510
511     shfc->predictShells = (getenv("GMX_NOPREDICT") == nullptr);
512     shfc->requireInit   = false;
513     if (!shfc->predictShells)
514     {
515         if (fplog)
516         {
517             fprintf(fplog, "\nWill never predict shell positions\n");
518         }
519     }
520     else
521     {
522         shfc->requireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != nullptr);
523         if (shfc->requireInit && fplog)
524         {
525             fprintf(fplog, "\nWill always initiate shell positions\n");
526         }
527     }
528
529     if (shfc->predictShells)
530     {
531         if (shfc->bInterCG)
532         {
533             if (fplog)
534             {
535                 fprintf(fplog,
536                         "\nNOTE: in the current version shell prediction during the crun is "
537                         "disabled\n\n");
538             }
539             /* Prediction improves performance, so we should implement either:
540              * 1. communication for the atoms needed for prediction
541              * 2. prediction using the velocities of shells; currently the
542              *    shell velocities are zeroed, it's a bit tricky to keep
543              *    track of the shell displacements and thus the velocity.
544              */
545             shfc->predictShells = false;
546         }
547     }
548
549     /* shfc->x is used as a coordinate buffer for the sim_util's `do_force` function, and
550      * when using PME it must be pinned. */
551     if (usingPmeOnGpu)
552     {
553         for (i = 0; i < 2; i++)
554         {
555             changePinningPolicy(&shfc->x[i], gmx::PinningPolicy::PinnedIfSupported);
556         }
557     }
558
559     return shfc;
560 }
561
562 void gmx::make_local_shells(const t_commrec* cr, const t_mdatoms* md, gmx_shellfc_t* shfc)
563 {
564     int           a0, a1;
565     gmx_domdec_t* dd = nullptr;
566
567     if (DOMAINDECOMP(cr))
568     {
569         dd = cr->dd;
570         a0 = 0;
571         a1 = dd_numHomeAtoms(*dd);
572     }
573     else
574     {
575         /* Single node: we need all shells, copy them */
576         shfc->shells = shfc->shell_gl;
577
578         return;
579     }
580
581     ArrayRef<const int> ind = shfc->shell_index_gl;
582
583     std::vector<t_shell>& shells = shfc->shells;
584     shells.clear();
585     for (int i = a0; i < a1; i++)
586     {
587         if (md->ptype[i] == eptShell)
588         {
589             if (dd)
590             {
591                 shells.push_back(shfc->shell_gl[ind[dd->globalAtomIndices[i]]]);
592             }
593             else
594             {
595                 shells.push_back(shfc->shell_gl[ind[i]]);
596             }
597             t_shell& shell = shells.back();
598
599             /* With inter-cg shells we can no do shell prediction,
600              * so we do not need the nuclei numbers.
601              */
602             if (!shfc->bInterCG)
603             {
604                 shell.nucl1 = i + shell.nucl1 - shell.shellIndex;
605                 if (shell.nnucl > 1)
606                 {
607                     shell.nucl2 = i + shell.nucl2 - shell.shellIndex;
608                 }
609                 if (shell.nnucl > 2)
610                 {
611                     shell.nucl3 = i + shell.nucl3 - shell.shellIndex;
612                 }
613             }
614             shell.shellIndex = i;
615         }
616     }
617 }
618
619 static void do_1pos(rvec xnew, const rvec xold, const rvec f, real step)
620 {
621     real xo, yo, zo;
622     real dx, dy, dz;
623
624     xo = xold[XX];
625     yo = xold[YY];
626     zo = xold[ZZ];
627
628     dx = f[XX] * step;
629     dy = f[YY] * step;
630     dz = f[ZZ] * step;
631
632     xnew[XX] = xo + dx;
633     xnew[YY] = yo + dy;
634     xnew[ZZ] = zo + dz;
635 }
636
637 static void do_1pos3(rvec xnew, const rvec xold, const rvec f, const rvec step)
638 {
639     real xo, yo, zo;
640     real dx, dy, dz;
641
642     xo = xold[XX];
643     yo = xold[YY];
644     zo = xold[ZZ];
645
646     dx = f[XX] * step[XX];
647     dy = f[YY] * step[YY];
648     dz = f[ZZ] * step[ZZ];
649
650     xnew[XX] = xo + dx;
651     xnew[YY] = yo + dy;
652     xnew[ZZ] = zo + dz;
653 }
654
655 static void directional_sd(ArrayRef<const RVec> xold,
656                            ArrayRef<RVec>       xnew,
657                            ArrayRef<const RVec> acc_dir,
658                            int                  homenr,
659                            real                 step)
660 {
661     const rvec* xo = as_rvec_array(xold.data());
662     rvec*       xn = as_rvec_array(xnew.data());
663
664     for (int i = 0; i < homenr; i++)
665     {
666         do_1pos(xn[i], xo[i], acc_dir[i], step);
667     }
668 }
669
670 static void shell_pos_sd(ArrayRef<const RVec> xcur,
671                          ArrayRef<RVec>       xnew,
672                          ArrayRef<const RVec> f,
673                          ArrayRef<t_shell>    shells,
674                          int                  count)
675 {
676     const real step_scale_min = 0.8, step_scale_increment = 0.2, step_scale_max = 1.2,
677                step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
678     int        d;
679     real       dx, df, k_est;
680     const real zero = 0;
681 #ifdef PRINT_STEP
682     real step_min, step_max;
683
684     step_min = 1e30;
685     step_max = 0;
686 #endif
687     for (t_shell& shell : shells)
688     {
689         const int ind = shell.shellIndex;
690         if (count == 1)
691         {
692             for (d = 0; d < DIM; d++)
693             {
694                 shell.step[d] = shell.k_1;
695 #ifdef PRINT_STEP
696                 step_min = std::min(step_min, shell.step[d]);
697                 step_max = std::max(step_max, shell.step[d]);
698 #endif
699             }
700         }
701         else
702         {
703             for (d = 0; d < DIM; d++)
704             {
705                 dx = xcur[ind][d] - shell.xold[d];
706                 df = f[ind][d] - shell.fold[d];
707                 /* -dx/df gets used to generate an interpolated value, but would
708                  * cause a NaN if df were binary-equal to zero. Values close to
709                  * zero won't cause problems (because of the min() and max()), so
710                  * just testing for binary inequality is OK. */
711                 if (zero != df)
712                 {
713                     k_est = -dx / df;
714                     /* Scale the step size by a factor interpolated from
715                      * step_scale_min to step_scale_max, as k_est goes from 0 to
716                      * step_scale_multiple * shell.step[d] */
717                     shell.step[d] = step_scale_min * shell.step[d]
718                                     + step_scale_increment
719                                               * std::min(step_scale_multiple * shell.step[d],
720                                                          std::max(k_est, zero));
721                 }
722                 else
723                 {
724                     /* Here 0 == df */
725                     if (gmx_numzero(dx)) /* 0 == dx */
726                     {
727                         /* Likely this will never happen, but if it does just
728                          * don't scale the step. */
729                     }
730                     else /* 0 != dx */
731                     {
732                         shell.step[d] *= step_scale_max;
733                     }
734                 }
735 #ifdef PRINT_STEP
736                 step_min = std::min(step_min, shell.step[d]);
737                 step_max = std::max(step_max, shell.step[d]);
738 #endif
739             }
740         }
741         copy_rvec(xcur[ind], shell.xold);
742         copy_rvec(f[ind], shell.fold);
743
744         do_1pos3(xnew[ind], xcur[ind], f[ind], shell.step);
745
746         if (gmx_debug_at)
747         {
748             fprintf(debug, "shell = %d\n", ind);
749             pr_rvec(debug, 0, "fshell", f[ind], DIM, TRUE);
750             pr_rvec(debug, 0, "xold", xcur[ind], DIM, TRUE);
751             pr_rvec(debug, 0, "step", shell.step, DIM, TRUE);
752             pr_rvec(debug, 0, "xnew", xnew[ind], DIM, TRUE);
753         }
754     }
755 #ifdef PRINT_STEP
756     printf("step %.3e %.3e\n", step_min, step_max);
757 #endif
758 }
759
760 static void decrease_step_size(ArrayRef<t_shell> shells)
761 {
762     for (t_shell& shell : shells)
763     {
764         svmul(0.8, shell.step, shell.step);
765     }
766 }
767
768 static void print_epot(FILE* fp, int64_t mdstep, int count, real epot, real df, int ndir, real sf_dir)
769 {
770     char buf[22];
771
772     fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e", gmx_step_str(mdstep, buf), count, epot, df);
773     if (ndir)
774     {
775         fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir / ndir));
776     }
777     else
778     {
779         fprintf(fp, "\n");
780     }
781 }
782
783
784 static real rms_force(const t_commrec*        cr,
785                       ArrayRef<const RVec>    force,
786                       ArrayRef<const t_shell> shells,
787                       int                     ndir,
788                       real*                   sf_dir,
789                       real*                   Epot)
790 {
791     double      buf[4];
792     const rvec* f = as_rvec_array(force.data());
793
794     buf[0] = *sf_dir;
795     for (const t_shell& shell : shells)
796     {
797         buf[0] += norm2(f[shell.shellIndex]);
798     }
799     int ntot = shells.ssize();
800
801     if (PAR(cr))
802     {
803         buf[1] = ntot;
804         buf[2] = *sf_dir;
805         buf[3] = *Epot;
806         gmx_sumd(4, buf, cr);
807         ntot    = gmx::roundToInt(buf[1]);
808         *sf_dir = buf[2];
809         *Epot   = buf[3];
810     }
811     ntot += ndir;
812
813     return (ntot ? std::sqrt(buf[0] / ntot) : 0);
814 }
815
816 static void dump_shells(FILE* fp, ArrayRef<RVec> f, real ftol, ArrayRef<const t_shell> shells)
817 {
818     real ft2, ff2;
819
820     ft2 = gmx::square(ftol);
821
822     for (const t_shell& shell : shells)
823     {
824         const int ind = shell.shellIndex;
825         ff2           = iprod(f[ind], f[ind]);
826         if (ff2 > ft2)
827         {
828             fprintf(fp,
829                     "SHELL %5d, force %10.5f  %10.5f  %10.5f, |f| %10.5f\n",
830                     ind,
831                     f[ind][XX],
832                     f[ind][YY],
833                     f[ind][ZZ],
834                     std::sqrt(ff2));
835         }
836     }
837 }
838
839 static void init_adir(gmx_shellfc_t*            shfc,
840                       gmx::Constraints*         constr,
841                       const t_inputrec*         ir,
842                       const t_commrec*          cr,
843                       int                       dd_ac1,
844                       int64_t                   step,
845                       const t_mdatoms*          md,
846                       int                       end,
847                       ArrayRefWithPadding<RVec> xOld,
848                       ArrayRef<RVec>            x_init,
849                       ArrayRefWithPadding<RVec> xCurrent,
850                       ArrayRef<RVec>            f,
851                       ArrayRef<RVec>            acc_dir,
852                       const matrix              box,
853                       ArrayRef<const real>      lambda,
854                       real*                     dvdlambda)
855 {
856     double          dt, w_dt;
857     int             n, d;
858     unsigned short* ptype;
859
860     if (DOMAINDECOMP(cr))
861     {
862         n = dd_ac1;
863     }
864     else
865     {
866         n = end;
867     }
868     shfc->adir_xnold.resizeWithPadding(n);
869     shfc->adir_xnew.resizeWithPadding(n);
870     rvec* xnold = as_rvec_array(shfc->adir_xnold.data());
871     rvec* xnew  = as_rvec_array(shfc->adir_xnew.data());
872     rvec* x_old = as_rvec_array(xOld.paddedArrayRef().data());
873     rvec* x     = as_rvec_array(xCurrent.paddedArrayRef().data());
874
875     ptype = md->ptype;
876
877     dt = ir->delta_t;
878
879     /* Does NOT work with freeze groups (yet) */
880     for (n = 0; n < end; n++)
881     {
882         w_dt = md->invmass[n] * dt;
883
884         for (d = 0; d < DIM; d++)
885         {
886             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
887             {
888                 xnold[n][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
889                 xnew[n][d]  = 2 * x[n][d] - x_old[n][d] + f[n][d] * w_dt * dt;
890             }
891             else
892             {
893                 xnold[n][d] = x[n][d];
894                 xnew[n][d]  = x[n][d];
895             }
896         }
897     }
898     bool needsLogging  = false;
899     bool computeEnergy = false;
900     bool computeVirial = false;
901     constr->apply(needsLogging,
902                   computeEnergy,
903                   step,
904                   0,
905                   1.0,
906                   xCurrent,
907                   shfc->adir_xnold.arrayRefWithPadding(),
908                   {},
909                   box,
910                   lambda[efptBONDED],
911                   &(dvdlambda[efptBONDED]),
912                   {},
913                   computeVirial,
914                   nullptr,
915                   gmx::ConstraintVariable::Positions);
916     constr->apply(needsLogging,
917                   computeEnergy,
918                   step,
919                   0,
920                   1.0,
921                   xCurrent,
922                   shfc->adir_xnew.arrayRefWithPadding(),
923                   {},
924                   box,
925                   lambda[efptBONDED],
926                   &(dvdlambda[efptBONDED]),
927                   {},
928                   computeVirial,
929                   nullptr,
930                   gmx::ConstraintVariable::Positions);
931
932     for (n = 0; n < end; n++)
933     {
934         for (d = 0; d < DIM; d++)
935         {
936             xnew[n][d] = -(2 * x[n][d] - xnold[n][d] - xnew[n][d]) / gmx::square(dt)
937                          - f[n][d] * md->invmass[n];
938         }
939         clear_rvec(acc_dir[n]);
940     }
941
942     /* Project the acceleration on the old bond directions */
943     constr->apply(needsLogging,
944                   computeEnergy,
945                   step,
946                   0,
947                   1.0,
948                   xOld,
949                   shfc->adir_xnew.arrayRefWithPadding(),
950                   acc_dir,
951                   box,
952                   lambda[efptBONDED],
953                   &(dvdlambda[efptBONDED]),
954                   {},
955                   computeVirial,
956                   nullptr,
957                   gmx::ConstraintVariable::Deriv_FlexCon);
958 }
959
960 void relax_shell_flexcon(FILE*                         fplog,
961                          const t_commrec*              cr,
962                          const gmx_multisim_t*         ms,
963                          gmx_bool                      bVerbose,
964                          gmx_enfrot*                   enforcedRotation,
965                          int64_t                       mdstep,
966                          const t_inputrec*             inputrec,
967                          gmx::ImdSession*              imdSession,
968                          pull_t*                       pull_work,
969                          gmx_bool                      bDoNS,
970                          int                           force_flags,
971                          const gmx_localtop_t*         top,
972                          gmx::Constraints*             constr,
973                          gmx_enerdata_t*               enerd,
974                          int                           natoms,
975                          ArrayRefWithPadding<RVec>     xPadded,
976                          ArrayRefWithPadding<RVec>     vPadded,
977                          const matrix                  box,
978                          ArrayRef<real>                lambda,
979                          const history_t*              hist,
980                          gmx::ForceBuffersView*        f,
981                          tensor                        force_vir,
982                          const t_mdatoms*              md,
983                          t_nrnb*                       nrnb,
984                          gmx_wallcycle_t               wcycle,
985                          gmx_shellfc_t*                shfc,
986                          t_forcerec*                   fr,
987                          gmx::MdrunScheduleWorkload*   runScheduleWork,
988                          double                        t,
989                          rvec                          mu_tot,
990                          gmx::VirtualSitesHandler*     vsite,
991                          const DDBalanceRegionHandler& ddBalanceRegionHandler)
992 {
993     real Epot[2], df[2];
994     real sf_dir, invdt;
995     real dum = 0;
996     char sbuf[22];
997     int  nat, dd_ac0, dd_ac1 = 0, i;
998     int  homenr = md->homenr, end = homenr;
999     int  d, Min = 0, count = 0;
1000 #define Try (1 - Min) /* At start Try = 1 */
1001
1002     const bool        bCont        = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1003     const bool        bInit        = (mdstep == inputrec->init_step) || shfc->requireInit;
1004     const real        ftol         = inputrec->em_tol;
1005     const int         number_steps = inputrec->niter;
1006     ArrayRef<t_shell> shells       = shfc->shells;
1007     const int         nflexcon     = shfc->nflexcon;
1008
1009     if (DOMAINDECOMP(cr))
1010     {
1011         nat = dd_natoms_vsite(cr->dd);
1012         if (nflexcon > 0)
1013         {
1014             dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1015             nat = std::max(nat, dd_ac1);
1016         }
1017     }
1018     else
1019     {
1020         nat = natoms;
1021     }
1022
1023     for (i = 0; (i < 2); i++)
1024     {
1025         shfc->x[i].resizeWithPadding(nat);
1026         shfc->f[i].resizeWithPadding(nat);
1027     }
1028
1029     /* Create views that we can swap for trail and minimum for positions and forces */
1030     ArrayRefWithPadding<RVec> posWithPadding[2];
1031     ArrayRefWithPadding<RVec> forceWithPadding[2];
1032     ArrayRef<RVec>            pos[2];
1033     ArrayRef<RVec>            force[2];
1034     for (i = 0; (i < 2); i++)
1035     {
1036         posWithPadding[i]   = shfc->x[i].arrayRefWithPadding();
1037         pos[i]              = posWithPadding[i].paddedArrayRef();
1038         forceWithPadding[i] = shfc->f[i].arrayRefWithPadding();
1039         force[i]            = forceWithPadding[i].paddedArrayRef();
1040     }
1041
1042     ArrayRef<RVec> x = xPadded.unpaddedArrayRef();
1043     ArrayRef<RVec> v = vPadded.unpaddedArrayRef();
1044
1045     if (bDoNS && inputrec->pbcType != PbcType::No && !DOMAINDECOMP(cr))
1046     {
1047         /* This is the only time where the coordinates are used
1048          * before do_force is called, which normally puts all
1049          * charge groups in the box.
1050          */
1051         put_atoms_in_box_omp(
1052                 fr->pbcType, box, x.subArray(0, md->homenr), gmx_omp_nthreads_get(emntDefault));
1053     }
1054
1055     if (nflexcon)
1056     {
1057         shfc->acc_dir.resize(nat);
1058         shfc->x_old.resizeWithPadding(nat);
1059         ArrayRef<RVec> x_old = shfc->x_old.arrayRefWithPadding().unpaddedArrayRef();
1060         for (i = 0; i < homenr; i++)
1061         {
1062             for (d = 0; d < DIM; d++)
1063             {
1064                 x_old[i][d] = x[i][d] - v[i][d] * inputrec->delta_t;
1065             }
1066         }
1067     }
1068
1069     /* Do a prediction of the shell positions, when appropriate.
1070      * Without velocities (EM, NM, BD) we only do initial prediction.
1071      */
1072     if (shfc->predictShells && !bCont && (EI_STATE_VELOCITY(inputrec->eI) || bInit))
1073     {
1074         predict_shells(fplog, x, v, inputrec->delta_t, shells, md->massT, nullptr, bInit);
1075     }
1076
1077     /* Calculate the forces first time around */
1078     if (gmx_debug_at)
1079     {
1080         pr_rvecs(debug, 0, "x b4 do_force", as_rvec_array(x.data()), homenr);
1081     }
1082     int                   shellfc_flags = force_flags | (bVerbose ? GMX_FORCE_ENERGY : 0);
1083     gmx::ForceBuffersView forceViewInit = gmx::ForceBuffersView(forceWithPadding[Min], {}, false);
1084     do_force(fplog,
1085              cr,
1086              ms,
1087              inputrec,
1088              nullptr,
1089              enforcedRotation,
1090              imdSession,
1091              pull_work,
1092              mdstep,
1093              nrnb,
1094              wcycle,
1095              top,
1096              box,
1097              xPadded,
1098              hist,
1099              &forceViewInit,
1100              force_vir,
1101              md,
1102              enerd,
1103              lambda,
1104              fr,
1105              runScheduleWork,
1106              vsite,
1107              mu_tot,
1108              t,
1109              nullptr,
1110              (bDoNS ? GMX_FORCE_NS : 0) | shellfc_flags,
1111              ddBalanceRegionHandler);
1112
1113     sf_dir = 0;
1114     if (nflexcon)
1115     {
1116         init_adir(shfc,
1117                   constr,
1118                   inputrec,
1119                   cr,
1120                   dd_ac1,
1121                   mdstep,
1122                   md,
1123                   end,
1124                   shfc->x_old.arrayRefWithPadding(),
1125                   x,
1126                   xPadded,
1127                   force[Min],
1128                   shfc->acc_dir,
1129                   box,
1130                   lambda,
1131                   &dum);
1132
1133         for (i = 0; i < end; i++)
1134         {
1135             sf_dir += md->massT[i] * norm2(shfc->acc_dir[i]);
1136         }
1137     }
1138     accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals);
1139     Epot[Min] = enerd->term[F_EPOT];
1140
1141     df[Min] = rms_force(cr, forceWithPadding[Min].paddedArrayRef(), shells, nflexcon, &sf_dir, &Epot[Min]);
1142     df[Try] = 0;
1143     if (debug)
1144     {
1145         fprintf(debug, "df = %g  %g\n", df[Min], df[Try]);
1146     }
1147
1148     if (gmx_debug_at)
1149     {
1150         pr_rvecs(debug, 0, "force0", as_rvec_array(force[Min].data()), md->nr);
1151     }
1152
1153     if (!shells.empty() || nflexcon > 0)
1154     {
1155         /* Copy x to pos[Min] & pos[Try]: during minimization only the
1156          * shell positions are updated, therefore the other particles must
1157          * be set here, in advance.
1158          */
1159         std::copy(xPadded.paddedArrayRef().begin(),
1160                   xPadded.paddedArrayRef().end(),
1161                   posWithPadding[Min].paddedArrayRef().begin());
1162         std::copy(xPadded.paddedArrayRef().begin(),
1163                   xPadded.paddedArrayRef().end(),
1164                   posWithPadding[Try].paddedArrayRef().begin());
1165     }
1166
1167     if (bVerbose && MASTER(cr))
1168     {
1169         print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1170     }
1171
1172     if (debug)
1173     {
1174         fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1175         fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1176         fprintf(debug, "%17s: %14.10e\n", interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1177         fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1178     }
1179
1180     /* First check whether we should do shells, or whether the force is
1181      * low enough even without minimization.
1182      */
1183     bool bConverged = (df[Min] < ftol);
1184
1185     for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1186     {
1187         if (vsite)
1188         {
1189             vsite->construct(pos[Min], inputrec->delta_t, v, box);
1190         }
1191
1192         if (nflexcon)
1193         {
1194             init_adir(shfc,
1195                       constr,
1196                       inputrec,
1197                       cr,
1198                       dd_ac1,
1199                       mdstep,
1200                       md,
1201                       end,
1202                       shfc->x_old.arrayRefWithPadding(),
1203                       x,
1204                       posWithPadding[Min],
1205                       force[Min],
1206                       shfc->acc_dir,
1207                       box,
1208                       lambda,
1209                       &dum);
1210
1211             directional_sd(pos[Min], pos[Try], shfc->acc_dir, end, fr->fc_stepsize);
1212         }
1213
1214         /* New positions, Steepest descent */
1215         shell_pos_sd(pos[Min], pos[Try], force[Min], shells, count);
1216
1217         if (gmx_debug_at)
1218         {
1219             pr_rvecs(debug, 0, "RELAX: pos[Min]  ", as_rvec_array(pos[Min].data()), homenr);
1220             pr_rvecs(debug, 0, "RELAX: pos[Try]  ", as_rvec_array(pos[Try].data()), homenr);
1221         }
1222         /* Try the new positions */
1223         gmx::ForceBuffersView forceViewTry = gmx::ForceBuffersView(forceWithPadding[Try], {}, false);
1224         do_force(fplog,
1225                  cr,
1226                  ms,
1227                  inputrec,
1228                  nullptr,
1229                  enforcedRotation,
1230                  imdSession,
1231                  pull_work,
1232                  1,
1233                  nrnb,
1234                  wcycle,
1235                  top,
1236                  box,
1237                  posWithPadding[Try],
1238                  hist,
1239                  &forceViewTry,
1240                  force_vir,
1241                  md,
1242                  enerd,
1243                  lambda,
1244                  fr,
1245                  runScheduleWork,
1246                  vsite,
1247                  mu_tot,
1248                  t,
1249                  nullptr,
1250                  shellfc_flags,
1251                  ddBalanceRegionHandler);
1252         accumulatePotentialEnergies(enerd, lambda, inputrec->fepvals);
1253         if (gmx_debug_at)
1254         {
1255             pr_rvecs(debug, 0, "RELAX: force[Min]", as_rvec_array(force[Min].data()), homenr);
1256             pr_rvecs(debug, 0, "RELAX: force[Try]", as_rvec_array(force[Try].data()), homenr);
1257         }
1258         sf_dir = 0;
1259         if (nflexcon)
1260         {
1261             init_adir(shfc,
1262                       constr,
1263                       inputrec,
1264                       cr,
1265                       dd_ac1,
1266                       mdstep,
1267                       md,
1268                       end,
1269                       shfc->x_old.arrayRefWithPadding(),
1270                       x,
1271                       posWithPadding[Try],
1272                       force[Try],
1273                       shfc->acc_dir,
1274                       box,
1275                       lambda,
1276                       &dum);
1277
1278             ArrayRef<const RVec> acc_dir = shfc->acc_dir;
1279             for (i = 0; i < end; i++)
1280             {
1281                 sf_dir += md->massT[i] * norm2(acc_dir[i]);
1282             }
1283         }
1284
1285         Epot[Try] = enerd->term[F_EPOT];
1286
1287         df[Try] = rms_force(cr, force[Try], shells, nflexcon, &sf_dir, &Epot[Try]);
1288
1289         if (debug)
1290         {
1291             fprintf(debug, "df = %g  %g\n", df[Min], df[Try]);
1292         }
1293
1294         if (debug)
1295         {
1296             if (gmx_debug_at)
1297             {
1298                 pr_rvecs(debug, 0, "F na do_force", as_rvec_array(force[Try].data()), homenr);
1299             }
1300             if (gmx_debug_at)
1301             {
1302                 fprintf(debug, "SHELL ITER %d\n", count);
1303                 dump_shells(debug, force[Try], ftol, shells);
1304             }
1305         }
1306
1307         if (bVerbose && MASTER(cr))
1308         {
1309             print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1310         }
1311
1312         bConverged = (df[Try] < ftol);
1313
1314         if ((df[Try] < df[Min]))
1315         {
1316             if (debug)
1317             {
1318                 fprintf(debug, "Swapping Min and Try\n");
1319             }
1320             if (nflexcon)
1321             {
1322                 /* Correct the velocities for the flexible constraints */
1323                 invdt = 1 / inputrec->delta_t;
1324                 for (i = 0; i < end; i++)
1325                 {
1326                     for (d = 0; d < DIM; d++)
1327                     {
1328                         v[i][d] += (pos[Try][i][d] - pos[Min][i][d]) * invdt;
1329                     }
1330                 }
1331             }
1332             Min = Try;
1333         }
1334         else
1335         {
1336             decrease_step_size(shells);
1337         }
1338     }
1339     shfc->numForceEvaluations += count;
1340     if (bConverged)
1341     {
1342         shfc->numConvergedIterations++;
1343     }
1344     if (MASTER(cr) && !(bConverged))
1345     {
1346         /* Note that the energies and virial are incorrect when not converged */
1347         if (fplog)
1348         {
1349             fprintf(fplog,
1350                     "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1351                     gmx_step_str(mdstep, sbuf),
1352                     number_steps,
1353                     df[Min]);
1354         }
1355         fprintf(stderr,
1356                 "step %s: EM did not converge in %d iterations, RMS force %6.2e\n",
1357                 gmx_step_str(mdstep, sbuf),
1358                 number_steps,
1359                 df[Min]);
1360     }
1361
1362     /* Copy back the coordinates and the forces */
1363     std::copy(pos[Min].begin(), pos[Min].end(), x.data());
1364     std::copy(force[Min].begin(), force[Min].end(), f->force().begin());
1365 }
1366
1367 void done_shellfc(FILE* fplog, gmx_shellfc_t* shfc, int64_t numSteps)
1368 {
1369     if (shfc && fplog && numSteps > 0)
1370     {
1371         double numStepsAsDouble = static_cast<double>(numSteps);
1372         fprintf(fplog,
1373                 "Fraction of iterations that converged:           %.2f %%\n",
1374                 (shfc->numConvergedIterations * 100.0) / numStepsAsDouble);
1375         fprintf(fplog,
1376                 "Average number of force evaluations per MD step: %.2f\n\n",
1377                 shfc->numForceEvaluations / numStepsAsDouble);
1378     }
1379
1380     delete shfc;
1381 }