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