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