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