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