Fixes for clang-tidy-9
[alexxy/gromacs.git] / src / gromacs / mdrun / minimize.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-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017 The GROMACS development team.
7  * Copyright (c) 2018,2019,2020, 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 /*! \internal \file
39  *
40  * \brief This file defines integrators for energy minimization
41  *
42  * \author Berk Hess <hess@kth.se>
43  * \author Erik Lindahl <erik@kth.se>
44  * \ingroup module_mdrun
45  */
46 #include "gmxpre.h"
47
48 #include "config.h"
49
50 #include <cmath>
51 #include <cstring>
52 #include <ctime>
53
54 #include <algorithm>
55 #include <vector>
56
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/collect.h"
59 #include "gromacs/domdec/dlbtiming.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/mdsetup.h"
63 #include "gromacs/domdec/partition.h"
64 #include "gromacs/ewald/pme_pp.h"
65 #include "gromacs/fileio/confio.h"
66 #include "gromacs/fileio/mtxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/imd/imd.h"
70 #include "gromacs/linearalgebra/sparsematrix.h"
71 #include "gromacs/listed_forces/manage_threading.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/mdlib/constr.h"
75 #include "gromacs/mdlib/dispersioncorrection.h"
76 #include "gromacs/mdlib/ebin.h"
77 #include "gromacs/mdlib/enerdata_utils.h"
78 #include "gromacs/mdlib/energyoutput.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/force_flags.h"
81 #include "gromacs/mdlib/forcerec.h"
82 #include "gromacs/mdlib/gmx_omp_nthreads.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/mdatoms.h"
85 #include "gromacs/mdlib/stat.h"
86 #include "gromacs/mdlib/tgroup.h"
87 #include "gromacs/mdlib/trajectory_writing.h"
88 #include "gromacs/mdlib/update.h"
89 #include "gromacs/mdlib/vsite.h"
90 #include "gromacs/mdrunutility/handlerestart.h"
91 #include "gromacs/mdrunutility/printtime.h"
92 #include "gromacs/mdtypes/commrec.h"
93 #include "gromacs/mdtypes/forcerec.h"
94 #include "gromacs/mdtypes/inputrec.h"
95 #include "gromacs/mdtypes/interaction_const.h"
96 #include "gromacs/mdtypes/md_enums.h"
97 #include "gromacs/mdtypes/mdatom.h"
98 #include "gromacs/mdtypes/mdrunoptions.h"
99 #include "gromacs/mdtypes/state.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/timing/wallcycle.h"
102 #include "gromacs/timing/walltime_accounting.h"
103 #include "gromacs/topology/mtop_util.h"
104 #include "gromacs/topology/topology.h"
105 #include "gromacs/utility/cstringutil.h"
106 #include "gromacs/utility/exceptions.h"
107 #include "gromacs/utility/fatalerror.h"
108 #include "gromacs/utility/logger.h"
109 #include "gromacs/utility/smalloc.h"
110
111 #include "legacysimulator.h"
112 #include "shellfc.h"
113
114 using gmx::ArrayRef;
115 using gmx::MdrunScheduleWorkload;
116 using gmx::RVec;
117
118 //! Utility structure for manipulating states during EM
119 typedef struct
120 {
121     //! Copy of the global state
122     t_state s;
123     //! Force array
124     PaddedHostVector<gmx::RVec> f;
125     //! Potential energy
126     real epot;
127     //! Norm of the force
128     real fnorm;
129     //! Maximum force
130     real fmax;
131     //! Direction
132     int a_fmax;
133 } em_state_t;
134
135 //! Print the EM starting conditions
136 static void print_em_start(FILE*                     fplog,
137                            const t_commrec*          cr,
138                            gmx_walltime_accounting_t walltime_accounting,
139                            gmx_wallcycle_t           wcycle,
140                            const char*               name)
141 {
142     walltime_accounting_start_time(walltime_accounting);
143     wallcycle_start(wcycle, ewcRUN);
144     print_start(fplog, cr, walltime_accounting, name);
145 }
146
147 //! Stop counting time for EM
148 static void em_time_end(gmx_walltime_accounting_t walltime_accounting, gmx_wallcycle_t wcycle)
149 {
150     wallcycle_stop(wcycle, ewcRUN);
151
152     walltime_accounting_end_time(walltime_accounting);
153 }
154
155 //! Printing a log file and console header
156 static void sp_header(FILE* out, const char* minimizer, real ftol, int nsteps)
157 {
158     fprintf(out, "\n");
159     fprintf(out, "%s:\n", minimizer);
160     fprintf(out, "   Tolerance (Fmax)   = %12.5e\n", ftol);
161     fprintf(out, "   Number of steps    = %12d\n", nsteps);
162 }
163
164 //! Print warning message
165 static void warn_step(FILE* fp, real ftol, real fmax, gmx_bool bLastStep, gmx_bool bConstrain)
166 {
167     constexpr bool realIsDouble = GMX_DOUBLE;
168     char           buffer[2048];
169
170     if (!std::isfinite(fmax))
171     {
172         sprintf(buffer,
173                 "\nEnergy minimization has stopped because the force "
174                 "on at least one atom is not finite. This usually means "
175                 "atoms are overlapping. Modify the input coordinates to "
176                 "remove atom overlap or use soft-core potentials with "
177                 "the free energy code to avoid infinite forces.\n%s",
178                 !realIsDouble ? "You could also be lucky that switching to double precision "
179                                 "is sufficient to obtain finite forces.\n"
180                               : "");
181     }
182     else if (bLastStep)
183     {
184         sprintf(buffer,
185                 "\nEnergy minimization reached the maximum number "
186                 "of steps before the forces reached the requested "
187                 "precision Fmax < %g.\n",
188                 ftol);
189     }
190     else
191     {
192         sprintf(buffer,
193                 "\nEnergy minimization has stopped, but the forces have "
194                 "not converged to the requested precision Fmax < %g (which "
195                 "may not be possible for your system). It stopped "
196                 "because the algorithm tried to make a new step whose size "
197                 "was too small, or there was no change in the energy since "
198                 "last step. Either way, we regard the minimization as "
199                 "converged to within the available machine precision, "
200                 "given your starting configuration and EM parameters.\n%s%s",
201                 ftol,
202                 !realIsDouble ? "\nDouble precision normally gives you higher accuracy, but "
203                                 "this is often not needed for preparing to run molecular "
204                                 "dynamics.\n"
205                               : "",
206                 bConstrain ? "You might need to increase your constraint accuracy, or turn\n"
207                              "off constraints altogether (set constraints = none in mdp file)\n"
208                            : "");
209     }
210
211     fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
212     fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
213 }
214
215 //! Print message about convergence of the EM
216 static void print_converged(FILE*             fp,
217                             const char*       alg,
218                             real              ftol,
219                             int64_t           count,
220                             gmx_bool          bDone,
221                             int64_t           nsteps,
222                             const em_state_t* ems,
223                             double            sqrtNumAtoms)
224 {
225     char buf[STEPSTRSIZE];
226
227     if (bDone)
228     {
229         fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n", alg, ftol, gmx_step_str(count, buf));
230     }
231     else if (count < nsteps)
232     {
233         fprintf(fp,
234                 "\n%s converged to machine precision in %s steps,\n"
235                 "but did not reach the requested Fmax < %g.\n",
236                 alg, gmx_step_str(count, buf), ftol);
237     }
238     else
239     {
240         fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n", alg, ftol,
241                 gmx_step_str(count, buf));
242     }
243
244 #if GMX_DOUBLE
245     fprintf(fp, "Potential Energy  = %21.14e\n", ems->epot);
246     fprintf(fp, "Maximum force     = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
247     fprintf(fp, "Norm of force     = %21.14e\n", ems->fnorm / sqrtNumAtoms);
248 #else
249     fprintf(fp, "Potential Energy  = %14.7e\n", ems->epot);
250     fprintf(fp, "Maximum force     = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
251     fprintf(fp, "Norm of force     = %14.7e\n", ems->fnorm / sqrtNumAtoms);
252 #endif
253 }
254
255 //! Compute the norm and max of the force array in parallel
256 static void get_f_norm_max(const t_commrec* cr,
257                            t_grpopts*       opts,
258                            t_mdatoms*       mdatoms,
259                            const rvec*      f,
260                            real*            fnorm,
261                            real*            fmax,
262                            int*             a_fmax)
263 {
264     double fnorm2, *sum;
265     real   fmax2, fam;
266     int    la_max, a_max, start, end, i, m, gf;
267
268     /* This routine finds the largest force and returns it.
269      * On parallel machines the global max is taken.
270      */
271     fnorm2 = 0;
272     fmax2  = 0;
273     la_max = -1;
274     start  = 0;
275     end    = mdatoms->homenr;
276     if (mdatoms->cFREEZE)
277     {
278         for (i = start; i < end; i++)
279         {
280             gf  = mdatoms->cFREEZE[i];
281             fam = 0;
282             for (m = 0; m < DIM; m++)
283             {
284                 if (!opts->nFreeze[gf][m])
285                 {
286                     fam += gmx::square(f[i][m]);
287                 }
288             }
289             fnorm2 += fam;
290             if (fam > fmax2)
291             {
292                 fmax2  = fam;
293                 la_max = i;
294             }
295         }
296     }
297     else
298     {
299         for (i = start; i < end; i++)
300         {
301             fam = norm2(f[i]);
302             fnorm2 += fam;
303             if (fam > fmax2)
304             {
305                 fmax2  = fam;
306                 la_max = i;
307             }
308         }
309     }
310
311     if (la_max >= 0 && DOMAINDECOMP(cr))
312     {
313         a_max = cr->dd->globalAtomIndices[la_max];
314     }
315     else
316     {
317         a_max = la_max;
318     }
319     if (PAR(cr))
320     {
321         snew(sum, 2 * cr->nnodes + 1);
322         sum[2 * cr->nodeid]     = fmax2;
323         sum[2 * cr->nodeid + 1] = a_max;
324         sum[2 * cr->nnodes]     = fnorm2;
325         gmx_sumd(2 * cr->nnodes + 1, sum, cr);
326         fnorm2 = sum[2 * cr->nnodes];
327         /* Determine the global maximum */
328         for (i = 0; i < cr->nnodes; i++)
329         {
330             if (sum[2 * i] > fmax2)
331             {
332                 fmax2 = sum[2 * i];
333                 a_max = gmx::roundToInt(sum[2 * i + 1]);
334             }
335         }
336         sfree(sum);
337     }
338
339     if (fnorm)
340     {
341         *fnorm = sqrt(fnorm2);
342     }
343     if (fmax)
344     {
345         *fmax = sqrt(fmax2);
346     }
347     if (a_fmax)
348     {
349         *a_fmax = a_max;
350     }
351 }
352
353 //! Compute the norm of the force
354 static void get_state_f_norm_max(const t_commrec* cr, t_grpopts* opts, t_mdatoms* mdatoms, em_state_t* ems)
355 {
356     get_f_norm_max(cr, opts, mdatoms, ems->f.rvec_array(), &ems->fnorm, &ems->fmax, &ems->a_fmax);
357 }
358
359 //! Initialize the energy minimization
360 static void init_em(FILE*                fplog,
361                     const gmx::MDLogger& mdlog,
362                     const char*          title,
363                     const t_commrec*     cr,
364                     t_inputrec*          ir,
365                     gmx::ImdSession*     imdSession,
366                     pull_t*              pull_work,
367                     t_state*             state_global,
368                     const gmx_mtop_t*    top_global,
369                     em_state_t*          ems,
370                     gmx_localtop_t*      top,
371                     t_nrnb*              nrnb,
372                     t_forcerec*          fr,
373                     gmx::MDAtoms*        mdAtoms,
374                     gmx_global_stat_t*   gstat,
375                     gmx_vsite_t*         vsite,
376                     gmx::Constraints*    constr,
377                     gmx_shellfc_t**      shellfc)
378 {
379     real dvdl_constr;
380
381     if (fplog)
382     {
383         fprintf(fplog, "Initiating %s\n", title);
384     }
385
386     if (MASTER(cr))
387     {
388         state_global->ngtc = 0;
389     }
390     initialize_lambdas(fplog, *ir, MASTER(cr), &(state_global->fep_state), state_global->lambda, nullptr);
391
392     if (ir->eI == eiNM)
393     {
394         GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
395
396         *shellfc = init_shell_flexcon(stdout, top_global, constr ? constr->numFlexibleConstraints() : 0,
397                                       ir->nstcalcenergy, DOMAINDECOMP(cr));
398     }
399     else
400     {
401         GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
402                    "This else currently only handles energy minimizers, consider if your algorithm "
403                    "needs shell/flexible-constraint support");
404
405         /* With energy minimization, shells and flexible constraints are
406          * automatically minimized when treated like normal DOFS.
407          */
408         if (shellfc != nullptr)
409         {
410             *shellfc = nullptr;
411         }
412     }
413
414     auto mdatoms = mdAtoms->mdatoms();
415     if (DOMAINDECOMP(cr))
416     {
417         dd_init_local_state(cr->dd, state_global, &ems->s);
418
419         /* Distribute the charge groups over the nodes from the master node */
420         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
421                             imdSession, pull_work, &ems->s, &ems->f, mdAtoms, top, fr, vsite,
422                             constr, nrnb, nullptr, FALSE);
423         dd_store_state(cr->dd, &ems->s);
424     }
425     else
426     {
427         state_change_natoms(state_global, state_global->natoms);
428         /* Just copy the state */
429         ems->s = *state_global;
430         state_change_natoms(&ems->s, ems->s.natoms);
431
432         mdAlgorithmsSetupAtomData(cr, ir, *top_global, top, fr, &ems->f, mdAtoms, constr, vsite,
433                                   shellfc ? *shellfc : nullptr);
434
435         if (vsite)
436         {
437             set_vsite_top(vsite, top, mdatoms);
438         }
439     }
440
441     update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
442
443     if (constr)
444     {
445         // TODO how should this cross-module support dependency be managed?
446         if (ir->eConstrAlg == econtSHAKE && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
447         {
448             gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
449                       econstr_names[econtSHAKE], econstr_names[econtLINCS]);
450         }
451
452         if (!ir->bContinuation)
453         {
454             /* Constrain the starting coordinates */
455             dvdl_constr = 0;
456             constr->apply(TRUE, TRUE, -1, 0, 1.0, ems->s.x.arrayRefWithPadding(),
457                           ems->s.x.arrayRefWithPadding(), ArrayRef<RVec>(), ems->s.box,
458                           ems->s.lambda[efptFEP], &dvdl_constr, gmx::ArrayRefWithPadding<RVec>(),
459                           nullptr, gmx::ConstraintVariable::Positions);
460         }
461     }
462
463     if (PAR(cr))
464     {
465         *gstat = global_stat_init(ir);
466     }
467     else
468     {
469         *gstat = nullptr;
470     }
471
472     calc_shifts(ems->s.box, fr->shift_vec);
473 }
474
475 //! Finalize the minimization
476 static void finish_em(const t_commrec*          cr,
477                       gmx_mdoutf_t              outf,
478                       gmx_walltime_accounting_t walltime_accounting,
479                       gmx_wallcycle_t           wcycle)
480 {
481     if (!thisRankHasDuty(cr, DUTY_PME))
482     {
483         /* Tell the PME only node to finish */
484         gmx_pme_send_finish(cr);
485     }
486
487     done_mdoutf(outf);
488
489     em_time_end(walltime_accounting, wcycle);
490 }
491
492 //! Swap two different EM states during minimization
493 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
494 {
495     em_state_t* tmp;
496
497     tmp   = *ems1;
498     *ems1 = *ems2;
499     *ems2 = tmp;
500 }
501
502 //! Save the EM trajectory
503 static void write_em_traj(FILE*               fplog,
504                           const t_commrec*    cr,
505                           gmx_mdoutf_t        outf,
506                           gmx_bool            bX,
507                           gmx_bool            bF,
508                           const char*         confout,
509                           const gmx_mtop_t*   top_global,
510                           t_inputrec*         ir,
511                           int64_t             step,
512                           em_state_t*         state,
513                           t_state*            state_global,
514                           ObservablesHistory* observablesHistory)
515 {
516     int mdof_flags = 0;
517
518     if (bX)
519     {
520         mdof_flags |= MDOF_X;
521     }
522     if (bF)
523     {
524         mdof_flags |= MDOF_F;
525     }
526
527     /* If we want IMD output, set appropriate MDOF flag */
528     if (ir->bIMD)
529     {
530         mdof_flags |= MDOF_IMD;
531     }
532
533     mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step,
534                                      static_cast<double>(step), &state->s, state_global,
535                                      observablesHistory, state->f);
536
537     if (confout != nullptr)
538     {
539         if (DOMAINDECOMP(cr))
540         {
541             /* If bX=true, x was collected to state_global in the call above */
542             if (!bX)
543             {
544                 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
545                 dd_collect_vec(cr->dd, &state->s, state->s.x, globalXRef);
546             }
547         }
548         else
549         {
550             /* Copy the local state pointer */
551             state_global = &state->s;
552         }
553
554         if (MASTER(cr))
555         {
556             if (ir->pbcType != PbcType::No && !ir->bPeriodicMols && DOMAINDECOMP(cr))
557             {
558                 /* Make molecules whole only for confout writing */
559                 do_pbc_mtop(ir->pbcType, state->s.box, top_global, state_global->x.rvec_array());
560             }
561
562             write_sto_conf_mtop(confout, *top_global->name, top_global,
563                                 state_global->x.rvec_array(), nullptr, ir->pbcType, state->s.box);
564         }
565     }
566 }
567
568 //! \brief Do one minimization step
569 //
570 // \returns true when the step succeeded, false when a constraint error occurred
571 static bool do_em_step(const t_commrec*                   cr,
572                        t_inputrec*                        ir,
573                        t_mdatoms*                         md,
574                        em_state_t*                        ems1,
575                        real                               a,
576                        const PaddedHostVector<gmx::RVec>* force,
577                        em_state_t*                        ems2,
578                        gmx::Constraints*                  constr,
579                        int64_t                            count)
580
581 {
582     t_state *s1, *s2;
583     int      start, end;
584     real     dvdl_constr;
585     int nthreads gmx_unused;
586
587     bool validStep = true;
588
589     s1 = &ems1->s;
590     s2 = &ems2->s;
591
592     if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
593     {
594         gmx_incons("state mismatch in do_em_step");
595     }
596
597     s2->flags = s1->flags;
598
599     if (s2->natoms != s1->natoms)
600     {
601         state_change_natoms(s2, s1->natoms);
602         ems2->f.resizeWithPadding(s2->natoms);
603     }
604     if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
605     {
606         s2->cg_gl.resize(s1->cg_gl.size());
607     }
608
609     copy_mat(s1->box, s2->box);
610     /* Copy free energy state */
611     s2->lambda = s1->lambda;
612     copy_mat(s1->box, s2->box);
613
614     start = 0;
615     end   = md->homenr;
616
617     nthreads = gmx_omp_nthreads_get(emntUpdate);
618 #pragma omp parallel num_threads(nthreads)
619     {
620         const rvec* x1 = s1->x.rvec_array();
621         rvec*       x2 = s2->x.rvec_array();
622         const rvec* f  = force->rvec_array();
623
624         int gf = 0;
625 #pragma omp for schedule(static) nowait
626         for (int i = start; i < end; i++)
627         {
628             try
629             {
630                 if (md->cFREEZE)
631                 {
632                     gf = md->cFREEZE[i];
633                 }
634                 for (int m = 0; m < DIM; m++)
635                 {
636                     if (ir->opts.nFreeze[gf][m])
637                     {
638                         x2[i][m] = x1[i][m];
639                     }
640                     else
641                     {
642                         x2[i][m] = x1[i][m] + a * f[i][m];
643                     }
644                 }
645             }
646             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
647         }
648
649         if (s2->flags & (1 << estCGP))
650         {
651             /* Copy the CG p vector */
652             const rvec* p1 = s1->cg_p.rvec_array();
653             rvec*       p2 = s2->cg_p.rvec_array();
654 #pragma omp for schedule(static) nowait
655             for (int i = start; i < end; i++)
656             {
657                 // Trivial OpenMP block that does not throw
658                 copy_rvec(p1[i], p2[i]);
659             }
660         }
661
662         if (DOMAINDECOMP(cr))
663         {
664             /* OpenMP does not supported unsigned loop variables */
665 #pragma omp for schedule(static) nowait
666             for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
667             {
668                 s2->cg_gl[i] = s1->cg_gl[i];
669             }
670         }
671     }
672
673     if (DOMAINDECOMP(cr))
674     {
675         s2->ddp_count       = s1->ddp_count;
676         s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
677     }
678
679     if (constr)
680     {
681         dvdl_constr = 0;
682         validStep   = constr->apply(
683                 TRUE, TRUE, count, 0, 1.0, s1->x.arrayRefWithPadding(), s2->x.arrayRefWithPadding(),
684                 ArrayRef<RVec>(), s2->box, s2->lambda[efptBONDED], &dvdl_constr,
685                 gmx::ArrayRefWithPadding<RVec>(), nullptr, gmx::ConstraintVariable::Positions);
686
687         if (cr->nnodes > 1)
688         {
689             /* This global reduction will affect performance at high
690              * parallelization, but we can not really avoid it.
691              * But usually EM is not run at high parallelization.
692              */
693             int reductionBuffer = static_cast<int>(!validStep);
694             gmx_sumi(1, &reductionBuffer, cr);
695             validStep = (reductionBuffer == 0);
696         }
697
698         // We should move this check to the different minimizers
699         if (!validStep && ir->eI != eiSteep)
700         {
701             gmx_fatal(FARGS,
702                       "The coordinates could not be constrained. Minimizer '%s' can not handle "
703                       "constraint failures, use minimizer '%s' before using '%s'.",
704                       EI(ir->eI), EI(eiSteep), EI(ir->eI));
705         }
706     }
707
708     return validStep;
709 }
710
711 //! Prepare EM for using domain decomposition parallellization
712 static void em_dd_partition_system(FILE*                fplog,
713                                    const gmx::MDLogger& mdlog,
714                                    int                  step,
715                                    const t_commrec*     cr,
716                                    const gmx_mtop_t*    top_global,
717                                    t_inputrec*          ir,
718                                    gmx::ImdSession*     imdSession,
719                                    pull_t*              pull_work,
720                                    em_state_t*          ems,
721                                    gmx_localtop_t*      top,
722                                    gmx::MDAtoms*        mdAtoms,
723                                    t_forcerec*          fr,
724                                    gmx_vsite_t*         vsite,
725                                    gmx::Constraints*    constr,
726                                    t_nrnb*              nrnb,
727                                    gmx_wallcycle_t      wcycle)
728 {
729     /* Repartition the domain decomposition */
730     dd_partition_system(fplog, mdlog, step, cr, FALSE, 1, nullptr, *top_global, ir, imdSession, pull_work,
731                         &ems->s, &ems->f, mdAtoms, top, fr, vsite, constr, nrnb, wcycle, FALSE);
732     dd_store_state(cr->dd, &ems->s);
733 }
734
735 namespace
736 {
737
738 /*! \brief Class to handle the work of setting and doing an energy evaluation.
739  *
740  * This class is a mere aggregate of parameters to pass to evaluate an
741  * energy, so that future changes to names and types of them consume
742  * less time when refactoring other code.
743  *
744  * Aggregate initialization is used, for which the chief risk is that
745  * if a member is added at the end and not all initializer lists are
746  * updated, then the member will be value initialized, which will
747  * typically mean initialization to zero.
748  *
749  * Use a braced initializer list to construct one of these. */
750 class EnergyEvaluator
751 {
752 public:
753     /*! \brief Evaluates an energy on the state in \c ems.
754      *
755      * \todo In practice, the same objects mu_tot, vir, and pres
756      * are always passed to this function, so we would rather have
757      * them as data members. However, their C-array types are
758      * unsuited for aggregate initialization. When the types
759      * improve, the call signature of this method can be reduced.
760      */
761     void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst);
762     //! Handles logging (deprecated).
763     FILE* fplog;
764     //! Handles logging.
765     const gmx::MDLogger& mdlog;
766     //! Handles communication.
767     const t_commrec* cr;
768     //! Coordinates multi-simulations.
769     const gmx_multisim_t* ms;
770     //! Holds the simulation topology.
771     const gmx_mtop_t* top_global;
772     //! Holds the domain topology.
773     gmx_localtop_t* top;
774     //! User input options.
775     t_inputrec* inputrec;
776     //! The Interactive Molecular Dynamics session.
777     gmx::ImdSession* imdSession;
778     //! The pull work object.
779     pull_t* pull_work;
780     //! Manages flop accounting.
781     t_nrnb* nrnb;
782     //! Manages wall cycle accounting.
783     gmx_wallcycle_t wcycle;
784     //! Coordinates global reduction.
785     gmx_global_stat_t gstat;
786     //! Handles virtual sites.
787     gmx_vsite_t* vsite;
788     //! Handles constraints.
789     gmx::Constraints* constr;
790     //! Handles strange things.
791     t_fcdata* fcd;
792     //! Per-atom data for this domain.
793     gmx::MDAtoms* mdAtoms;
794     //! Handles how to calculate the forces.
795     t_forcerec* fr;
796     //! Schedule of force-calculation work each step for this task.
797     MdrunScheduleWorkload* runScheduleWork;
798     //! Stores the computed energies.
799     gmx_enerdata_t* enerd;
800 };
801
802 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst)
803 {
804     real     t;
805     gmx_bool bNS;
806     tensor   force_vir, shake_vir, ekin;
807     real     dvdl_constr;
808     real     terminate = 0;
809
810     /* Set the time to the initial time, the time does not change during EM */
811     t = inputrec->init_t;
812
813     if (bFirst || (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
814     {
815         /* This is the first state or an old state used before the last ns */
816         bNS = TRUE;
817     }
818     else
819     {
820         bNS = FALSE;
821         if (inputrec->nstlist > 0)
822         {
823             bNS = TRUE;
824         }
825     }
826
827     if (vsite)
828     {
829         construct_vsites(vsite, ems->s.x.rvec_array(), 1, nullptr, top->idef.iparams, top->idef.il,
830                          fr->pbcType, fr->bMolPBC, cr, ems->s.box);
831     }
832
833     if (DOMAINDECOMP(cr) && bNS)
834     {
835         /* Repartition the domain decomposition */
836         em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work,
837                                ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
838     }
839
840     /* Calc force & energy on new trial position  */
841     /* do_force always puts the charge groups in the box and shifts again
842      * We do not unshift, so molecules are always whole in congrad.c
843      */
844     do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession, pull_work, count, nrnb, wcycle,
845              top, ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist,
846              ems->f.arrayRefWithPadding(), force_vir, mdAtoms->mdatoms(), enerd, fcd, ems->s.lambda,
847              fr, runScheduleWork, vsite, mu_tot, t, nullptr,
848              GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY
849                      | (bNS ? GMX_FORCE_NS : 0),
850              DDBalanceRegionHandler(cr));
851
852     /* Clear the unused shake virial and pressure */
853     clear_mat(shake_vir);
854     clear_mat(pres);
855
856     /* Communicate stuff when parallel */
857     if (PAR(cr) && inputrec->eI != eiNM)
858     {
859         wallcycle_start(wcycle, ewcMoveE);
860
861         global_stat(gstat, cr, enerd, force_vir, shake_vir, inputrec, nullptr, nullptr, nullptr, 1,
862                     &terminate, nullptr, FALSE, CGLO_ENERGY | CGLO_PRESSURE | CGLO_CONSTRAINT);
863
864         wallcycle_stop(wcycle, ewcMoveE);
865     }
866
867     if (fr->dispersionCorrection)
868     {
869         /* Calculate long range corrections to pressure and energy */
870         const DispersionCorrection::Correction correction =
871                 fr->dispersionCorrection->calculate(ems->s.box, ems->s.lambda[efptVDW]);
872
873         enerd->term[F_DISPCORR] = correction.energy;
874         enerd->term[F_EPOT] += correction.energy;
875         enerd->term[F_PRES] += correction.pressure;
876         enerd->term[F_DVDL] += correction.dvdl;
877     }
878     else
879     {
880         enerd->term[F_DISPCORR] = 0;
881     }
882
883     ems->epot = enerd->term[F_EPOT];
884
885     if (constr)
886     {
887         /* Project out the constraint components of the force */
888         dvdl_constr = 0;
889         auto f      = ems->f.arrayRefWithPadding();
890         constr->apply(FALSE, FALSE, count, 0, 1.0, ems->s.x.arrayRefWithPadding(), f,
891                       f.unpaddedArrayRef(), ems->s.box, ems->s.lambda[efptBONDED], &dvdl_constr,
892                       gmx::ArrayRefWithPadding<RVec>(), &shake_vir, gmx::ConstraintVariable::ForceDispl);
893         enerd->term[F_DVDL_CONSTR] += dvdl_constr;
894         m_add(force_vir, shake_vir, vir);
895     }
896     else
897     {
898         copy_mat(force_vir, vir);
899     }
900
901     clear_mat(ekin);
902     enerd->term[F_PRES] = calc_pres(fr->pbcType, inputrec->nwall, ems->s.box, ekin, vir, pres);
903
904     sum_dhdl(enerd, ems->s.lambda, *inputrec->fepvals);
905
906     if (EI_ENERGY_MINIMIZATION(inputrec->eI))
907     {
908         get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
909     }
910 }
911
912 } // namespace
913
914 //! Parallel utility summing energies and forces
915 static double reorder_partsum(const t_commrec*  cr,
916                               t_grpopts*        opts,
917                               const gmx_mtop_t* top_global,
918                               em_state_t*       s_min,
919                               em_state_t*       s_b)
920 {
921     if (debug)
922     {
923         fprintf(debug, "Doing reorder_partsum\n");
924     }
925
926     const rvec* fm = s_min->f.rvec_array();
927     const rvec* fb = s_b->f.rvec_array();
928
929     /* Collect fm in a global vector fmg.
930      * This conflicts with the spirit of domain decomposition,
931      * but to fully optimize this a much more complicated algorithm is required.
932      */
933     const int natoms = top_global->natoms;
934     rvec*     fmg;
935     snew(fmg, natoms);
936
937     gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
938     int                      i          = 0;
939     for (int a : indicesMin)
940     {
941         copy_rvec(fm[i], fmg[a]);
942         i++;
943     }
944     gmx_sum(top_global->natoms * 3, fmg[0], cr);
945
946     /* Now we will determine the part of the sum for the cgs in state s_b */
947     gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
948
949     double partsum                        = 0;
950     i                                     = 0;
951     int                                gf = 0;
952     gmx::ArrayRef<const unsigned char> grpnrFREEZE =
953             top_global->groups.groupNumbers[SimulationAtomGroupType::Freeze];
954     for (int a : indicesB)
955     {
956         if (!grpnrFREEZE.empty())
957         {
958             gf = grpnrFREEZE[i];
959         }
960         for (int m = 0; m < DIM; m++)
961         {
962             if (!opts->nFreeze[gf][m])
963             {
964                 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
965             }
966         }
967         i++;
968     }
969
970     sfree(fmg);
971
972     return partsum;
973 }
974
975 //! Print some stuff, like beta, whatever that means.
976 static real pr_beta(const t_commrec*  cr,
977                     t_grpopts*        opts,
978                     t_mdatoms*        mdatoms,
979                     const gmx_mtop_t* top_global,
980                     em_state_t*       s_min,
981                     em_state_t*       s_b)
982 {
983     double sum;
984
985     /* This is just the classical Polak-Ribiere calculation of beta;
986      * it looks a bit complicated since we take freeze groups into account,
987      * and might have to sum it in parallel runs.
988      */
989
990     if (!DOMAINDECOMP(cr)
991         || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
992     {
993         const rvec* fm = s_min->f.rvec_array();
994         const rvec* fb = s_b->f.rvec_array();
995         sum            = 0;
996         int gf         = 0;
997         /* This part of code can be incorrect with DD,
998          * since the atom ordering in s_b and s_min might differ.
999          */
1000         for (int i = 0; i < mdatoms->homenr; i++)
1001         {
1002             if (mdatoms->cFREEZE)
1003             {
1004                 gf = mdatoms->cFREEZE[i];
1005             }
1006             for (int m = 0; m < DIM; m++)
1007             {
1008                 if (!opts->nFreeze[gf][m])
1009                 {
1010                     sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1011                 }
1012             }
1013         }
1014     }
1015     else
1016     {
1017         /* We need to reorder cgs while summing */
1018         sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1019     }
1020     if (PAR(cr))
1021     {
1022         gmx_sumd(1, &sum, cr);
1023     }
1024
1025     return sum / gmx::square(s_min->fnorm);
1026 }
1027
1028 namespace gmx
1029 {
1030
1031 void LegacySimulator::do_cg()
1032 {
1033     const char* CG = "Polak-Ribiere Conjugate Gradients";
1034
1035     gmx_localtop_t    top(top_global->ffparams);
1036     gmx_global_stat_t gstat;
1037     double            tmp, minstep;
1038     real              stepsize;
1039     real              a, b, c, beta = 0.0;
1040     real              epot_repl = 0;
1041     real              pnorm;
1042     gmx_bool          converged, foundlower;
1043     rvec              mu_tot = { 0 };
1044     gmx_bool          do_log = FALSE, do_ene = FALSE, do_x, do_f;
1045     tensor            vir, pres;
1046     int               number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1047     int               m, step, nminstep;
1048     auto              mdatoms = mdAtoms->mdatoms();
1049
1050     GMX_LOG(mdlog.info)
1051             .asParagraph()
1052             .appendText(
1053                     "Note that activating conjugate gradient energy minimization via the "
1054                     "integrator .mdp option and the command gmx mdrun may "
1055                     "be available in a different form in a future version of GROMACS, "
1056                     "e.g. gmx minimize and an .mdp option.");
1057
1058     step = 0;
1059
1060     if (MASTER(cr))
1061     {
1062         // In CG, the state is extended with a search direction
1063         state_global->flags |= (1 << estCGP);
1064
1065         // Ensure the extra per-atom state array gets allocated
1066         state_change_natoms(state_global, state_global->natoms);
1067
1068         // Initialize the search direction to zero
1069         for (RVec& cg_p : state_global->cg_p)
1070         {
1071             cg_p = { 0, 0, 0 };
1072         }
1073     }
1074
1075     /* Create 4 states on the stack and extract pointers that we will swap */
1076     em_state_t  s0{}, s1{}, s2{}, s3{};
1077     em_state_t* s_min = &s0;
1078     em_state_t* s_a   = &s1;
1079     em_state_t* s_b   = &s2;
1080     em_state_t* s_c   = &s3;
1081
1082     /* Init em and store the local state in s_min */
1083     init_em(fplog, mdlog, CG, cr, inputrec, imdSession, pull_work, state_global, top_global, s_min,
1084             &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
1085     const bool        simulationsShareState = false;
1086     gmx_mdoutf*       outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
1087                                    mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1088                                    StartingBehavior::NewSimulation, simulationsShareState, ms);
1089     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
1090                                    false, StartingBehavior::NewSimulation, mdModulesNotifier);
1091
1092     /* Print to log file */
1093     print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1094
1095     /* Max number of steps */
1096     number_steps = inputrec->nsteps;
1097
1098     if (MASTER(cr))
1099     {
1100         sp_header(stderr, CG, inputrec->em_tol, number_steps);
1101     }
1102     if (fplog)
1103     {
1104         sp_header(fplog, CG, inputrec->em_tol, number_steps);
1105     }
1106
1107     EnergyEvaluator energyEvaluator{
1108         fplog, mdlog,  cr,    ms,    top_global, &top, inputrec, imdSession, pull_work,
1109         nrnb,  wcycle, gstat, vsite, constr,     fcd,  mdAtoms,  fr,         runScheduleWork,
1110         enerd
1111     };
1112     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1113     /* do_force always puts the charge groups in the box and shifts again
1114      * We do not unshift, so molecules are always whole in congrad.c
1115      */
1116     energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1117
1118     if (MASTER(cr))
1119     {
1120         /* Copy stuff to the energy bin for easy printing etc. */
1121         matrix nullBox = {};
1122         energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1123                                          enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1124                                          nullptr, vir, pres, nullptr, mu_tot, constr);
1125
1126         EnergyOutput::printHeader(fplog, step, step);
1127         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step,
1128                                            step, fcd, nullptr);
1129     }
1130
1131     /* Estimate/guess the initial stepsize */
1132     stepsize = inputrec->em_stepsize / s_min->fnorm;
1133
1134     if (MASTER(cr))
1135     {
1136         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1137         fprintf(stderr, "   F-max             = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1138         fprintf(stderr, "   F-Norm            = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1139         fprintf(stderr, "\n");
1140         /* and copy to the log file too... */
1141         fprintf(fplog, "   F-max             = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1142         fprintf(fplog, "   F-Norm            = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1143         fprintf(fplog, "\n");
1144     }
1145     /* Start the loop over CG steps.
1146      * Each successful step is counted, and we continue until
1147      * we either converge or reach the max number of steps.
1148      */
1149     converged = FALSE;
1150     for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1151     {
1152
1153         /* start taking steps in a new direction
1154          * First time we enter the routine, beta=0, and the direction is
1155          * simply the negative gradient.
1156          */
1157
1158         /* Calculate the new direction in p, and the gradient in this direction, gpa */
1159         rvec*       pm  = s_min->s.cg_p.rvec_array();
1160         const rvec* sfm = s_min->f.rvec_array();
1161         double      gpa = 0;
1162         int         gf  = 0;
1163         for (int i = 0; i < mdatoms->homenr; i++)
1164         {
1165             if (mdatoms->cFREEZE)
1166             {
1167                 gf = mdatoms->cFREEZE[i];
1168             }
1169             for (m = 0; m < DIM; m++)
1170             {
1171                 if (!inputrec->opts.nFreeze[gf][m])
1172                 {
1173                     pm[i][m] = sfm[i][m] + beta * pm[i][m];
1174                     gpa -= pm[i][m] * sfm[i][m];
1175                     /* f is negative gradient, thus the sign */
1176                 }
1177                 else
1178                 {
1179                     pm[i][m] = 0;
1180                 }
1181             }
1182         }
1183
1184         /* Sum the gradient along the line across CPUs */
1185         if (PAR(cr))
1186         {
1187             gmx_sumd(1, &gpa, cr);
1188         }
1189
1190         /* Calculate the norm of the search vector */
1191         get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1192
1193         /* Just in case stepsize reaches zero due to numerical precision... */
1194         if (stepsize <= 0)
1195         {
1196             stepsize = inputrec->em_stepsize / pnorm;
1197         }
1198
1199         /*
1200          * Double check the value of the derivative in the search direction.
1201          * If it is positive it must be due to the old information in the
1202          * CG formula, so just remove that and start over with beta=0.
1203          * This corresponds to a steepest descent step.
1204          */
1205         if (gpa > 0)
1206         {
1207             beta = 0;
1208             step--;   /* Don't count this step since we are restarting */
1209             continue; /* Go back to the beginning of the big for-loop */
1210         }
1211
1212         /* Calculate minimum allowed stepsize, before the average (norm)
1213          * relative change in coordinate is smaller than precision
1214          */
1215         minstep      = 0;
1216         auto s_min_x = makeArrayRef(s_min->s.x);
1217         for (int i = 0; i < mdatoms->homenr; i++)
1218         {
1219             for (m = 0; m < DIM; m++)
1220             {
1221                 tmp = fabs(s_min_x[i][m]);
1222                 if (tmp < 1.0)
1223                 {
1224                     tmp = 1.0;
1225                 }
1226                 tmp = pm[i][m] / tmp;
1227                 minstep += tmp * tmp;
1228             }
1229         }
1230         /* Add up from all CPUs */
1231         if (PAR(cr))
1232         {
1233             gmx_sumd(1, &minstep, cr);
1234         }
1235
1236         minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global->natoms));
1237
1238         if (stepsize < minstep)
1239         {
1240             converged = TRUE;
1241             break;
1242         }
1243
1244         /* Write coordinates if necessary */
1245         do_x = do_per_step(step, inputrec->nstxout);
1246         do_f = do_per_step(step, inputrec->nstfout);
1247
1248         write_em_traj(fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min,
1249                       state_global, observablesHistory);
1250
1251         /* Take a step downhill.
1252          * In theory, we should minimize the function along this direction.
1253          * That is quite possible, but it turns out to take 5-10 function evaluations
1254          * for each line. However, we dont really need to find the exact minimum -
1255          * it is much better to start a new CG step in a modified direction as soon
1256          * as we are close to it. This will save a lot of energy evaluations.
1257          *
1258          * In practice, we just try to take a single step.
1259          * If it worked (i.e. lowered the energy), we increase the stepsize but
1260          * the continue straight to the next CG step without trying to find any minimum.
1261          * If it didn't work (higher energy), there must be a minimum somewhere between
1262          * the old position and the new one.
1263          *
1264          * Due to the finite numerical accuracy, it turns out that it is a good idea
1265          * to even accept a SMALL increase in energy, if the derivative is still downhill.
1266          * This leads to lower final energies in the tests I've done. / Erik
1267          */
1268         s_a->epot = s_min->epot;
1269         a         = 0.0;
1270         c         = a + stepsize; /* reference position along line is zero */
1271
1272         if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1273         {
1274             em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec, imdSession,
1275                                    pull_work, s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
1276         }
1277
1278         /* Take a trial step (new coords in s_c) */
1279         do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c, constr, -1);
1280
1281         neval++;
1282         /* Calculate energy for the trial step */
1283         energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1284
1285         /* Calc derivative along line */
1286         const rvec* pc  = s_c->s.cg_p.rvec_array();
1287         const rvec* sfc = s_c->f.rvec_array();
1288         double      gpc = 0;
1289         for (int i = 0; i < mdatoms->homenr; i++)
1290         {
1291             for (m = 0; m < DIM; m++)
1292             {
1293                 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1294             }
1295         }
1296         /* Sum the gradient along the line across CPUs */
1297         if (PAR(cr))
1298         {
1299             gmx_sumd(1, &gpc, cr);
1300         }
1301
1302         /* This is the max amount of increase in energy we tolerate */
1303         tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1304
1305         /* Accept the step if the energy is lower, or if it is not significantly higher
1306          * and the line derivative is still negative.
1307          */
1308         if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1309         {
1310             foundlower = TRUE;
1311             /* Great, we found a better energy. Increase step for next iteration
1312              * if we are still going down, decrease it otherwise
1313              */
1314             if (gpc < 0)
1315             {
1316                 stepsize *= 1.618034; /* The golden section */
1317             }
1318             else
1319             {
1320                 stepsize *= 0.618034; /* 1/golden section */
1321             }
1322         }
1323         else
1324         {
1325             /* New energy is the same or higher. We will have to do some work
1326              * to find a smaller value in the interval. Take smaller step next time!
1327              */
1328             foundlower = FALSE;
1329             stepsize *= 0.618034;
1330         }
1331
1332
1333         /* OK, if we didn't find a lower value we will have to locate one now - there must
1334          * be one in the interval [a=0,c].
1335          * The same thing is valid here, though: Don't spend dozens of iterations to find
1336          * the line minimum. We try to interpolate based on the derivative at the endpoints,
1337          * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1338          *
1339          * I also have a safeguard for potentially really pathological functions so we never
1340          * take more than 20 steps before we give up ...
1341          *
1342          * If we already found a lower value we just skip this step and continue to the update.
1343          */
1344         double gpb;
1345         if (!foundlower)
1346         {
1347             nminstep = 0;
1348
1349             do
1350             {
1351                 /* Select a new trial point.
1352                  * If the derivatives at points a & c have different sign we interpolate to zero,
1353                  * otherwise just do a bisection.
1354                  */
1355                 if (gpa < 0 && gpc > 0)
1356                 {
1357                     b = a + gpa * (a - c) / (gpc - gpa);
1358                 }
1359                 else
1360                 {
1361                     b = 0.5 * (a + c);
1362                 }
1363
1364                 /* safeguard if interpolation close to machine accuracy causes errors:
1365                  * never go outside the interval
1366                  */
1367                 if (b <= a || b >= c)
1368                 {
1369                     b = 0.5 * (a + c);
1370                 }
1371
1372                 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1373                 {
1374                     /* Reload the old state */
1375                     em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec, imdSession, pull_work,
1376                                            s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
1377                 }
1378
1379                 /* Take a trial step to this new point - new coords in s_b */
1380                 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b, constr, -1);
1381
1382                 neval++;
1383                 /* Calculate energy for the trial step */
1384                 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1385
1386                 /* p does not change within a step, but since the domain decomposition
1387                  * might change, we have to use cg_p of s_b here.
1388                  */
1389                 const rvec* pb  = s_b->s.cg_p.rvec_array();
1390                 const rvec* sfb = s_b->f.rvec_array();
1391                 gpb             = 0;
1392                 for (int i = 0; i < mdatoms->homenr; i++)
1393                 {
1394                     for (m = 0; m < DIM; m++)
1395                     {
1396                         gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1397                     }
1398                 }
1399                 /* Sum the gradient along the line across CPUs */
1400                 if (PAR(cr))
1401                 {
1402                     gmx_sumd(1, &gpb, cr);
1403                 }
1404
1405                 if (debug)
1406                 {
1407                     fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot,
1408                             s_c->epot, gpb);
1409                 }
1410
1411                 epot_repl = s_b->epot;
1412
1413                 /* Keep one of the intervals based on the value of the derivative at the new point */
1414                 if (gpb > 0)
1415                 {
1416                     /* Replace c endpoint with b */
1417                     swap_em_state(&s_b, &s_c);
1418                     c   = b;
1419                     gpc = gpb;
1420                 }
1421                 else
1422                 {
1423                     /* Replace a endpoint with b */
1424                     swap_em_state(&s_b, &s_a);
1425                     a   = b;
1426                     gpa = gpb;
1427                 }
1428
1429                 /*
1430                  * Stop search as soon as we find a value smaller than the endpoints.
1431                  * Never run more than 20 steps, no matter what.
1432                  */
1433                 nminstep++;
1434             } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1435
1436             if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1437             {
1438                 /* OK. We couldn't find a significantly lower energy.
1439                  * If beta==0 this was steepest descent, and then we give up.
1440                  * If not, set beta=0 and restart with steepest descent before quitting.
1441                  */
1442                 if (beta == 0.0)
1443                 {
1444                     /* Converged */
1445                     converged = TRUE;
1446                     break;
1447                 }
1448                 else
1449                 {
1450                     /* Reset memory before giving up */
1451                     beta = 0.0;
1452                     continue;
1453                 }
1454             }
1455
1456             /* Select min energy state of A & C, put the best in B.
1457              */
1458             if (s_c->epot < s_a->epot)
1459             {
1460                 if (debug)
1461                 {
1462                     fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot,
1463                             s_a->epot);
1464                 }
1465                 swap_em_state(&s_b, &s_c);
1466                 gpb = gpc;
1467             }
1468             else
1469             {
1470                 if (debug)
1471                 {
1472                     fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot,
1473                             s_c->epot);
1474                 }
1475                 swap_em_state(&s_b, &s_a);
1476                 gpb = gpa;
1477             }
1478         }
1479         else
1480         {
1481             if (debug)
1482             {
1483                 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1484             }
1485             swap_em_state(&s_b, &s_c);
1486             gpb = gpc;
1487         }
1488
1489         /* new search direction */
1490         /* beta = 0 means forget all memory and restart with steepest descents. */
1491         if (nstcg && ((step % nstcg) == 0))
1492         {
1493             beta = 0.0;
1494         }
1495         else
1496         {
1497             /* s_min->fnorm cannot be zero, because then we would have converged
1498              * and broken out.
1499              */
1500
1501             /* Polak-Ribiere update.
1502              * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1503              */
1504             beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1505         }
1506         /* Limit beta to prevent oscillations */
1507         if (fabs(beta) > 5.0)
1508         {
1509             beta = 0.0;
1510         }
1511
1512
1513         /* update positions */
1514         swap_em_state(&s_min, &s_b);
1515         gpa = gpb;
1516
1517         /* Print it if necessary */
1518         if (MASTER(cr))
1519         {
1520             if (mdrunOptions.verbose)
1521             {
1522                 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1523                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step,
1524                         s_min->epot, s_min->fnorm / sqrtNumAtoms, s_min->fmax, s_min->a_fmax + 1);
1525                 fflush(stderr);
1526             }
1527             /* Store the new (lower) energies */
1528             matrix nullBox = {};
1529             energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1530                                              enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1531                                              nullptr, vir, pres, nullptr, mu_tot, constr);
1532
1533             do_log = do_per_step(step, inputrec->nstlog);
1534             do_ene = do_per_step(step, inputrec->nstenergy);
1535
1536             imdSession->fillEnergyRecord(step, TRUE);
1537
1538             if (do_log)
1539             {
1540                 EnergyOutput::printHeader(fplog, step, step);
1541             }
1542             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1543                                                do_log ? fplog : nullptr, step, step, fcd, nullptr);
1544         }
1545
1546         /* Send energies and positions to the IMD client if bIMD is TRUE. */
1547         if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0))
1548         {
1549             imdSession->sendPositionsAndEnergies();
1550         }
1551
1552         /* Stop when the maximum force lies below tolerance.
1553          * If we have reached machine precision, converged is already set to true.
1554          */
1555         converged = converged || (s_min->fmax < inputrec->em_tol);
1556
1557     } /* End of the loop */
1558
1559     if (converged)
1560     {
1561         step--; /* we never took that last step in this case */
1562     }
1563     if (s_min->fmax > inputrec->em_tol)
1564     {
1565         if (MASTER(cr))
1566         {
1567             warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1568         }
1569         converged = FALSE;
1570     }
1571
1572     if (MASTER(cr))
1573     {
1574         /* If we printed energy and/or logfile last step (which was the last step)
1575          * we don't have to do it again, but otherwise print the final values.
1576          */
1577         if (!do_log)
1578         {
1579             /* Write final value to log since we didn't do anything the last step */
1580             EnergyOutput::printHeader(fplog, step, step);
1581         }
1582         if (!do_ene || !do_log)
1583         {
1584             /* Write final energy file entries */
1585             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1586                                                !do_log ? fplog : nullptr, step, step, fcd, nullptr);
1587         }
1588     }
1589
1590     /* Print some stuff... */
1591     if (MASTER(cr))
1592     {
1593         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1594     }
1595
1596     /* IMPORTANT!
1597      * For accurate normal mode calculation it is imperative that we
1598      * store the last conformation into the full precision binary trajectory.
1599      *
1600      * However, we should only do it if we did NOT already write this step
1601      * above (which we did if do_x or do_f was true).
1602      */
1603     /* Note that with 0 < nstfout != nstxout we can end up with two frames
1604      * in the trajectory with the same step number.
1605      */
1606     do_x = !do_per_step(step, inputrec->nstxout);
1607     do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1608
1609     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec,
1610                   step, s_min, state_global, observablesHistory);
1611
1612
1613     if (MASTER(cr))
1614     {
1615         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1616         print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1617         print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1618
1619         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1620     }
1621
1622     finish_em(cr, outf, walltime_accounting, wcycle);
1623
1624     /* To print the actual number of steps we needed somewhere */
1625     walltime_accounting_set_nsteps_done(walltime_accounting, step);
1626 }
1627
1628
1629 void LegacySimulator::do_lbfgs()
1630 {
1631     static const char* LBFGS = "Low-Memory BFGS Minimizer";
1632     em_state_t         ems;
1633     gmx_localtop_t     top(top_global->ffparams);
1634     gmx_global_stat_t  gstat;
1635     int                ncorr, nmaxcorr, point, cp, neval, nminstep;
1636     double             stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1637     real *             rho, *alpha, *p, *s, **dx, **dg;
1638     real               a, b, c, maxdelta, delta;
1639     real               diag, Epot0;
1640     real               dgdx, dgdg, sq, yr, beta;
1641     gmx_bool           converged;
1642     rvec               mu_tot = { 0 };
1643     gmx_bool           do_log, do_ene, do_x, do_f, foundlower, *frozen;
1644     tensor             vir, pres;
1645     int                start, end, number_steps;
1646     int                i, k, m, n, gf, step;
1647     int                mdof_flags;
1648     auto               mdatoms = mdAtoms->mdatoms();
1649
1650     GMX_LOG(mdlog.info)
1651             .asParagraph()
1652             .appendText(
1653                     "Note that activating L-BFGS energy minimization via the "
1654                     "integrator .mdp option and the command gmx mdrun may "
1655                     "be available in a different form in a future version of GROMACS, "
1656                     "e.g. gmx minimize and an .mdp option.");
1657
1658     if (PAR(cr))
1659     {
1660         gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1661     }
1662
1663     if (nullptr != constr)
1664     {
1665         gmx_fatal(
1666                 FARGS,
1667                 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1668                 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1669     }
1670
1671     n        = 3 * state_global->natoms;
1672     nmaxcorr = inputrec->nbfgscorr;
1673
1674     snew(frozen, n);
1675
1676     snew(p, n);
1677     snew(rho, nmaxcorr);
1678     snew(alpha, nmaxcorr);
1679
1680     snew(dx, nmaxcorr);
1681     for (i = 0; i < nmaxcorr; i++)
1682     {
1683         snew(dx[i], n);
1684     }
1685
1686     snew(dg, nmaxcorr);
1687     for (i = 0; i < nmaxcorr; i++)
1688     {
1689         snew(dg[i], n);
1690     }
1691
1692     step  = 0;
1693     neval = 0;
1694
1695     /* Init em */
1696     init_em(fplog, mdlog, LBFGS, cr, inputrec, imdSession, pull_work, state_global, top_global,
1697             &ems, &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
1698     const bool        simulationsShareState = false;
1699     gmx_mdoutf*       outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
1700                                    mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1701                                    StartingBehavior::NewSimulation, simulationsShareState, ms);
1702     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
1703                                    false, StartingBehavior::NewSimulation, mdModulesNotifier);
1704
1705     start = 0;
1706     end   = mdatoms->homenr;
1707
1708     /* We need 4 working states */
1709     em_state_t  s0{}, s1{}, s2{}, s3{};
1710     em_state_t* sa   = &s0;
1711     em_state_t* sb   = &s1;
1712     em_state_t* sc   = &s2;
1713     em_state_t* last = &s3;
1714     /* Initialize by copying the state from ems (we could skip x and f here) */
1715     *sa = ems;
1716     *sb = ems;
1717     *sc = ems;
1718
1719     /* Print to log file */
1720     print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1721
1722     do_log = do_ene = do_x = do_f = TRUE;
1723
1724     /* Max number of steps */
1725     number_steps = inputrec->nsteps;
1726
1727     /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1728     gf = 0;
1729     for (i = start; i < end; i++)
1730     {
1731         if (mdatoms->cFREEZE)
1732         {
1733             gf = mdatoms->cFREEZE[i];
1734         }
1735         for (m = 0; m < DIM; m++)
1736         {
1737             frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
1738         }
1739     }
1740     if (MASTER(cr))
1741     {
1742         sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1743     }
1744     if (fplog)
1745     {
1746         sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1747     }
1748
1749     if (vsite)
1750     {
1751         construct_vsites(vsite, state_global->x.rvec_array(), 1, nullptr, top.idef.iparams,
1752                          top.idef.il, fr->pbcType, fr->bMolPBC, cr, state_global->box);
1753     }
1754
1755     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1756     /* do_force always puts the charge groups in the box and shifts again
1757      * We do not unshift, so molecules are always whole
1758      */
1759     neval++;
1760     EnergyEvaluator energyEvaluator{
1761         fplog, mdlog,  cr,    ms,    top_global, &top, inputrec, imdSession, pull_work,
1762         nrnb,  wcycle, gstat, vsite, constr,     fcd,  mdAtoms,  fr,         runScheduleWork,
1763         enerd
1764     };
1765     energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1766
1767     if (MASTER(cr))
1768     {
1769         /* Copy stuff to the energy bin for easy printing etc. */
1770         matrix nullBox = {};
1771         energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1772                                          enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1773                                          nullptr, vir, pres, nullptr, mu_tot, constr);
1774
1775         EnergyOutput::printHeader(fplog, step, step);
1776         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step,
1777                                            step, fcd, nullptr);
1778     }
1779
1780     /* Set the initial step.
1781      * since it will be multiplied by the non-normalized search direction
1782      * vector (force vector the first time), we scale it by the
1783      * norm of the force.
1784      */
1785
1786     if (MASTER(cr))
1787     {
1788         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1789         fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1790         fprintf(stderr, "   F-max             = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1791         fprintf(stderr, "   F-Norm            = %12.5e\n", ems.fnorm / sqrtNumAtoms);
1792         fprintf(stderr, "\n");
1793         /* and copy to the log file too... */
1794         fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1795         fprintf(fplog, "   F-max             = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1796         fprintf(fplog, "   F-Norm            = %12.5e\n", ems.fnorm / sqrtNumAtoms);
1797         fprintf(fplog, "\n");
1798     }
1799
1800     // Point is an index to the memory of search directions, where 0 is the first one.
1801     point = 0;
1802
1803     // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1804     real* fInit = static_cast<real*>(ems.f.rvec_array()[0]);
1805     for (i = 0; i < n; i++)
1806     {
1807         if (!frozen[i])
1808         {
1809             dx[point][i] = fInit[i]; /* Initial search direction */
1810         }
1811         else
1812         {
1813             dx[point][i] = 0;
1814         }
1815     }
1816
1817     // Stepsize will be modified during the search, and actually it is not critical
1818     // (the main efficiency in the algorithm comes from changing directions), but
1819     // we still need an initial value, so estimate it as the inverse of the norm
1820     // so we take small steps where the potential fluctuates a lot.
1821     stepsize = 1.0 / ems.fnorm;
1822
1823     /* Start the loop over BFGS steps.
1824      * Each successful step is counted, and we continue until
1825      * we either converge or reach the max number of steps.
1826      */
1827
1828     ncorr = 0;
1829
1830     /* Set the gradient from the force */
1831     converged = FALSE;
1832     for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1833     {
1834
1835         /* Write coordinates if necessary */
1836         do_x = do_per_step(step, inputrec->nstxout);
1837         do_f = do_per_step(step, inputrec->nstfout);
1838
1839         mdof_flags = 0;
1840         if (do_x)
1841         {
1842             mdof_flags |= MDOF_X;
1843         }
1844
1845         if (do_f)
1846         {
1847             mdof_flags |= MDOF_F;
1848         }
1849
1850         if (inputrec->bIMD)
1851         {
1852             mdof_flags |= MDOF_IMD;
1853         }
1854
1855         mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step,
1856                                          static_cast<real>(step), &ems.s, state_global,
1857                                          observablesHistory, ems.f);
1858
1859         /* Do the linesearching in the direction dx[point][0..(n-1)] */
1860
1861         /* make s a pointer to current search direction - point=0 first time we get here */
1862         s = dx[point];
1863
1864         real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
1865         real* ff = static_cast<real*>(ems.f.rvec_array()[0]);
1866
1867         // calculate line gradient in position A
1868         for (gpa = 0, i = 0; i < n; i++)
1869         {
1870             gpa -= s[i] * ff[i];
1871         }
1872
1873         /* Calculate minimum allowed stepsize along the line, before the average (norm)
1874          * relative change in coordinate is smaller than precision
1875          */
1876         for (minstep = 0, i = 0; i < n; i++)
1877         {
1878             tmp = fabs(xx[i]);
1879             if (tmp < 1.0)
1880             {
1881                 tmp = 1.0;
1882             }
1883             tmp = s[i] / tmp;
1884             minstep += tmp * tmp;
1885         }
1886         minstep = GMX_REAL_EPS / sqrt(minstep / n);
1887
1888         if (stepsize < minstep)
1889         {
1890             converged = TRUE;
1891             break;
1892         }
1893
1894         // Before taking any steps along the line, store the old position
1895         *last       = ems;
1896         real* lastx = static_cast<real*>(last->s.x.data()[0]);
1897         real* lastf = static_cast<real*>(last->f.data()[0]);
1898         Epot0       = ems.epot;
1899
1900         *sa = ems;
1901
1902         /* Take a step downhill.
1903          * In theory, we should find the actual minimum of the function in this
1904          * direction, somewhere along the line.
1905          * That is quite possible, but it turns out to take 5-10 function evaluations
1906          * for each line. However, we dont really need to find the exact minimum -
1907          * it is much better to start a new BFGS step in a modified direction as soon
1908          * as we are close to it. This will save a lot of energy evaluations.
1909          *
1910          * In practice, we just try to take a single step.
1911          * If it worked (i.e. lowered the energy), we increase the stepsize but
1912          * continue straight to the next BFGS step without trying to find any minimum,
1913          * i.e. we change the search direction too. If the line was smooth, it is
1914          * likely we are in a smooth region, and then it makes sense to take longer
1915          * steps in the modified search direction too.
1916          *
1917          * If it didn't work (higher energy), there must be a minimum somewhere between
1918          * the old position and the new one. Then we need to start by finding a lower
1919          * value before we change search direction. Since the energy was apparently
1920          * quite rough, we need to decrease the step size.
1921          *
1922          * Due to the finite numerical accuracy, it turns out that it is a good idea
1923          * to accept a SMALL increase in energy, if the derivative is still downhill.
1924          * This leads to lower final energies in the tests I've done. / Erik
1925          */
1926
1927         // State "A" is the first position along the line.
1928         // reference position along line is initially zero
1929         a = 0.0;
1930
1931         // Check stepsize first. We do not allow displacements
1932         // larger than emstep.
1933         //
1934         do
1935         {
1936             // Pick a new position C by adding stepsize to A.
1937             c = a + stepsize;
1938
1939             // Calculate what the largest change in any individual coordinate
1940             // would be (translation along line * gradient along line)
1941             maxdelta = 0;
1942             for (i = 0; i < n; i++)
1943             {
1944                 delta = c * s[i];
1945                 if (delta > maxdelta)
1946                 {
1947                     maxdelta = delta;
1948                 }
1949             }
1950             // If any displacement is larger than the stepsize limit, reduce the step
1951             if (maxdelta > inputrec->em_stepsize)
1952             {
1953                 stepsize *= 0.1;
1954             }
1955         } while (maxdelta > inputrec->em_stepsize);
1956
1957         // Take a trial step and move the coordinate array xc[] to position C
1958         real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
1959         for (i = 0; i < n; i++)
1960         {
1961             xc[i] = lastx[i] + c * s[i];
1962         }
1963
1964         neval++;
1965         // Calculate energy for the trial step in position C
1966         energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
1967
1968         // Calc line gradient in position C
1969         real* fc = static_cast<real*>(sc->f.rvec_array()[0]);
1970         for (gpc = 0, i = 0; i < n; i++)
1971         {
1972             gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
1973         }
1974         /* Sum the gradient along the line across CPUs */
1975         if (PAR(cr))
1976         {
1977             gmx_sumd(1, &gpc, cr);
1978         }
1979
1980         // This is the max amount of increase in energy we tolerate.
1981         // By allowing VERY small changes (close to numerical precision) we
1982         // frequently find even better (lower) final energies.
1983         tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
1984
1985         // Accept the step if the energy is lower in the new position C (compared to A),
1986         // or if it is not significantly higher and the line derivative is still negative.
1987         foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
1988         // If true, great, we found a better energy. We no longer try to alter the
1989         // stepsize, but simply accept this new better position. The we select a new
1990         // search direction instead, which will be much more efficient than continuing
1991         // to take smaller steps along a line. Set fnorm based on the new C position,
1992         // which will be used to update the stepsize to 1/fnorm further down.
1993
1994         // If false, the energy is NOT lower in point C, i.e. it will be the same
1995         // or higher than in point A. In this case it is pointless to move to point C,
1996         // so we will have to do more iterations along the same line to find a smaller
1997         // value in the interval [A=0.0,C].
1998         // Here, A is still 0.0, but that will change when we do a search in the interval
1999         // [0.0,C] below. That search we will do by interpolation or bisection rather
2000         // than with the stepsize, so no need to modify it. For the next search direction
2001         // it will be reset to 1/fnorm anyway.
2002
2003         if (!foundlower)
2004         {
2005             // OK, if we didn't find a lower value we will have to locate one now - there must
2006             // be one in the interval [a,c].
2007             // The same thing is valid here, though: Don't spend dozens of iterations to find
2008             // the line minimum. We try to interpolate based on the derivative at the endpoints,
2009             // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2010             // I also have a safeguard for potentially really pathological functions so we never
2011             // take more than 20 steps before we give up.
2012             // If we already found a lower value we just skip this step and continue to the update.
2013             real fnorm = 0;
2014             nminstep   = 0;
2015             do
2016             {
2017                 // Select a new trial point B in the interval [A,C].
2018                 // If the derivatives at points a & c have different sign we interpolate to zero,
2019                 // otherwise just do a bisection since there might be multiple minima/maxima
2020                 // inside the interval.
2021                 if (gpa < 0 && gpc > 0)
2022                 {
2023                     b = a + gpa * (a - c) / (gpc - gpa);
2024                 }
2025                 else
2026                 {
2027                     b = 0.5 * (a + c);
2028                 }
2029
2030                 /* safeguard if interpolation close to machine accuracy causes errors:
2031                  * never go outside the interval
2032                  */
2033                 if (b <= a || b >= c)
2034                 {
2035                     b = 0.5 * (a + c);
2036                 }
2037
2038                 // Take a trial step to point B
2039                 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2040                 for (i = 0; i < n; i++)
2041                 {
2042                     xb[i] = lastx[i] + b * s[i];
2043                 }
2044
2045                 neval++;
2046                 // Calculate energy for the trial step in point B
2047                 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2048                 fnorm = sb->fnorm;
2049
2050                 // Calculate gradient in point B
2051                 real* fb = static_cast<real*>(sb->f.rvec_array()[0]);
2052                 for (gpb = 0, i = 0; i < n; i++)
2053                 {
2054                     gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2055                 }
2056                 /* Sum the gradient along the line across CPUs */
2057                 if (PAR(cr))
2058                 {
2059                     gmx_sumd(1, &gpb, cr);
2060                 }
2061
2062                 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2063                 // at the new point B, and rename the endpoints of this new interval A and C.
2064                 if (gpb > 0)
2065                 {
2066                     /* Replace c endpoint with b */
2067                     c = b;
2068                     /* copy state b to c */
2069                     *sc = *sb;
2070                 }
2071                 else
2072                 {
2073                     /* Replace a endpoint with b */
2074                     a = b;
2075                     /* copy state b to a */
2076                     *sa = *sb;
2077                 }
2078
2079                 /*
2080                  * Stop search as soon as we find a value smaller than the endpoints,
2081                  * or if the tolerance is below machine precision.
2082                  * Never run more than 20 steps, no matter what.
2083                  */
2084                 nminstep++;
2085             } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2086
2087             if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2088             {
2089                 /* OK. We couldn't find a significantly lower energy.
2090                  * If ncorr==0 this was steepest descent, and then we give up.
2091                  * If not, reset memory to restart as steepest descent before quitting.
2092                  */
2093                 if (ncorr == 0)
2094                 {
2095                     /* Converged */
2096                     converged = TRUE;
2097                     break;
2098                 }
2099                 else
2100                 {
2101                     /* Reset memory */
2102                     ncorr = 0;
2103                     /* Search in gradient direction */
2104                     for (i = 0; i < n; i++)
2105                     {
2106                         dx[point][i] = ff[i];
2107                     }
2108                     /* Reset stepsize */
2109                     stepsize = 1.0 / fnorm;
2110                     continue;
2111                 }
2112             }
2113
2114             /* Select min energy state of A & C, put the best in xx/ff/Epot
2115              */
2116             if (sc->epot < sa->epot)
2117             {
2118                 /* Use state C */
2119                 ems        = *sc;
2120                 step_taken = c;
2121             }
2122             else
2123             {
2124                 /* Use state A */
2125                 ems        = *sa;
2126                 step_taken = a;
2127             }
2128         }
2129         else
2130         {
2131             /* found lower */
2132             /* Use state C */
2133             ems        = *sc;
2134             step_taken = c;
2135         }
2136
2137         /* Update the memory information, and calculate a new
2138          * approximation of the inverse hessian
2139          */
2140
2141         /* Have new data in Epot, xx, ff */
2142         if (ncorr < nmaxcorr)
2143         {
2144             ncorr++;
2145         }
2146
2147         for (i = 0; i < n; i++)
2148         {
2149             dg[point][i] = lastf[i] - ff[i];
2150             dx[point][i] *= step_taken;
2151         }
2152
2153         dgdg = 0;
2154         dgdx = 0;
2155         for (i = 0; i < n; i++)
2156         {
2157             dgdg += dg[point][i] * dg[point][i];
2158             dgdx += dg[point][i] * dx[point][i];
2159         }
2160
2161         diag = dgdx / dgdg;
2162
2163         rho[point] = 1.0 / dgdx;
2164         point++;
2165
2166         if (point >= nmaxcorr)
2167         {
2168             point = 0;
2169         }
2170
2171         /* Update */
2172         for (i = 0; i < n; i++)
2173         {
2174             p[i] = ff[i];
2175         }
2176
2177         cp = point;
2178
2179         /* Recursive update. First go back over the memory points */
2180         for (k = 0; k < ncorr; k++)
2181         {
2182             cp--;
2183             if (cp < 0)
2184             {
2185                 cp = ncorr - 1;
2186             }
2187
2188             sq = 0;
2189             for (i = 0; i < n; i++)
2190             {
2191                 sq += dx[cp][i] * p[i];
2192             }
2193
2194             alpha[cp] = rho[cp] * sq;
2195
2196             for (i = 0; i < n; i++)
2197             {
2198                 p[i] -= alpha[cp] * dg[cp][i];
2199             }
2200         }
2201
2202         for (i = 0; i < n; i++)
2203         {
2204             p[i] *= diag;
2205         }
2206
2207         /* And then go forward again */
2208         for (k = 0; k < ncorr; k++)
2209         {
2210             yr = 0;
2211             for (i = 0; i < n; i++)
2212             {
2213                 yr += p[i] * dg[cp][i];
2214             }
2215
2216             beta = rho[cp] * yr;
2217             beta = alpha[cp] - beta;
2218
2219             for (i = 0; i < n; i++)
2220             {
2221                 p[i] += beta * dx[cp][i];
2222             }
2223
2224             cp++;
2225             if (cp >= ncorr)
2226             {
2227                 cp = 0;
2228             }
2229         }
2230
2231         for (i = 0; i < n; i++)
2232         {
2233             if (!frozen[i])
2234             {
2235                 dx[point][i] = p[i];
2236             }
2237             else
2238             {
2239                 dx[point][i] = 0;
2240             }
2241         }
2242
2243         /* Print it if necessary */
2244         if (MASTER(cr))
2245         {
2246             if (mdrunOptions.verbose)
2247             {
2248                 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2249                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step,
2250                         ems.epot, ems.fnorm / sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2251                 fflush(stderr);
2252             }
2253             /* Store the new (lower) energies */
2254             matrix nullBox = {};
2255             energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
2256                                              enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
2257                                              nullptr, vir, pres, nullptr, mu_tot, constr);
2258
2259             do_log = do_per_step(step, inputrec->nstlog);
2260             do_ene = do_per_step(step, inputrec->nstenergy);
2261
2262             imdSession->fillEnergyRecord(step, TRUE);
2263
2264             if (do_log)
2265             {
2266                 EnergyOutput::printHeader(fplog, step, step);
2267             }
2268             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2269                                                do_log ? fplog : nullptr, step, step, fcd, nullptr);
2270         }
2271
2272         /* Send x and E to IMD client, if bIMD is TRUE. */
2273         if (imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0) && MASTER(cr))
2274         {
2275             imdSession->sendPositionsAndEnergies();
2276         }
2277
2278         // Reset stepsize in we are doing more iterations
2279         stepsize = 1.0;
2280
2281         /* Stop when the maximum force lies below tolerance.
2282          * If we have reached machine precision, converged is already set to true.
2283          */
2284         converged = converged || (ems.fmax < inputrec->em_tol);
2285
2286     } /* End of the loop */
2287
2288     if (converged)
2289     {
2290         step--; /* we never took that last step in this case */
2291     }
2292     if (ems.fmax > inputrec->em_tol)
2293     {
2294         if (MASTER(cr))
2295         {
2296             warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2297         }
2298         converged = FALSE;
2299     }
2300
2301     /* If we printed energy and/or logfile last step (which was the last step)
2302      * we don't have to do it again, but otherwise print the final values.
2303      */
2304     if (!do_log) /* Write final value to log since we didn't do anythin last step */
2305     {
2306         EnergyOutput::printHeader(fplog, step, step);
2307     }
2308     if (!do_ene || !do_log) /* Write final energy file entries */
2309     {
2310         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2311                                            !do_log ? fplog : nullptr, step, step, fcd, nullptr);
2312     }
2313
2314     /* Print some stuff... */
2315     if (MASTER(cr))
2316     {
2317         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2318     }
2319
2320     /* IMPORTANT!
2321      * For accurate normal mode calculation it is imperative that we
2322      * store the last conformation into the full precision binary trajectory.
2323      *
2324      * However, we should only do it if we did NOT already write this step
2325      * above (which we did if do_x or do_f was true).
2326      */
2327     do_x = !do_per_step(step, inputrec->nstxout);
2328     do_f = !do_per_step(step, inputrec->nstfout);
2329     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec,
2330                   step, &ems, state_global, observablesHistory);
2331
2332     if (MASTER(cr))
2333     {
2334         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2335         print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2336         print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2337
2338         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2339     }
2340
2341     finish_em(cr, outf, walltime_accounting, wcycle);
2342
2343     /* To print the actual number of steps we needed somewhere */
2344     walltime_accounting_set_nsteps_done(walltime_accounting, step);
2345 }
2346
2347 void LegacySimulator::do_steep()
2348 {
2349     const char*       SD = "Steepest Descents";
2350     gmx_localtop_t    top(top_global->ffparams);
2351     gmx_global_stat_t gstat;
2352     real              stepsize;
2353     real              ustep;
2354     gmx_bool          bDone, bAbort, do_x, do_f;
2355     tensor            vir, pres;
2356     rvec              mu_tot = { 0 };
2357     int               nsteps;
2358     int               count          = 0;
2359     int               steps_accepted = 0;
2360     auto              mdatoms        = mdAtoms->mdatoms();
2361
2362     GMX_LOG(mdlog.info)
2363             .asParagraph()
2364             .appendText(
2365                     "Note that activating steepest-descent energy minimization via the "
2366                     "integrator .mdp option and the command gmx mdrun may "
2367                     "be available in a different form in a future version of GROMACS, "
2368                     "e.g. gmx minimize and an .mdp option.");
2369
2370     /* Create 2 states on the stack and extract pointers that we will swap */
2371     em_state_t  s0{}, s1{};
2372     em_state_t* s_min = &s0;
2373     em_state_t* s_try = &s1;
2374
2375     /* Init em and store the local state in s_try */
2376     init_em(fplog, mdlog, SD, cr, inputrec, imdSession, pull_work, state_global, top_global, s_try,
2377             &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
2378     const bool        simulationsShareState = false;
2379     gmx_mdoutf*       outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
2380                                    mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2381                                    StartingBehavior::NewSimulation, simulationsShareState, ms);
2382     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
2383                                    false, StartingBehavior::NewSimulation, mdModulesNotifier);
2384
2385     /* Print to log file  */
2386     print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2387
2388     /* Set variables for stepsize (in nm). This is the largest
2389      * step that we are going to make in any direction.
2390      */
2391     ustep    = inputrec->em_stepsize;
2392     stepsize = 0;
2393
2394     /* Max number of steps  */
2395     nsteps = inputrec->nsteps;
2396
2397     if (MASTER(cr))
2398     {
2399         /* Print to the screen  */
2400         sp_header(stderr, SD, inputrec->em_tol, nsteps);
2401     }
2402     if (fplog)
2403     {
2404         sp_header(fplog, SD, inputrec->em_tol, nsteps);
2405     }
2406     EnergyEvaluator energyEvaluator{
2407         fplog, mdlog,  cr,    ms,    top_global, &top, inputrec, imdSession, pull_work,
2408         nrnb,  wcycle, gstat, vsite, constr,     fcd,  mdAtoms,  fr,         runScheduleWork,
2409         enerd
2410     };
2411
2412     /**** HERE STARTS THE LOOP ****
2413      * count is the counter for the number of steps
2414      * bDone will be TRUE when the minimization has converged
2415      * bAbort will be TRUE when nsteps steps have been performed or when
2416      * the stepsize becomes smaller than is reasonable for machine precision
2417      */
2418     count  = 0;
2419     bDone  = FALSE;
2420     bAbort = FALSE;
2421     while (!bDone && !bAbort)
2422     {
2423         bAbort = (nsteps >= 0) && (count == nsteps);
2424
2425         /* set new coordinates, except for first step */
2426         bool validStep = true;
2427         if (count > 0)
2428         {
2429             validStep = do_em_step(cr, inputrec, mdatoms, s_min, stepsize, &s_min->f, s_try, constr, count);
2430         }
2431
2432         if (validStep)
2433         {
2434             energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2435         }
2436         else
2437         {
2438             // Signal constraint error during stepping with energy=inf
2439             s_try->epot = std::numeric_limits<real>::infinity();
2440         }
2441
2442         if (MASTER(cr))
2443         {
2444             EnergyOutput::printHeader(fplog, count, count);
2445         }
2446
2447         if (count == 0)
2448         {
2449             s_min->epot = s_try->epot;
2450         }
2451
2452         /* Print it if necessary  */
2453         if (MASTER(cr))
2454         {
2455             if (mdrunOptions.verbose)
2456             {
2457                 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2458                         count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax + 1,
2459                         ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2460                 fflush(stderr);
2461             }
2462
2463             if ((count == 0) || (s_try->epot < s_min->epot))
2464             {
2465                 /* Store the new (lower) energies  */
2466                 matrix nullBox = {};
2467                 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(count), mdatoms->tmass,
2468                                                  enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
2469                                                  nullptr, vir, pres, nullptr, mu_tot, constr);
2470
2471                 imdSession->fillEnergyRecord(count, TRUE);
2472
2473                 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2474                 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2475                 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or,
2476                                                    fplog, count, count, fcd, nullptr);
2477                 fflush(fplog);
2478             }
2479         }
2480
2481         /* Now if the new energy is smaller than the previous...
2482          * or if this is the first step!
2483          * or if we did random steps!
2484          */
2485
2486         if ((count == 0) || (s_try->epot < s_min->epot))
2487         {
2488             steps_accepted++;
2489
2490             /* Test whether the convergence criterion is met...  */
2491             bDone = (s_try->fmax < inputrec->em_tol);
2492
2493             /* Copy the arrays for force, positions and energy  */
2494             /* The 'Min' array always holds the coords and forces of the minimal
2495                sampled energy  */
2496             swap_em_state(&s_min, &s_try);
2497             if (count > 0)
2498             {
2499                 ustep *= 1.2;
2500             }
2501
2502             /* Write to trn, if necessary */
2503             do_x = do_per_step(steps_accepted, inputrec->nstxout);
2504             do_f = do_per_step(steps_accepted, inputrec->nstfout);
2505             write_em_traj(fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min,
2506                           state_global, observablesHistory);
2507         }
2508         else
2509         {
2510             /* If energy is not smaller make the step smaller...  */
2511             ustep *= 0.5;
2512
2513             if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2514             {
2515                 /* Reload the old state */
2516                 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
2517                                        pull_work, s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
2518             }
2519         }
2520
2521         // If the force is very small after finishing minimization,
2522         // we risk dividing by zero when calculating the step size.
2523         // So we check first if the minimization has stopped before
2524         // trying to obtain a new step size.
2525         if (!bDone)
2526         {
2527             /* Determine new step  */
2528             stepsize = ustep / s_min->fmax;
2529         }
2530
2531         /* Check if stepsize is too small, with 1 nm as a characteristic length */
2532 #if GMX_DOUBLE
2533         if (count == nsteps || ustep < 1e-12)
2534 #else
2535         if (count == nsteps || ustep < 1e-6)
2536 #endif
2537         {
2538             if (MASTER(cr))
2539             {
2540                 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
2541             }
2542             bAbort = TRUE;
2543         }
2544
2545         /* Send IMD energies and positions, if bIMD is TRUE. */
2546         if (imdSession->run(count, TRUE, state_global->box,
2547                             MASTER(cr) ? state_global->x.rvec_array() : nullptr, 0)
2548             && MASTER(cr))
2549         {
2550             imdSession->sendPositionsAndEnergies();
2551         }
2552
2553         count++;
2554     } /* End of the loop  */
2555
2556     /* Print some data...  */
2557     if (MASTER(cr))
2558     {
2559         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2560     }
2561     write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2562                   top_global, inputrec, count, s_min, state_global, observablesHistory);
2563
2564     if (MASTER(cr))
2565     {
2566         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2567
2568         print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2569         print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2570     }
2571
2572     finish_em(cr, outf, walltime_accounting, wcycle);
2573
2574     /* To print the actual number of steps we needed somewhere */
2575     inputrec->nsteps = count;
2576
2577     walltime_accounting_set_nsteps_done(walltime_accounting, count);
2578 }
2579
2580 void LegacySimulator::do_nm()
2581 {
2582     const char*         NM = "Normal Mode Analysis";
2583     int                 nnodes;
2584     gmx_localtop_t      top(top_global->ffparams);
2585     gmx_global_stat_t   gstat;
2586     tensor              vir, pres;
2587     rvec                mu_tot = { 0 };
2588     rvec*               dfdx;
2589     gmx_bool            bSparse; /* use sparse matrix storage format */
2590     size_t              sz;
2591     gmx_sparsematrix_t* sparse_matrix = nullptr;
2592     real*               full_matrix   = nullptr;
2593
2594     /* added with respect to mdrun */
2595     int  row, col;
2596     real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
2597     real x_min;
2598     bool bIsMaster = MASTER(cr);
2599     auto mdatoms   = mdAtoms->mdatoms();
2600
2601     GMX_LOG(mdlog.info)
2602             .asParagraph()
2603             .appendText(
2604                     "Note that activating normal-mode analysis via the integrator "
2605                     ".mdp option and the command gmx mdrun may "
2606                     "be available in a different form in a future version of GROMACS, "
2607                     "e.g. gmx normal-modes.");
2608
2609     if (constr != nullptr)
2610     {
2611         gmx_fatal(
2612                 FARGS,
2613                 "Constraints present with Normal Mode Analysis, this combination is not supported");
2614     }
2615
2616     gmx_shellfc_t* shellfc;
2617
2618     em_state_t state_work{};
2619
2620     /* Init em and store the local state in state_minimum */
2621     init_em(fplog, mdlog, NM, cr, inputrec, imdSession, pull_work, state_global, top_global,
2622             &state_work, &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, &shellfc);
2623     const bool  simulationsShareState = false;
2624     gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
2625                                    mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2626                                    StartingBehavior::NewSimulation, simulationsShareState, ms);
2627
2628     std::vector<int>       atom_index = get_atom_index(top_global);
2629     std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
2630     snew(dfdx, atom_index.size());
2631
2632 #if !GMX_DOUBLE
2633     if (bIsMaster)
2634     {
2635         fprintf(stderr,
2636                 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2637                 "      which MIGHT not be accurate enough for normal mode analysis.\n"
2638                 "      GROMACS now uses sparse matrix storage, so the memory requirements\n"
2639                 "      are fairly modest even if you recompile in double precision.\n\n");
2640     }
2641 #endif
2642
2643     /* Check if we can/should use sparse storage format.
2644      *
2645      * Sparse format is only useful when the Hessian itself is sparse, which it
2646      * will be when we use a cutoff.
2647      * For small systems (n<1000) it is easier to always use full matrix format, though.
2648      */
2649     if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2650     {
2651         GMX_LOG(mdlog.warning)
2652                 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2653         bSparse = FALSE;
2654     }
2655     else if (atom_index.size() < 1000)
2656     {
2657         GMX_LOG(mdlog.warning)
2658                 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2659                                      atom_index.size());
2660         bSparse = FALSE;
2661     }
2662     else
2663     {
2664         GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2665         bSparse = TRUE;
2666     }
2667
2668     /* Number of dimensions, based on real atoms, that is not vsites or shell */
2669     sz = DIM * atom_index.size();
2670
2671     fprintf(stderr, "Allocating Hessian memory...\n\n");
2672
2673     if (bSparse)
2674     {
2675         sparse_matrix                       = gmx_sparsematrix_init(sz);
2676         sparse_matrix->compressed_symmetric = TRUE;
2677     }
2678     else
2679     {
2680         snew(full_matrix, sz * sz);
2681     }
2682
2683     /* Write start time and temperature */
2684     print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2685
2686     /* fudge nr of steps to nr of atoms */
2687     inputrec->nsteps = atom_index.size() * 2;
2688
2689     if (bIsMaster)
2690     {
2691         fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2692                 *(top_global->name), inputrec->nsteps);
2693     }
2694
2695     nnodes = cr->nnodes;
2696
2697     /* Make evaluate_energy do a single node force calculation */
2698     cr->nnodes = 1;
2699     EnergyEvaluator energyEvaluator{
2700         fplog, mdlog,  cr,    ms,    top_global, &top, inputrec, imdSession, pull_work,
2701         nrnb,  wcycle, gstat, vsite, constr,     fcd,  mdAtoms,  fr,         runScheduleWork,
2702         enerd
2703     };
2704     energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2705     cr->nnodes = nnodes;
2706
2707     /* if forces are not small, warn user */
2708     get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2709
2710     GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2711     if (state_work.fmax > 1.0e-3)
2712     {
2713         GMX_LOG(mdlog.warning)
2714                 .appendText(
2715                         "The force is probably not small enough to "
2716                         "ensure that you are at a minimum.\n"
2717                         "Be aware that negative eigenvalues may occur\n"
2718                         "when the resulting matrix is diagonalized.");
2719     }
2720
2721     /***********************************************************
2722      *
2723      *      Loop over all pairs in matrix
2724      *
2725      *      do_force called twice. Once with positive and
2726      *      once with negative displacement
2727      *
2728      ************************************************************/
2729
2730     /* Steps are divided one by one over the nodes */
2731     bool bNS          = true;
2732     auto state_work_x = makeArrayRef(state_work.s.x);
2733     auto state_work_f = makeArrayRef(state_work.f);
2734     for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
2735     {
2736         size_t atom = atom_index[aid];
2737         for (size_t d = 0; d < DIM; d++)
2738         {
2739             int64_t step        = 0;
2740             int     force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2741             double  t           = 0;
2742
2743             x_min = state_work_x[atom][d];
2744
2745             for (unsigned int dx = 0; (dx < 2); dx++)
2746             {
2747                 if (dx == 0)
2748                 {
2749                     state_work_x[atom][d] = x_min - der_range;
2750                 }
2751                 else
2752                 {
2753                     state_work_x[atom][d] = x_min + der_range;
2754                 }
2755
2756                 /* Make evaluate_energy do a single node force calculation */
2757                 cr->nnodes = 1;
2758                 if (shellfc)
2759                 {
2760                     /* Now is the time to relax the shells */
2761                     relax_shell_flexcon(
2762                             fplog, cr, ms, mdrunOptions.verbose, nullptr, step, inputrec, imdSession,
2763                             pull_work, bNS, force_flags, &top, constr, enerd, fcd, state_work.s.natoms,
2764                             state_work.s.x.arrayRefWithPadding(), state_work.s.v.arrayRefWithPadding(),
2765                             state_work.s.box, state_work.s.lambda, &state_work.s.hist,
2766                             state_work.f.arrayRefWithPadding(), vir, mdatoms, nrnb, wcycle, shellfc,
2767                             fr, runScheduleWork, t, mu_tot, vsite, DDBalanceRegionHandler(nullptr));
2768                     bNS = false;
2769                     step++;
2770                 }
2771                 else
2772                 {
2773                     energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE);
2774                 }
2775
2776                 cr->nnodes = nnodes;
2777
2778                 if (dx == 0)
2779                 {
2780                     std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(),
2781                               fneg.begin());
2782                 }
2783             }
2784
2785             /* x is restored to original */
2786             state_work_x[atom][d] = x_min;
2787
2788             for (size_t j = 0; j < atom_index.size(); j++)
2789             {
2790                 for (size_t k = 0; (k < DIM); k++)
2791                 {
2792                     dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
2793                 }
2794             }
2795
2796             if (!bIsMaster)
2797             {
2798 #if GMX_MPI
2799 #    define mpi_type GMX_MPI_REAL
2800                 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid,
2801                          cr->mpi_comm_mygroup);
2802 #endif
2803             }
2804             else
2805             {
2806                 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
2807                 {
2808                     if (node > 0)
2809                     {
2810 #if GMX_MPI
2811                         MPI_Status stat;
2812                         MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node,
2813                                  cr->mpi_comm_mygroup, &stat);
2814 #    undef mpi_type
2815 #endif
2816                     }
2817
2818                     row = (aid + node) * DIM + d;
2819
2820                     for (size_t j = 0; j < atom_index.size(); j++)
2821                     {
2822                         for (size_t k = 0; k < DIM; k++)
2823                         {
2824                             col = j * DIM + k;
2825
2826                             if (bSparse)
2827                             {
2828                                 if (col >= row && dfdx[j][k] != 0.0)
2829                                 {
2830                                     gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
2831                                 }
2832                             }
2833                             else
2834                             {
2835                                 full_matrix[row * sz + col] = dfdx[j][k];
2836                             }
2837                         }
2838                     }
2839                 }
2840             }
2841
2842             if (mdrunOptions.verbose && fplog)
2843             {
2844                 fflush(fplog);
2845             }
2846         }
2847         /* write progress */
2848         if (bIsMaster && mdrunOptions.verbose)
2849         {
2850             fprintf(stderr, "\rFinished step %d out of %td",
2851                     std::min<int>(atom + nnodes, atom_index.size()), ssize(atom_index));
2852             fflush(stderr);
2853         }
2854     }
2855
2856     if (bIsMaster)
2857     {
2858         fprintf(stderr, "\n\nWriting Hessian...\n");
2859         gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2860     }
2861
2862     finish_em(cr, outf, walltime_accounting, wcycle);
2863
2864     walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);
2865 }
2866
2867 } // namespace gmx