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