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