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