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