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