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