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