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