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