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