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